Home > fvcom_prepro > interp_NEMO2FVCOM.m

interp_NEMO2FVCOM

PURPOSE ^

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

SYNOPSIS ^

function Mobj = interp_NEMO2FVCOM(Mobj, nemo, start_date, varlist)

DESCRIPTION ^

 Use an FVCOM restart file to seed a model run with spatially varying
 versions of otherwise constant variables.

 function interp_NEMO2FVCOM(Mobj, nemo, start_date, varlist)

 DESCRIPTION:
    FVCOM does not yet support spatially varying 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. For this, we
    need variables interpolated onto the unstructured grid. The
    interpolated data can then be written out with write_FVCOM_restart.m.

 INPUT:
   Mobj        = MATLAB mesh structure which must contain:
                   - Mobj.lon, Mobj.lat and/or Mobj.lonc, Mobj.latc - node
                   or element centre coordinates (long/lat). Dependent on
                   the variable being interpolated (u and v need lonc,
                   latc whilst temperature and salinity need lat and lon).
                   - Mobj.siglayz - sigma layer depths for all model
                   nodes.
                   - Mobj.ts_times - Modified Julian Day times for the
                   model run.
   nemo        = Struct output from some other function. Must include
                 fields:
                   - nemo.lon, nemo.lat - rectangular arrays.
                   - nemo.Depth - NEMO depth levels.
   start_date  = Gregorian start date array (YYYY, MM, DD, hh, mm, ss).
   varlist     = cell array of fields to use from the NEMO struct.

 OUTPUT:
   Mobj.restart = struct whose field names are the variables which have
   been interpolated (e.g. Mobj.restart.temperature for NEMO
   temperature).

 EXAMPLE USAGE
   interp_NEMO2FVCOM(Mobj, nemo, [2006, 01, 01, 00, 00, 00], ...
       {'lon', 'lat', 'time', 'temperature', 'salinity'})

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

 Revision history
   2018-05-23 First version based on interp_HYCOMS2FVCOM.m.

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

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function Mobj = interp_NEMO2FVCOM(Mobj, nemo, start_date, varlist)
0002 % Use an FVCOM restart file to seed a model run with spatially varying
0003 % versions of otherwise constant variables.
0004 %
0005 % function interp_NEMO2FVCOM(Mobj, nemo, start_date, varlist)
0006 %
0007 % DESCRIPTION:
0008 %    FVCOM does not yet support spatially varying initial conditions. To
0009 %    avoid having to run a model for a long time in order for temperature
0010 %    and salinity to settle within the model from the atmospheric and
0011 %    boundary forcing, we can use a restart file to cheat. For this, we
0012 %    need variables interpolated onto the unstructured grid. The
0013 %    interpolated data can then be written out with write_FVCOM_restart.m.
0014 %
0015 % INPUT:
0016 %   Mobj        = MATLAB mesh structure which must contain:
0017 %                   - Mobj.lon, Mobj.lat and/or Mobj.lonc, Mobj.latc - node
0018 %                   or element centre coordinates (long/lat). Dependent on
0019 %                   the variable being interpolated (u and v need lonc,
0020 %                   latc whilst temperature and salinity need lat and lon).
0021 %                   - Mobj.siglayz - sigma layer depths for all model
0022 %                   nodes.
0023 %                   - Mobj.ts_times - Modified Julian Day times for the
0024 %                   model run.
0025 %   nemo        = Struct output from some other function. Must include
0026 %                 fields:
0027 %                   - nemo.lon, nemo.lat - rectangular arrays.
0028 %                   - nemo.Depth - NEMO depth levels.
0029 %   start_date  = Gregorian start date array (YYYY, MM, DD, hh, mm, ss).
0030 %   varlist     = cell array of fields to use from the NEMO struct.
0031 %
0032 % OUTPUT:
0033 %   Mobj.restart = struct whose field names are the variables which have
0034 %   been interpolated (e.g. Mobj.restart.temperature for NEMO
0035 %   temperature).
0036 %
0037 % EXAMPLE USAGE
0038 %   interp_NEMO2FVCOM(Mobj, nemo, [2006, 01, 01, 00, 00, 00], ...
0039 %       {'lon', 'lat', 'time', 'temperature', 'salinity'})
0040 %
0041 % Author(s):
0042 %   Ricardo Torres (Plymouth Marine Laboratory)
0043 %   Pierre Cazenave (Plymouth Marine Laboratory)
0044 %
0045 % Revision history
0046 %   2018-05-23 First version based on interp_HYCOMS2FVCOM.m.
0047 %
0048 %==========================================================================
0049 
0050 subname = 'interp_NEMO2FVCOM';
0051 
0052 global ftbverbose
0053 if ftbverbose
0054     fprintf('\nbegin : %s\n', subname)
0055 end
0056 
0057 if nargin == 0
0058     error('Not enough input arguments. See HELP %s', subname)
0059 end
0060 
0061 % Run jobs on multiple workers if we have that functionality. Not sure if
0062 % it's necessary, but check we have the Parallel Toolbox first.
0063 wasOpened = false;
0064 if license('test', 'Distrib_Computing_Toolbox')
0065     % We have the Parallel Computing Toolbox, so launch a bunch of workers.
0066     % New version for MATLAB 2014a (I think) onwards.
0067     if isempty(gcp('nocreate'))
0068         pool = parpool('local');
0069         wasOpened = true;
0070     end
0071 end
0072 
0073 % Given our input time (in start_date), find the nearest time index for
0074 % the NEMO data.
0075 stime = greg2mjulian(start_date(1), start_date(2), ...
0076     start_date(3), start_date(4), ...
0077     start_date(5), start_date(6));
0078 [~, tidx] = min(abs(nemo.time - stime));
0079 
0080 for vv = 1:length(varlist)
0081 
0082     currvar = varlist{vv};
0083 
0084     switch currvar
0085         case {'lon', 'lat', 'longitude', 'latitude', 't_1', 'time', 'Depth', 'MT'}
0086             continue
0087 
0088         otherwise
0089             %--------------------------------------------------------------
0090             % Interpolate the regularly gridded data onto the FVCOM grid
0091             % (vertical grid first).
0092             %--------------------------------------------------------------
0093 
0094             tic
0095 
0096             % Set up all the constants which are based on the model and
0097             % data grids.
0098             [~, fz] = size(Mobj.siglayz);
0099             if strcmpi(currvar, 'u') || strcmpi(currvar, 'v')
0100                 plon = nemo.lon(:);
0101                 plat = nemo.lat(:);
0102                 flon = Mobj.lonc;
0103                 flat = Mobj.latc;
0104                 fn = numel(flon);
0105                 fvtemp = nan(fn, fz);
0106             else
0107                 plon = nemo.lon(:);
0108                 plat = nemo.lat(:);
0109                 flon = Mobj.lon;
0110                 flat = Mobj.lat;
0111                 fn = numel(flon);
0112                 fvtemp = nan(fn, fz);
0113             end
0114 
0115             if ftbverbose
0116                 fprintf('%s : interpolate %s ... ', subname, currvar)
0117             end
0118 
0119             % Calculate resolution of parent coarse model
0120             hdx = max([max(diff(plon(:, 1))), max(diff(plat(1, :)))]);
0121 
0122             xlimits = [min(flon) - 3 * hdx, max(flon) + 3 * hdx];
0123             ylimits = [min(flat) - 3 * hdx, max(flat) + 3 * hdx];
0124             xidxgood = find(nemo.lon<xlimits(2) & nemo.lon>xlimits(1) & nemo.lat<ylimits(2) & nemo.lat>ylimits(1));
0125 
0126             hdepth = nemo.Depth;
0127 
0128             % Do not interpolate 2D variables in depth.
0129             if ~ismatrix(nemo.(currvar))
0130                 hyinterp.(currvar) = grid_vert_interp(Mobj, ...
0131                     nemo.lon, nemo.lat, ...
0132                     squeeze(nemo.(currvar)(:, :, :, tidx)), ...
0133                     hdepth, isnan(nemo.land_mask), 'extrapolate', ...
0134                     [flon, flat]);
0135 
0136                 %----------------------------------------------------------
0137                 % Now we have vertically interpolated data, we can
0138                 % interpolate each sigma layer onto the FVCOM unstructured
0139                 % grid ready to write out to NetCDF. We'll use the
0140                 % triangular interpolation in MATLAB with the natural
0141                 % method (gives pretty good results, at least
0142                 % qualitatively).
0143                 %----------------------------------------------------------
0144 
0145                 if ftbverbose
0146                     fprintf('horizontally... ')
0147                 end
0148 
0149                 hytempz = hyinterp.(currvar);
0150 
0151                 parfor zi = 1:fz
0152                     % Get the current depth layer's data and mask out the
0153                     % NaN values.
0154                     hytempzcurrent = hytempz(:, :, zi);
0155 
0156                     % Strip out NaNs so we can extrapolate with
0157                     % TriScatteredInterp.
0158                     nanmask = ~isnan(hytempzcurrent);
0159                     plonclean = plon(nanmask);
0160                     platclean = plat(nanmask);
0161                     hytempzclean = hytempzcurrent(nanmask);
0162                     % Set up the interpolation object and interpolate the
0163                     % current variable to the FVCOM unstructured grid.
0164                     ft = scatteredInterpolant(plonclean, platclean, ...
0165                         hytempzclean, 'natural', 'none');
0166                     fvtemp(:, zi) = ft(flon, flat);
0167                 end
0168             else
0169                 % Remove points outside our area of interest and
0170                 % interpolate this non-depth-resolved data.
0171                 datatmp = nemo.(currvar)(xidxgood);
0172                 plon = nemo.lon(xidxgood);
0173                 plat = nemo.lat(xidxgood);
0174                 % Strip out NaNs so we can extrapolate with
0175                 % TriScatteredInterp.
0176                 nanmask = ~isnan(datatmp);
0177                 plonclean = plon(nanmask);
0178                 platclean = plat(nanmask);
0179                 hytempzclean = datatmp(nanmask);
0180 
0181                 ft = scatteredInterpolant(plonclean, platclean, ...
0182                     hytempzclean, 'natural', 'none');
0183                 fvtemp = ft(flon, flat);
0184             end
0185 
0186             % Unfortunately, TriScatteredInterp won't extrapolate,
0187             % returning instead NaNs outside the original data's extents.
0188             % So, for each NaN position, find the nearest non-NaN value and
0189             % use that instead. The order in which the NaN-nodes are found
0190             % will determine the spatial pattern of the extrapolation.
0191 
0192             % We can assume that all layers will have NaNs in the same
0193             % place (horizontally), so just use the surface layer (1) for
0194             % the identification of NaNs. Also store the finite values so
0195             % we can find the nearest real value to the current NaN node
0196             % and use its temperature and salinity values.
0197             fvidx = 1:fn;
0198             fvnanidx = fvidx(isnan(fvtemp(:, 1)));
0199             fvfinidx = fvidx(~isnan(fvtemp(:, 1)));
0200             if isfield(Mobj,'dist')
0201                 [~,idx2coast]=sort(Mobj.dist(fvnanidx),'descend');
0202                 fvnanidx = fvnanidx(idx2coast);
0203             end
0204             for ni = 1:length(fvnanidx)
0205                 % Find the nearest non-nan data (temp, salinity, u or v)
0206                 % value.
0207                 xx = flon(fvnanidx(ni));
0208                 yy = flat(fvnanidx(ni));
0209 
0210                 [~, di] = min(sqrt((flon(fvfinidx) - xx).^2 + (flat(fvfinidx) - yy).^2));
0211                 fvtemp(fvnanidx(ni), :) = fvtemp(fvfinidx(di), :);
0212                 fvfinidx = fvidx(~isnan(fvtemp(:, 1)));
0213 
0214             end
0215 
0216             clear plon plat flon flat ptempz
0217 
0218             Mobj.restart.(currvar) = fvtemp;
0219 
0220             te = toc;
0221 
0222             if ftbverbose
0223                 fprintf('done. (elapsed time = %.2f seconds)\n', te)
0224             end
0225 
0226     end
0227 end
0228 
0229 % Close the MATLAB pool if we opened it.
0230 if wasOpened
0231     pool.delete
0232 end
0233 
0234 if ftbverbose
0235     fprintf('end   : %s\n', subname)
0236 end
0237 
0238 %% Debugging figure
0239 
0240 % close all
0241 %
0242 % ri = 85; % column index
0243 % ci = 95; % row index
0244 %
0245 % [~, idx] = min(sqrt((Mobj.lon - lon(ri, ci)).^2 + (Mobj.lat - lat(ri, ci)).^2));
0246 %
0247 % % Vertical profiles
0248 % figure
0249 % clf
0250 %
0251 % % The top row shows the temperature/salinity values as plotted against
0252 % % index (i.e. position in the array). I had thought POLCOMS stored its data
0253 % % from seabed to sea surface, but looking at the NetCDF files in ncview,
0254 % % it's clear that the data are in fact stored surface to seabed (like
0255 % % FVCOM). As such, the two plots in the top row should be upside down (i.e.
0256 % % surface values at the bottom of the plot). The two lower rows should have
0257 % % three lines which all match: the raw POLCOMS data, the POLCOMS data for
0258 % % the current time step (i.e. that in 'temperature' and 'salinity') and the
0259 % % interpolated FVCOM data against depth.
0260 % %
0261 % % Also worth noting, the pc.* have the rows and columns flipped, so
0262 % % (ci, ri) in pc.* and (ri, ci) in 'temperature', 'salinity' and
0263 % % 'depth'. Needless to say, the two lines in the lower plots should
0264 % % overlap.
0265 %
0266 % subplot(2,2,1)
0267 % plot(squeeze(pc.ETWD(ci, ri, :, tidx)), 1:size(depth, 3), 'rx:')
0268 % xlabel('Temperature (^{\circ}C)')
0269 % ylabel('Array index')
0270 % title('Array Temperature')
0271 %
0272 % subplot(2,2,2)
0273 % plot(squeeze(pc.x1XD(ci, ri, :, tidx)), 1:size(depth, 3), 'rx:')
0274 % xlabel('Salinity')
0275 % ylabel('Array index')
0276 % title('Array Salinity')
0277 %
0278 % subplot(2,2,3)
0279 % % Although POLCOMS stores its temperature values from seabed to surface,
0280 % % the depths are stored surface to seabed. Nice. Flip the
0281 % % temperature/salinity data accordingly. The lines here should match one
0282 % % another.
0283 % plot(squeeze(pc.ETWD(ci, ri, :, tidx)), squeeze(pc.depth(ci, ri, :, tidx)), 'rx-')
0284 % hold on
0285 % plot(squeeze(temperature(ri, ci, :)), squeeze(depth(ri, ci, :)), '.:')
0286 % % Add the equivalent interpolated FVCOM data point
0287 % plot(fvtemp(idx, :), Mobj.siglayz(idx, :), 'g.-')
0288 % xlabel('Temperature (^{\circ}C)')
0289 % ylabel('Depth (m)')
0290 % title('Depth Temperature')
0291 % legend('pc', 'temp', 'fvcom', 'location', 'north')
0292 % legend('boxoff')
0293 %
0294 % subplot(2,2,4)
0295 % % Although POLCOMS stores its temperature values from seabed to surface,
0296 % % the depths are stored surface to seabed. Nice. Flip the
0297 % % temperature/salinity data accordingly. The lines here should match one
0298 % % another.
0299 % plot(squeeze(pc.x1XD(ci, ri, :, tidx)), squeeze(pc.depth(ci, ri, :, tidx)), 'rx-')
0300 % hold on
0301 % plot(squeeze(salinity(ri, ci, :)), squeeze(depth(ri, ci, :)), '.:')
0302 % % Add the equivalent interpolated FVCOM data point
0303 % plot(fvsalt(idx, :), Mobj.siglayz(idx, :), 'g.-')
0304 % xlabel('Salinity')
0305 % ylabel('Depth (m)')
0306 % title('Depth Salinity')
0307 % legend('pc', 'salt', 'fvcom', 'location', 'north')
0308 % legend('boxoff')
0309 %
0310 % %% Plot the sample location
0311 % figure
0312 % dx = mean(diff(pc.lon));
0313 % dy = mean(diff(pc.lat));
0314 % z = depth(:, :, end); % water depth (bottom layer depth)
0315 % z(mask) = 0; % clear out nonsense values
0316 % pcolor(lon - (dx / 2), lat - (dy / 2), z)
0317 % shading flat
0318 % axis('equal', 'tight')
0319 % daspect([1.5, 1, 1])
0320 % colorbar
0321 % caxis([-150, 0])
0322 % hold on
0323 % plot(lon(ri, ci), lat(ri, ci), 'ko', 'MarkerFaceColor', 'w')
0324

Generated on Wed 20-Feb-2019 16:06:01 by m2html © 2005