


Duplicate an FVCOM restart file, replacing variable values with those
specified in the struct data.
function write_FVCOM_restart(fv_restart, out_restart, indata)
DESCRIPTION:
Use an existing FVCOM restart file as a template, export all existing
data except for variables whose names match the data in the struct
'data'.
INPUT:
fv_restart = full path to an existing FVCOM restart file.
out_restart = full path to the restart file to be created.
indata = struct whose field names are the variable names to be
replaced.
OPTIONAL INPUT (keyword-value pairs):
'out_date' = [optional] reset the restart file times to this date
([YYYY, MM, DD, HH, MM, SS]). This can be a single date only or a time
series whose length matches the provided restart data in indata. If new
data are being provided, the data time series length must match the
out_date series. If a only a single time is given, the output file will
include three time steps to bracket the specified time by 30 minutes
each way (to allow FVCOM some wiggle room when loading the data). The
data will be replicated over the three time steps.
OUTPUT:
FVCOM restart file.
EXAMPLE USAGE:
Replace temperature and salinity but leave the times as is:
indata.temp = interpolated_temp;
indata.salinity = interpolated_salinity;
write_FVCOM_restart('/tmp/fvcom_restart.nc', ...
'/tmp/fvcom_restart_interp.nc', indata)
Replace temperature only and change the restart times to bracketed
values:
indata.temp = interpolated_temp;
write_FVCOM_restart('/tmp/fvcom_restart.nc', ...
'/tmp/fvcom_restart_interp.nc', indata, 'out_date', ...
[2003, 05, 25, 13, 34, 07])
Replace all temperatures with a single value leaving times as they are:
indata.temp = 13;
write_FVCOM_restart('/tmp/fvcom_restart.nc', ...
'/tmp/fvcom_restart_interp.nc', indata)
Replace temperature with a constant at all positions and times and all
times with a new array of times (for a restart file with five time
steps):
newtime = datevec(datenum([2004, 06, 20, 11, 09, 44]) + ...
[-0.5, -0.25, 0, 0.25, 0.5]);
indata.temp = 13;
write_FVCOM_restart('/tmp/fvcom_restart.nc', ...
'/tmp/fvcom_restart_interp.nc', indata, 'out_time', newtime)
Author(s):
Pierre Cazenave (Plymouth Marine Laboratory)
Karen Amoudry (National Oceanography Centre, Liverpool)
Revision history:
2013-02-08 First version.
2013-02-15 Fix bug wherein only the last field in the new data would
only be added to the output netCDF file.
2013-03-13 Make the time rewriting optional and not just commented out.
2014-01-23 Add functionality to specify length of time series in output
file (Karen Amoudry).
2014-01-31 Add functionality to replace user-specified variables in the
output file with user-specified constant values (Karen Amoudry).
2014-02-04 Incorporate Karen's functionality (see above) but with the
ability to retain the existing behaviour (where a new start time is
still optional). User-specified constants are also supported but
instead of specifying a new input argument, if a single scalar value is
given in the input struct but the output is non-scalar (i.e. an array),
then that scalar is tiled to the size of the expected output array.
2014-02-13 Fix output when only a single time step is in the input
template restart file (fv_restart). Single time step files have the
indata put in the time step (previously the ramping from model to
indata values meant that the output contained only the model data,
which makes no sense.
2014-07-08 Fix the creation of the Itime2 variable (milliseconds since
midnight). Also update the help to better describe how the optional
'out_date' argument works.
2017-08-01 Add netCDF4 compression to save on file size.
==========================================================================

0001 function write_FVCOM_restart(fv_restart, out_restart, indata, varargin) 0002 % Duplicate an FVCOM restart file, replacing variable values with those 0003 % specified in the struct data. 0004 % 0005 % function write_FVCOM_restart(fv_restart, out_restart, indata) 0006 % 0007 % DESCRIPTION: 0008 % Use an existing FVCOM restart file as a template, export all existing 0009 % data except for variables whose names match the data in the struct 0010 % 'data'. 0011 % 0012 % INPUT: 0013 % fv_restart = full path to an existing FVCOM restart file. 0014 % out_restart = full path to the restart file to be created. 0015 % indata = struct whose field names are the variable names to be 0016 % replaced. 0017 % OPTIONAL INPUT (keyword-value pairs): 0018 % 'out_date' = [optional] reset the restart file times to this date 0019 % ([YYYY, MM, DD, HH, MM, SS]). This can be a single date only or a time 0020 % series whose length matches the provided restart data in indata. If new 0021 % data are being provided, the data time series length must match the 0022 % out_date series. If a only a single time is given, the output file will 0023 % include three time steps to bracket the specified time by 30 minutes 0024 % each way (to allow FVCOM some wiggle room when loading the data). The 0025 % data will be replicated over the three time steps. 0026 % 0027 % OUTPUT: 0028 % FVCOM restart file. 0029 % 0030 % EXAMPLE USAGE: 0031 % 0032 % Replace temperature and salinity but leave the times as is: 0033 % indata.temp = interpolated_temp; 0034 % indata.salinity = interpolated_salinity; 0035 % write_FVCOM_restart('/tmp/fvcom_restart.nc', ... 0036 % '/tmp/fvcom_restart_interp.nc', indata) 0037 % 0038 % Replace temperature only and change the restart times to bracketed 0039 % values: 0040 % indata.temp = interpolated_temp; 0041 % write_FVCOM_restart('/tmp/fvcom_restart.nc', ... 0042 % '/tmp/fvcom_restart_interp.nc', indata, 'out_date', ... 0043 % [2003, 05, 25, 13, 34, 07]) 0044 % 0045 % Replace all temperatures with a single value leaving times as they are: 0046 % indata.temp = 13; 0047 % write_FVCOM_restart('/tmp/fvcom_restart.nc', ... 0048 % '/tmp/fvcom_restart_interp.nc', indata) 0049 % 0050 % Replace temperature with a constant at all positions and times and all 0051 % times with a new array of times (for a restart file with five time 0052 % steps): 0053 % newtime = datevec(datenum([2004, 06, 20, 11, 09, 44]) + ... 0054 % [-0.5, -0.25, 0, 0.25, 0.5]); 0055 % indata.temp = 13; 0056 % write_FVCOM_restart('/tmp/fvcom_restart.nc', ... 0057 % '/tmp/fvcom_restart_interp.nc', indata, 'out_time', newtime) 0058 % 0059 % Author(s): 0060 % Pierre Cazenave (Plymouth Marine Laboratory) 0061 % Karen Amoudry (National Oceanography Centre, Liverpool) 0062 % 0063 % Revision history: 0064 % 2013-02-08 First version. 0065 % 2013-02-15 Fix bug wherein only the last field in the new data would 0066 % only be added to the output netCDF file. 0067 % 2013-03-13 Make the time rewriting optional and not just commented out. 0068 % 2014-01-23 Add functionality to specify length of time series in output 0069 % file (Karen Amoudry). 0070 % 2014-01-31 Add functionality to replace user-specified variables in the 0071 % output file with user-specified constant values (Karen Amoudry). 0072 % 2014-02-04 Incorporate Karen's functionality (see above) but with the 0073 % ability to retain the existing behaviour (where a new start time is 0074 % still optional). User-specified constants are also supported but 0075 % instead of specifying a new input argument, if a single scalar value is 0076 % given in the input struct but the output is non-scalar (i.e. an array), 0077 % then that scalar is tiled to the size of the expected output array. 0078 % 2014-02-13 Fix output when only a single time step is in the input 0079 % template restart file (fv_restart). Single time step files have the 0080 % indata put in the time step (previously the ramping from model to 0081 % indata values meant that the output contained only the model data, 0082 % which makes no sense. 0083 % 2014-07-08 Fix the creation of the Itime2 variable (milliseconds since 0084 % midnight). Also update the help to better describe how the optional 0085 % 'out_date' argument works. 0086 % 2017-08-01 Add netCDF4 compression to save on file size. 0087 % 0088 %========================================================================== 0089 0090 [~, subname] = fileparts(mfilename('fullpath')); 0091 0092 global ftbverbose 0093 if ftbverbose 0094 fprintf('\nbegin : %s \n', subname) 0095 end 0096 0097 % Default to not bracketing data. If we've got a single time input, then 0098 % this will be set to true and the data will be bracketed around that time. 0099 % If we've got a new input time array, the new data will be ramped to the 0100 % new input data values and bracketed is left false. If we have no new 0101 % times, then the new data will also be ramped up. 0102 bracketed = false; 0103 0104 if nargin > 3 0105 assert(rem(length(varargin), 2) == 0, 'Incorrect keyword-value arguments.') 0106 for aa = 1:2:length(varargin) 0107 key = varargin{aa}; 0108 val = varargin{aa + 1}; 0109 switch key 0110 case 'out_date' 0111 out_date = val; 0112 if isscalar(val(:, 1)) 0113 % Bracket the date by 30 minutes either way. 0114 bracketed = true; 0115 tOffset = 30; 0116 new_date = datevec([... 0117 datenum(... 0118 val(1), val(2), val(3), val(4), val(5) - tOffset, val(6)); ... 0119 datenum(... 0120 val(1), val(2), val(3), val(4), val(5), val(6));... 0121 datenum(... 0122 val(1), val(2), val(3), val(4), val(5) + tOffset, val(6))... 0123 ]); 0124 else 0125 % Assume we've been given a time series. We need to 0126 % make sure (later) that this is the same length as the 0127 % input file times. 0128 new_date = val; 0129 end 0130 otherwise 0131 warning('Unrecognised input argument keyword ''%s'' and value ''%s''.', key, varargin{v + 1}) 0132 end 0133 end 0134 end 0135 0136 % Get the fieldnames which must match the variable names to be replaced 0137 % (case sensitive). 0138 fnames = fieldnames(indata); 0139 nf = length(fnames); 0140 0141 if exist(out_restart, 'file') 0142 delete(out_restart) 0143 end 0144 nc = netcdf.open(fv_restart, 'NOWRITE'); 0145 ncout = netcdf.create(out_restart, 'NETCDF4'); 0146 0147 [numdims, numvars, numglobatts, unlimdimID] = netcdf.inq(nc); 0148 0149 % Define the dimensions for all variables. 0150 dimid = nan(numdims, 1); 0151 dimnames = cell(numdims, 1); 0152 dimlengths = nan(numdims, 1); 0153 for ii = 1:numdims 0154 [dimname, dimlen] = netcdf.inqDim(nc, ii - 1); 0155 % If we've been asked to rewrite the times, then set the length 0156 % of the time dimension to three (bracketed by two time steps each 0157 % way). 0158 if exist('out_date', 'var') && ii == unlimdimID + 1 0159 if length(out_date(:, 1)) == dimlen 0160 % We're replacing times with those in new_date 0161 out_date = new_date; 0162 elseif isvector(out_date) 0163 % We have bracketed dates 0164 out_date = new_date; 0165 end 0166 dimlen = length(out_date(:, 1)); 0167 end 0168 0169 if ii ~= unlimdimID + 1 % netCDF indices start at zero 0170 dimid(ii) = netcdf.defDim(ncout, dimname, dimlen); 0171 else 0172 dimid(ii) = netcdf.defDim(ncout, dimname, netcdf.getConstant('NC_UNLIMITED')); 0173 end 0174 dimnames{ii} = dimname; 0175 dimlengths(ii) = dimlen; 0176 end 0177 0178 % Now define the variables and attributes. 0179 for ii = 1:numvars 0180 0181 % Find name of the current variable. 0182 [varname, xtype, varDimIDs, varAtts] = netcdf.inqVar(nc, ii - 1); 0183 0184 % Create the variables. 0185 varid = netcdf.defVar(ncout, varname, xtype, varDimIDs); 0186 % Enable compression. 0187 netcdf.defVarDeflate(ncout, varid, true, true, 7); 0188 0189 % Get each attribute and add it to the current variable. 0190 for j = 1:varAtts 0191 0192 attname = netcdf.inqAttName(nc, varid, j - 1); 0193 attval = netcdf.getAtt(nc, varid, attname); 0194 if startsWith(attname ,{'_FillValue'}) 0195 netcdf.defVarFill(ncout, varid,false,attval); 0196 else 0197 netcdf.putAtt(ncout, varid, attname, attval); 0198 end 0199 end 0200 end 0201 0202 % Do the global attributes 0203 for ii = 1:numglobatts 0204 0205 % Find the current global attribute's name and value. 0206 gattname = netcdf.inqAttName(nc, netcdf.getConstant('NC_GLOBAL'), ii - 1); 0207 gattval = netcdf.getAtt(nc, netcdf.getConstant('NC_GLOBAL'), gattname); 0208 0209 % Put that back into the output netCDF file. 0210 netcdf.putAtt(ncout, netcdf.getConstant('NC_GLOBAL'), gattname, gattval); 0211 end 0212 0213 netcdf.endDef(ncout); 0214 0215 % Get the existing data and output to the new netCDF file, except for 0216 % variables which match the fieldnames in the data struct. 0217 for ii = 1:numvars 0218 0219 [varname, ~, varDimIDs, ~] = netcdf.inqVar(nc, ii - 1); 0220 varid = netcdf.inqVarID(nc, varname); 0221 0222 if ftbverbose 0223 fprintf('\tvariable %s... ', varname) 0224 end 0225 0226 % Get the size of the data and the dimension names. 0227 currDimsNames = dimnames(varDimIDs + 1); 0228 currDimsLengths = dimlengths(varDimIDs + 1); 0229 0230 % Find whether we've got an unlimited dimension in this data. 0231 wasUnlimited = -1; 0232 for jj = varDimIDs 0233 if numel(unlimdimID) > 1 0234 error('Do not currently support multiple unlimited dimensions.') 0235 end 0236 if strcmpi(dimnames(jj + 1), dimnames(unlimdimID + 1)) 0237 wasUnlimited = jj; 0238 end 0239 end 0240 0241 % Since the restart file has a number of time values, we'll ramp up the 0242 % replacement data from the existing start condition to the actual 0243 % value over the time steps. So, we need to know how many time steps we 0244 % actually have. 0245 0246 % Get the dimension data ready for the replacement arrays. 0247 tIdx = strncmp(dimnames(unlimdimID + 1), currDimsNames, length(dimnames{unlimdimID + 1})); 0248 nt = currDimsLengths(tIdx); 0249 % Not sure about the hardcoded strings below... 0250 sIdx = strncmp('siglay', currDimsNames, length(dimnames{unlimdimID + 1})); 0251 nIdx = strncmp('node', currDimsNames, length(dimnames{unlimdimID + 1})); 0252 ns = currDimsLengths(sIdx); 0253 nd = currDimsLengths(nIdx); 0254 if isempty(nd) 0255 % We've got data on the elements (i.e. u and v) 0256 nIdx = strncmp('nele', currDimsNames, length(dimnames{unlimdimID + 1})); 0257 nd = currDimsLengths(nIdx); 0258 end 0259 0260 % Iterate through the field names to see if we're on one of the 0261 % variables to be replaced. 0262 0263 % Set variable so we know if we've already written this variable to the 0264 % output file. 0265 writtenAlready = 0; 0266 for vv = 1:nf 0267 if strcmp(varname, fnames{vv}) && writtenAlready == 0 0268 if ftbverbose 0269 fprintf('NEW DATA... ') 0270 end 0271 0272 % Only repeat the values when we've bracketed the data, 0273 % otherwise scale up from the first time step. 0274 if bracketed 0275 % Use the new input data repeated the required number of 0276 % times. 0277 if wasUnlimited < 0 % no time, easy peasy. 0278 sfvdata = indata.(fnames{vv}); 0279 else 0280 % Damn. We've got to filter based on the shape of the 0281 % input data first (do we have a scalar input we need 0282 % to repeat to the shape of the output?) and then on 0283 % the basis of the expected output (is it 1D, 2D or 0284 % 3D?). 0285 if isscalar(indata.(fnames{vv})) 0286 if isempty(ns) && isempty(nd) 0287 sfvdata = indata.(fnames{vv}); 0288 elseif isempty(ns) || isempty(nd) 0289 if ftbverbose 0290 fprintf('tiling input scalar to non-scalar array... ') 0291 end 0292 sfvdata = repmat(indata.(fnames{vv}), [nd, nt]); 0293 elseif ~isempty(ns) && ~isempty(nd) 0294 if ftbverbose 0295 fprintf('tiling input scalar to non-scalar array... ') 0296 end 0297 sfvdata = repmat(indata.(fnames{vv}), [nd, ns, nt]); 0298 else 0299 error('Hmmm... FVCOM doesn''t have 4D arrays...') 0300 end 0301 else 0302 % We don't have scalar input so repeat the last 0303 % time in the given input nt times. 0304 if isempty(ns) && isempty(nd) 0305 sfvdata = indata.(fnames{vv}); 0306 elseif isempty(ns) || isempty(nd) 0307 sfvdata = repmat(indata.(fnames{vv})(:, :, end), [1, nt]); 0308 elseif ~isempty(ns) && ~isempty(nd) 0309 sfvdata = repmat(indata.(fnames{vv})(:, :, end), [1, 1, nt]); 0310 else 0311 error('Hmmm... FVCOM doesn''t have 4D arrays...') 0312 end 0313 end 0314 if ftbverbose 0315 fprintf('clipping in time... ') 0316 end 0317 end 0318 else 0319 % To make the scaling go from the initial value to the 0320 % supplied data value, we need to scale the difference 0321 % between the end members by the scaling factor at each 0322 % time and add to the current time's value. 0323 data = netcdf.getVar(nc, varid); 0324 if ~isscalar(data) && ~isempty(nt) 0325 if isscalar(indata.(fnames{vv})) && ftbverbose 0326 fprintf('tiling input scalar to non-scalar array... ') 0327 end 0328 sfvdata = nan(nd, ns, nt); 0329 % Scale from 0 to 1 if we have more than one time step, 0330 % otherwise just dump the data we've been given as 0331 % input. 0332 if nt == 1 0333 % Use the observed data. 0334 sfvdata = indata.(fnames{vv}); 0335 else 0336 ss = 0:1 / (nt - 1):1; % scale from 0 to 1. 0337 % Use the first modelled time step. 0338 startdata = squeeze(data(:, :, 1)); 0339 td = indata.(fnames{vv}) - startdata; 0340 for tt = 1:nt 0341 if tt == 1 0342 sfvdata(:, :, 1) = startdata; 0343 else 0344 sfvdata(:, :, tt) = startdata + (ss(tt) .* td); 0345 end 0346 end 0347 if ftbverbose 0348 fprintf('ramping data in time... ') 0349 end 0350 end 0351 elseif ~isscalar(data) && isempty(nt) 0352 sfvdata = indata.(fnames{vv}); 0353 else 0354 sfvdata = data; 0355 end 0356 end 0357 % Replace the values with the scaled interpolated values, 0358 % checking for unlimited dimensions as we go. 0359 if wasUnlimited < 0 0360 netcdf.putVar(ncout, varid, sfvdata) 0361 else 0362 netcdf.putVar(ncout, varid, zeros(length(currDimsLengths), 1), currDimsLengths, sfvdata) 0363 end 0364 0365 writtenAlready = 1; 0366 0367 % We might also want to replace the time. If so, supply a fourth 0368 % argument (out_date) to replace the times in the existing 0369 % restart file with an arbitrary time period. 0370 elseif strcmpi(varname, 'time') && writtenAlready == 0 && exist('out_date', 'var') 0371 if ftbverbose 0372 fprintf('NEW DATA... ') 0373 end 0374 tmp_time = greg2mjulian(out_date(:, 1), out_date(:, 2), out_date(:, 3), out_date(:, 4), out_date(:, 5), out_date(:, 6)); 0375 netcdf.putVar(ncout, varid, tmp_time) 0376 0377 writtenAlready = 1; 0378 0379 elseif strcmpi(varname, 'Times') && writtenAlready == 0 && exist('out_date', 'var') 0380 if ftbverbose 0381 fprintf('NEW DATA... ') 0382 end 0383 tmp_time = []; 0384 for i = 1:nt 0385 tmp_time = [tmp_time, sprintf('%-026s', datestr(datenum(out_date(i, :)), 'yyyy-mm-dd HH:MM:SS.FFF'))]; 0386 end 0387 netcdf.putVar(ncout, varid, tmp_time) 0388 0389 writtenAlready = 1; 0390 0391 elseif strcmpi(varname, 'Itime') && writtenAlready == 0 && exist('out_date', 'var') 0392 if ftbverbose 0393 fprintf('NEW DATA... ') 0394 end 0395 tmp_time = greg2mjulian(out_date(:, 1), out_date(:, 2), out_date(:, 3), out_date(:, 4), out_date(:, 5), out_date(:, 6)); 0396 netcdf.putVar(ncout, varid, floor(tmp_time)) 0397 0398 writtenAlready = 1; 0399 0400 elseif strcmpi(varname, 'Itime2') && writtenAlready == 0 && exist('out_date', 'var') 0401 if ftbverbose 0402 fprintf('NEW DATA... ') 0403 end 0404 tmp_time_mjd = greg2mjulian(out_date(:, 1), out_date(:, 2), out_date(:, 3), out_date(:, 4), out_date(:, 5), out_date(:, 6)); 0405 tmp_time = round((mod(tmp_time_mjd, 1) * 24 * 3600 * 1000) / (3600 * 1000)) * (3600 * 1000); 0406 netcdf.putVar(ncout, varid, floor(tmp_time)) 0407 0408 writtenAlready = 1; 0409 0410 end 0411 end 0412 0413 % If writtenAlready is zero, we haven't had one of the variables we're 0414 % replacing, so just dump the existing data. 0415 if writtenAlready == 0 0416 if ftbverbose 0417 fprintf('existing data... ') 0418 end 0419 0420 % Grab the data. 0421 data = netcdf.getVar(nc, varid); 0422 0423 % We need to check if the dimension is unlimited, and use a 0424 % start and end with netcdf.putVar if it is. This is largely 0425 % because MATLAB can't handle unlimited dimensions in the same 0426 % way as it does finite dimensions. 0427 if wasUnlimited < 0 0428 % We can just dump the entire data without specifying over 0429 % what indices. 0430 netcdf.putVar(ncout, varid, data); 0431 else 0432 % If we're clipping in time, use only the last value repeated 0433 % three times for the bracketed times. Otherwise, we just dump 0434 % it as is. FVCOM currently has no 4D variables: since the two 0435 % spatial dimensions are collapsed into one, the maximum number 0436 % of dimensions is 3: time, depth and position. 0437 if any(ismember(currDimsNames, 'time')) && exist('out_date', 'var') 0438 if ftbverbose 0439 fprintf('clipping in time... ') 0440 end 0441 if isvector(data) % 1D data 0442 data = data(end - (nt - 1):end); 0443 elseif ismatrix(data) % 2D data 0444 data = repmat(data(:, end), [1, nt]); 0445 else % 3D data 0446 data = repmat(data(:, :, end), [1, 1, nt]); 0447 end 0448 end 0449 0450 % Use the dimension length we extracted above to output the 0451 % data with the valid unlimited dimension format. 0452 netcdf.putVar(ncout, varid, zeros(length(currDimsLengths), 1), currDimsLengths, data); 0453 end 0454 end 0455 0456 if ftbverbose 0457 fprintf('done.\n') 0458 end 0459 end 0460 0461 netcdf.close(nc) 0462 netcdf.close(ncout) 0463 0464 if ftbverbose 0465 fprintf('end : %s \n', subname) 0466 end