0001 function write_nesting_bdy_file(Mobj, MJD, OutFile, MyTitle)
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 global ftbverbose
0041 report = false;
0042 if(ftbverbose); report = true; end
0043 subname = 'write_nesting_bdy_file';
0044 if(report); fprintf('\n'); end
0045 if(report); fprintf(['begin : ' subname '\n']); end
0046
0047
0048
0049
0050
0051 tmpObcNodes = Mobj.nnodesID';
0052
0053
0054 ObcNodes = tmpObcNodes(tmpObcNodes~=0)';
0055
0056
0057 ObcNodes = Mobj.nnodesID;
0058 nElems = Mobj.nElems;
0059
0060
0061
0062
0063 nTimes = numel(MJD);
0064 if(report); fprintf('Number of time steps %d\n',nTimes); end
0065
0066 nObcs = numel(ObcNodes);
0067 if(report); fprintf('Number of Open Boundary Nodes %d\n',nObcs); end
0068
0069 [chk1, chk2] = size(Mobj.surfaceElevation);
0070 if nObcs ~= chk1 || nTimes ~= chk2
0071 fprintf('surface Elevation dimensions do not match time series and number of boundary nodes.\n')
0072 fprintf('surface Elevation nodes and time sizes: (%d, %d)\n', chk1, chk2)
0073 fprintf('Boundary nodes size: %d\n', nObcs)
0074 fprintf('Times size: %d\n', nTimes)
0075 error('Input data sizes do not match. Check and try again.');
0076 end
0077
0078
0079
0080
0081
0082 nc=netcdf.create(OutFile,'clobber');
0083
0084
0085 netcdf.putAtt(nc,netcdf.getConstant('NC_GLOBAL'),'title',MyTitle)
0086 netcdf.putAtt(nc,netcdf.getConstant('NC_GLOBAL'),'type','FVCOM TIME SERIES NESTING FORCING FILE')
0087 netcdf.putAtt(nc,netcdf.getConstant('NC_GLOBAL'),'history', sprintf('File created using %s from the MATLAB fvcom-toolbox', subname))
0088 netcdf.putAtt(nc,netcdf.getConstant('NC_GLOBAL'),'Conventions','CF-1.0')
0089 netcdf.putAtt(nc,netcdf.getConstant('NC_GLOBAL'),'CoordinateSystem','Spherical')
0090
0091
0092
0093
0094
0095
0096
0097
0098 nlevels = Mobj.nsiglev;
0099 nlayers = nlevels - 1;
0100 nele_dimid = netcdf.defDim(nc,'nele',nElems);
0101 nobc_dimid = netcdf.defDim(nc,'node',nObcs);
0102 siglay_dimid= netcdf.defDim(nc,'siglay',nlayers);
0103 siglev_dimid= netcdf.defDim(nc,'siglev',nlevels);
0104 three_dimid = netcdf.defDim(nc,'three',3);
0105 time_dimid = netcdf.defDim(nc,'time',netcdf.getConstant('NC_UNLIMITED'));
0106 one_dimid = netcdf.defDim(nc,'one',1);
0107 date_str_len_dimid=netcdf.defDim(nc,'DateStrLen',26);
0108
0109
0110
0111
0112 disp('writing nprocs attributes');
0113 nprocs_varid=netcdf.defVar(nc,'nprocs','NC_INT',one_dimid);
0114 netcdf.putAtt(nc,nprocs_varid,'long_name','number of processors');
0115
0116 disp('writing npartition attributes');
0117 npart_varid=netcdf.defVar(nc,'partition','NC_INT',nele_dimid);
0118 netcdf.putAtt(nc,npart_varid,'long_name','partition');
0119
0120 disp('writing x attributes');
0121 x_varid=netcdf.defVar(nc,'x','NC_FLOAT',nobc_dimid);
0122 netcdf.putAtt(nc,x_varid,'long_name','nodal x-coordinate');
0123 netcdf.putAtt(nc,x_varid,'units','meters');
0124
0125 disp('writing y attributes');
0126 y_varid=netcdf.defVar(nc,'y','NC_FLOAT',nobc_dimid);
0127 netcdf.putAtt(nc,y_varid,'long_name','nodal y-coordinate');
0128 netcdf.putAtt(nc,y_varid,'units','meters');
0129
0130 disp('writing lon attributes');
0131 lon_varid=netcdf.defVar(nc,'lon','NC_FLOAT',nobc_dimid);
0132 netcdf.putAtt(nc,lon_varid,'long_name','nodal longitude');
0133 netcdf.putAtt(nc,lon_varid,'standard_name','longitude');
0134 netcdf.putAtt(nc,lon_varid,'units','degrees_east');
0135
0136 disp('writing lat attributes');
0137 lat_varid=netcdf.defVar(nc,'lat','NC_FLOAT',nobc_dimid);
0138 netcdf.putAtt(nc,lat_varid,'long_name','nodal latitude');
0139 netcdf.putAtt(nc,lat_varid,'standard_name','latitude');
0140 netcdf.putAtt(nc,lat_varid,'units','degrees_north');
0141
0142 disp('writing xc attributes');
0143 xc_varid=netcdf.defVar(nc,'xc','NC_FLOAT',nele_dimid);
0144 netcdf.putAtt(nc,xc_varid,'long_name','zonal x-coordinate');
0145 netcdf.putAtt(nc,xc_varid,'units','meters');
0146
0147 disp('writing yc attributes');
0148 yc_varid=netcdf.defVar(nc,'yc','NC_FLOAT',nele_dimid);
0149 netcdf.putAtt(nc,yc_varid,'long_name','zonal y-coordinate');
0150 netcdf.putAtt(nc,yc_varid,'units','meters');
0151
0152 disp('writing lonc attributes');
0153 lonc_varid=netcdf.defVar(nc,'lonc','NC_FLOAT',nele_dimid);
0154 netcdf.putAtt(nc,lonc_varid,'long_name','zonal longitude');
0155 netcdf.putAtt(nc,lonc_varid,'standard_name','longitude');
0156 netcdf.putAtt(nc,lonc_varid,'units','degrees_east');
0157
0158 disp('writing latc attributes');
0159 latc_varid=netcdf.defVar(nc,'latc','NC_FLOAT',nele_dimid);
0160 netcdf.putAtt(nc,latc_varid,'long_name','zonal latitude');
0161 netcdf.putAtt(nc,latc_varid,'standard_name','latitude');
0162 netcdf.putAtt(nc,latc_varid,'units','degrees_north');
0163
0164 disp('writing siglay attributes');
0165 siglay_varid=netcdf.defVar(nc,'siglay','NC_FLOAT',[nobc_dimid, siglay_dimid ]);
0166 netcdf.putAtt(nc,siglay_varid,'long_name','Sigma Layers');
0167 netcdf.putAtt(nc,siglay_varid,'standard_name','ocean_sigma/general_coordinate');
0168 netcdf.putAtt(nc,siglay_varid,'positive','up');
0169 netcdf.putAtt(nc,siglay_varid,'valid_min','-1.f');
0170 netcdf.putAtt(nc,siglay_varid,'valid_max','0.f');
0171 netcdf.putAtt(nc,siglay_varid,'formula_terms','sigma: siglay eta: zeta depth: h');
0172
0173 disp('writing siglev attributes');
0174 siglev_varid=netcdf.defVar(nc,'siglev','NC_FLOAT',[nobc_dimid,siglev_dimid ]);
0175 netcdf.putAtt(nc,siglev_varid,'long_name','Sigma Levels');
0176 netcdf.putAtt(nc,siglev_varid,'standard_name','ocean_sigma/general_coordinate');
0177 netcdf.putAtt(nc,siglev_varid,'positive','up');
0178 netcdf.putAtt(nc,siglev_varid,'valid_min','-1.f');
0179 netcdf.putAtt(nc,siglev_varid,'valid_max','0.f');
0180 netcdf.putAtt(nc,siglev_varid,'formula_terms','sigma: siglay eta: zeta depth: h');
0181
0182 disp('writing h attributes');
0183 h_varid=netcdf.defVar(nc,'h','NC_FLOAT',nobc_dimid);
0184 netcdf.putAtt(nc,h_varid,'long_name','Bathymetry');
0185 netcdf.putAtt(nc,h_varid,'standard_name','sea_floor_depth_below_geoid');
0186 netcdf.putAtt(nc,h_varid,'units','m');
0187 netcdf.putAtt(nc,h_varid,'positive','down');
0188 netcdf.putAtt(nc,h_varid,'grid','Bathymetry_Mesh');
0189 netcdf.putAtt(nc,h_varid,'coordinates','lat lon');
0190 netcdf.putAtt(nc,h_varid,'type','data');
0191
0192 disp('writing nv attributes');
0193
0194 nv_varid=netcdf.defVar(nc,'nv','NC_INT',[nele_dimid, three_dimid]);
0195 netcdf.putAtt(nc,nv_varid,'long_name','nodes surrounding element');
0196
0197
0198 disp('writing obc_nodes attributes');
0199 nobc_varid=netcdf.defVar(nc,'obc_nodes','NC_INT',nobc_dimid);
0200 netcdf.putAtt(nc,nobc_varid,'long_name','Open Boundary Node Number');
0201 netcdf.putAtt(nc,nobc_varid,'grid','obc_grid');
0202
0203 disp('writing iint attributes');
0204 iint_varid=netcdf.defVar(nc,'iint','NC_INT',time_dimid);
0205 netcdf.putAtt(nc,iint_varid,'long_name','internal mode iteration number');
0206
0207 disp('writing time attributes');
0208 time_varid=netcdf.defVar(nc,'time','NC_FLOAT',time_dimid);
0209 netcdf.putAtt(nc,time_varid,'long_name','time');
0210 netcdf.putAtt(nc,time_varid,'units','days since 1858-11-17 00:00:00');
0211 netcdf.putAtt(nc,time_varid,'format','modified julian day (MJD)');
0212 netcdf.putAtt(nc,time_varid,'time_zone','UTC');
0213
0214 disp('writing Itime attributes');
0215 itime_varid=netcdf.defVar(nc,'Itime','NC_INT',time_dimid);
0216 netcdf.putAtt(nc,itime_varid,'units','days since 1858-11-17 00:00:00');
0217 netcdf.putAtt(nc,itime_varid,'format','modified julian day (MJD)');
0218 netcdf.putAtt(nc,itime_varid,'time_zone','UTC');
0219
0220 disp('writing Itime2 attributes');
0221 itime2_varid=netcdf.defVar(nc,'Itime2','NC_INT',time_dimid);
0222 netcdf.putAtt(nc,itime2_varid,'units','msec since 00:00:00');
0223 netcdf.putAtt(nc,itime2_varid,'time_zone','UTC');
0224
0225 disp('writing Times attributes');
0226 Times_varid=netcdf.defVar(nc,'Times','NC_CHAR',[date_str_len_dimid, time_dimid]);
0227 netcdf.putAtt(nc,Times_varid,'time_zone','UTC');
0228
0229 disp('writing salinity attributes');
0230 salinity_varid=netcdf.defVar(nc,'salinity','NC_FLOAT',[nobc_dimid, siglay_dimid, time_dimid]);
0231 netcdf.putAtt(nc,salinity_varid,'long_name','salinity');
0232 netcdf.putAtt(nc,salinity_varid,'standard_name','sea_water_salinity');
0233 netcdf.putAtt(nc,salinity_varid,'units','1e-3');
0234 netcdf.putAtt(nc,salinity_varid,'grid','fvcom_grid');
0235 netcdf.putAtt(nc,salinity_varid,'coordinates','lat lon');
0236 netcdf.putAtt(nc,salinity_varid,'type','data');
0237
0238 disp('writing temperature attributes');
0239 temperature_varid=netcdf.defVar(nc,'temp','NC_FLOAT',[nobc_dimid,siglay_dimid, time_dimid]);
0240 netcdf.putAtt(nc,temperature_varid,'long_name','temperature');
0241 netcdf.putAtt(nc,temperature_varid,'standard_name','sea_water_temperature');
0242 netcdf.putAtt(nc,temperature_varid,'units','degrees_C');
0243 netcdf.putAtt(nc,temperature_varid,'grid','fvcom_grid');
0244 netcdf.putAtt(nc,temperature_varid,'coordinates','lat lon');
0245 netcdf.putAtt(nc,temperature_varid,'type','data');
0246
0247 disp('writing u velocity attributes');
0248 u_varid=netcdf.defVar(nc,'u','NC_FLOAT',[nele_dimid,siglay_dimid, time_dimid]);
0249 netcdf.putAtt(nc,u_varid,'long_name','Eastward Water Velocity');
0250 netcdf.putAtt(nc,u_varid,'units','meters s-1');
0251 netcdf.putAtt(nc,u_varid,'grid','fvcom_grid');
0252 netcdf.putAtt(nc,u_varid,'type','data');
0253
0254 disp('writing v velocity attributes');
0255 v_varid=netcdf.defVar(nc,'v','NC_FLOAT',[nele_dimid,siglay_dimid, time_dimid]);
0256 netcdf.putAtt(nc,v_varid,'long_name','Northward Water Velocity');
0257 netcdf.putAtt(nc,v_varid,'units','meters s-1');
0258 netcdf.putAtt(nc,v_varid,'grid','fvcom_grid');
0259 netcdf.putAtt(nc,v_varid,'type','data');
0260
0261 disp('writing vertical velocity attributes');
0262 omega_varid=netcdf.defVar(nc,'omega','NC_FLOAT',[nobc_dimid,siglev_dimid, time_dimid]);
0263 netcdf.putAtt(nc,omega_varid,'long_name','hydro static vertical velocity');
0264 netcdf.putAtt(nc,omega_varid,'units','meters s-1');
0265 netcdf.putAtt(nc,omega_varid,'grid','fvcom_grid');
0266 netcdf.putAtt(nc,omega_varid,'coordinates','lat lon');
0267 netcdf.putAtt(nc,omega_varid,'type','data');
0268
0269 disp('writing vertical velocity attributes');
0270 hyw_varid=netcdf.defVar(nc,'hyw','NC_FLOAT',[nobc_dimid,siglev_dimid, time_dimid]);
0271 netcdf.putAtt(nc,hyw_varid,'long_name','hydro static vertical velocity');
0272 netcdf.putAtt(nc,hyw_varid,'units','meters s-1');
0273 netcdf.putAtt(nc,hyw_varid,'grid','fvcom_grid');
0274 netcdf.putAtt(nc,hyw_varid,'coordinates','lat lon');
0275 netcdf.putAtt(nc,hyw_varid,'type','data');
0276
0277 disp('writing depth averaged u velocity attributes');
0278 ua_varid=netcdf.defVar(nc,'ua','NC_FLOAT',[nele_dimid, time_dimid]);
0279 netcdf.putAtt(nc,ua_varid,'long_name','Vertically Averaged x-velocity');
0280 netcdf.putAtt(nc,ua_varid,'units','meters s-1');
0281 netcdf.putAtt(nc,ua_varid,'grid','fvcom_grid');
0282 netcdf.putAtt(nc,ua_varid,'type','data');
0283
0284 disp('writing depth averaged v velocity attributes');
0285 va_varid=netcdf.defVar(nc,'va','NC_FLOAT',[nele_dimid, time_dimid]);
0286 netcdf.putAtt(nc,va_varid,'long_name','Vertically Averaged y-velocity');
0287 netcdf.putAtt(nc,va_varid,'units','meters s-1');
0288 netcdf.putAtt(nc,va_varid,'grid','fvcom_grid');
0289 netcdf.putAtt(nc,va_varid,'type','data');
0290
0291 disp('writing water level attributes');
0292 zeta_varid=netcdf.defVar(nc,'zeta','NC_FLOAT',[nobc_dimid, time_dimid]);
0293 netcdf.putAtt(nc,zeta_varid,'long_name','Water Surface Elevation');
0294 netcdf.putAtt(nc,zeta_varid,'units','meters');
0295 netcdf.putAtt(nc,zeta_varid,'positive','up');
0296 netcdf.putAtt(nc,zeta_varid,'standard_name','sea_surface_elevation');
0297 netcdf.putAtt(nc,zeta_varid,'grid','SSH_Mesh');
0298 netcdf.putAtt(nc,zeta_varid,'coordinates','lat lon');
0299 netcdf.putAtt(nc,zeta_varid,'type','data');
0300
0301
0302
0303 if (Mobj.NEST_TYPE=='TYPE3')
0304
0305 disp('writing node weighting attributes');
0306 node_weight_varid=netcdf.defVar(nc,'weight_node','NC_FLOAT',[nobc_dimid time_dimid]);
0307 netcdf.putAtt(nc,node_weight_varid,'long_name','nesting node weighting');
0308 netcdf.putAtt(nc,node_weight_varid,'units','none');
0309 netcdf.putAtt(nc,node_weight_varid,'grid','fvcom_grid');
0310 netcdf.putAtt(nc,node_weight_varid,'coordinates','lat lon');
0311 netcdf.putAtt(nc,node_weight_varid,'type','data');
0312
0313 disp('writing element weighting attributes');
0314 elem_weight_varid=netcdf.defVar(nc,'weight_cell','NC_FLOAT',[nele_dimid time_dimid]);
0315 netcdf.putAtt(nc,elem_weight_varid,'long_name','nesting element weighting');
0316 netcdf.putAtt(nc,elem_weight_varid,'units','none');
0317 netcdf.putAtt(nc,elem_weight_varid,'grid','fvcom_grid');
0318 netcdf.putAtt(nc,node_weight_varid,'coordinates','latc lonc');
0319 netcdf.putAtt(nc,elem_weight_varid,'type','data');
0320
0321
0322 end
0323
0324
0325 netcdf.endDef(nc);
0326
0327
0328
0329
0330 disp('writing nprocs - not used??');
0331 nprocs(1)= 3;
0332 netcdf.putVar(nc,nprocs_varid,nprocs);
0333
0334 disp('writing partition - not used??');
0335 partition = ones(1,nElems);
0336 netcdf.putVar(nc,npart_varid,partition);
0337
0338 disp('writing co-ordinate data at nodal points: x, y, lon, lat');
0339 netcdf.putVar(nc,x_varid,Mobj.x)
0340 netcdf.putVar(nc,y_varid,Mobj.y)
0341 netcdf.putVar(nc,lon_varid,Mobj.lon)
0342 netcdf.putVar(nc,lat_varid,Mobj.lat)
0343
0344 disp('writing co-ordinate data at element centres: xc, yc, lonc, latc');
0345 netcdf.putVar(nc,xc_varid,Mobj.xc)
0346 netcdf.putVar(nc,yc_varid,Mobj.yc)
0347 netcdf.putVar(nc,lonc_varid,Mobj.lonc)
0348 netcdf.putVar(nc,latc_varid,Mobj.latc)
0349
0350 disp('writing vertical schematization data: sigma layers, sigma levels');
0351 inc = 1./real(nlayers);
0352 siglev = 0:-inc:-1;
0353 for i=1:nlevels
0354 Msiglev(i,1:nObcs) = siglev(i);
0355 end;
0356 for i=1:nlayers
0357 Msiglay(i,1:nObcs) = mean(siglev(i:i+1));
0358 end;
0359
0360 Msiglay = -1*ones(nlayers,nObcs);
0361 tlayer = 1/nlayers;
0362 temp = 0.5*tlayer;
0363 for i=1:nlayers
0364 Msiglay(i,:)= temp*Msiglay(i,:);
0365 temp = temp + tlayer;
0366 end
0367
0368 Msiglev = -1*ones(nlevels,nObcs);
0369 temp = 0.0;
0370 for i=1:nlevels
0371 Msiglev(i,:)= temp*Msiglev(i,:);
0372 temp = temp + tlayer;
0373 end
0374
0375 netcdf.putVar(nc,siglay_varid,Msiglay)
0376 netcdf.putVar(nc,siglev_varid,Msiglev)
0377
0378 disp('writing h_varid data');
0379 netcdf.putVar(nc,h_varid,Mobj.h);
0380
0381 disp('writing nv_varid data');
0382 Mobj.nv = Mobj.nv';
0383 netcdf.putVar(nc,nv_varid,Mobj.nv);
0384
0385
0386 disp('writing nobc_varid data');
0387 netcdf.putVar(nc,nobc_varid,ObcNodes);
0388
0389 disp('writing iint_varid data');
0390 netcdf.putVar(nc,iint_varid,0,nTimes,1:nTimes);
0391
0392 disp('writing time_varid data');
0393 netcdf.putVar(nc,time_varid,0,nTimes,MJD);
0394
0395 disp('writing itime_varid data');
0396 netcdf.putVar(nc,itime_varid,floor(MJD));
0397
0398 disp('writing itime2_varid data');
0399 netcdf.putVar(nc,itime2_varid,0,nTimes,mod(MJD,1)*24*3600*1000);
0400
0401 nStringOut = char();
0402 for i=1:nTimes
0403 [nYr, nMon, nDay, nHour, nMin, nSec] = mjulian2greg(MJD(i));
0404 if strcmp(sprintf('%02i', nSec), '60')
0405
0406
0407
0408
0409 if mod(nMin + 1, 60) == 0
0410
0411 nHour = mod(nHour + 1, 24);
0412 end
0413 nMin = mod(nMin + 1, 60);
0414 nSec = 0;
0415 end
0416 nDate = [nYr, nMon, nDay, nHour, nMin, nSec];
0417 nStringOut = [nStringOut, sprintf('%04i/%02i/%02i %02i:%02i:%02i ',nDate)];
0418 end
0419
0420 disp('writing Times_varid data');
0421 netcdf.putVar(nc,Times_varid,nStringOut);
0422
0423
0424
0425
0426 disp('writing salinity_varid data');
0427 netcdf.putVar(nc,salinity_varid,Mobj.salinity);
0428
0429 disp('writing temperature_varid data');
0430 netcdf.putVar(nc,temperature_varid,Mobj.temperature);
0431
0432 disp('writing u_varid data');
0433 netcdf.putVar(nc,u_varid,Mobj.u);
0434
0435 disp('writing v_varid data');
0436 netcdf.putVar(nc,v_varid,Mobj.v);
0437
0438 disp('writing ua_varid data');
0439 netcdf.putVar(nc,ua_varid,Mobj.daUvel);
0440
0441 disp('writing va_varid data');
0442 netcdf.putVar(nc,va_varid,Mobj.daVvel);
0443
0444 disp('writing omega_varid data');
0445 netcdf.putVar(nc,omega_varid,Mobj.w);
0446
0447 disp('writing hyw_varid data');
0448 netcdf.putVar(nc,hyw_varid,Mobj.w);
0449
0450 disp('writing zeta_varid data');
0451 netcdf.putVar(nc,zeta_varid,Mobj.surfaceElevation);
0452
0453
0454 if (Mobj.NEST_TYPE=='TYPE3')
0455
0456 disp('writing node_weight_varid data');
0457 netcdf.putVar(nc,node_weight_varid,Mobj.weight_node_val);
0458
0459 disp('writing element weighting attributes');
0460 netcdf.putVar(nc,elem_weight_varid,Mobj.weight_elem_val);
0461
0462
0463
0464 end
0465
0466
0467
0468 netcdf.close(nc);
0469
0470 if(report); fprintf(['end : ' subname '\n']); end;
0471
0472
0473 end