


Use an FVCOM restart file to seed a model run with spatially varying
versions of otherwise constant variables.
function interp_HYCOM2FVCOM(Mobj, hycom, start_date, 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. For this, we need temperature and salinity
(potentially other variables too) 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.
hycom = Struct output by get_HYCOM_forcing. Must include fields:
- hycom.lon, hycom.lat - rectangular arrays.
- hycom.Depth - HYCOM depth levels.
The data fields (specified in varlist) must be shaped
(x, y, z, time).
start_date = Gregorian start date array (YYYY, MM, DD, hh, mm, ss).
varlist = cell array of fields to use from the HYCOM struct.
OUTPUT:
Mobj.restart = struct whose field names are the variables which have
been interpolated (e.g. Mobj.restart.temperature for HYCOM
temperature).
EXAMPLE USAGE
interp_HYCOM2FVCOM(Mobj, hycom, [2006, 01, 01, 00, 00, 00], ...
{'lon', 'lat', 'time', 'temperature', 'salinity'})
Author(s):
Pierre Cazenave (Plymouth Marine Laboratory)
Revision history
2013-09-10 First version based on interp_POLCOMS2FVCOM.m.
2013-12-10 Fix the identification of the time index in the HYCOM data
(use hycom.time instead of Mobj.ts_times). Also ignore a field name of
'MT' if supplied in varlist.
2014-04-28 Fix bug when interpolating velocity data due to incorrectly
sized preallocated array. Make the way the data to be interpolated onto
the element centres (i.e. the velocity data), as opposed to the element
nodes, more understandable. This is done by having a single block which
deals with setting up the variables dependent on the number of points
onto which to interpolate. Also update the parallel pool code to use
the new parpool function instead of matlabpool in anticipation of the
latter's eventual removal from MATLAB.
2016-03-18 Clarify help on the shape of the required input data arrays.
2016-03-29 Update the interpolation to extrapolate by removing NaNs
from the input data and setting the 'nearest' flag on
scatteredInterpolant.
==========================================================================

0001 function Mobj = interp_HYCOM2FVCOM(Mobj, hycom, 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_HYCOM2FVCOM(Mobj, hycom, start_date, varlist) 0006 % 0007 % DESCRIPTION: 0008 % FVCOM does not yet support spatially varying temperature and salinity 0009 % inputs as initial conditions. To avoid having to run a model for a 0010 % long time in order for temperature and salinity to settle within the 0011 % model from the atmospheric and boundary forcing, we can use a restart 0012 % file to cheat. For this, we need temperature and salinity 0013 % (potentially other variables too) interpolated onto the unstructured 0014 % grid. The interpolated data can then be written out with 0015 % write_FVCOM_restart.m. 0016 % 0017 % INPUT: 0018 % Mobj = MATLAB mesh structure which must contain: 0019 % - Mobj.lon, Mobj.lat and/or Mobj.lonc, Mobj.latc - node 0020 % or element centre coordinates (long/lat). Dependent on 0021 % the variable being interpolated (u and v need lonc, 0022 % latc whilst temperature and salinity need lat and lon). 0023 % - Mobj.siglayz - sigma layer depths for all model 0024 % nodes. 0025 % - Mobj.ts_times - Modified Julian Day times for the 0026 % model run. 0027 % hycom = Struct output by get_HYCOM_forcing. Must include fields: 0028 % - hycom.lon, hycom.lat - rectangular arrays. 0029 % - hycom.Depth - HYCOM depth levels. 0030 % The data fields (specified in varlist) must be shaped 0031 % (x, y, z, time). 0032 % start_date = Gregorian start date array (YYYY, MM, DD, hh, mm, ss). 0033 % varlist = cell array of fields to use from the HYCOM struct. 0034 % 0035 % OUTPUT: 0036 % Mobj.restart = struct whose field names are the variables which have 0037 % been interpolated (e.g. Mobj.restart.temperature for HYCOM 0038 % temperature). 0039 % 0040 % EXAMPLE USAGE 0041 % interp_HYCOM2FVCOM(Mobj, hycom, [2006, 01, 01, 00, 00, 00], ... 0042 % {'lon', 'lat', 'time', 'temperature', 'salinity'}) 0043 % 0044 % Author(s): 0045 % Pierre Cazenave (Plymouth Marine Laboratory) 0046 % 0047 % Revision history 0048 % 2013-09-10 First version based on interp_POLCOMS2FVCOM.m. 0049 % 2013-12-10 Fix the identification of the time index in the HYCOM data 0050 % (use hycom.time instead of Mobj.ts_times). Also ignore a field name of 0051 % 'MT' if supplied in varlist. 0052 % 2014-04-28 Fix bug when interpolating velocity data due to incorrectly 0053 % sized preallocated array. Make the way the data to be interpolated onto 0054 % the element centres (i.e. the velocity data), as opposed to the element 0055 % nodes, more understandable. This is done by having a single block which 0056 % deals with setting up the variables dependent on the number of points 0057 % onto which to interpolate. Also update the parallel pool code to use 0058 % the new parpool function instead of matlabpool in anticipation of the 0059 % latter's eventual removal from MATLAB. 0060 % 2016-03-18 Clarify help on the shape of the required input data arrays. 0061 % 2016-03-29 Update the interpolation to extrapolate by removing NaNs 0062 % from the input data and setting the 'nearest' flag on 0063 % scatteredInterpolant. 0064 % 0065 %========================================================================== 0066 0067 subname = 'interp_HYCOM2FVCOM'; 0068 0069 global ftbverbose; 0070 if ftbverbose 0071 fprintf('\nbegin : %s\n', subname) 0072 end 0073 0074 if nargin == 0 0075 error('Not enough input arguments. See HELP %s', subname) 0076 end 0077 0078 % Run jobs on multiple workers if we have that functionality. Not sure if 0079 % it's necessary, but check we have the Parallel Toolbox first. 0080 wasOpened = false; 0081 if license('test', 'Distrib_Computing_Toolbox') 0082 % We have the Parallel Computing Toolbox, so launch a bunch of workers. 0083 % New version for MATLAB 2014a (I think) onwards. 0084 if isempty(gcp('nocreate')) 0085 pool = parpool('local'); 0086 wasOpened = true; 0087 end 0088 end 0089 0090 % Given our input time (in start_date), find the nearest time index for 0091 % the HYCOM data. 0092 stime = greg2mjulian(start_date(1), start_date(2), ... 0093 start_date(3), start_date(4), ... 0094 start_date(5), start_date(6)); 0095 [~, tidx] = min(abs(hycom.time - stime)); 0096 0097 for vv = 1:length(varlist) 0098 0099 currvar = varlist{vv}; 0100 0101 switch currvar 0102 case {'lon', 'lat', 'longitude', 'latitude', 't_1', 'time', 'Depth', 'MT'} 0103 continue 0104 0105 otherwise 0106 0107 %-------------------------------------------------------------- 0108 % Interpolate the regularly gridded data onto the FVCOM grid 0109 % (vertical grid first). 0110 %-------------------------------------------------------------- 0111 0112 % Set up all the constants which are based on the model and 0113 % data grids. 0114 [~, fz] = size(Mobj.siglayz); 0115 if strcmpi(currvar, 'u') || strcmpi(currvar, 'v') 0116 plon = hycom.lon(:); 0117 plat = hycom.lat(:); 0118 flon = Mobj.lonc; 0119 flat = Mobj.latc; 0120 fn = numel(flon); 0121 fvtemp = nan(fn, fz); 0122 else 0123 plon = hycom.lon(:); 0124 plat = hycom.lat(:); 0125 flon = Mobj.lon; 0126 flat = Mobj.lat; 0127 fn = numel(flon); 0128 fvtemp = nan(fn, fz); 0129 end 0130 0131 if ftbverbose 0132 fprintf('%s : interpolate %s vertically... ', subname, currvar) 0133 end 0134 0135 [hx, hy, ~, ~] = size(hycom.(currvar).data); 0136 if iscolumn(hycom.Depth.data) 0137 hdepth = permute(repmat(-hycom.Depth.data, [1, hx, hy]), [2, 3, 1]); 0138 else 0139 hdepth = -hycom.Depth.data; 0140 end 0141 0142 % Make a land mask of the HYCOM domain (based on the surface 0143 % layer from the first time step). 0144 landmask = hycom.(currvar).data(:, :, 1, 1) > 1.26e29; 0145 0146 % Create a mask of all the layers which are deeper than the 0147 % water depth. 0148 zmask = hycom.(currvar).data(:, :, :, 1) > 1.26e29; 0149 0150 % Set invalid depths to NaN. 0151 hdepth(zmask) = nan; 0152 0153 ftbverbose = false; 0154 hyinterp.(currvar) = grid_vert_interp(Mobj, ... 0155 hycom.lon, hycom.lat, ... 0156 squeeze(hycom.(currvar).data(:, :, :, tidx)), ... 0157 hdepth, landmask, 'extrapolate', [flon, flat]); 0158 ftbverbose = true; 0159 0160 %-------------------------------------------------------------- 0161 % Now we have vertically interpolated data, we can interpolate 0162 % each sigma layer onto the FVCOM unstructured grid ready to 0163 % write out to NetCDF. We'll use the triangular interpolation 0164 % in MATLAB with the natural method (gives pretty good results, 0165 % at least qualitatively). 0166 %-------------------------------------------------------------- 0167 0168 if ftbverbose 0169 fprintf('horizontally... ') 0170 end 0171 0172 hytempz = hyinterp.(currvar); 0173 0174 tic 0175 parfor zi = 1:fz 0176 % Get the current depth layer's data and mask out the NaN 0177 % values. 0178 hytempzcurrent = hytempz(:, :, zi); 0179 0180 % Strip out NaNs so we can extrapolate with TriScatteredInterp. 0181 nanmask = ~isnan(hytempzcurrent); 0182 plonclean = plon(nanmask); 0183 platclean = plat(nanmask); 0184 hytempzclean = hytempzcurrent(nanmask); 0185 % Set up the interpolation object and interpolate the 0186 % current variable to the FVCOM unstructured grid. 0187 ft = scatteredInterpolant(plonclean, platclean, ... 0188 hytempzclean, 'natural', 'nearest'); 0189 fvtemp(:, zi) = ft(flon, flat); 0190 end 0191 0192 % Unfortunately, TriScatteredInterp won't extrapolate, 0193 % returning instead NaNs outside the original data's extents. 0194 % So, for each NaN position, find the nearest non-NaN value and 0195 % use that instead. The order in which the NaN-nodes are found 0196 % will determine the spatial pattern of the extrapolation. 0197 0198 % We can assume that all layers will have NaNs in the same 0199 % place (horizontally), so just use the surface layer (1) for 0200 % the identification of NaNs. Also store the finite values so 0201 % we can find the nearest real value to the current NaN node 0202 % and use its temperature and salinity values. 0203 fvidx = 1:fn; 0204 fvnanidx = fvidx(isnan(fvtemp(:, 1))); 0205 fvfinidx = fvidx(~isnan(fvtemp(:, 1))); 0206 0207 for ni = 1:length(fvnanidx) 0208 % Find the nearest non-nan data (temp, salinity, u or v) 0209 % value. 0210 xx = flon(fvnanidx(ni)); 0211 yy = flat(fvnanidx(ni)); 0212 [~, di] = min(sqrt((flon(fvfinidx) - xx).^2 + (flat(fvfinidx) - yy).^2)); 0213 fvtemp(fvnanidx(ni), :) = fvtemp(fvfinidx(di), :); 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.*.data have the rows and columns flipped, so 0262 % % (ci, ri) in pc.*.data 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.data(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.data(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.data(ci, ri, :, tidx)), squeeze(pc.depth.data(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.data(ci, ri, :, tidx)), squeeze(pc.depth.data(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.data)); 0313 % dy = mean(diff(pc.lat.data)); 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