Home > fvcom_prepro > replace_FVCOM_restart_vars.m

replace_FVCOM_restart_vars

PURPOSE ^

Use an FVCOM restart file to seed a model run with spatially varying

SYNOPSIS ^

function replace_FVCOM_restart_vars(Mobj, polcoms_ts, polcoms_z, start_date, fv_restart, varlist)

DESCRIPTION ^

 Use an FVCOM restart file to seed a model run with spatially varying
 versions of otherwise constant variables (temperature and salinity only
 for the time being).

 function replace_FVCOM_restart_vars(Mobj, polcoms_ts, polcoms_z, ...
    start_date, fv_restart, varlist)

 DESCRIPTION:
    FVCOM does not yet support spatially varying temperature and salinity
    inputs as initial conditions. To avoid having to run a model for a
    long time in order for temperature and salinity to settle within the
    model from the atmospheric and boundary forcing, we can use a restart
    file to cheat. If we run a model for a week or so (until the
    hydrodynamics has stabilised, we can use the restart file generated
    from that run as the basis for a new run, except we replace the
    currently computed temperature and salinity and replace them with data
    interpolated from another source (in this case, POLCOMS). 

 INPUT:
   Mobj        = MATLAB mesh structure which must contain:
                   - Mobj.siglayz - sigma layer depths for all model
                   nodes.
                   - Mobj.lon, Mobj.lat - node coordinates (long/lat).
   polcoms_ts  = POLCOMS NetCDF file in which 4D variables of temperature
   and salinity (called 'ETW' and 'x1X') exist. Their shape should be (y,
   x, sigma, time).
   polcoms_z   = POLCOMS NetCDF file in which 4D variables of bathymetry
   and sigma layer thickness can be found. They should be called 'depth'
   and 'pdepth' respectively.
   start_date  = Gregorian start date array (YYYY, MM, DD, hh, mm, ss).
   fv_restart  = full path to the FVCOM restart file.
   varlist     = cell array of variables to extract from the POLCOMS
   NetCDF files (polcoms_ts only, pdepth and depth will be extracted from
   polcoms_z).
 
 OUTPUT:
   NetCDF file in which the temperature and salinity variables have been
   replaced with the POLCOMS versions. The file name is the input restart
   file name appended _polcoms.nc.

 EXAMPLE USAGE
   replace_FVCOM_restart_vars(Mobj, '/tmp/pc_ts.nc', '/tmp/pc_z.nc', ...
       '2006-01-01 00:00:00', '/tmp/fvcom_restart.nc', ...
       {'lon', 'lat', 'ETW', 'x1X', 'time'})

 Author(s):
   Pierre Cazenave (Plymouth Marine Laboratory)

 Revision history
   2013-01-28 First version (based on initial_tsobc.m).

==========================================================================

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function replace_FVCOM_restart_vars(Mobj, polcoms_ts, polcoms_z, start_date, fv_restart, varlist)
0002 % Use an FVCOM restart file to seed a model run with spatially varying
0003 % versions of otherwise constant variables (temperature and salinity only
0004 % for the time being).
0005 %
0006 % function replace_FVCOM_restart_vars(Mobj, polcoms_ts, polcoms_z, ...
0007 %    start_date, fv_restart, varlist)
0008 %
0009 % DESCRIPTION:
0010 %    FVCOM does not yet support spatially varying temperature and salinity
0011 %    inputs as initial conditions. To avoid having to run a model for a
0012 %    long time in order for temperature and salinity to settle within the
0013 %    model from the atmospheric and boundary forcing, we can use a restart
0014 %    file to cheat. If we run a model for a week or so (until the
0015 %    hydrodynamics has stabilised, we can use the restart file generated
0016 %    from that run as the basis for a new run, except we replace the
0017 %    currently computed temperature and salinity and replace them with data
0018 %    interpolated from another source (in this case, POLCOMS).
0019 %
0020 % INPUT:
0021 %   Mobj        = MATLAB mesh structure which must contain:
0022 %                   - Mobj.siglayz - sigma layer depths for all model
0023 %                   nodes.
0024 %                   - Mobj.lon, Mobj.lat - node coordinates (long/lat).
0025 %   polcoms_ts  = POLCOMS NetCDF file in which 4D variables of temperature
0026 %   and salinity (called 'ETW' and 'x1X') exist. Their shape should be (y,
0027 %   x, sigma, time).
0028 %   polcoms_z   = POLCOMS NetCDF file in which 4D variables of bathymetry
0029 %   and sigma layer thickness can be found. They should be called 'depth'
0030 %   and 'pdepth' respectively.
0031 %   start_date  = Gregorian start date array (YYYY, MM, DD, hh, mm, ss).
0032 %   fv_restart  = full path to the FVCOM restart file.
0033 %   varlist     = cell array of variables to extract from the POLCOMS
0034 %   NetCDF files (polcoms_ts only, pdepth and depth will be extracted from
0035 %   polcoms_z).
0036 %
0037 % OUTPUT:
0038 %   NetCDF file in which the temperature and salinity variables have been
0039 %   replaced with the POLCOMS versions. The file name is the input restart
0040 %   file name appended _polcoms.nc.
0041 %
0042 % EXAMPLE USAGE
0043 %   replace_FVCOM_restart_vars(Mobj, '/tmp/pc_ts.nc', '/tmp/pc_z.nc', ...
0044 %       '2006-01-01 00:00:00', '/tmp/fvcom_restart.nc', ...
0045 %       {'lon', 'lat', 'ETW', 'x1X', 'time'})
0046 %
0047 % Author(s):
0048 %   Pierre Cazenave (Plymouth Marine Laboratory)
0049 %
0050 % Revision history
0051 %   2013-01-28 First version (based on initial_tsobc.m).
0052 %
0053 %==========================================================================
0054 
0055 subname = 'replace_FVCOM_restart_vars';
0056 
0057 global ftbverbose;
0058 if ftbverbose
0059     fprintf('\n')
0060     fprintf(['begin : ' subname '\n'])
0061 end
0062 
0063 %--------------------------------------------------------------------------
0064 % Extract the POLCOMS data specified in varlist
0065 %--------------------------------------------------------------------------
0066 
0067 if ftbverbose
0068     fprintf('%s : read POLCOMS data... ', subname)
0069 end
0070 
0071 nc = netcdf.open(polcoms_ts, 'NOWRITE');
0072 
0073 for var=1:numel(varlist)
0074     
0075     getVar = varlist{var};
0076     varid_pc = netcdf.inqVarID(nc, getVar);
0077     
0078     data = netcdf.getVar(nc, varid_pc, 'single');
0079     pc.(getVar).data = double(data);
0080     % Try to get some units (important for the calculation of MJD).
0081     try
0082         units = netcdf.getAtt(nc,varid_pc,'units');
0083     catch err
0084         units = [];
0085     end
0086     pc.(getVar).units = units;
0087 end
0088 clear data getVar varid_pc var
0089 
0090 netcdf.close(nc)
0091 
0092 % Extract the bathymetry ('pdepth' is cell thickness, 'depth' is cell
0093 % centre depth).
0094 nc = netcdf.open(polcoms_z, 'NOWRITE');
0095 varid_pc = netcdf.inqVarID(nc, 'depth'); 
0096 pc.depth.data = double(netcdf.getVar(nc, varid_pc, 'single'));
0097 try
0098     pc.depth.units = netcdf.getAtt(nc, varid_pc, 'units');
0099 catch err
0100     pc.depth.units = [];
0101 end
0102 varid_pc = netcdf.inqVarID(nc, 'pdepth'); 
0103 pc.pdepth.data = double(netcdf.getVar(nc, varid_pc, 'single'));
0104 try
0105     pc.pdepth.units = netcdf.getAtt(nc, varid_pc, 'units');
0106 catch err
0107     pc.pdepth.units = [];
0108 end
0109 netcdf.close(nc)
0110 
0111 % Data format:
0112 %
0113 %   pc.ETW.data and pc.x1X.data are y, x, sigma, time
0114 %
0115 [ny, nx, ~, ~] = size(pc.ETW.data);
0116 
0117 % Number of sigma layers.
0118 [fn, fz] = size(Mobj.siglayz);
0119 
0120 % Make rectangular arrays for the nearest point lookup.
0121 [lon, lat] = meshgrid(pc.lon.data, pc.lat.data);
0122 
0123 % Convert the POLCOMS times to Modified Julian Day (this is a very ugly).
0124 pc.time.yyyymmdd = strtrim(regexp(pc.time.units, 'since', 'split'));
0125 pc.time.strtime = regexp(pc.time.yyyymmdd{end}, '-', 'split');
0126 % This new version of the time has the year in a weird format (yr.#). We
0127 % thus need to split it again to get the decimal year (post-2000 only?).
0128 pc.time.strtimeyr = regexp(pc.time.strtime, '\.', 'split');
0129 pc.time.yyyymmdd = str2double([pc.time.strtimeyr{1}(2), pc.time.strtime{2:3}]);
0130 pc.time.yyyymmdd(1) = pc.time.yyyymmdd(1) + 2000; % add full year.
0131 Mobj.ts_times = greg2mjulian(pc.time.yyyymmdd(1), pc.time.yyyymmdd(2), pc.time.yyyymmdd(3), 0, 0, 0) + (pc.time.data / 3600 / 24);
0132 
0133 % Given our intput time (in start_date), find the nearest time
0134 % index for the POLCOMS data.
0135 stime = greg2mjulian(start_date(1), start_date(2), ...
0136     start_date(3), start_date(4), ...
0137     start_date(5), start_date(6));
0138 [~, tidx] = min(abs(Mobj.ts_times - stime));
0139 
0140 if ftbverbose
0141     fprintf('done.\n') 
0142 end
0143 
0144 %--------------------------------------------------------------------------
0145 % Interpolate the regularly gridded POLCOMS data onto the FVCOM grid
0146 % (vertical grid first).
0147 %--------------------------------------------------------------------------
0148 
0149 if ftbverbose
0150     fprintf('%s : interpolate POLCOMS onto FVCOM''s vertical grid... ', subname)
0151 end
0152 
0153 % Preallocate the output arrays
0154 pctempz = nan(nx, ny, fz);
0155 pcsalz = nan(nx, ny, fz);
0156 
0157 for xi = 1:nx
0158     % For the parallel loop, get all the y and z dimension data for the
0159     % current x position (temperature, salinity and depth).
0160     % N.B. The POLCOMS data is stored y, x, z, t (first two dimensions the
0161     % wrong way around).
0162     xtemp = squeeze(pc.ETW.data(:, xi, :, tidx));
0163     xsalt = squeeze(pc.x1X.data(:, xi, :, tidx));
0164     xdepth = squeeze(pc.depth.data(:, xi, :, tidx));
0165     
0166     % Preallocate the arrays for the inner loop
0167     ytemp = nan(ny, fz);
0168     ysalt = nan(ny, fz);
0169     for yi = 1:ny
0170         % First things first, check we're within the POLCOMS domain
0171         % (assuming anything beyond a maximum depth of 20km is outside the
0172         % domain).
0173         if xdepth(yi, end) < -20000 % use deepest depth
0174             continue
0175         end
0176 
0177         % Find the nearest sigma layer z values from the unstructured grid.
0178         [md, mi] = min(sqrt((Mobj.lon - lon(xi, yi)).^2 + (Mobj.lat - lat(xi, yi)).^2));
0179 
0180         % Skip data point if the closest FVCOM node is more than 10 minutes
0181         % away.
0182         if md > 10 / 60
0183             continue
0184         else
0185             % Use the FVCOM node's sigma depths to interpolate this POLCOMS
0186             % position's temperature and salinity data.
0187             
0188             % Get the FVCOM depths closest to this POLCOMS grid position.
0189             tfz = Mobj.siglayz(mi, :);
0190             % Now get the POLCOMS depth at this node for the time index
0191             % identified above.
0192             tpz = xdepth(yi, :); 
0193             ytemp(yi, :) = interp1(tpz, xtemp(yi, :), tfz, 'linear', 'extrap');
0194             
0195             ysalt(yi, :) = interp1(tpz, xsalt(yi, :), tfz, 'linear', 'extrap');
0196         end
0197     end
0198     pctempz(xi, :, :) = ytemp;
0199     pcsalz(xi, :, :) = ysalt;
0200 end
0201 % Tidy up the namespace a bit.
0202 clear ytemp ysalt tfz tpz md mi xtemp xsalt xdepth xi yi zi
0203 
0204 if ftbverbose
0205     fprintf('done.\n') 
0206 end
0207 
0208 %--------------------------------------------------------------------------
0209 % Now we have vertically interpolated POLCOMS data, we can interpolate each
0210 % sigma layer onto the FVCOM unstructured grid ready to write out to
0211 % NetCDF. We'll use the triangular interpolation in MATLAB with the natural
0212 % method (gives pretty good results, at least qualitatively).
0213 %--------------------------------------------------------------------------
0214 
0215 if ftbverbose
0216     fprintf('%s : interpolate POLCOMS onto FVCOM''s horizontal grid... ', subname)
0217 end
0218 
0219 fvtemp = nan(fn, fz);
0220 fvsalt = nan(fn, fz);
0221 for zi = 1:fz
0222     % Set up the interpolation object.
0223     ft = TriScatteredInterp(lon(:), lat(:), reshape(pctempz(:, :, zi), [], 1), 'natural');
0224     fs = TriScatteredInterp(lon(:), lat(:), reshape(pcsalz(:, :, zi), [], 1), 'natural');
0225     % Interpolate temperature and salinity onto the unstructured grid.
0226     fvtemp(:, zi) = ft(Mobj.lon, Mobj.lat);
0227     fvsalt(:, zi) = fs(Mobj.lon, Mobj.lat);
0228 end
0229 
0230 % Unfortunately, TriScatteredInterp won't extrapolate, returning instead
0231 % NaNs outside the original data's extents. So, for each NaN position, find
0232 % the nearest non-NaN value and use that instead. The order in which the
0233 % NaN-nodes are found will determine the spatial pattern of the
0234 % extrapolation.
0235 
0236 % We can assume that all layers will have NaNs in the same place
0237 % (horizontally), so just use the surface layer (1) for the identification
0238 % of NaNs. Also store the finite values so we can find the nearest real
0239 % value to the current NaN node and use its temperature and salinity
0240 % values.
0241 fvidx = 1:fn;
0242 fvnanidx = fvidx(isnan(fvtemp(:, 1)));
0243 fvfinidx = fvidx(~isnan(fvtemp(:, 1)));
0244 for ni = 1:length(fvnanidx);
0245     % Current position
0246     xx = Mobj.lon(fvnanidx(ni));
0247     yy = Mobj.lat(fvnanidx(ni));
0248     % Find the nearest non-nan temperature and salinity value.
0249     [~, di] = min(sqrt((Mobj.lon(fvfinidx) - xx).^2 + (Mobj.lat(fvfinidx) - yy).^2));
0250     % Replace the temperature and salinity values at all depths at the
0251     % current NaN position with the closest non-nan value.
0252     fvtemp(fvnanidx(ni), :) = fvtemp(fvfinidx(di), :);
0253     fvsalt(fvnanidx(ni), :) = fvsalt(fvfinidx(di), :);
0254 end
0255 
0256 if ftbverbose
0257     fprintf('done.\n') 
0258 end
0259 
0260 %--------------------------------------------------------------------------
0261 % Get the restart data and replace with the interpolated data.
0262 %--------------------------------------------------------------------------
0263 
0264 if ftbverbose
0265     fprintf('%s : export interpolated data to NetCDF:\n', subname)
0266 end
0267 
0268 nc = netcdf.open(fv_restart, 'NOWRITE');
0269 [pp, nn, ee] = fileparts(fv_restart);
0270 ncout = netcdf.create(fullfile(pp, [nn, '_polcoms', ee]), 'clobber');
0271 % ncout = netcdf.create(fullfile('/tmp/pica/', [nn, '_polcoms', ee]), 'clobber');
0272 
0273 [numdims, numvars, numglobatts, unlimdimID] = netcdf.inq(nc);
0274 
0275 % Define the dimensions for all variables.
0276 dimid = nan(numdims, 1);
0277 dimnames = cell(numdims, 1);
0278 dimlengths = nan(numdims, 1);
0279 for ii = 1:numdims
0280     [dimname, dimlen] = netcdf.inqDim(nc, ii - 1);
0281     if ii ~= unlimdimID + 1 % NetCDF indices start at zero
0282         dimid(ii) = netcdf.defDim(ncout, dimname, dimlen);
0283     else
0284         dimid(ii) = netcdf.defDim(ncout, dimname, netcdf.getConstant('NC_UNLIMITED'));
0285     end
0286     dimnames{ii} = dimname;
0287     dimlengths(ii) = dimlen;
0288 end
0289 
0290 % Now define the variables and attributes.
0291 for ii = 1:numvars
0292 
0293     % Find name of the current variable.
0294     [varname, xtype, varDimIDs, varAtts] = netcdf.inqVar(nc, ii - 1);
0295 
0296     % Create the variables.
0297     varid = netcdf.defVar(ncout, varname, xtype, varDimIDs);
0298 
0299     % Get each attribute and add it to the current variable.
0300     for j = 1:varAtts
0301 
0302         attname = netcdf.inqAttName(nc, varid, j - 1);
0303         attval = netcdf.getAtt(nc, varid, attname);
0304 
0305         netcdf.putAtt(ncout, varid, attname, attval);
0306     end
0307 end
0308 
0309 % Do the global attributes
0310 for ii = 1:numglobatts
0311     
0312     % Find the current global attribute's name and value.
0313     gattname = netcdf.inqAttName(nc, netcdf.getConstant('NC_GLOBAL'), ii - 1);
0314     gattval = netcdf.getAtt(nc, netcdf.getConstant('NC_GLOBAL'), gattname);
0315     
0316     % Put that back into the output NetCDF file.
0317     netcdf.putAtt(ncout, netcdf.getConstant('NC_GLOBAL'), gattname, gattval);
0318 end
0319 
0320 netcdf.endDef(ncout);
0321 
0322 % Get the existing data and output to the new NetCDF file, except for
0323 % temperature and salinity, where we replace the values with the
0324 % interpolated POLCOMS data.
0325 for ii = 1:numvars
0326 
0327     [varname, ~, varDimIDs, ~] = netcdf.inqVar(nc, ii - 1);
0328     varid = netcdf.inqVarID(nc, varname);
0329 
0330     if ftbverbose
0331         fprintf('\tvariable %s... ', varname)
0332     end
0333 
0334     % We need the data irrespective of whether we're replacing it or not,
0335     % so grab it outside the if statement below.
0336     data = netcdf.getVar(nc, varid);
0337 
0338     % Get the size of the data and the dimension names.
0339     currDimsNames = dimnames(varDimIDs + 1);
0340     currDimsLengths = dimlengths(varDimIDs + 1);
0341 
0342     % Find whether we've got an unlimited dimension in this data.
0343     wasUnlimited = -1;
0344     for jj = varDimIDs
0345         if numel(unlimdimID) > 1
0346             error('Do not currently support multiple unlimited dimensions.')
0347         end
0348         if strcmpi(dimnames(jj + 1), dimnames(unlimdimID + 1))
0349             wasUnlimited = jj;
0350         end
0351     end
0352 
0353     % Since the restart file has a number of time values, we'll ramp up
0354     % the temperature from some constant to the actual value over the
0355     % time steps. So, we need to know how many time steps we actually
0356     % have.
0357 
0358     % Get the dimension data ready for the temperature and salinity arrays.
0359     tIdx = strncmp(dimnames(unlimdimID + 1), currDimsNames, length(dimnames{unlimdimID + 1}));
0360     % Not sure about the hardcoded strings below...
0361     sIdx = strncmp('siglay', currDimsNames, length(dimnames{unlimdimID + 1}));
0362     nIdx = strncmp('node', currDimsNames, length(dimnames{unlimdimID + 1}));
0363     nt = currDimsLengths(tIdx);
0364     ns = currDimsLengths(sIdx);
0365     nd = currDimsLengths(nIdx);
0366 
0367     if strcmpi(varname, 'temp') || strcmpi(varname, 'salinity')
0368         % To make the scaling go from the initial value to the POLCOMS value,
0369         % we need to take the scale the difference between the end members by
0370         % the scaling factor at each time and add to the current time's value.
0371         sfvdata = nan(nd, ns, nt);
0372         ss = 0:1 / (nt - 1):1; % scale from 0 to 1.
0373         startdata = squeeze(data(:, :, 1)); % use the first modelled time step
0374         for tt = 1:nt;
0375             if tt == 1
0376                 sfvdata(:, :, 1) = startdata;
0377             else
0378                 td = fvtemp - startdata;
0379                 sfvdata(:, :, tt) = startdata + (ss(tt) .* td);
0380             end
0381         end
0382 
0383         % Replace the values with the scaled interpolated values.
0384         netcdf.putVar(ncout, varid, sfvdata)
0385     else
0386         % We need to check if the dimension is unlimited, and use a start
0387         % and end with netcdf.putVar if it is. This is largely because
0388         % MATLAB can't handle unlimited dimensions in the same way as it
0389         % does finite dimensions.
0390         if wasUnlimited < 0
0391             % We can just dump the entire data without specifying over what
0392             % indices.
0393             netcdf.putVar(ncout, varid, data);
0394         else
0395             % Use the dimension length we extracted above to output the
0396             % data with the valid unlimited dimension format.
0397             netcdf.putVar(ncout, varid, zeros(length(currDimsLengths), 1), currDimsLengths, data);
0398         end
0399     end
0400 
0401     if ftbverbose
0402         fprintf('done.\n')
0403     end
0404 end
0405 
0406 netcdf.close(nc)
0407 netcdf.close(ncout)
0408 
0409 if ftbverbose
0410     fprintf('%s : export complete.\n', subname)
0411 end
0412 
0413 if ftbverbose
0414     fprintf(['end   : ' subname '\n'])
0415 end

Generated on Mon 04-Feb-2013 14:22:28 by m2html © 2005