-
Notifications
You must be signed in to change notification settings - Fork 7
/
Copy pathhgpt2.m
271 lines (241 loc) · 9.32 KB
/
hgpt2.m
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
function [P, T, RH, Tm, ZHD, ZWD, PWV] = hgpt2(dt, x0, y0, z0, z0_type)
% hgpt2.m
%
% This routine determines the surface pressure (P, hPa), surface air temperature (T, K), relative humidity (RH, %)
% weighed mean temperature (Tm, K), zenith hydrostatic delay (ZHD, m), zenith wet delay (ZWD, m), and precipitable water vapor (PWV, m)
% from binary coefficient files
% As available from:
% https://github.com/pjmateus/hgpt_model (release v2.0)
% press_grid.bin; temp_grid.bin; tm_grid.bin; and rh_grid.bin
%
% It is admitted that the binary files with the coefficients are in the same directory as this script.
% In alternative you can define the "coeffiles" variable
%
% The epoch can be an array of size 1, and in this case is the Modified Julian Date (MJD)
% or can be an array of size 6, with the Gregorian Calendar in the following format (year, month, day, hour, min, sec)
%
% All parameters are bilinear interpolated to the input ellipsoidal longitude and latitude
%
% Reference for HGPT and HGPT2:
% Mateus, P.; Catalão, J.; Mendes, V.B.; Nico, G. An ERA5-Based Hourly Global Pressure and Temperature (HGPT) Model.
% Remote Sens. 2020, 12, 1098. https://doi.org/10.3390/rs12071098
% HGPT2: an ERA5-based global model to estimate relative humidity (Remote Sensing, MDPI)
%
% INPUT:
% dt : if size(dt)=1 => modified julian date
% if size(dt)=6 => year, month, day, hour, min, sec
% x0 : ellipsoidal longitude (degrees)
% y0 : ellipsoidal latitude (degrees)
% z0 : height (m)
% z0_type : ‘orth’ for orthometric height or ‘elli’ for ellipsoidal height
%
% OUTPUT:
% P : surface pressure valid at (x0, y0, z0), in hPa
% T : surface air temperature valid at (x0, y0, z0), in Kelvins
% RH : relative humidity valid at (x0, y0, z0), in %
% Tm : weighed mean temperature valid at (x0, y0, z0), in Kelvins
% ZHD : zenith hydrostatic delay, valid at (x0, y0, z0), in meters
% ZWD : zenith wet delay, valid at (x0, y0, z0), in meters
% PWV : precipitable water vapor, valid at (x0, y0, z0), in meters
%
%--------------------------------------------------------------------------
% Example:
% y0 = 38.5519;
% x0 = -9.0147;
% z0 = 25;
% dt = 58119.5; or dt = [2018, 1, 1, 12, 0, 0]
% [P, T, RH, Tm, ZHD, ZWD, PWV] = hgpt2(dt, x0, y0, z0, 'orth')
% [1020.262,
% 286.589,
% 74.908,
% 277.236,
% 2.324,
% 0.117,
% 0.018]
%--------------------------------------------------------------------------
% written by Pedro Mateus (2021/05/15)
% Instituto Dom Luiz (IDL), Faculdade de Ciências, Universidade de Lisboa, 1749-016 Lisboa, Portugal
%
% Location of coefficient files
coeffiles = '';
% Constants
row = 721;
col = 1440;
p1 = 365.250;
p2 = 182.625;
p3 = 91.3125;
% Geographic coordinates ( equal to ERA5 )
lon = linspace(-180, 179.75, col);
lat = linspace(-90, 90, row);
% Input datetime format
if length(dt) == 6
mjd = juliandate( datetime(dt(1), dt(2), dt(3), dt(4), dt(5), dt(6)) ) - 2400000.5;
hour = dt(4);
elseif length(dt) == 1
mjd = dt;
mjdconv = datevec(mjd+678942);
hour = mjdconv(4);
else
error('Use 1) Modified Julian Date (MJD) or 2) Gregorian date (y,m,d,HH,MM,SS).')
end
% Finding indexes for bilinear interpolation
% x-location
[~, indx] = sort( sqrt( (lon-x0).^2 ) );
ix1 = indx( 1 ); ix2 = indx( 2 );
x1 = lon( ix1 ); x2 = lon( ix2 );
x = [ix1, ix1, ix2, ix2];
% y-location
[~, indy] = sort( sqrt( (lat-y0).^2 ) );
jy1 = indy( 1 ); jy2 = indy( 2 );
y1 = lat( jy1 ); y2 = lat( jy2 );
y = [jy1, jy2, jy1, jy2];
% xy-distances (weights)
dx1y1= 1/sqrt( (x1 - x0)^2 + (y1 - y0)^2 );
dx1y2= 1/sqrt( (x1 - x0)^2 + (y2 - y0)^2 );
dx2y1= 1/sqrt( (x2 - x0)^2 + (y1 - y0)^2 );
dx2y2= 1/sqrt( (x2 - x0)^2 + (y2 - y0)^2 );
dxy = [dx1y1, dx1y2, dx2y1, dx2y2];
if sum(isinf(dxy)) > 0
% Exact point grid
dxy = [1,0,0,0];
end
% ******************************************************************
% Open and read the surface air temperature coefficients file
% ******************************************************************
if ~isempty(coeffiles)
if isunix && ~strcmp(coeffiles(end),'/'), coeffiles(end+1) = '/'; end
if ~isunix && ~strcmp(coeffiles(end),'\'), coeffiles(end+1) = '\'; end
end
[fid, errmsg] = fopen([coeffiles,'temp_grid.bin'], 'r');
if fid == -1
error(errmsg)
else
fseek(fid, (row*col*26)*hour, -1);
a0 = fread(fid, [row,col], 'single');
b0 = fread(fid, [row,col], 'single');
a1 = fread(fid, [row,col], 'single');
f1 = (fread(fid, [row,col], 'int16'))./10000;
a2 = fread(fid, [row,col], 'single');
f2 = (fread(fid, [row,col], 'int16'))./10000;
a3 = fread(fid, [row,col], 'single');
f3 = (fread(fid, [row,col], 'int16'))./10000;
fclose(fid);
end
% Surface air temperature model
fun_t = @(a0, b0, a1, f1, a2, f2, a3, f3) a0 + b0*(mjd - 51178) + a1*cos(2*pi*(mjd - 51178)/p1+f1) + ...
a2*cos(2*pi*(mjd - 51178)/p2+f2) + a3*cos(2*pi*(mjd - 51178)/p3+f3);
% Applying the bilinear interpolation
for j = 1 : length(x)
tij(j) = fun_t(a0(y(j),x(j)), b0(y(j),x(j)), a1(y(j),x(j)), f1(y(j),x(j)), ...
a2(y(j),x(j)), f2(y(j),x(j)), a3(y(j),x(j)), f3(y(j),x(j)));
end
T = sum(tij.*dxy)/sum(dxy);
% ******************************************************************
% Open and read the surface pressure coefficients file
% ******************************************************************
[fid, errmsg] = fopen([coeffiles,'press_grid.bin'], 'r');
if fid == -1
error(errmsg)
else
fseek(fid, (row*col*20)*hour, -1);
a0 = fread(fid, [row,col], 'single');
b0 = fread(fid, [row,col], 'single');
a1 = fread(fid, [row,col], 'single');
f1 = (fread(fid, [row,col], 'int16'))./10000;
a2 = fread(fid, [row,col], 'single');
f2 = (fread(fid, [row,col], 'int16'))./10000;
fclose(fid);
end
% Surface pressure model
fun_p = @(a0, b0, a1, f1, a2, f2) a0 + b0*(mjd - 51178) + a1*cos(2*pi*(mjd - 51178)/p1+f1) + ...
a2*cos(2*pi*(mjd - 51178)/p2+f2);
% Applying the bilinear interpolation
for j = 1 : length(x)
pij(j) = fun_p(a0(y(j),x(j)), b0(y(j),x(j)), a1(y(j),x(j)), f1(y(j),x(j)), a2(y(j),x(j)), f2(y(j),x(j)));
end
P = sum(pij.*dxy)/sum(dxy);
% ******************************************************************
% Open and read the surface relative humidity coefficients file
% ******************************************************************
[fid, errmsg] = fopen([coeffiles,'rh_grid.bin'], 'r');
if fid == -1
error(errmsg)
else
fseek(fid, (row*col*22)*hour, -1);
a0 = fread(fid, [row,col], 'single');
a1 = fread(fid, [row,col], 'single');
f1 = (fread(fid, [row,col], 'int16'))./10000;
a2 = fread(fid, [row,col], 'single');
f2 = (fread(fid, [row,col], 'int16'))./10000;
a3 = fread(fid, [row,col], 'single');
f3 = (fread(fid, [row,col], 'int16'))./10000;
fclose(fid);
end
% Surface relative humidity model
fun_rh = @(a0, a1, f1, a2, f2, a3, f3) a0 + a1*cos(2*pi*(mjd - 51178)/p1+f1) + ...
a2*cos(2*pi*(mjd - 51178)/p2+f2) + a3*cos(2*pi*(mjd - 51178)/p3+f3);
% Applying the bilinear interpolation
for j = 1 : length(x)
rhij(j) = fun_rh(a0(y(j),x(j)), a1(y(j),x(j)), f1(y(j),x(j)), ...
a2(y(j),x(j)), f2(y(j),x(j)), a3(y(j),x(j)), f3(y(j),x(j)));
end
RH = sum(rhij.*dxy)/sum(dxy);
% ******************************************************************
% Open and read the Tm coefficients and undulation file
% ******************************************************************
[fid, errmsg] = fopen([coeffiles,'tm_grid.bin'], 'r');
if fid == -1
error(errmsg)
else
a = fread(fid, [row,col], 'single');
b = fread(fid, [row,col], 'single');
orography = fread(fid, [row,col], 'single');
undu = fread(fid, [row,col], 'single');
fclose(fid);
end
% Applying the bilinear interpolation
for j = 1 : length(x)
aij(j) = a(y(j),x(j));
bij(j) = b(y(j),x(j));
hij(j) = orography(y(j),x(j));
nij(j) = undu(y(j),x(j));
end
a = sum(aij.*dxy)/sum(dxy);
b = sum(bij.*dxy)/sum(dxy);
geo_height = sum(hij.*dxy)/sum(dxy);
N = sum(nij.*dxy)/sum(dxy);
% Zenith hydrostatic delay (ZHD), Saastamoinen model
if strcmp(z0_type, 'orth')
H_orth = z0;
elseif strcmp(z0_type, 'elli')
H_orth = z0 - N;
else
error('Use 1) <<orth>> for Orthometric height or 2) <<elli>> for Ellipsoidal height (in m).')
end
% Correction to P, T, and RH (see Guochanf Xu, GPS Theory, Algorithms and Applications, 3nd Edition, page 82)
dh = H_orth - geo_height;
% Temperature lapse rate by latitude
% rate = (6.25 - 2*sin( deg2rad(y0) )^4)/1000;
rate = (10.3 + 0.03182*(T-273.15) - 0.00436*P)/1000;
% Altimetric corrections
P = (P*100 * ( 1 - (rate*dh)/T )^5.6004)/100;
T = T - rate * ( dh );
% Ensure that we eliminate any noise outside this range
% Gamit has an unidentified error when assuming rh = 100% (met-files)
if RH >= 100, RH = 99.9; end
if RH < 0, RH = 0.0; end
% Call external function
es = es_wexler(T, P);
e = es * (RH/100);
% Weight mean temperature, Tm
Tm = a + b*T;
% ZHD using the Saastamoinen Model (see Saastamoinen, 1973)
ZHD = (0.0022768 * P)/(1 - 0.0026*cos(2*deg2rad(y0))-0.00000028*H_orth);
% ZWD using the Saastamoinen Tropospheric Model (see Saastamoinen, 1972, 1973)
ZWD = 0.002277 * ((1255/T + 0.05) * e);
% PWV
% Dimensionless constant
const = 10^8 / (461525 * (22.9744 + 375463/Tm));
PWV = ZWD * const;
return