1 |
% |
2 |
% $Id: assemble_SH.m,v 1.2 2006/01/31 03:44:54 edhill Exp $ |
3 |
% |
4 |
% Ed Hill |
5 |
% |
6 |
% Generate the APE SH01 -- SH30 fields: |
7 |
% |
8 |
% For all time steps: |
9 |
% For all variables: |
10 |
% For all facets: |
11 |
% compute the "derived" fields: |
12 |
% sum: |
13 |
% Compute the averages from the sums: |
14 |
% Rotate vector components as necessary: |
15 |
% Convert units as necessary: |
16 |
% Regrid onto lat-lon: |
17 |
% Write netCDF output with all attributes: |
18 |
% |
19 |
|
20 |
%====================================================================== |
21 |
% |
22 |
% Define the connections between the APE "SH" variables and the |
23 |
% MITgcm diagnostics output |
24 |
% |
25 |
nfacets = 6; |
26 |
months.i_start = 7; |
27 |
months.i_end = 42; |
28 |
bname = 'data/SHave.t%03d.nc'; |
29 |
gname = 'data/grid.t%03d.nc'; |
30 |
r_con_name = 'scrip_regrid/rmp_CS32_to_LL128x64_conserv.nc'; |
31 |
r_dwt_name = 'scrip_regrid/rmp_CS32_to_LL128x64_distwgt.nc'; |
32 |
oname = { 'SH_fields.nc' }; |
33 |
|
34 |
flags.regrid_type_for_vec = 2; |
35 |
flags.debug = 0; |
36 |
|
37 |
% Get the following from: |
38 |
% ! (cd ../../ape_data_specs/ ; ./sh_parse.sh) |
39 |
vars = {}; |
40 |
vars{1} = {'SH01','sh_sw_toai','RADSWT','toa_incoming_shortwave_flux','W m-2','-','1'}; |
41 |
vars{2} = {'SH02','sh_sw_toar','OSR','toa_outgoing_shortwave_flux','W m-2','-','-1'}; |
42 |
vars{3} = {'SH04','sh_lw_toa','OLR','toa_outgoing_longwave_flux','W m-2','-','1'}; |
43 |
vars{4} = {'SH07','sh_cld_frac','CLDFRC','cloud_area_fraction','1','-','1'}; |
44 |
vars{5} = {'SH08','sh_cldw','---','atmosphere_cloud_condensed_water_co','kg m-2','-','1'}; |
45 |
vars{6} = {'SH09','sh_cldi','---','atmosphere_cloud_ice_content','kg m-2','-','1'}; |
46 |
vars{7} = {'SH11','sh_cppn','PRECON','convective_precipitation_flux','kg m-2 s-1','-','0.000011574074'}; |
47 |
vars{8} = {'SH12','sh_dppn','SH12','large_scale_precipitation_flux','kg m-2 s-1','-','0.000011574074'}; |
48 |
vars{9} = {'SH13','sh_evap','EVAP','evaporation_flux','kg m-2 s-1','-','0.000011574074'}; |
49 |
vars{10} = {'SH15','sh_sswi','SH15','surface_downwelling_shortwave_flux','W m-2','-','1'}; |
50 |
vars{11} = {'SH16','sh_sswr','SH16','surface_upwelling_shortwave_flux','W m-2','-','-1'}; |
51 |
vars{12} = {'SH18','sh_slwd','LWGDOWN','surface_downwelling_longwave_flux','W m-2','-','-1'}; |
52 |
vars{13} = {'SH19','sh_slwu','LWGUP','surface_upwelling_longwave_flux','W m-2','-','1'}; |
53 |
vars{14} = {'SH21','sh_slh','EFLUX','surface_upward_latent_heat_flux','W m-2','-','1'}; |
54 |
vars{15} = {'SH22','sh_ssh','HFLUX','surface_upward_sensible_heat_flux','W m-2','-','1'}; |
55 |
vars{16} = {'SH24','sh_ts2m','T2M','surface_air_temperature','K','-','1'}; |
56 |
vars{17} = {'SH25','sh_q2m','Q2M','specific_humidity','1','-','0.001'}; |
57 |
vars{18} = {'SH26','sh_tauu','UFLUX','surface_downward_eastward_stress','Pa','u','-1'}; |
58 |
vars{19} = {'SH27','sh_tauv','VFLUX','surface_downward_northward_stress','Pa','v','-1'}; |
59 |
vars{20} = {'SH28','sh_u10m','U10M','eastward_wind','m s-1','u','1'}; |
60 |
vars{21} = {'SH29','sh_v10m','V10M','northward_wind','m s-1','v','1'}; |
61 |
vars{22} = {'SH30','sh_ps','PS','surface_air_pressure','Pa','-','100'}; |
62 |
|
63 |
|
64 |
|
65 |
%====================================================================== |
66 |
% |
67 |
% Open the input files and get the (minimum) number of time steps |
68 |
% available across all the files |
69 |
% |
70 |
disp('Finding available time steps') |
71 |
|
72 |
ncl = {}; |
73 |
nTmin = 0; |
74 |
ifacet = 1; |
75 |
for ifacet = 1:nfacets |
76 |
fname = sprintf(bname,ifacet); |
77 |
nc = netcdf(fname, 'nowrite'); |
78 |
ncl{ifacet} = nc; |
79 |
T = nc{'T'}(:); |
80 |
nT = length(T); |
81 |
if (ifacet == 1) |
82 |
nTmin = nT; |
83 |
else |
84 |
nTmin = min(nTmin,nT); |
85 |
end |
86 |
nc = []; |
87 |
end |
88 |
disp(sprintf(' total iterations found = %d',nTmin)); |
89 |
|
90 |
|
91 |
%====================================================================== |
92 |
% |
93 |
% Compute the derived fields and the temporal sums of all fields on |
94 |
% the original MITgcm grid |
95 |
% |
96 |
disp('Computing Sums') |
97 |
|
98 |
dat = []; |
99 |
iz = 1; |
100 |
it = 1; |
101 |
iv = 1; |
102 |
ifacet = 1; |
103 |
for it = [ 1 ] % months.i_start:min(nTmin,months.i_end) |
104 |
disp(sprintf(' iteration = %d',it)); |
105 |
for iv = 1:length(vars) |
106 |
for ifacet = 1:nfacets |
107 |
|
108 |
% cull unknown variables |
109 |
if ( vars{iv}{3}(1) == '-' ) |
110 |
continue |
111 |
end |
112 |
% disp(sprintf('%d, %d',iv,ifacet)); |
113 |
|
114 |
% skip over computed quantities |
115 |
if strcmp( vars{iv}{1} , vars{iv}{3} ) |
116 |
continue |
117 |
end |
118 |
|
119 |
% WARN about missing variables |
120 |
nc = ncl{ifacet}; |
121 |
if prod(size(nc{vars{iv}{3}})) == 0 |
122 |
if it == months.i_start |
123 |
str = 'var "%s" does not exist in file "%s"'; |
124 |
disp(sprintf([' Warning: ' str], vars{iv}{3}, name(nc))); |
125 |
end |
126 |
continue |
127 |
end |
128 |
|
129 |
if it == 1 |
130 |
dat(iv,ifacet).n = 1; |
131 |
dat(iv,ifacet).a = squeeze( nc{ vars{iv}{3} }(it,iz,:,:) ); |
132 |
else |
133 |
dat(iv,ifacet).n = dat(iv,ifacet).n + 1; |
134 |
dat(iv,ifacet).a = dat(iv,ifacet).a + ... |
135 |
squeeze( nc{ vars{iv}{3} }(it,iz,:,:) ); |
136 |
end |
137 |
|
138 |
end |
139 |
end |
140 |
|
141 |
% fill in the derived fields |
142 |
for iv = 1:length(vars) |
143 |
for ifacet = 1:nfacets |
144 |
tmp = []; |
145 |
nc = ncl{ifacet}; |
146 |
|
147 |
if strncmp(vars{iv}{1},'SH12',4) |
148 |
% SH12 = PREACC - PRECON |
149 |
tmp = squeeze( nc{ 'PREACC' }(it,iz,:,:) ) ... |
150 |
- squeeze( nc{ 'PRECON' }(it,iz,:,:) ); |
151 |
end |
152 |
|
153 |
if strncmp(vars{iv}{1},'SH15',4) || strncmp(vars{iv}{1},'SH16',4) |
154 |
albedo = ( squeeze( nc{ 'ALBNIRDF' }(it,iz,:,:)) ... |
155 |
+ squeeze( nc{ 'ALBVISDF' }(it,iz,:,:)) ) ./ 2.0; |
156 |
if vars{iv}{1} == 'SH15' |
157 |
% SH15 = RADSWG/(1 - ALBEDO) |
158 |
tmp = squeeze( nc{ 'RADSWG' }(it,iz,:,:) ) ... |
159 |
./ (1.0 - albedo); |
160 |
end |
161 |
if strncmp(vars{iv}{1},'SH16',4) |
162 |
% SH16 = RADSWG*(1/(1 - ALBEDO) - 1) |
163 |
tmp = squeeze( nc{ 'RADSWG' }(it,iz,:,:) ) ... |
164 |
.* (1.0./(1.0 - albedo) - 1.0); |
165 |
end |
166 |
end |
167 |
|
168 |
if prod(size(tmp)) > 0 |
169 |
if it == 1 |
170 |
dat(iv,ifacet).n = 1; |
171 |
dat(iv,ifacet).a = tmp; |
172 |
else |
173 |
dat(iv,ifacet).n = dat(iv,ifacet).n + 1; |
174 |
dat(iv,ifacet).a = dat(iv,ifacet).a + tmp; |
175 |
end |
176 |
end |
177 |
|
178 |
end |
179 |
end |
180 |
end |
181 |
clear tmp |
182 |
|
183 |
% close the input files |
184 |
for ifacet = 1:nfacets |
185 |
nc = close( ncl{ifacet} ); |
186 |
end |
187 |
clear ncl; |
188 |
|
189 |
%====================================================================== |
190 |
% |
191 |
% Compute the averages from the temporal sums |
192 |
% |
193 |
disp('Computing Averages from Sums') |
194 |
|
195 |
ave = []; |
196 |
for iv = 1:length(vars) |
197 |
ind = length(ave) + 1; |
198 |
for ifacet = 1:nfacets |
199 |
if length(dat(iv,ifacet).n) == 0 |
200 |
continue |
201 |
end |
202 |
ave(ind).v(:,:,ifacet) = dat(iv,ifacet).a ./ dat(iv,ifacet).n; |
203 |
ave(ind).ivar = iv; |
204 |
end |
205 |
end |
206 |
|
207 |
clear dat; |
208 |
|
209 |
%====================================================================== |
210 |
% |
211 |
% Rotate vector components & convert units |
212 |
% |
213 |
clear ig; |
214 |
gvars = { 'XC','YC','dxF','dyF','rA','XG','YG','dxV', ... |
215 |
'dyU','rAz','dxC','dyC','rAw','rAs','dxG','dyG' }; |
216 |
for iface = 1:6 |
217 |
fname = sprintf(gname, iface); |
218 |
nc = netcdf(fname, 'nowrite'); |
219 |
if iface == 1 |
220 |
ig.ne = length(nc('X')); |
221 |
end |
222 |
for jj = 1:length(gvars) |
223 |
comm = [ 'ig.' gvars{jj} '(:,:,iface) = nc{''' ... |
224 |
gvars{jj} '''}(:);']; |
225 |
eval(comm); |
226 |
end |
227 |
nc = close(nc); |
228 |
end |
229 |
|
230 |
ne = ig.ne; |
231 |
n1 = ne - 1; |
232 |
dux = zeros(size(ig.XC)); |
233 |
duy = zeros(size(ig.XC)); |
234 |
dvx = zeros(size(ig.XC)); |
235 |
dvy = zeros(size(ig.XC)); |
236 |
% dux(:,:,:) = diff(ig.XG(:,1:ne,:),1,1); |
237 |
% dvx(:,:,:) = diff(ig.XG(1:ne,:,:),1,2); |
238 |
% duy(:,:,:) = diff(ig.YG(:,1:ne,:),1,1); |
239 |
% dvy(:,:,:) = diff(ig.YG(1:ne,:,:),1,2); |
240 |
% |
241 |
% Note the first two dimensions below are permuted relative to the |
242 |
% ordering in the *.mitgrid files (commented-out above). |
243 |
dux(:,:,:) = diff(ig.XG(1:ne,:,:),1,2); |
244 |
dvx(:,:,:) = diff(ig.XG(:,1:ne,:),1,1); |
245 |
duy(:,:,:) = diff(ig.YG(1:ne,:,:),1,2); |
246 |
dvy(:,:,:) = diff(ig.YG(:,1:ne,:),1,1); |
247 |
dux = dux + 360*double(dux < 180); |
248 |
dux = dux - 360*double(dux > 180); % squeeze([ min(min(dux)) max(max(dux)) ]) |
249 |
duy = duy + 360*double(duy < 180); |
250 |
duy = duy - 360*double(duy > 180); % squeeze([ min(min(duy)) max(max(duy)) ]) |
251 |
dvx = dvx + 360*double(dvx < 180); |
252 |
dvx = dvx - 360*double(dvx > 180); % squeeze([ min(min(dvx)) max(max(dvx)) ]) |
253 |
dvy = dvy + 360*double(dvy < 180); |
254 |
dvy = dvy - 360*double(dvy > 180); % squeeze([ min(min(dvy)) max(max(dvy)) ]) |
255 |
dut = sqrt(dux.^2 + duy.^2); % squeeze([ min(min(dut)) max(max(dut)) ]) |
256 |
dvt = sqrt(dvx.^2 + dvy.^2); % squeeze([ min(min(dvt)) max(max(dvt)) ]) |
257 |
ig.llux = dux ./ dut; % squeeze([ min(min(ig.llux)) max(max(ig.llux)) ]) |
258 |
ig.lluy = duy ./ dut; % squeeze([ min(min(ig.lluy)) max(max(ig.lluy)) ]) |
259 |
ig.llvx = dvx ./ dvt; % squeeze([ min(min(ig.llvx)) max(max(ig.llvx)) ]) |
260 |
ig.llvy = dvy ./ dvt; % squeeze([ min(min(ig.llvy)) max(max(ig.llvy)) ]) |
261 |
clear dux duy dvx dvy dut dvt ne n1 |
262 |
|
263 |
if flags.regrid_type_for_vec ~= 2 |
264 |
% Vector rotation |
265 |
for ia = 1:length(ave) |
266 |
if strcmp(vars{ ave(ia).ivar }{6}, 'u') |
267 |
% llu = u .* llux + v .* llvx; |
268 |
% llv = u .* lluy + v .* llvy; |
269 |
llu = ave(ia).v .* ig.llux + ave(ia+1).v .* ig.llvx; |
270 |
llv = ave(ia).v .* ig.lluy + ave(ia+1).v .* ig.llvy; |
271 |
ave(ia).v = llu; |
272 |
ave(ia+1).v = llv; |
273 |
end |
274 |
end |
275 |
end |
276 |
clear llu llv |
277 |
|
278 |
% Units conversion |
279 |
for ia = 1:length(ave) |
280 |
if strcmp(vars{ ave(ia).ivar }{7}, '-') ~= 1 |
281 |
eval(['fac = ' vars{ ave(ia).ivar }{7} ';']); |
282 |
if fac ~= 1 |
283 |
ave(ia).v = fac .* ave(ia).v; |
284 |
end |
285 |
end |
286 |
end |
287 |
|
288 |
%====================================================================== |
289 |
% |
290 |
% Regrid |
291 |
% |
292 |
% JMC suggested an evenly spaced Lat-Lon at: 128x64 |
293 |
% |
294 |
disp('Regridding') |
295 |
|
296 |
og = []; |
297 |
og.nlat = 64; |
298 |
og.nlon = 128; |
299 |
og.latcell = 180/og.nlat; |
300 |
og.loncell = 360/og.nlon; |
301 |
og.lat = linspace(-90+og.latcell/2, 90-og.latcell/2, og.nlat); |
302 |
og.lon = linspace(0+og.loncell/2, 360-og.loncell/2, og.nlon); |
303 |
og.latbnd(:,1) = og.lat - og.latcell/2.0; |
304 |
og.latbnd(:,2) = og.lat + og.latcell/2.0; |
305 |
og.lonbnd(:,1) = og.lon - og.loncell/2.0; |
306 |
og.lonbnd(:,2) = og.lon + og.loncell/2.0; |
307 |
|
308 |
rgvar = { 'src_address', 'dst_address', 'remap_matrix' }; |
309 |
rtype = { { 'con', r_con_name}, {'dwt', r_dwt_name } }; |
310 |
for kk = 1:length(rtype) |
311 |
nc = netcdf(rtype{kk}{2}, 'nowrite'); |
312 |
for jj = 1:length(rgvar) |
313 |
comm = sprintf('rgrid.%s.%s = nc{''%s''}(:);',... |
314 |
rtype{kk}{1},rgvar{jj},rgvar{jj}); |
315 |
eval(comm); |
316 |
% disp(comm); |
317 |
end |
318 |
nc = close(nc); |
319 |
end |
320 |
|
321 |
if ( flags.regrid_type_for_vec > 0 ) |
322 |
% >> x=rdmds('XC'); |
323 |
% >> y=rdmds('YC'); |
324 |
% >> t=rdmds('Ttave.0000513360'); |
325 |
% >> xi=-179:2:180;yi=-89:2:90; |
326 |
% >> del=cube2latlon_preprocess(x,y,xi,yi); |
327 |
% >> ti=cube2latlon_fast(del,t); |
328 |
aa_regrid.x = rdmds('aa_regrid/XC'); |
329 |
aa_regrid.y = rdmds('aa_regrid/YC'); |
330 |
% u = rdmds('aa_regrid/U.0000000000'); |
331 |
% v = rdmds('aa_regrid/V.0000000000'); |
332 |
% del = cube2latlon_preprocess(LON,LAT,xc,yc); |
333 |
% |
334 |
% aa_regrid.del = ... |
335 |
% cube2latlon_preprocess(aa_regrid.x,aa_regrid.y, ... |
336 |
% og.lon,og.lat); |
337 |
end |
338 |
|
339 |
% save zzz |
340 |
% zzz |
341 |
|
342 |
ii = 1; |
343 |
for ii = 1:length(ave) |
344 |
|
345 |
dcp = zeros(prod(size( ave(ii).v )),1); |
346 |
nn = 0; |
347 |
for kk = 1:size( ave(ii).v, 3 ) |
348 |
for jk = 1:size( ave(ii).v, 2 ) |
349 |
for ik = 1:size( ave(ii).v, 1 ) |
350 |
nn = nn + 1; |
351 |
dcp(nn) = ave(ii).v(ik,jk,kk); |
352 |
end |
353 |
end |
354 |
end |
355 |
lld = zeros(og.nlat*og.nlon,1); |
356 |
|
357 |
if ( ( flags.regrid_type_for_vec > 0 ) ... |
358 |
&& ( vars{ave(ii).ivar}{6} == 'u' ... |
359 |
|| vars{ave(ii).ivar}{6} == 'v' ) ) |
360 |
|
361 |
if flags.regrid_type_for_vec == 1 |
362 |
|
363 |
% use distributed weights for vectors |
364 |
for jj = 1:length(rgrid.dwt.src_address) |
365 |
lld(rgrid.dwt.dst_address(jj)) = ... |
366 |
lld(rgrid.dwt.dst_address(jj)) ... |
367 |
+ ( dcp(rgrid.dwt.src_address(jj)) ... |
368 |
* rgrid.dwt.remap_matrix(jj,1) ); |
369 |
end |
370 |
ave(ii).r = reshape(lld',og.nlat,og.nlon); |
371 |
|
372 |
elseif ( flags.regrid_type_for_vec == 2 ... |
373 |
&& vars{ave(ii).ivar}{6} == 'u' ) |
374 |
|
375 |
% use Alistair's uvcube2latlon() for vectors |
376 |
nsiz = size(ave(ii).v); |
377 |
by6u = zeros([ nsiz(1)*nsiz(3) nsiz(2) ]); |
378 |
by6v = zeros([ nsiz(1)*nsiz(3) nsiz(2) ]); |
379 |
for kk = 1:size( ave(ii).v, 3 ) |
380 |
is = 1 + (kk-1)*nsiz(1); |
381 |
ie = is + nsiz(1) - 1; |
382 |
by6u(is:ie,:) = ave(ii ).v(:,:,kk)'; |
383 |
by6v(is:ie,:) = ave(ii+1).v(:,:,kk)'; |
384 |
end |
385 |
% xi=-179:2:180; yi=-89:2:90; |
386 |
xi = og.lon; |
387 |
ind = find(xi > 179.5); |
388 |
xi(ind) = xi(ind) - 360; |
389 |
yi = og.lat; |
390 |
[ul,vl] = uvacube2latlon(aa_regrid.x,aa_regrid.y, ... |
391 |
by6u,by6v, xi,yi); |
392 |
ave(ii ).r = ul'; |
393 |
ave(ii+1).r = vl'; |
394 |
|
395 |
if flags.debug > 0 |
396 |
figure(1), subplot(1,1,1) |
397 |
subplot(2,1,1), surf(ul'), axis equal, view(2), shading flat |
398 |
subplot(2,1,2), surf(vl'), axis equal, view(2), shading flat |
399 |
end |
400 |
|
401 |
elseif ( vars{ave(ii).ivar}{6} == 'v' ) |
402 |
|
403 |
disp(sprintf( ' %s was computed along with %s', ... |
404 |
vars{ave(ii).ivar}{2}, vars{ave(ii-1).ivar}{2} )); |
405 |
|
406 |
elseif ( vars{ave(ii).ivar}{6} == 'u' ) |
407 |
|
408 |
disp('Error: please specify regrid type for vectors'); |
409 |
disp(' using "flags.regrid_type_for_vec"'); |
410 |
|
411 |
end |
412 |
else |
413 |
% conservative for scalars |
414 |
for jj = 1:length(rgrid.con.src_address) |
415 |
lld(rgrid.con.dst_address(jj)) = ... |
416 |
lld(rgrid.con.dst_address(jj)) ... |
417 |
+ ( dcp(rgrid.con.src_address(jj)) ... |
418 |
* rgrid.con.remap_matrix(jj,1) ); |
419 |
end |
420 |
ave(ii).r = reshape(lld',og.nlat,og.nlon); |
421 |
end |
422 |
|
423 |
if flags.debug > 0 |
424 |
figure(2) |
425 |
surf(og.lon,og.lat,ll128x64), view(2), shading flat |
426 |
end |
427 |
|
428 |
end |
429 |
|
430 |
|
431 |
%====================================================================== |
432 |
% |
433 |
% Write netCDF output |
434 |
% |
435 |
disp('Writing netCDF output for SH') |
436 |
|
437 |
nc = netcdf(oname{1}, 'clobber'); |
438 |
nc.title = [ 'Aqua Planet: ' ... |
439 |
'Single-Level 2-D Means from "Example" Experiment' ]; |
440 |
nc.institution = 'MIT Dept. of EAPS, Cambridge, MA, USA'; |
441 |
nc.source = [ 'MITgcm: ' ]; |
442 |
nc.Conventions = 'CF-1.0'; |
443 |
nc.history = [ 'Original data produced: ' '2002/08/20' ]; |
444 |
|
445 |
nc('lon') = og.nlon; |
446 |
nc('lat') = og.nlat; |
447 |
nc('bnd') = 2; |
448 |
|
449 |
nc{'lon'} = ncdouble('lon'); |
450 |
nc{'lon'}.standard_name = 'longitude'; |
451 |
nc{'lon'}.units = 'degrees_east'; |
452 |
nc{'lon'}.bounds = 'lon_bnds'; |
453 |
nc{'lon'}(:) = og.lon; |
454 |
|
455 |
nc{'lat'} = ncdouble('lat'); |
456 |
nc{'lat'}.standard_name = 'latitude'; |
457 |
nc{'lat'}.units = 'degrees_north'; |
458 |
nc{'lat'}.bounds = 'lat_bnds'; |
459 |
nc{'lat'}(:) = og.lat; |
460 |
|
461 |
nc{'lon_bnds'} = ncdouble('lon','bnd'); |
462 |
nc{'lon_bnds'}.ape_name = 'longitude cell bounds'; |
463 |
nc{'lon_bnds'}.units = 'degrees_east'; |
464 |
nc{'lon_bnds'}(:) = og.lonbnd; |
465 |
|
466 |
nc{'lat_bnds'} = ncdouble('lat','bnd'); |
467 |
nc{'lat_bnds'}.ape_name = 'latitude cell bounds'; |
468 |
nc{'lat_bnds'}.units = 'degrees_north'; |
469 |
nc{'lat_bnds'}(:) = og.latbnd; |
470 |
|
471 |
for ii = 1:length(ave) |
472 |
|
473 |
iv = ave(ii).ivar; |
474 |
nc{ vars{iv}{2} } = ncfloat( 'lat', 'lon' ); |
475 |
nc{ vars{iv}{2} }.ape_name = vars{iv}{4}; |
476 |
nc{ vars{iv}{2} }.units = vars{iv}{5}; |
477 |
nc{ vars{iv}{2} }.FillValue_ = 1.0e20; |
478 |
|
479 |
% deal with missing values |
480 |
ind = find(isnan( ave(ii).r )); |
481 |
ave(ii).r(ind) = 1.0e20; |
482 |
|
483 |
nc{ vars{iv}{2} }(:) = ave(ii).r; |
484 |
|
485 |
end |
486 |
|
487 |
nc = close(nc); |
488 |
|