0001 function write_FVCOM_forcing(Mobj, fileprefix, data, infos, fver, varargin)
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053
0054
0055
0056
0057
0058
0059
0060
0061
0062
0063
0064
0065
0066
0067
0068
0069
0070
0071
0072
0073
0074
0075
0076
0077
0078
0079
0080
0081
0082
0083
0084
0085
0086
0087
0088
0089
0090
0091
0092
0093
0094
0095
0096
0097
0098
0099
0100
0101
0102
0103
0104
0105
0106
0107
0108
0109
0110
0111
0112
0113
0114
0115
0116
0117
0118
0119
0120
0121
0122
0123
0124
0125
0126
0127
0128
0129
0130
0131 multi_out = false;
0132 if (nargin < 4 || nargin > 5) && isempty(varargin)
0133 error('Incorrect number of arguments')
0134 elseif nargin == 5
0135 if strcmpi(fver, '3.1.0')
0136 multi_out = true;
0137 end
0138 end
0139
0140 subname = 'write_FVCOM_forcing';
0141
0142 global ftbverbose;
0143 if ftbverbose
0144 fprintf('\nbegin : %s \n', subname)
0145 end
0146
0147
0148 strtime = true;
0149 inttime = false;
0150 floattime = false;
0151 for vv = 1:2:length(varargin)
0152 switch varargin{vv}
0153 case 'strtime'
0154 strtime = true;
0155 case 'inttime'
0156 inttime = true;
0157 case 'floattime'
0158 floattime = true;
0159 end
0160 end
0161
0162 tri = Mobj.tri;
0163 nNodes = Mobj.nVerts;
0164 nElems = Mobj.nElems;
0165 ntimes = numel(data.time);
0166
0167 if strcmpi(Mobj.nativeCoords, 'cartesian')
0168 x = Mobj.x;
0169 y = Mobj.y;
0170 else
0171 x = Mobj.lon;
0172 y = Mobj.lat;
0173 end
0174
0175 coordString = sprintf('FVCOM %s coordinates', Mobj.nativeCoords);
0176
0177
0178 xc = nodes2elems(x, Mobj);
0179 yc = nodes2elems(y, Mobj);
0180
0181
0182
0183
0184
0185 if multi_out
0186 suffixes = {'_wnd', '_hfx', '_evap', '_air_press'};
0187 else
0188 suffixes = {'_wnd'};
0189 end
0190
0191
0192
0193 nshf = 0;
0194
0195 for i=1:length(suffixes)
0196 nc = netcdf.create([fileprefix, suffixes{i}, '.nc'], 'clobber');
0197
0198
0199 netcdf.putAtt(nc,netcdf.getConstant('NC_GLOBAL'),'title','FVCOM Forcing File')
0200 netcdf.putAtt(nc,netcdf.getConstant('NC_GLOBAL'),'institution','Plymouth Marine Laboratory')
0201 netcdf.putAtt(nc,netcdf.getConstant('NC_GLOBAL'),'source','FVCOM grid (unstructured) surface forcing')
0202 netcdf.putAtt(nc,netcdf.getConstant('NC_GLOBAL'),'history', sprintf('File created with %s from the MATLAB fvcom-toolbox', subname))
0203 netcdf.putAtt(nc,netcdf.getConstant('NC_GLOBAL'),'references','http://fvcom.smast.umassd.edu, http://codfish.smast.umassd.edu')
0204 netcdf.putAtt(nc,netcdf.getConstant('NC_GLOBAL'),'Conventions','CF-1.0')
0205
0206 netcdf.putAtt(nc,netcdf.getConstant('NC_GLOBAL'),'CoordinateSystem',Mobj.nativeCoords)
0207 netcdf.putAtt(nc,netcdf.getConstant('NC_GLOBAL'),'CoordinateProjection','init=epsg:4326')
0208
0209
0210 nele_dimid=netcdf.defDim(nc,'nele',nElems);
0211 node_dimid=netcdf.defDim(nc,'node',nNodes);
0212 three_dimid=netcdf.defDim(nc,'three',3);
0213 time_dimid=netcdf.defDim(nc,'time',netcdf.getConstant('NC_UNLIMITED'));
0214 datestrlen_dimid=netcdf.defDim(nc,'DateStrLen',26);
0215
0216
0217 x_varid=netcdf.defVar(nc,'x','NC_FLOAT',node_dimid);
0218 netcdf.putAtt(nc,x_varid,'long_name','nodal x-coordinate');
0219 netcdf.putAtt(nc,x_varid,'units','meters');
0220
0221 y_varid=netcdf.defVar(nc,'y','NC_FLOAT',node_dimid);
0222 netcdf.putAtt(nc,y_varid,'long_name','nodal y-coordinate');
0223 netcdf.putAtt(nc,y_varid,'units','meters');
0224
0225 xc_varid=netcdf.defVar(nc,'xc','NC_FLOAT',nele_dimid);
0226 netcdf.putAtt(nc,xc_varid,'long_name','zonal x-coordinate');
0227 netcdf.putAtt(nc,xc_varid,'units','meters');
0228
0229 yc_varid=netcdf.defVar(nc,'yc','NC_FLOAT',nele_dimid);
0230 netcdf.putAtt(nc,yc_varid,'long_name','zonal y-coordinate');
0231 netcdf.putAtt(nc,yc_varid,'units','meters');
0232
0233 nv_varid=netcdf.defVar(nc,'nv','NC_INT',[nele_dimid, three_dimid]);
0234 netcdf.putAtt(nc,nv_varid,'long_name','nodes surrounding element');
0235
0236
0237 if floattime
0238 time_varid=netcdf.defVar(nc,'time','NC_FLOAT',time_dimid);
0239 netcdf.putAtt(nc,time_varid,'long_name','time');
0240 netcdf.putAtt(nc,time_varid,'units','days since 1858-11-17 00:00:00');
0241 netcdf.putAtt(nc,time_varid,'format','modified julian day (MJD)');
0242 netcdf.putAtt(nc,time_varid,'time_zone','UTC');
0243 end
0244
0245 if inttime
0246 itime_varid=netcdf.defVar(nc,'Itime','NC_INT',time_dimid);
0247 netcdf.putAtt(nc,itime_varid,'units','days since 1858-11-17 00:00:00');
0248 netcdf.putAtt(nc,itime_varid,'format','modified julian day (MJD)');
0249 netcdf.putAtt(nc,itime_varid,'time_zone','UTC');
0250
0251 itime2_varid=netcdf.defVar(nc,'Itime2','NC_INT',time_dimid);
0252 netcdf.putAtt(nc,itime2_varid,'units','msec since 00:00:00');
0253 netcdf.putAtt(nc,itime2_varid,'time_zone','UTC');
0254 netcdf.putAtt(nc,itime2_varid,'long_name','time');
0255 end
0256
0257 if strtime
0258 times_varid=netcdf.defVar(nc,'Times','NC_CHAR',[datestrlen_dimid,time_dimid]);
0259 netcdf.putAtt(nc,times_varid,'long_name','Calendar Date');
0260 netcdf.putAtt(nc,times_varid,'format','String: Calendar Time');
0261 netcdf.putAtt(nc,times_varid,'time_zone','UTC');
0262 end
0263
0264
0265
0266 fnames = fieldnames(data);
0267 used_varids = cell(0);
0268 used_fnames = cell(0);
0269 used_dims = cell(0);
0270
0271 for vv=1:length(fnames)
0272
0273
0274
0275 switch fnames{vv}
0276 case {'uwnd', 'u10'}
0277 if strcmpi(suffixes{i}, '_wnd') || ~multi_out
0278
0279
0280
0281 u10_varid=netcdf.defVar(nc,'U10','NC_FLOAT',[nele_dimid, time_dimid]);
0282 netcdf.putAtt(nc,u10_varid,'long_name','Eastward Wind Speed');
0283 netcdf.putAtt(nc,u10_varid,'units','m/s');
0284 netcdf.putAtt(nc,u10_varid,'grid','fvcom_grid');
0285 netcdf.putAtt(nc,u10_varid,'coordinates',coordString);
0286 netcdf.putAtt(nc,u10_varid,'type','data');
0287
0288 v10_varid=netcdf.defVar(nc,'V10','NC_FLOAT',[nele_dimid, time_dimid]);
0289 netcdf.putAtt(nc,v10_varid,'long_name','Northward Wind Speed');
0290 netcdf.putAtt(nc,v10_varid,'units','m/s');
0291 netcdf.putAtt(nc,v10_varid,'grid','fvcom_grid');
0292 netcdf.putAtt(nc,v10_varid,'coordinates',coordString);
0293 netcdf.putAtt(nc,v10_varid,'type','data');
0294
0295 uwind_varid=netcdf.defVar(nc,'uwind_speed','NC_FLOAT',[nele_dimid, time_dimid]);
0296 netcdf.putAtt(nc,uwind_varid,'long_name','Eastward Wind Speed');
0297 netcdf.putAtt(nc,uwind_varid,'standard_name','Wind Speed');
0298 netcdf.putAtt(nc,uwind_varid,'units','m/s');
0299 netcdf.putAtt(nc,uwind_varid,'grid','fvcom_grid');
0300 netcdf.putAtt(nc,uwind_varid,'type','data');
0301
0302 vwind_varid=netcdf.defVar(nc,'vwind_speed','NC_FLOAT',[nele_dimid, time_dimid]);
0303 netcdf.putAtt(nc,vwind_varid,'long_name','Northward Wind Speed');
0304 netcdf.putAtt(nc,vwind_varid,'standard_name','Wind Speed');
0305 netcdf.putAtt(nc,vwind_varid,'units','m/s');
0306 netcdf.putAtt(nc,vwind_varid,'grid','fvcom_grid');
0307 netcdf.putAtt(nc,vwind_varid,'type','data');
0308
0309
0310
0311 used_varids = [used_varids, {'u10_varid', 'v10_varid', 'uwind_varid', 'vwind_varid'}];
0312
0313
0314 if isfield(data, 'uwnd') && ~isfield(data, 'u10')
0315
0316 used_fnames = [used_fnames, {'uwnd', 'vwnd', 'uwnd', 'vwnd'}];
0317 elseif isfield(data, 'u10') && ~isfield(data, 'uwnd')
0318
0319 used_fnames = [used_fnames, {'u10', 'v10', 'u10', 'v10'}];
0320 elseif isfield(data, 'u10') && isfield(data, 'uwnd')
0321 error('Supply only one set of wind fields: ''uwnd'' and ''vwnd'' or ''u10'' and ''v10''.')
0322 else
0323 error('Unrecognised wind field names.')
0324 end
0325 used_dims = [used_dims, {'nElems', 'nElems', 'nElems', 'nElems'}];
0326 end
0327
0328 case {'vwnd', 'v10'}
0329
0330
0331 true;
0332
0333 case {'slp', 'pres'}
0334 if strcmpi(suffixes{i}, '_air_press') || ~multi_out
0335
0336 slp_varid=netcdf.defVar(nc,'air_pressure','NC_FLOAT',[node_dimid, time_dimid]);
0337 netcdf.putAtt(nc,slp_varid,'long_name','Surface air pressure');
0338 netcdf.putAtt(nc,slp_varid,'units','Pa');
0339 netcdf.putAtt(nc,slp_varid,'grid','fvcom_grid');
0340 netcdf.putAtt(nc,slp_varid,'coordinates',coordString);
0341 netcdf.putAtt(nc,slp_varid,'type','data');
0342
0343 used_varids = [used_varids, 'slp_varid'];
0344 used_fnames = [used_fnames, fnames{vv}];
0345 used_dims = [used_dims, 'nNodes'];
0346 end
0347
0348 case 'lcc'
0349 if strcmpi(suffixes{i}, '_low_cloud_cover') || ~multi_out
0350
0351 slp_varid=netcdf.defVar(nc,'low_cloud_cover','NC_FLOAT',[node_dimid, time_dimid]);
0352 netcdf.putAtt(nc,slp_varid,'long_name','Low cloud cover');
0353 netcdf.putAtt(nc,slp_varid,'units','Fraction');
0354 netcdf.putAtt(nc,slp_varid,'grid','fvcom_grid');
0355 netcdf.putAtt(nc,slp_varid,'coordinates',coordString);
0356 netcdf.putAtt(nc,slp_varid,'type','data');
0357
0358 used_varids = [used_varids, 'lcc_varid'];
0359 used_fnames = [used_fnames, fnames{vv}];
0360 used_dims = [used_dims, 'nNodes'];
0361 end
0362
0363 case 'shum'
0364 if strcmpi(suffixes{i}, '_specific_humidity') || ~multi_out
0365
0366 slp_varid=netcdf.defVar(nc,'specific_humidity','NC_FLOAT',[node_dimid, time_dimid]);
0367 netcdf.putAtt(nc,slp_varid,'long_name','Specific humidity');
0368 netcdf.putAtt(nc,slp_varid,'units','Kg kg^-1');
0369 netcdf.putAtt(nc,slp_varid,'grid','fvcom_grid');
0370 netcdf.putAtt(nc,slp_varid,'coordinates',coordString);
0371 netcdf.putAtt(nc,slp_varid,'type','data');
0372
0373 used_varids = [used_varids, 'shum_varid'];
0374 used_fnames = [used_fnames, fnames{vv}];
0375 used_dims = [used_dims, 'nNodes'];
0376 end
0377
0378 case {'Et', 'evap'}
0379 if strcmpi(suffixes{i}, '_evap') || ~multi_out
0380
0381 pevpr_varid=netcdf.defVar(nc,'evap','NC_FLOAT',[node_dimid, time_dimid]);
0382 netcdf.putAtt(nc,pevpr_varid,'long_name','Evaporation');
0383 netcdf.putAtt(nc,pevpr_varid,'description','Evaporation, ocean lose water is negative');
0384 netcdf.putAtt(nc,pevpr_varid,'units','m s-1');
0385 netcdf.putAtt(nc,pevpr_varid,'grid','fvcom_grid');
0386 netcdf.putAtt(nc,pevpr_varid,'coordinates',coordString);
0387 netcdf.putAtt(nc,pevpr_varid,'type','data');
0388
0389 used_varids = [used_varids, 'pevpr_varid'];
0390 used_fnames = [used_fnames, fnames{vv}];
0391 used_dims = [used_dims, 'nNodes'];
0392 end
0393
0394 case {'prate', 'P_E'}
0395 if strcmpi(suffixes{i}, '_evap') || ~multi_out
0396
0397 prate_varid=netcdf.defVar(nc,'precip','NC_FLOAT',[node_dimid, time_dimid]);
0398 netcdf.putAtt(nc,prate_varid,'long_name','Precipitation');
0399 netcdf.putAtt(nc,prate_varid,'description','Precipitation, ocean lose water is negative');
0400 netcdf.putAtt(nc,prate_varid,'units','m s-1');
0401 netcdf.putAtt(nc,prate_varid,'grid','fvcom_grid');
0402 netcdf.putAtt(nc,prate_varid,'coordinates',coordString);
0403 netcdf.putAtt(nc,prate_varid,'type','data');
0404
0405 used_varids = [used_varids, 'prate_varid'];
0406 used_fnames = [used_fnames, fnames{vv}];
0407 used_dims = [used_dims, 'nNodes'];
0408 end
0409
0410 case 'nswrs'
0411 if strcmpi(suffixes{i}, '_hfx') || ~multi_out
0412
0413 nswrs_varid=netcdf.defVar(nc,'short_wave','NC_FLOAT',[node_dimid, time_dimid]);
0414 netcdf.putAtt(nc,nswrs_varid,'long_name','Short Wave Radiation');
0415 netcdf.putAtt(nc,nswrs_varid,'units','W m-2');
0416 netcdf.putAtt(nc,nswrs_varid,'grid','fvcom_grid');
0417 netcdf.putAtt(nc,nswrs_varid,'coordinates',coordString);
0418 netcdf.putAtt(nc,nswrs_varid,'type','data');
0419
0420 used_varids = [used_varids, 'nswrs_varid'];
0421 used_fnames = [used_fnames, fnames{vv}];
0422 used_dims = [used_dims, 'nNodes'];
0423 end
0424
0425 case 'nlwrs'
0426 if strcmpi(suffixes{i}, '_hfx') || ~multi_out
0427
0428 nlwrs_varid=netcdf.defVar(nc,'long_wave','NC_FLOAT',[node_dimid, time_dimid]);
0429 netcdf.putAtt(nc,nlwrs_varid,'long_name','Long Wave Radiation');
0430 netcdf.putAtt(nc,nlwrs_varid,'units','W m-2');
0431 netcdf.putAtt(nc,nlwrs_varid,'grid','fvcom_grid');
0432 netcdf.putAtt(nc,nlwrs_varid,'coordinates',coordString);
0433 netcdf.putAtt(nc,nlwrs_varid,'type','data');
0434
0435 used_varids = [used_varids, 'nlwrs_varid'];
0436 used_fnames = [used_fnames, fnames{vv}];
0437 used_dims = [used_dims, 'nNodes'];
0438 end
0439
0440 case 'air'
0441 if strcmpi(suffixes{i}, '_hfx') || ~multi_out
0442
0443 airt_varid = netcdf.defVar(nc, 'air_temperature', 'NC_FLOAT', [node_dimid, time_dimid]);
0444 netcdf.putAtt(nc, airt_varid, 'long_name', 'Surface air temperature');
0445 netcdf.putAtt(nc, airt_varid, 'units', 'Celsius Degree');
0446 netcdf.putAtt(nc, airt_varid, 'grid', 'fvcom_grid');
0447 netcdf.putAtt(nc, airt_varid, 'coordinates', coordString);
0448 netcdf.putAtt(nc, airt_varid, 'type', 'data');
0449 used_varids = [used_varids, 'airt_varid'];
0450 used_fnames = [used_fnames, fnames{vv}];
0451 used_dims = [used_dims, 'nNodes'];
0452 end
0453
0454 case 'rhum'
0455 if strcmpi(suffixes{i}, '_hfx') || ~multi_out
0456
0457
0458 rhum_varid = netcdf.defVar(nc, 'relative_humidity', 'NC_FLOAT', [node_dimid, time_dimid]);
0459 netcdf.putAtt(nc, rhum_varid, 'long_name', 'surface air relative humidity');
0460 netcdf.putAtt(nc, rhum_varid, 'units', 'percentage');
0461 netcdf.putAtt(nc, rhum_varid, 'grid', 'fvcom_grid');
0462 netcdf.putAtt(nc, rhum_varid, 'coordinates', coordString);
0463 netcdf.putAtt(nc, rhum_varid, 'type', 'data');
0464
0465 used_varids = [used_varids, 'rhum_varid'];
0466 used_fnames = [used_fnames, fnames{vv}];
0467 used_dims = [used_dims, 'nNodes'];
0468 end
0469
0470 case 'dswrf'
0471 if strcmpi(suffixes{i}, '_hfx') || ~multi_out
0472
0473
0474 dswrf_varid = netcdf.defVar(nc, 'short_wave', 'NC_FLOAT', [node_dimid, time_dimid]);
0475 netcdf.putAtt(nc, dswrf_varid, 'long_name', 'Downward solar shortwave radiation flux');
0476 netcdf.putAtt(nc, dswrf_varid, 'units', 'Watts meter-2');
0477 netcdf.putAtt(nc, dswrf_varid, 'grid', 'fvcom_grid');
0478 netcdf.putAtt(nc, dswrf_varid, 'coordinates', coordString);
0479 netcdf.putAtt(nc, dswrf_varid, 'type', 'data');
0480
0481 used_varids = [used_varids, 'dswrf_varid'];
0482 used_fnames = [used_fnames, fnames{vv}];
0483 used_dims = [used_dims, 'nNodes'];
0484 end
0485
0486 case 'dlwrf'
0487 if strcmpi(suffixes{i}, '_hfx') || ~multi_out
0488
0489
0490 dlwrf_varid = netcdf.defVar(nc, 'long_wave', 'NC_FLOAT', [node_dimid, time_dimid]);
0491 netcdf.putAtt(nc, dlwrf_varid, 'long_name', 'Downward solar longwave radiation flux');
0492 netcdf.putAtt(nc, dlwrf_varid, 'units', 'Watts meter-2');
0493 netcdf.putAtt(nc, dlwrf_varid, 'grid', 'fvcom_grid');
0494 netcdf.putAtt(nc, dlwrf_varid, 'coordinates', coordString);
0495 netcdf.putAtt(nc, dlwrf_varid, 'type', 'data');
0496
0497 used_varids = [used_varids, 'dlwrf_varid'];
0498 used_fnames = [used_fnames, fnames{vv}];
0499 used_dims = [used_dims, 'nNodes'];
0500 end
0501
0502 case {'shtfl', 'lhtfl', 'nshf'}
0503
0504
0505
0506
0507
0508
0509
0510
0511
0512
0513 if strcmpi(fnames{vv}, 'nshf')
0514 nshf = 1;
0515 end
0516 try
0517
0518
0519
0520
0521 if strcmpi(suffixes{i}, '_hfx') || ~multi_out
0522
0523 nhf_varid=netcdf.defVar(nc,'net_heat_flux','NC_FLOAT',[node_dimid, time_dimid]);
0524 netcdf.putAtt(nc,nhf_varid,'long_name','Surface Net Heat Flux');
0525 netcdf.putAtt(nc,nhf_varid,'units','W m-2');
0526 netcdf.putAtt(nc,nhf_varid,'grid','fvcom_grid');
0527 netcdf.putAtt(nc,nhf_varid,'coordinates',coordString);
0528 netcdf.putAtt(nc,nhf_varid,'type','data');
0529 end
0530 end
0531
0532 if ~multi_out
0533
0534
0535 used_varids = [used_varids, 'nhf_varid'];
0536 used_fnames = [used_fnames, fnames{vv}];
0537 used_dims = [used_dims, 'nNodes'];
0538 end
0539
0540 case {'time', 'lon', 'lat', 'x', 'y'}
0541 continue
0542
0543 otherwise
0544 if ftbverbose
0545 warning('Unknown or possibly unused input data type: %s', fnames{vv})
0546 end
0547 end
0548 end
0549
0550
0551 netcdf.endDef(nc);
0552
0553
0554 netcdf.putVar(nc, nv_varid, tri);
0555 netcdf.putVar(nc,x_varid,x);
0556 netcdf.putVar(nc,y_varid,y);
0557 netcdf.putVar(nc,xc_varid,xc);
0558 netcdf.putVar(nc,yc_varid,yc);
0559
0560 if floattime
0561 netcdf.putVar(nc,time_varid,0,ntimes,data.time);
0562 end
0563 if inttime
0564 netcdf.putVar(nc,itime_varid,0,ntimes,floor(data.time));
0565
0566
0567
0568 netcdf.putVar(nc,itime2_varid,0,ntimes,round((mod(data.time,1)*24*3600*1000)/(3600*1000))*(3600*1000));
0569 end
0570
0571 if strtime
0572 nStringOut = char();
0573 [nYr, nMon, nDay, nHour, nMin, nSec] = mjulian2greg(data.time);
0574 for i=1:ntimes
0575 nDate = [nYr(i), nMon(i), nDay(i), nHour(i), nMin(i), nSec(i)];
0576 nStringOut = [nStringOut, sprintf('%04i/%02i/%02i %02i:%02i:%09.6f', nDate)];
0577 end
0578 netcdf.putVar(nc,times_varid,[0, 0], [26, ntimes], nStringOut);
0579 end
0580
0581
0582
0583
0584 hf_done = 0;
0585 wnd_done = 0;
0586 for ff = 1:length(used_fnames)
0587 if ftbverbose
0588 fprintf('write : %s... ', used_fnames{ff})
0589 end
0590 if strcmpi(used_fnames{ff}, 'shtfl') || strcmpi(used_fnames{ff}, 'lhtfl') || strcmpi(used_fnames{ff}, 'nlwrs') || strcmpi(used_fnames{ff}, 'nswrs')
0591
0592 hf_done = hf_done + 1;
0593
0594 if hf_done == 4 && nshf == 0
0595 if ftbverbose
0596 fprintf('combining heat flux ... ')
0597 end
0598
0599
0600
0601
0602 hf = data.nlwrs.node + data.nswrs.node - ...
0603 data.shtfl.node - data.lhtfl.node;
0604 netcdf.putVar(nc,nhf_varid,[0,0],[nNodes,ntimes],hf)
0605 elseif strcmpi(used_fnames{ff}, 'nswrs') || strcmpi(used_fnames{ff}, 'nlwrs')
0606
0607
0608
0609 if strcmpi(used_dims{ff}, 'nNodes')
0610 eval(['netcdf.putVar(nc,',used_varids{ff},',[0,0],[',used_dims{ff},',ntimes],data.',used_fnames{ff},'.node);'])
0611 else
0612 eval(['netcdf.putVar(nc,',used_varids{ff},',[0,0],[',used_dims{ff},',ntimes],data.',used_fnames{ff},'.data);'])
0613 end
0614 else
0615
0616
0617
0618
0619 end
0620 elseif strcmpi(used_fnames{ff}, 'nshf') && nshf == 1
0621 if ftbverbose
0622 fprintf('existing combined heat flux ... ')
0623 end
0624
0625
0626
0627
0628 hf_done = 4;
0629 netcdf.putVar(nc, nhf_varid, [0, 0], [nNodes, ntimes], data.nshf.node)
0630 else
0631
0632
0633 if strcmpi(used_dims{ff}, 'nNodes')
0634 eval(['netcdf.putVar(nc,',used_varids{ff},',[0,0],[',used_dims{ff},',ntimes],data.',used_fnames{ff},'.node);'])
0635 else
0636 try
0637 eval(['netcdf.putVar(nc,',used_varids{ff},',[0,0],[',used_dims{ff},',ntimes],data.',used_fnames{ff},'.data);'])
0638 wnd_done = wnd_done + 1;
0639 catch err
0640 fprintf('%s', err.message)
0641 end
0642 end
0643 end
0644 if ftbverbose
0645 fprintf('done.\n')
0646 end
0647 end
0648 if hf_done < 4 && nshf == 1
0649
0650
0651
0652 warning('Did not have all the required heat flux parameters for HEATING_ON. Need ''shtfl'', ''lhtfl'', ''nlwrs'' and ''nwsrs''.')
0653 end
0654
0655 if wnd_done < 2;
0656 warning('No wind data was provided (or one component was missing). Expected fields u10 and v10 or uwnd and vwnd.')
0657 end
0658
0659
0660 netcdf.close(nc);
0661 end
0662
0663 if ftbverbose
0664 fprintf('end : %s \n', subname)
0665 end