


Creates an FVCOM nesting file.
function modify_FVCOM_nested_forcing(conf, ncfile, nesttype)
DESCRIPTION:
Modifies nested file to change weights or number of levels and writes a
new Nested file suitable for types different from 3. you can only do
one or the other. Changing weights and reducing the number of levels is
not permitted in the current implementation. You can always cut the
number of levels and then call it again to change the values of the
weights.
INPUT:
conf = with either new_weights field or conf.nest.levels and
method
ncfile2read = full path to the nesting file to be created.
nesttype = [optional] nesting type (defaults to 1 = direct nesting).
OUTPUT:
FVCOM nesting file.
EXAMPLE USAGE:
conf.method = one of weights, levels or positions
conf.new_weight_cell = weights see manual for explanation [0-1]. It
only requires nlevel values, not spatial or temporal dimension is
expected.
conf.new_weight_node = weights see manual for explanation [0-1]. It
only requires nlevel values, not spatial or temporal dimension is
expected.
conf.nest.levels = [1]; The value reflect the inner most level to keep.
1 will keep the most external boundary, 4 will keep 1-4 levels.
conf.nest.latlon_nodes_file = '/users/modellers/rito/Models/FVCOM/fvcom-projects/rosa/run/input/aqua_v16/tapas_nodes_v0latlon_v16_nest.dat'
conf.nest.latlon_elems_file = '/users/modellers/rito/Models/FVCOM/fvcom-projects/rosa/run/input/aqua_v16/tapas_elems_v0latlon_v16_nest.dat'
modify_FVCOM_nested_forcing(conf, '/tmp/fvcom_nested.nc', 1)
Author(s):
Ricardo Torres (Plymouth Marine Laboratory)
Pierre Cazenave (Plymouth Marine Laboratory)
Revision history:
2017-01-17 First version based on Pierre's write_FVCOM_nested_forcing.m
script. Compress the time series data to save space. Requires
netCDF4 in FVCOM.
==========================================================================

0001 function modify_FVCOM_nested_forcing(conf, ncfile2read, nesttype) 0002 % Creates an FVCOM nesting file. 0003 % 0004 % function modify_FVCOM_nested_forcing(conf, ncfile, nesttype) 0005 % 0006 % DESCRIPTION: 0007 % Modifies nested file to change weights or number of levels and writes a 0008 % new Nested file suitable for types different from 3. you can only do 0009 % one or the other. Changing weights and reducing the number of levels is 0010 % not permitted in the current implementation. You can always cut the 0011 % number of levels and then call it again to change the values of the 0012 % weights. 0013 % 0014 % 0015 % INPUT: 0016 % conf = with either new_weights field or conf.nest.levels and 0017 % method 0018 % ncfile2read = full path to the nesting file to be created. 0019 % nesttype = [optional] nesting type (defaults to 1 = direct nesting). 0020 % 0021 % OUTPUT: 0022 % FVCOM nesting file. 0023 % 0024 % EXAMPLE USAGE: 0025 % conf.method = one of weights, levels or positions 0026 % conf.new_weight_cell = weights see manual for explanation [0-1]. It 0027 % only requires nlevel values, not spatial or temporal dimension is 0028 % expected. 0029 % conf.new_weight_node = weights see manual for explanation [0-1]. It 0030 % only requires nlevel values, not spatial or temporal dimension is 0031 % expected. 0032 % conf.nest.levels = [1]; The value reflect the inner most level to keep. 0033 % 1 will keep the most external boundary, 4 will keep 1-4 levels. 0034 % conf.nest.latlon_nodes_file = '/users/modellers/rito/Models/FVCOM/fvcom-projects/rosa/run/input/aqua_v16/tapas_nodes_v0latlon_v16_nest.dat' 0035 % conf.nest.latlon_elems_file = '/users/modellers/rito/Models/FVCOM/fvcom-projects/rosa/run/input/aqua_v16/tapas_elems_v0latlon_v16_nest.dat' 0036 % 0037 % modify_FVCOM_nested_forcing(conf, '/tmp/fvcom_nested.nc', 1) 0038 % 0039 % Author(s): 0040 % Ricardo Torres (Plymouth Marine Laboratory) 0041 % Pierre Cazenave (Plymouth Marine Laboratory) 0042 % 0043 % Revision history: 0044 % 2017-01-17 First version based on Pierre's write_FVCOM_nested_forcing.m 0045 % script. Compress the time series data to save space. Requires 0046 % netCDF4 in FVCOM. 0047 % 0048 %========================================================================== 0049 0050 % The following variables are read from the file: 0051 % 0052 % lon, lat: Grid node positions [node] 0053 % lonc, latc: Grid element positions [nele] 0054 % h: Grid node depths [node] 0055 % hc: Grid element depth [nele] 0056 % nv: Triangulation table [nele, 3] 0057 % zeta: Sea surface elevation [node, time] 0058 % ua: Vertically averaged x velocity [node, time] 0059 % va: Vertically averaged y velocity [nele, time] 0060 % u: Eastward Water Velocity [nele, siglay, time] 0061 % v: Northward Water Velocity [nele, siglay, time] 0062 % temp: Temperature [node, siglay, time] 0063 % salinity: Salinity [node, siglay, time] 0064 % hyw: Hydrostatic vertical velocity [node, siglev, time] 0065 % weight_cell: Weighting for elements [nele] 0066 % weight_node: Weighting for nodes [node] 0067 % Itime: Days since 1858-11-17 00:00:00 [time] 0068 % Itime2: msec since 00:00:00 [time] 0069 % 0070 % We include these optional ones for humans: 0071 % time: Modified Julian Day [time] 0072 % Times: Gregorian dates [time, datestrlen] 0073 % ... and all the ERSEM variables if present... 0074 [~, subname] = fileparts(mfilename('fullpath')); 0075 0076 global ftbverbose 0077 if ftbverbose 0078 fprintf('\nbegin : %s\n', subname) 0079 end 0080 if nargin == 2 0081 nesttype = 1; 0082 elseif nargin < 2 || nargin > 3 0083 error(['Incorrect input arguments. Supply netCDF file path, ', ... 0084 'conf struct and optionally the nesting type (1, 2 or 3).']) 0085 end 0086 0087 % Can't use CLOBBER and NETCDF4 at the same time (the bitwise or didn't 0088 % work). Fall back to a horrible delete and then create instead. 0089 if exist(ncfile2read, 'file') 0090 nc2read = netcdf.open(ncfile2read,'NOWRITE'); 0091 [PATHSTR,NAME,EXT] = fileparts(ncfile2read); 0092 if ftbverbose 0093 fprintf('\nProcessing file : %s\n', ncfile2read) 0094 end 0095 0096 % new file to hold the modified nesting data 0097 ncfile = fullfile(PATHSTR,[NAME,'modified',EXT]); 0098 else 0099 sprintf('I am stoping. The nesting file %s doesn''t exist',ncfile2read); 0100 error(['Aborting the script ',subname]) 0101 end 0102 % List of data we need. 0103 % work with existing variables in netcdf 0104 [numdims,numvars,numglobalatts,unlimdimid] = netcdf.inq(nc2read); 0105 required = cell(1,numvars); 0106 for vv=0:numvars-1 0107 [required{vv+1},xtype,dimids,natts] = netcdf.inqVar(nc2read,vv); 0108 end 0109 % 0110 % required = {'time', 'x', 'y', 'lon', 'lat', 'xc', 'yc', 'lonc', 'latc', ... 0111 % 'nv','zeta', 'h', 'h_center', 'u', 'v', 'ua', 'va', 'temp', 'salinity', 'hyw', ... 0112 % 'weight_cell', 'weight_node', 'siglay', 'siglay_center', 'siglev', 'siglev_center'}; 0113 % for f = required 0114 % nest.(f{1})=[]; 0115 % end 0116 % % Read all variables from the existing file 0117 % if ~isfield (conf,'nest') 0118 % error('Missing conf.nest, aborting... ') 0119 % end 0120 % if ~isfield (conf.nest,'levels') 0121 % error('Missing conf.nest.levels, aborting... ') 0122 % end 0123 0124 for f = required 0125 try 0126 varid = netcdf.inqVarID(nc2read,f{1}); 0127 nest.(f{1}) = netcdf.getVar(nc2read,varid,'double'); 0128 disp(['Extracting variable ',f{1}, ' from file']) 0129 catch 0130 warning( 'Missing %s variable from nesting file', f{1}); 0131 % remove field from list 0132 required(startsWith(required,f{1}))=[]; 0133 end 0134 end 0135 [~, nsiglay, ~] = size(nest.u); 0136 nsiglev = nsiglay + 1; 0137 [~, ~] = size(nest.zeta); 0138 0139 % Change weights if requested 0140 switch lower(conf.method) 0141 case {'weights'} 0142 if isfield(conf,'new_weight_cell') && isfield(conf,'new_weight_node') 0143 varid = netcdf.inqVarID(nc2read,'weight_cell'); 0144 nest.weight_cell = netcdf.getVar(nc2read,varid,'double'); 0145 weight_cell = unique(conf.new_weight_cell); % these are sorted 0146 Nweight_cell = unique(nest.weight_cell); % these are sorted 0147 if length(weight_cell) ~= length(Nweight_cell) 0148 error('new_weight_cell doesn''t have the same number of UNIQUE levels as the nesting file') 0149 end 0150 for rr=1:length(Nweight_cell) 0151 nest.weight_cell(nest.weight_cell==Nweight_cell(rr)) = weight_cell(rr); 0152 end 0153 % now do the Node's weights 0154 varid = netcdf.inqVarID(nc2read,'weight_node'); 0155 nest.weight_node = netcdf.getVar(nc2read,varid,'double'); 0156 weight_node = unique(conf.new_weight_node); % these are sorted 0157 Nweight_node = unique(nest.weight_node); % these are sorted 0158 if length(weight_node) ~= length(Nweight_node) 0159 error('new_weight_node doesn''t have the same number of UNIQUE levels as the nesting file') 0160 end 0161 for rr=1:length(Nweight_node) 0162 nest.weight_node(nest.weight_node==Nweight_node(rr)) = weight_node(rr); 0163 end 0164 sprintf('The weights have been changed in the original nesting file %s ',ncfile2read); 0165 disp('Returning to the main script...') 0166 0167 end 0168 case {'levels'} 0169 0170 if isfield(conf.nest,'levels') 0171 % Split nodes into the different levels 0172 Nweight_node = unique(nest.weight_node); % these are sorted increasing 0173 Nweight_cell = unique(nest.weight_cell); % these are sorted increasing 0174 nlevs = unique(conf.nest.levels); 0175 if length(Nweight_node)==conf.nest.levels 0176 warning(['conf.nest.levels have been provided but it has the same number of levels',... 0177 'as the original file. Nothing has changed unless you have provided new weights']) 0178 else 0179 % Decide how many levels to keep 0180 warning(['Chopping levels from nesting file. Assuming unique weight values for',... 0181 'each nesting level. If there are repeated values, this function will not work']) 0182 levels2keepN = Nweight_node(end:-1:end-nlevs); % these are now decreasing in order 0183 levels2keepC = Nweight_cell(end:-1:end-nlevs+1); % these are now decreasing in order 0184 % remove nodes and cells from each variable 0185 nodeid = netcdf.inqDimID(nc2read,'node'); 0186 neleid = netcdf.inqDimID(nc2read,'nele'); 0187 dropidxnode = nest.weight_node(:,1)<min(levels2keepN); 0188 dropidxele = nest.weight_cell(:,1)<min(levels2keepC); 0189 for f = required 0190 varid = netcdf.inqVarID(nc2read,f{1}); 0191 [~,~,dimids,~] = netcdf.inqVar(nc2read,varid); 0192 % check the variable has a node or element dimension 0193 if any(dimids==nodeid) 0194 switch ndims (nest.(f{1})) 0195 case 1 0196 nest.(f{1})(dropidxnode)=[]; 0197 case 2 0198 nest.(f{1})(dropidxnode,:)=[]; 0199 0200 case 3 0201 nest.(f{1})(dropidxnode,:,:)=[]; 0202 end 0203 elseif any(dimids==neleid) 0204 switch ndims (nest.(f{1})) 0205 case 1 0206 nest.(f{1})(dropidxele)=[]; 0207 case 2 0208 nest.(f{1})(dropidxele,:)=[]; 0209 0210 case 3 0211 nest.(f{1})(dropidxele,:,:)=[]; 0212 end 0213 0214 end 0215 end 0216 end 0217 end 0218 % if weight_nodes doesn't exist in nesting files and we want to only 0219 % extract a set a lat and lon positions read a list of nodes and delete the 0220 % ones that are not identified in the netcdf 0221 case {'positions'} 0222 if isfield(conf.nest,'latlon_nodes_file') 0223 % read nodes ids 0224 fid = fopen (conf.nest.latlon_nodes_file,'rt'); 0225 C = textscan(fid,"%*d %f %f",'Headerlines',1) 0226 fclose(fid) 0227 nodes.x = C{1};nodes.y=C{2}; 0228 fid = fopen (conf.nest.latlon_elems_file,'rt'); 0229 C = textscan(fid,"%*d %f %f",'Headerlines',1) 0230 fclose(fid) 0231 elems.x = C{1};elems.y=C{2}; 0232 % Find node numbers of nested boundary nodes in coarse mesh 0233 missmatch=nan(size(nodes.x)); 0234 for nn=1:length(nodes.x) 0235 % Find the nodes that are identical in coarse domain 0236 [missmatch(nn),keepidx.node(nn)]=min(abs(complex( nodes.x(nn)-nest.x,nodes.y(nn)-nest.y ) )); 0237 end 0238 warning('Maximum difference in node locations between find_nesting and FVCOM output is %f meters',max(missmatch)) 0239 for nn=1:length(elems.x) 0240 % Find the nodes that are identical in coarse domain 0241 [missmatch(nn),keepidx.elems(nn)]=min(abs(complex( elems.x(nn)-nest.xc,elems.y(nn)-nest.yc ) )); 0242 end 0243 warning('Maximum difference in element locations between find_nesting and FVCOM output is %f meters',max(missmatch)) 0244 else 0245 error('We are selecting nodes according to position but you haven''t provided a file to read them from') 0246 end 0247 nodeid = netcdf.inqDimID(nc2read,'node'); 0248 neleid = netcdf.inqDimID(nc2read,'nele'); 0249 dropidxnode = setdiff([1:length(nest.x)],keepidx.node); 0250 dropidxele = setdiff([1:length(nest.xc)],keepidx.elems); 0251 for f = required 0252 varid = netcdf.inqVarID(nc2read,f{1}); 0253 [~,~,dimids,~] = netcdf.inqVar(nc2read,varid); 0254 % check the variable has a node or element dimension 0255 if any(dimids==nodeid) 0256 switch ndims (nest.(f{1})) 0257 case 1 0258 nest.(f{1})(dropidxnode)=[]; 0259 case 2 0260 nest.(f{1})(dropidxnode,:)=[]; 0261 0262 case 3 0263 nest.(f{1})(dropidxnode,:,:)=[]; 0264 end 0265 elseif any(dimids==neleid) 0266 switch ndims (nest.(f{1})) 0267 case 1 0268 nest.(f{1})(dropidxele)=[]; 0269 case 2 0270 nest.(f{1})(dropidxele,:)=[]; 0271 0272 case 3 0273 nest.(f{1})(dropidxele,:,:)=[]; 0274 end 0275 0276 end 0277 end 0278 0279 0280 end 0281 % Now identify indexes to drop 0282 %% Can't use CLOBBER and NETCDF4 at the same time (the bitwise or didn't 0283 % work). Fall back to a horrible delete and then create instead. 0284 if exist(ncfile, 'file') 0285 delete(ncfile) 0286 end 0287 [elems, nsiglay, ntimes] = size(nest.u); 0288 nsiglev = nsiglay + 1; 0289 [nodes, ~] = size(nest.zeta); 0290 0291 % output the new nesting file 0292 nc = netcdf.create(ncfile, 'NETCDF4'); 0293 0294 % define global attributes 0295 netcdf.putAtt(nc, netcdf.getConstant('NC_GLOBAL'), 'type', ... 0296 'FVCOM nesting TIME SERIES FILE') 0297 netcdf.putAtt(nc, netcdf.getConstant('NC_GLOBAL'), 'title', ... 0298 sprintf('FVCOM nesting TYPE %d TIME SERIES data for open boundary', ... 0299 nesttype)) 0300 netcdf.putAtt(nc, netcdf.getConstant('NC_GLOBAL'), 'history', ... 0301 sprintf('File created using %s from the MATLAB fvcom-toolbox', subname)) 0302 netcdf.putAtt(nc, netcdf.getConstant('NC_GLOBAL'), 'filename', ncfile) 0303 netcdf.putAtt(nc, netcdf.getConstant('NC_GLOBAL'), 'Conventions', 'CF-1.0') 0304 0305 % define dimensions 0306 % These need to happen in the same order as in the FVCOM output nesting 0307 % file!! 0308 elem_dimid = netcdf.defDim(nc, 'nele', elems); 0309 node_dimid = netcdf.defDim(nc, 'node', nodes); 0310 siglay_dimid = netcdf.defDim(nc, 'siglay', nsiglay); 0311 siglev_dimid = netcdf.defDim(nc, 'siglev', nsiglev); 0312 three_dimid = netcdf.defDim(nc, 'three', 3); 0313 time_dimid = netcdf.defDim(nc, 'time', netcdf.getConstant('NC_UNLIMITED')); 0314 datestrlen_dimid = netcdf.defDim(nc, 'DateStrLen', 26); 0315 0316 % define Times variable manually as it is not output directly by FVCOM 0317 0318 Times_varid = netcdf.defVar(nc, 'Times' ,'NC_CHAR', ... 0319 [datestrlen_dimid, time_dimid]); 0320 netcdf.putAtt(nc, Times_varid, 'time_zone', 'UTC'); 0321 0322 0323 0324 % define variables 0325 name=cell(1,numvars); 0326 xtype = nan(numvars,1); 0327 dimids = cell(1,numvars); 0328 natts = nan(numvars,1); 0329 required2=required; 0330 for vv=1:numvars 0331 [name{vv},xtype(vv),dimids{vv},natts(vv)]=netcdf.inqVar(nc2read,netcdf.inqVarID(nc2read,required{vv})); 0332 if isempty(dimids{vv}) 0333 % drop from list... 0334 name(vv)=[];xtype(vv)=[];dimids(vv)=[];natts(vv)=[]; 0335 required2(vv)=[]; 0336 else 0337 eval([name{vv},'_varid = netcdf.defVar(nc, ''',name{vv},''',xtype(vv),dimids{vv});']) 0338 for aa=1:natts(vv) 0339 varid = netcdf.inqVarID(nc2read,name{vv}); 0340 eval(['att_name= netcdf.inqAttName(nc2read,varid,aa-1);']) 0341 att_value= netcdf.getAtt(nc2read,varid,att_name) 0342 if startsWith(att_name ,{'_FillValue'}) 0343 eval(['netcdf.defVarFill(nc, ',name{vv},'_varid,false,att_value);']) 0344 else 0345 eval(['netcdf.putAtt(nc, ',name{vv},'_varid, att_name, att_value);']) 0346 end 0347 end 0348 % enable compression on the big variables. These are the ones with 3 0349 % dimensions 0350 if length(dimids(vv))>=2 0351 eval(['netcdf.putAtt(nc, ',name{vv},'_varid, att_name, att_value);']) 0352 0353 eval(['netcdf.defVarDeflate(nc, ',name{vv},'_varid, true, true, 7);']) 0354 % netcdf.defVarDeflate(nc, u_varid, true, true, 7); 0355 % netcdf.defVarDeflate(nc, v_varid, true, true, 7); 0356 % netcdf.defVarDeflate(nc, temp_varid, true, true, 7); 0357 % netcdf.defVarDeflate(nc, salinity_varid, true, true, 7); 0358 % netcdf.defVarDeflate(nc, hyw_varid, true, true, 7); 0359 end 0360 end 0361 end 0362 % end definitions 0363 netcdf.endDef(nc); 0364 disp('Finished with definition section') 0365 % write time data 0366 nStringOut = char(); 0367 [nYr, nMon, nDay, nHour, nMin, nSec] = mjulian2greg(nest.Itime+(nest.Itime2/1000)./(24*60*60)); 0368 for i = 1:ntimes 0369 nDate = [nYr(i), nMon(i), nDay(i), nHour(i), nMin(i), nSec(i)]; 0370 nStringOut = [nStringOut, sprintf('%-026s', datestr(datenum(nDate), 'yyyy-mm-dd HH:MM:SS.FFF'))]; 0371 end 0372 % do time manually and remove from list of variables to process 0373 netcdf.putVar(nc, Itime_varid, 0, numel(nest.time), floor(nest.time)); 0374 required2(startsWith(required2,'itime'))=[]; 0375 netcdf.putVar(nc, Itime2_varid, 0, numel(nest.time), ... 0376 mod(nest.time, 1)*24*3600*1000); 0377 required2(startsWith(required2,'itime2'))=[]; 0378 0379 netcdf.putVar(nc, Times_varid, nStringOut); 0380 required2(startsWith(required2,'Times'))=[]; 0381 netcdf.putVar(nc, time_varid, nest.time); 0382 required2(startsWith(required2,'time'))=[]; 0383 % now loop through all variables again and put them in 0384 for vv=1:length(required2) 0385 [name{vv},xtype(vv),dimids{vv},natts(vv)]=netcdf.inqVar(nc2read,netcdf.inqVarID(nc2read,required2{vv})); 0386 eval(['netcdf.putVar(nc, ',required2{vv},'_varid, nest.(required2{vv}));']) 0387 disp(['Adding variable ',required2{vv},' to file']) 0388 0389 end 0390 % netcdf.putVar(nc, x_varid, nest.x); 0391 % netcdf.putVar(nc, y_varid, nest.y); 0392 % netcdf.putVar(nc, xc_varid, nest.xc); 0393 % netcdf.putVar(nc, yc_varid, nest.yc); 0394 % netcdf.putVar(nc, lon_varid, nest.lon); 0395 % netcdf.putVar(nc, lat_varid, nest.lat); 0396 % netcdf.putVar(nc, lonc_varid, nest.lonc); 0397 % netcdf.putVar(nc, latc_varid, nest.latc); 0398 % 0399 % % dump data 0400 % netcdf.putVar(nc, zeta_varid, nest.zeta); 0401 % netcdf.putVar(nc, ua_varid, nest.ua); 0402 % netcdf.putVar(nc, va_varid, nest.va); 0403 % netcdf.putVar(nc, u_varid, nest.u); 0404 % netcdf.putVar(nc, v_varid, nest.v); 0405 % netcdf.putVar(nc, temp_varid, nest.temp); 0406 % netcdf.putVar(nc, salinity_varid, nest.salinity); 0407 % netcdf.putVar(nc, hyw_varid, nest.hyw); 0408 % netcdf.putVar(nc, siglay_varid, nest.siglay); 0409 % netcdf.putVar(nc, siglayc_varid, nest.siglay_center); 0410 % netcdf.putVar(nc, siglev_varid, nest.siglev); 0411 % netcdf.putVar(nc, siglevc_varid, nest.siglev_center); 0412 % netcdf.putVar(nc, h_varid, nest.h); 0413 % netcdf.putVar(nc, hc_varid, nest.h_center); 0414 % if nesttype > 2 0415 % netcdf.putVar(nc, cweights_varid, nest.weight_cell); 0416 % netcdf.putVar(nc, nweights_varid, nest.weight_node); 0417 % end 0418 0419 % close file 0420 netcdf.close(nc) 0421 0422 if ftbverbose 0423 fprintf('end : %s\n', subname) 0424 end 0425 % 0426 % time_varid = netcdf.defVar(nc, 'time', 'NC_FLOAT', time_dimid); 0427 % netcdf.putAtt(nc, time_varid, 'long_name', 'time'); 0428 % netcdf.putAtt(nc, time_varid, 'units', 'days since 1858-11-17 00:00:00'); 0429 % netcdf.putAtt(nc, time_varid, 'format', 'modified julian day (MJD)'); 0430 % netcdf.putAtt(nc, time_varid, 'time_zone', 'UTC'); 0431 % 0432 % itime_varid = netcdf.defVar(nc, 'Itime', 'NC_INT', ... 0433 % time_dimid); 0434 % netcdf.putAtt(nc, itime_varid, 'units', 'days since 1858-11-17 00:00:00'); 0435 % netcdf.putAtt(nc, itime_varid, 'format', 'modified julian day (MJD)'); 0436 % netcdf.putAtt(nc, itime_varid, 'time_zone', 'UTC'); 0437 % 0438 % itime2_varid = netcdf.defVar(nc, 'Itime2', 'NC_INT', ... 0439 % time_dimid); 0440 % netcdf.putAtt(nc, itime2_varid, 'units', 'msec since 00:00:00'); 0441 % netcdf.putAtt(nc, itime2_varid, 'time_zone', 'UTC'); 0442 % 0443 % Times_varid = netcdf.defVar(nc, 'Times' ,'NC_CHAR', ... 0444 % [datestrlen_dimid, time_dimid]); 0445 % netcdf.putAtt(nc, Times_varid, 'time_zone', 'UTC'); 0446 % 0447 % x_varid = netcdf.defVar(nc, 'x', 'NC_FLOAT', ... 0448 % node_dimid); 0449 % netcdf.putAtt(nc, x_varid, 'units', 'meters'); 0450 % netcdf.putAtt(nc, x_varid, 'long_name', 'nodal x-coordinate'); 0451 % 0452 % y_varid = netcdf.defVar(nc, 'y', 'NC_FLOAT', ... 0453 % node_dimid); 0454 % netcdf.putAtt(nc, y_varid, 'units', 'meters'); 0455 % netcdf.putAtt(nc, y_varid, 'long_name', 'nodal y-coordinate'); 0456 % 0457 % xc_varid = netcdf.defVar(nc, 'xc', 'NC_FLOAT', ... 0458 % elem_dimid); 0459 % netcdf.putAtt(nc, xc_varid, 'units', 'meters'); 0460 % netcdf.putAtt(nc, xc_varid, 'long_name', 'zonal x-coordinate'); 0461 % 0462 % yc_varid = netcdf.defVar(nc, 'yc', 'NC_FLOAT', ... 0463 % elem_dimid); 0464 % netcdf.putAtt(nc, yc_varid, 'units', 'meters'); 0465 % netcdf.putAtt(nc, yc_varid, 'long_name', 'zonal y-coordinate'); 0466 % 0467 % lon_varid = netcdf.defVar(nc, 'lon', 'NC_FLOAT', ... 0468 % node_dimid); 0469 % netcdf.putAtt(nc, lon_varid, 'units', 'degrees_east'); 0470 % netcdf.putAtt(nc, lon_varid, 'standard_name', 'longitude'); 0471 % netcdf.putAtt(nc, lon_varid, 'long_name', 'nodal longitude'); 0472 % 0473 % lat_varid = netcdf.defVar(nc, 'lat', 'NC_FLOAT', ... 0474 % node_dimid); 0475 % netcdf.putAtt(nc, lat_varid, 'units', 'degrees_north'); 0476 % netcdf.putAtt(nc, lat_varid, 'standard_name', 'latitude'); 0477 % netcdf.putAtt(nc, lat_varid, 'long_name', 'nodal latitude'); 0478 % 0479 % lonc_varid = netcdf.defVar(nc, 'lonc', 'NC_FLOAT', ... 0480 % elem_dimid); 0481 % netcdf.putAtt(nc, lonc_varid, 'units', 'degrees_east'); 0482 % netcdf.putAtt(nc, lonc_varid, 'standard_name', 'longitude'); 0483 % netcdf.putAtt(nc, lonc_varid, 'long_name', 'zonal longitude'); 0484 % 0485 % latc_varid = netcdf.defVar(nc, 'latc', 'NC_FLOAT', ... 0486 % elem_dimid); 0487 % netcdf.putAtt(nc, latc_varid, 'units', 'degrees_north'); 0488 % netcdf.putAtt(nc, latc_varid, 'standard_name', 'latitude'); 0489 % netcdf.putAtt(nc, latc_varid, 'long_name', 'zonal latitude'); 0490 % 0491 % nv_varid = netcdf.defVar(nc, 'nv', 'NC_INT', ... 0492 % [elem_dimid, three_dimid]); 0493 % netcdf.putAtt(nc, xc_varid, 'units', 'no units'); 0494 % netcdf.putAtt(nc, xc_varid, 'long_name', 'elements nodes indices'); 0495 % 0496 % zeta_varid = netcdf.defVar(nc, 'zeta', 'NC_FLOAT', ... 0497 % [node_dimid, time_dimid]); 0498 % netcdf.putAtt(nc, zeta_varid, 'long_name', 'Water Surface Elevation'); 0499 % netcdf.putAtt(nc, zeta_varid, 'units', 'meters'); 0500 % netcdf.putAtt(nc, zeta_varid, 'positive', 'up'); 0501 % netcdf.putAtt(nc, zeta_varid, 'standard_name', ... 0502 % 'sea_surface_height_above_geoid'); 0503 % netcdf.putAtt(nc, zeta_varid, 'grid', 'Bathymetry_Mesh'); 0504 % netcdf.putAtt(nc, zeta_varid, 'coordinates', 'time lat lon'); 0505 % netcdf.putAtt(nc, zeta_varid, 'type', 'data'); 0506 % netcdf.putAtt(nc, zeta_varid, 'location', 'node'); 0507 % 0508 % ua_varid = netcdf.defVar(nc, 'ua', 'NC_FLOAT', ... 0509 % [elem_dimid, time_dimid]); 0510 % netcdf.putAtt(nc, ua_varid, 'long_name', 'Vertically Averaged x-velocity'); 0511 % netcdf.putAtt(nc, ua_varid, 'units', 'meters s-1'); 0512 % netcdf.putAtt(nc, ua_varid, 'grid', 'fvcom_grid'); 0513 % netcdf.putAtt(nc, ua_varid, 'type', 'data'); 0514 % 0515 % va_varid = netcdf.defVar(nc, 'va', 'NC_FLOAT', ... 0516 % [elem_dimid, time_dimid]); 0517 % netcdf.putAtt(nc, va_varid, 'long_name', 'Vertically Averaged y-velocity'); 0518 % netcdf.putAtt(nc, va_varid, 'units', 'meters s-1'); 0519 % netcdf.putAtt(nc, va_varid, 'grid', 'fvcom_grid'); 0520 % netcdf.putAtt(nc, va_varid, 'type', 'data'); 0521 % 0522 % u_varid = netcdf.defVar(nc, 'u', 'NC_FLOAT', ... 0523 % [elem_dimid, siglay_dimid, time_dimid]); 0524 % netcdf.putAtt(nc, u_varid, 'long_name', 'Eastward Water Velocity'); 0525 % netcdf.putAtt(nc, u_varid, 'units', 'meters s-1'); 0526 % netcdf.putAtt(nc, u_varid, 'standard_name', 'eastward_sea_water_velocity'); 0527 % netcdf.putAtt(nc, u_varid, 'grid', 'fvcom_grid'); 0528 % netcdf.putAtt(nc, u_varid, 'coordinates', 'time siglay latc lonc'); 0529 % netcdf.putAtt(nc, u_varid, 'type', 'data'); 0530 % netcdf.putAtt(nc, u_varid, 'location', 'face'); 0531 % 0532 % v_varid = netcdf.defVar(nc, 'v', 'NC_FLOAT', ... 0533 % [elem_dimid, siglay_dimid, time_dimid]); 0534 % netcdf.putAtt(nc, v_varid, 'long_name', 'Northward Water Velocity'); 0535 % netcdf.putAtt(nc, v_varid, 'units', 'meters s-1'); 0536 % netcdf.putAtt(nc, v_varid, 'standard_name', ... 0537 % 'Northward_sea_water_velocity'); 0538 % netcdf.putAtt(nc, v_varid, 'grid', 'fvcom_grid'); 0539 % netcdf.putAtt(nc, v_varid, 'coordinates', 'time siglay latc lonc'); 0540 % netcdf.putAtt(nc, v_varid, 'type', 'data'); 0541 % netcdf.putAtt(nc, v_varid, 'location', 'face'); 0542 % 0543 % temp_varid = netcdf.defVar(nc, 'temp', 'NC_FLOAT', ... 0544 % [node_dimid, siglay_dimid, time_dimid]); 0545 % netcdf.putAtt(nc, temp_varid, 'long_name', 'Temperature'); 0546 % netcdf.putAtt(nc, temp_varid, 'standard_name', 'sea_water_temperature'); 0547 % netcdf.putAtt(nc, temp_varid, 'units', 'degrees Celcius'); 0548 % netcdf.putAtt(nc, temp_varid, 'grid', 'fvcom_grid'); 0549 % netcdf.putAtt(nc, temp_varid, 'coordinates', 'time siglay lat lon'); 0550 % netcdf.putAtt(nc, temp_varid, 'type', 'data'); 0551 % netcdf.putAtt(nc, temp_varid, 'location', 'node'); 0552 % 0553 % salinity_varid = netcdf.defVar(nc, 'salinity', 'NC_FLOAT', ... 0554 % [node_dimid, siglay_dimid, time_dimid]); 0555 % netcdf.putAtt(nc, salinity_varid, 'long_name', 'Salinity'); 0556 % netcdf.putAtt(nc, salinity_varid, 'standard_name', 'sea_water_salinity'); 0557 % netcdf.putAtt(nc, salinity_varid, 'units', '1e-3'); 0558 % netcdf.putAtt(nc, salinity_varid, 'grid', 'fvcom_grid'); 0559 % netcdf.putAtt(nc, salinity_varid, 'coordinates', 'time siglay lat lon'); 0560 % netcdf.putAtt(nc, salinity_varid, 'type', 'data'); 0561 % netcdf.putAtt(nc, salinity_varid, 'location', 'node'); 0562 % 0563 % hyw_varid = netcdf.defVar(nc, 'hyw', 'NC_FLOAT', ... 0564 % [node_dimid, siglev_dimid, time_dimid]); 0565 % netcdf.putAtt(nc, hyw_varid, 'long_name', ... 0566 % 'hydro static vertical velocity'); 0567 % netcdf.putAtt(nc, hyw_varid, 'units', 'meters s-1'); 0568 % netcdf.putAtt(nc, hyw_varid, 'grid', 'fvcom_grid'); 0569 % netcdf.putAtt(nc, hyw_varid, 'type', 'data'); 0570 % netcdf.putAtt(nc, hyw_varid, 'coordinates', 'time siglay lat lon'); 0571 % 0572 % siglay_varid = netcdf.defVar(nc, 'siglay', 'NC_FLOAT', ... 0573 % [node_dimid, siglay_dimid]); 0574 % netcdf.putAtt(nc, siglay_varid, 'long_name', 'Sigma Layers'); 0575 % netcdf.putAtt(nc, siglay_varid, 'standard_name', 'ocean_sigma/general_coordinate'); 0576 % netcdf.putAtt(nc, siglay_varid, 'positive', 'up'); 0577 % netcdf.putAtt(nc, siglay_varid, 'valid_min', '-1.f'); 0578 % netcdf.putAtt(nc, siglay_varid, 'valid_max', '0.f'); 0579 % netcdf.putAtt(nc, siglay_varid, 'formula_terms', 'sigma: siglay eta: zeta depth: h'); 0580 % 0581 % siglayc_varid = netcdf.defVar(nc, 'siglay_center', 'NC_FLOAT', ... 0582 % [elem_dimid, siglay_dimid]); 0583 % netcdf.putAtt(nc, siglayc_varid, 'long_name', 'Sigma Layers'); 0584 % netcdf.putAtt(nc, siglayc_varid, 'standard_name', 'ocean_sigma/general_coordinate'); 0585 % netcdf.putAtt(nc, siglayc_varid, 'positive', 'up'); 0586 % netcdf.putAtt(nc, siglayc_varid, 'valid_min', '-1.f'); 0587 % netcdf.putAtt(nc, siglayc_varid, 'valid_max', '0.f'); 0588 % netcdf.putAtt(nc, siglayc_varid, 'formula_terms', 'sigma: siglay_center eta: zeta_center depth: h_center'); 0589 % 0590 % siglev_varid = netcdf.defVar(nc, 'siglev', 'NC_FLOAT', ... 0591 % [node_dimid, siglev_dimid]); 0592 % netcdf.putAtt(nc, siglev_varid, 'long_name', 'Sigma Levels'); 0593 % netcdf.putAtt(nc, siglev_varid, 'standard_name', 'ocean_sigma/general_coordinate'); 0594 % netcdf.putAtt(nc, siglev_varid, 'positive', 'up'); 0595 % netcdf.putAtt(nc, siglev_varid, 'valid_min', '-1.f'); 0596 % netcdf.putAtt(nc, siglev_varid, 'valid_max', '0.f'); 0597 % netcdf.putAtt(nc, siglev_varid, 'formula_terms', 'sigma:siglev eta: zeta depth: h'); 0598 % 0599 % siglevc_varid = netcdf.defVar(nc, 'siglev_center', 'NC_FLOAT', ... 0600 % [elem_dimid, siglev_dimid]); 0601 % netcdf.putAtt(nc, siglevc_varid, 'long_name', 'Sigma Layers'); 0602 % netcdf.putAtt(nc, siglevc_varid, 'standard_name', 'ocean_sigma/general_coordinate'); 0603 % netcdf.putAtt(nc, siglevc_varid, 'positive', 'up'); 0604 % netcdf.putAtt(nc, siglevc_varid, 'valid_min', '-1.f'); 0605 % netcdf.putAtt(nc, siglevc_varid, 'valid_max', '0.f'); 0606 % netcdf.putAtt(nc, siglevc_varid, 'formula_terms', 'sigma: siglev_center eta: zeta_center depth: h_center'); 0607 % 0608 % h_varid = netcdf.defVar(nc, 'h', 'NC_FLOAT', ... 0609 % node_dimid); 0610 % netcdf.putAtt(nc, h_varid, 'long_name', 'Bathymetry'); 0611 % netcdf.putAtt(nc, h_varid, 'standard_name', 'sea_floor_depth_below_geoid'); 0612 % netcdf.putAtt(nc, h_varid, 'units', 'm'); 0613 % netcdf.putAtt(nc, h_varid, 'positive', 'down'); 0614 % netcdf.putAtt(nc, h_varid, 'grid', 'Bathymetry_mesh'); 0615 % netcdf.putAtt(nc, h_varid, 'coordinates', 'x y'); 0616 % netcdf.putAtt(nc, h_varid, 'type', 'data'); 0617 % 0618 % hc_varid = netcdf.defVar(nc, 'h_center', 'NC_FLOAT', ... 0619 % elem_dimid); 0620 % netcdf.putAtt(nc, hc_varid, 'long_name', 'Bathymetry'); 0621 % netcdf.putAtt(nc, hc_varid, 'standard_name', 'sea_floor_depth_below_geoid'); 0622 % netcdf.putAtt(nc, hc_varid, 'units', 'm'); 0623 % netcdf.putAtt(nc, hc_varid, 'positive', 'down'); 0624 % netcdf.putAtt(nc, hc_varid, 'grid', 'grid1 grid3'); 0625 % netcdf.putAtt(nc, hc_varid, 'coordinates', 'latc lonc'); 0626 % netcdf.putAtt(nc, hc_varid, 'grid_location', 'center'); 0627 % 0628 % if nesttype > 2 0629 % cweights_varid = netcdf.defVar(nc, 'weight_cell', 'NC_FLOAT', ... 0630 % [elem_dimid, time_dimid]); 0631 % netcdf.putAtt(nc, cweights_varid, 'long_name', ... 0632 % 'Weights for elements in relaxation zone'); 0633 % netcdf.putAtt(nc, cweights_varid, 'units', 'no units'); 0634 % netcdf.putAtt(nc, cweights_varid, 'grid', 'fvcom_grid'); 0635 % netcdf.putAtt(nc, cweights_varid, 'type', 'data'); 0636 % 0637 % nweights_varid = netcdf.defVar(nc, 'weight_node', 'NC_FLOAT', ... 0638 % [node_dimid, time_dimid]); 0639 % netcdf.putAtt(nc, nweights_varid, 'long_name', ... 0640 % 'Weights for nodes in relaxation zone'); 0641 % netcdf.putAtt(nc, nweights_varid, 'units', 'no units'); 0642 % netcdf.putAtt(nc, nweights_varid, 'grid', 'fvcom_grid'); 0643 % netcdf.putAtt(nc, nweights_varid, 'type', 'data'); 0644 % end 0645