


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.
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.
==========================================================================

0001 function Mobj = interp_POLCOMS2FVCOMv1(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 % start_date = Gregorian start date array (YYYY, MM, DD, hh, mm, ss). 0031 % varlist = cell array of fields to use from the HYCOM struct. 0032 % 0033 % OUTPUT: 0034 % Mobj.restart = struct whose field names are the variables which have 0035 % been interpolated (e.g. Mobj.restart.temperature for HYCOM 0036 % temperature). 0037 % 0038 % EXAMPLE USAGE 0039 % interp_HYCOM2FVCOM(Mobj, hycom, [2006, 01, 01, 00, 00, 00], ... 0040 % {'lon', 'lat', 'time', 'temperature', 'salinity'}) 0041 % 0042 % Author(s): 0043 % Pierre Cazenave (Plymouth Marine Laboratory) 0044 % 0045 % Revision history 0046 % 2013-09-10 First version based on interp_POLCOMS2FVCOM.m. 0047 % 2013-12-10 Fix the identification of the time index in the HYCOM data 0048 % (use hycom.time instead of Mobj.ts_times). Also ignore a field name of 0049 % 'MT' if supplied in varlist. 0050 % 2014-04-28 Fix bug when interpolating velocity data due to incorrectly 0051 % sized preallocated array. Make the way the data to be interpolated onto 0052 % the element centres (i.e. the velocity data), as opposed to the element 0053 % nodes, more understandable. This is done by having a single block which 0054 % deals with setting up the variables dependent on the number of points 0055 % onto which to interpolate. Also update the parallel pool code to use 0056 % the new parpool function instead of matlabpool in anticipation of the 0057 % latter's eventual removal from MATLAB. 0058 % 0059 %========================================================================== 0060 0061 subname = 'interp_POLCOMS2FVCOMv1'; 0062 0063 global ftbverbose; 0064 if ftbverbose 0065 fprintf('\nbegin : %s\n', subname) 0066 end 0067 0068 if nargin == 0 0069 error('Not enough input arguments. See HELP %s', subname) 0070 end 0071 0072 % Run jobs on multiple workers if we have that functionality. Not sure if 0073 % it's necessary, but check we have the Parallel Toolbox first. 0074 wasOpened = false; 0075 if license('test', 'Distrib_Computing_Toolbox') 0076 % We have the Parallel Computing Toolbox, so launch a bunch of workers. 0077 % New version for MATLAB 2014a (I think) onwards. 0078 if isempty(gcp('nocreate')) 0079 pool = parpool('local'); 0080 wasOpened = true; 0081 end 0082 end 0083 0084 % Given our input time (in start_date), find the nearest time index for 0085 % the HYCOM data. 0086 stime = greg2mjulian(start_date(1), start_date(2), ... 0087 start_date(3), start_date(4), ... 0088 start_date(5), start_date(6)); 0089 [~, tidx] = min(abs(hycom.time - stime)); 0090 0091 for vv = 1:length(varlist); 0092 0093 currvar = varlist{vv}; 0094 0095 switch currvar 0096 case {'lon', 'lat', 'longitude', 'latitude', 't_1', 'time', 'Depth', 'MT'} 0097 continue 0098 0099 otherwise 0100 0101 %-------------------------------------------------------------- 0102 % Interpolate the regularly gridded data onto the FVCOM grid 0103 % (vertical grid first). 0104 %-------------------------------------------------------------- 0105 0106 % Set up all the constants which are based on the model and 0107 % data grids. 0108 [~, fz] = size(Mobj.siglayz); 0109 if strcmpi(currvar, 'u') || strcmpi(currvar, 'v') 0110 plon = hycom.lon(:); 0111 plat = hycom.lat(:); 0112 flon = Mobj.lonc; 0113 flat = Mobj.latc; 0114 fn = numel(flon); 0115 fvtemp = nan(fn, fz); 0116 else 0117 plon = hycom.lon(:); 0118 plat = hycom.lat(:); 0119 flon = Mobj.lon; 0120 flat = Mobj.lat; 0121 fn = numel(flon); 0122 fvtemp = nan(fn, fz); 0123 end 0124 0125 if ftbverbose 0126 fprintf('%s : interpolate %s vertically... ', subname, currvar) 0127 end 0128 0129 [hx, hy, ~, ~] = size(hycom.(currvar)); 0130 hdepth = hycom.Depth; 0131 % Make a land mask of the HYCOM domain (based on the surface 0132 % layer from the first time step). 0133 landmask = isnan(hycom.(currvar)(:, :, 1, 1)); 0134 0135 0136 ftbverbose = false; 0137 % flip data upside down 0138 hdepth = flip(hdepth,3); 0139 % CHECK HERE IF 2D OR 3D VARIABLE 0140 if (ndims(hycom.(currvar))>2) 0141 datain= squeeze(hycom.(currvar)(:, :, :, tidx)); 0142 datain = flip(datain,3); 0143 % DO NOT INTERPOLATE 2D VARIABLES!!! 0144 0145 hyinterp.(currvar) = grid_vert_interp(Mobj, ... 0146 hycom.lon, hycom.lat, ... 0147 datain, ... 0148 hdepth, landmask, 'extrapolate', [flon, flat]); 0149 ftbverbose = true; 0150 %-------------------------------------------------------------- 0151 % Now we have vertically interpolated data, we can interpolate 0152 % each sigma layer onto the FVCOM unstructured grid ready to 0153 % write out to NetCDF. We'll use the triangular interpolation 0154 % in MATLAB with the natural method (gives pretty good results, 0155 % at least qualitatively). 0156 %-------------------------------------------------------------- 0157 0158 if ftbverbose 0159 fprintf('horizontally... ') 0160 end 0161 0162 hytempz = hyinterp.(currvar); 0163 0164 tic 0165 parfor zi = 1:fz 0166 % Get the current depth layer's data and mask out the NaN 0167 % values. 0168 hytempzcurrent = hytempz(:, :, zi); 0169 0170 % Strip out NaNs so we can extrapolate with TriScatteredInterp. 0171 nanmask = ~isnan(hytempzcurrent); 0172 plonclean = plon(nanmask); 0173 platclean = plat(nanmask); 0174 hytempzclean = hytempzcurrent(nanmask); 0175 % Set up the interpolation object and interpolate the 0176 % current variable to the FVCOM unstructured grid. 0177 ft = scatteredInterpolant(plonclean, platclean, ... 0178 hytempzclean, 'natural', 'nearest'); 0179 fvtemp(:, zi) = ft(flon, flat); 0180 end 0181 else 0182 ft = scatteredInterpolant(plon, plat, reshape(hycom.(currvar), [], 1), 'nearest','linear'); 0183 fvtemp = ft(flon, flat); 0184 ftbverbose = true; 0185 end 0186 0187 % Unfortunately, TriScatteredInterp won't extrapolate, 0188 % returning instead NaNs outside the original data's extents. 0189 % So, for each NaN position, find the nearest non-NaN value and 0190 % use that instead. The order in which the NaN-nodes are found 0191 % will determine the spatial pattern of the extrapolation. 0192 0193 % We can assume that all layers will have NaNs in the same 0194 % place (horizontally), so just use the surface layer (1) for 0195 % the identification of NaNs. Also store the finite values so 0196 % we can find the nearest real value to the current NaN node 0197 % and use its temperature and salinity values. 0198 fvidx = 1:fn; 0199 fvnanidx = fvidx(isnan(fvtemp(:, 1))); 0200 fvfinidx = fvidx(~isnan(fvtemp(:, 1))); 0201 0202 for ni = 1:length(fvnanidx) 0203 % Find the nearest non-nan data (temp, salinity, u or v) 0204 % value. 0205 xx = flon(fvnanidx(ni)); 0206 yy = flat(fvnanidx(ni)); 0207 [~, di] = min(sqrt((flon(fvfinidx) - xx).^2 + (flat(fvfinidx) - yy).^2)); 0208 fvtemp(fvnanidx(ni), :) = fvtemp(fvfinidx(di), :); 0209 end 0210 0211 clear plon plat flon flat ptempz 0212 0213 Mobj.restart.(currvar) = fvtemp; 0214 0215 te = toc; 0216 0217 if ftbverbose 0218 fprintf('done. (elapsed time = %.2f seconds)\n', te) 0219 end 0220 0221 end 0222 end 0223 0224 % Close the MATLAB pool if we opened it. 0225 if wasOpened 0226 pool.delete 0227 end 0228 0229 if ftbverbose 0230 fprintf('end : %s\n', subname) 0231 end 0232 0233 %% Debugging figure 0234 0235 % close all 0236 % 0237 % ri = 85; % column index 0238 % ci = 95; % row index 0239 % 0240 % [~, idx] = min(sqrt((Mobj.lon - lon(ri, ci)).^2 + (Mobj.lat - lat(ri, ci)).^2)); 0241 % 0242 % % Vertical profiles 0243 % figure 0244 % clf 0245 % 0246 % % The top row shows the temperature/salinity values as plotted against 0247 % % index (i.e. position in the array). I had thought POLCOMS stored its data 0248 % % from seabed to sea surface, but looking at the NetCDF files in ncview, 0249 % % it's clear that the data are in fact stored surface to seabed (like 0250 % % FVCOM). As such, the two plots in the top row should be upside down (i.e. 0251 % % surface values at the bottom of the plot). The two lower rows should have 0252 % % three lines which all match: the raw POLCOMS data, the POLCOMS data for 0253 % % the current time step (i.e. that in 'temperature' and 'salinity') and the 0254 % % interpolated FVCOM data against depth. 0255 % % 0256 % % Also worth noting, the pc.* have the rows and columns flipped, so 0257 % % (ci, ri) in pc.* and (ri, ci) in 'temperature', 'salinity' and 0258 % % 'depth'. Needless to say, the two lines in the lower plots should 0259 % % overlap. 0260 % 0261 % subplot(2,2,1) 0262 % plot(squeeze(pc.ETWD(ci, ri, :, tidx)), 1:size(depth, 3), 'rx:') 0263 % xlabel('Temperature (^{\circ}C)') 0264 % ylabel('Array index') 0265 % title('Array Temperature') 0266 % 0267 % subplot(2,2,2) 0268 % plot(squeeze(pc.x1XD(ci, ri, :, tidx)), 1:size(depth, 3), 'rx:') 0269 % xlabel('Salinity') 0270 % ylabel('Array index') 0271 % title('Array Salinity') 0272 % 0273 % subplot(2,2,3) 0274 % % Although POLCOMS stores its temperature values from seabed to surface, 0275 % % the depths are stored surface to seabed. Nice. Flip the 0276 % % temperature/salinity data accordingly. The lines here should match one 0277 % % another. 0278 % plot(squeeze(pc.ETWD(ci, ri, :, tidx)), squeeze(pc.depth(ci, ri, :, tidx)), 'rx-') 0279 % hold on 0280 % plot(squeeze(temperature(ri, ci, :)), squeeze(depth(ri, ci, :)), '.:') 0281 % % Add the equivalent interpolated FVCOM data point 0282 % plot(fvtemp(idx, :), Mobj.siglayz(idx, :), 'g.-') 0283 % xlabel('Temperature (^{\circ}C)') 0284 % ylabel('Depth (m)') 0285 % title('Depth Temperature') 0286 % legend('pc', 'temp', 'fvcom', 'location', 'north') 0287 % legend('boxoff') 0288 % 0289 % subplot(2,2,4) 0290 % % Although POLCOMS stores its temperature values from seabed to surface, 0291 % % the depths are stored surface to seabed. Nice. Flip the 0292 % % temperature/salinity data accordingly. The lines here should match one 0293 % % another. 0294 % plot(squeeze(pc.x1XD(ci, ri, :, tidx)), squeeze(pc.depth(ci, ri, :, tidx)), 'rx-') 0295 % hold on 0296 % plot(squeeze(salinity(ri, ci, :)), squeeze(depth(ri, ci, :)), '.:') 0297 % % Add the equivalent interpolated FVCOM data point 0298 % plot(fvsalt(idx, :), Mobj.siglayz(idx, :), 'g.-') 0299 % xlabel('Salinity') 0300 % ylabel('Depth (m)') 0301 % title('Depth Salinity') 0302 % legend('pc', 'salt', 'fvcom', 'location', 'north') 0303 % legend('boxoff') 0304 % 0305 % %% Plot the sample location 0306 % figure 0307 % dx = mean(diff(pc.lon)); 0308 % dy = mean(diff(pc.lat)); 0309 % z = depth(:, :, end); % water depth (bottom layer depth) 0310 % z(mask) = 0; % clear out nonsense values 0311 % pcolor(lon - (dx / 2), lat - (dy / 2), z) 0312 % shading flat 0313 % axis('equal', 'tight') 0314 % daspect([1.5, 1, 1]) 0315 % colorbar 0316 % caxis([-150, 0]) 0317 % hold on 0318 % plot(lon(ri, ci), lat(ri, ci), 'ko', 'MarkerFaceColor', 'w') 0319