


Read an arbitrary 4D field (x, y, z, time) from the coarse struct and
interpolate onto the open boundaries in Mobj.
function Mobj = interp_coarse_to_obc(Mobj, coarse, varlist)
DESCRIPTION:
Interpolate 4D field values onto the FVCOM open boundaries at all
sigma levels.
INPUT:
Mobj = MATLAB mesh structure which must contain:
- Mobj.siglayz - sigma layer depths for all model nodes.
- Mobj.siglayzc - sigma layer depths for all model elements.
- Mobj.lon, Mobj.lat - node coordinates (lat/long).
- Mobj.lonc, Mobj.latc - element coordinates (lat/long).
- Mobj.read_obc_nodes - cell array of open boundary nodes.
- Mobj.read_obc_elems - cell array of open boundary
elements (only if using velocities - use find_nested_region
to get these indices).
- Mobj.h - water depths at nodes.
- Mobj.tri - triangulation table for the grid (nv in FVCOM
terms).
- Mobj.nObcNodes - number of nodes in each open boundary.
coarse = Struct with coarse data covering the model domain. Unless
varlist is specified (see below), all 4D fields will be
interpolated onto the open boundaries (1-3D data will be
ignored).
varlist = [optional] cell array of variable (field) names to use from
coarse. If omitted, all fields are processed.
OUTPUT:
Mobj = MATLAB structure with new fields whose names match those given
in coarse. The fields have sizes (sum(Mobj.nObcNodes), sigma, time).
The time dimension is determined based on the number of time steps in
coarse. The ts_time variable is just the input file times in Modified
Julian Day.
EXAMPLE USAGE
coarse = get_HYCOM_forcing(Mobj, [51500, 51531]);
Mobj = interp_coarse_to_obc(Mobj, coarse, {'u', 'v', 'temperature', 'salinity'})
Author(s):
Pierre Cazenave (Plymouth Marine Laboratory)
Revision history
2013-09-03 First version based on get_POLCOMS_tsobc.m.
2014-04-28 Update interp1 function to use pchip instead of csap as the
latter will be removed in a future version of MATLAB and the
innumerable warnings were doing my nut in. I checked the output using
the new interp1 call and it's identical to the old version. Also
update the parallel toolbox stuff for the same reason (future
removal).
2015-05-21 Remove the old parallel processing bits and replace with
the current versions.
2016-03-15 Add fallback interpolation to inverse distance weighted if
the triangular interpolation fails (which can happen if the points
supplied are all in a line, for example).
2017-01-27 Subset the coarse data (HYCOM, CMEMS etc.). This yields a
significant performance improvement (both in memory usage and time).
2017-02-16 Further performance improvement by only using the coarse
data in the vicinity of the open boundary nodes.
2017-10-12 Fix bug in indexing the open boundary positions which may
have caused the interpolation to fail as the identified positions
would be too far from the open boundary nodes.
==========================================================================

0001 function Mobj = interp_coarse_to_obc(Mobj, coarse, varlist) 0002 % Read an arbitrary 4D field (x, y, z, time) from the coarse struct and 0003 % interpolate onto the open boundaries in Mobj. 0004 % 0005 % function Mobj = interp_coarse_to_obc(Mobj, coarse, varlist) 0006 % 0007 % DESCRIPTION: 0008 % Interpolate 4D field values onto the FVCOM open boundaries at all 0009 % sigma levels. 0010 % 0011 % INPUT: 0012 % Mobj = MATLAB mesh structure which must contain: 0013 % - Mobj.siglayz - sigma layer depths for all model nodes. 0014 % - Mobj.siglayzc - sigma layer depths for all model elements. 0015 % - Mobj.lon, Mobj.lat - node coordinates (lat/long). 0016 % - Mobj.lonc, Mobj.latc - element coordinates (lat/long). 0017 % - Mobj.read_obc_nodes - cell array of open boundary nodes. 0018 % - Mobj.read_obc_elems - cell array of open boundary 0019 % elements (only if using velocities - use find_nested_region 0020 % to get these indices). 0021 % - Mobj.h - water depths at nodes. 0022 % - Mobj.tri - triangulation table for the grid (nv in FVCOM 0023 % terms). 0024 % - Mobj.nObcNodes - number of nodes in each open boundary. 0025 % coarse = Struct with coarse data covering the model domain. Unless 0026 % varlist is specified (see below), all 4D fields will be 0027 % interpolated onto the open boundaries (1-3D data will be 0028 % ignored). 0029 % varlist = [optional] cell array of variable (field) names to use from 0030 % coarse. If omitted, all fields are processed. 0031 % 0032 % OUTPUT: 0033 % Mobj = MATLAB structure with new fields whose names match those given 0034 % in coarse. The fields have sizes (sum(Mobj.nObcNodes), sigma, time). 0035 % The time dimension is determined based on the number of time steps in 0036 % coarse. The ts_time variable is just the input file times in Modified 0037 % Julian Day. 0038 % 0039 % EXAMPLE USAGE 0040 % coarse = get_HYCOM_forcing(Mobj, [51500, 51531]); 0041 % Mobj = interp_coarse_to_obc(Mobj, coarse, {'u', 'v', 'temperature', 'salinity'}) 0042 % 0043 % Author(s): 0044 % Pierre Cazenave (Plymouth Marine Laboratory) 0045 % 0046 % Revision history 0047 % 2013-09-03 First version based on get_POLCOMS_tsobc.m. 0048 % 2014-04-28 Update interp1 function to use pchip instead of csap as the 0049 % latter will be removed in a future version of MATLAB and the 0050 % innumerable warnings were doing my nut in. I checked the output using 0051 % the new interp1 call and it's identical to the old version. Also 0052 % update the parallel toolbox stuff for the same reason (future 0053 % removal). 0054 % 2015-05-21 Remove the old parallel processing bits and replace with 0055 % the current versions. 0056 % 2016-03-15 Add fallback interpolation to inverse distance weighted if 0057 % the triangular interpolation fails (which can happen if the points 0058 % supplied are all in a line, for example). 0059 % 2017-01-27 Subset the coarse data (HYCOM, CMEMS etc.). This yields a 0060 % significant performance improvement (both in memory usage and time). 0061 % 2017-02-16 Further performance improvement by only using the coarse 0062 % data in the vicinity of the open boundary nodes. 0063 % 2017-10-12 Fix bug in indexing the open boundary positions which may 0064 % have caused the interpolation to fail as the identified positions 0065 % would be too far from the open boundary nodes. 0066 % 0067 %========================================================================== 0068 0069 [~, subname] = fileparts(mfilename('fullpath')); 0070 0071 global ftbverbose 0072 if ftbverbose 0073 fprintf('\nbegin : %s\n', subname) 0074 end 0075 0076 if nargin == 2 0077 fields = fieldnames(coarse); 0078 else 0079 % We always want Depth because we need it to interpolate the vertical 0080 % component. 0081 fields = unique(['Depth', varlist], 'stable'); 0082 end 0083 0084 assert(isfield(coarse, 'Depth'), 'Require a depth vector to interpolate vertically.') 0085 0086 % Find the first 4D array and use it to get the number of vertical levels 0087 % and time steps. 0088 for ff = 1:length(fields) 0089 if isfield(coarse.(fields{ff}), 'data') && ndims(coarse.(fields{ff}).data) > 3 0090 [nx, ny, nz, nt] = size(coarse.(fields{ff}).data); 0091 break 0092 end 0093 end 0094 0095 % Use the existing rectangular arrays for the nearest point lookup. 0096 [lon, lat] = deal(coarse.lon, coarse.lat); 0097 0098 % Find the indices of all the coarse data which cover the open boundary 0099 % nodes given by some distance. This should massively increase performance 0100 % by chopping out the vast majority of the data. Find the 6 closest ones 0101 % for each node/element. 0102 0103 obc_coarse_idx = []; 0104 if isfield(Mobj, 'read_obc_elems') 0105 oElems = [Mobj.read_obc_elems{:}]; 0106 else 0107 oElems = []; 0108 end 0109 oNodes = [Mobj.read_obc_nodes{:}]; 0110 for obc_i = 1:length(oNodes) 0111 % Find the 20 closest ones. 0112 [~, h_idx] = sort(sqrt((lon(:) - Mobj.lon(oNodes(obc_i))).^2 + ... 0113 (lat(:) - Mobj.lat(oNodes(obc_i))).^2)); 0114 obc_coarse_idx = [obc_coarse_idx; h_idx(1:20)]; 0115 end 0116 for obc_i = 1:length(oElems) 0117 % Find the 20 closest ones. 0118 [~, h_idx] = sort(sqrt((lon(:) - Mobj.lonc(oElems(obc_i))).^2 + ... 0119 (lat(:) - Mobj.latc(oElems(obc_i))).^2)); 0120 obc_coarse_idx = [obc_coarse_idx; h_idx(1:20)]; 0121 end 0122 obc_coarse_idx = unique(obc_coarse_idx, 'stable'); 0123 [xidx, yidx] = ind2sub(size(lon), obc_coarse_idx); 0124 0125 % Calculate resolution of parent coarse model 0126 hdx = max([max(diff(lon(:,1))),max(diff(lat(1,:)))]); 0127 0128 % Number of sigma layers. 0129 fz = size(Mobj.siglayz, 2); 0130 0131 % Make a 3D array of the coarse depths and mask where we don't have data. 0132 % This can then be used in the interpolation instead of trying to deal with 0133 % this on the fly. 0134 hdepth = permute(repmat(coarse.Depth.data, [1, nx, ny]), [2, 3, 1]); 0135 mask = coarse.(fields{ff}).data(:, :, :, 1) > 1.26e29; 0136 % Used to use the landmask to check whether to extrapolate data onto land. 0137 % Realised this is no longer necessary, so comment this out for future 0138 % removal. 0139 % hdepth(mask) = nan; 0140 % landmask = coarse.(fields{ff}).data(:, :, 1, 1) > 1.26e29; 0141 0142 % Ignore depth layers in the coarse data which are below the maximum depth 0143 % of the open boundary nodes so we don't waste time trying to interpolate 0144 % those data. Add a buffer of one layer in case I'm an idiot. 0145 if ~isfield(Mobj, 'hc') 0146 Mobj.hc = nodes2elems(Mobj.h, Mobj); 0147 end 0148 max_obc_depth = max([max(Mobj.h([Mobj.read_obc_nodes{:}])), max(Mobj.hc([Mobj.read_obc_elems{:}]))]); 0149 [~, zidx] = min(abs(coarse.Depth.data - max_obc_depth)); 0150 zidx = min([zidx + 1, length(coarse.Depth.data)]); 0151 0152 if ftbverbose 0153 tic 0154 end 0155 % Only do warnings about removing values outside some ranges once per 0156 % variable. 0157 warned = true(4, 1); 0158 for v = 1:length(fields) 0159 0160 if ~(isfield(coarse.(fields{v}), 'data') && ndims(coarse.(fields{v}).data) > 3) 0161 continue 0162 end 0163 0164 % Supply FVCOM grid positions depending on whether we're working on the 0165 % velocity data (on elements) or anything else (on nodes). 0166 if any(strcmpi(fields{v}, {'u', 'v'})) 0167 fvlon = Mobj.lonc(oElems); 0168 fvlat = Mobj.latc(oElems); 0169 sigma = Mobj.siglayzc; 0170 grid_ids = oElems; 0171 % Number of boundary elements to process 0172 nf = length(oElems); 0173 0174 if ftbverbose 0175 fprintf('Variable %s on elements (%d positions)\n', fields{v}, length(fvlon)) 0176 end 0177 else 0178 % Make sure the nodes are listed in the same way as in 0179 % casename_obc.dat. 0180 fvlon = Mobj.lon(oNodes); 0181 fvlat = Mobj.lat(oNodes); 0182 sigma = Mobj.siglayz; 0183 grid_ids = oNodes; 0184 % Number of boundary nodes 0185 nf = length(oNodes); 0186 0187 if ftbverbose 0188 fprintf('Variable %s on nodes (%d positions)\n', fields{v}, length(fvlon)) 0189 end 0190 end 0191 0192 fvtemp = nan(nf, fz, nt); % FVCOM interpolated data 0193 0194 for t = 1:nt 0195 0196 if ftbverbose 0197 fprintf('%s : %i of %i %s timesteps... \n', subname, t, nt, fields{v}) 0198 end 0199 % Get the current 3D array of coarse results. Only load the data we 0200 % need in both the vertical and horizontal (i.e. clip to the 0201 % locations which cover the open boundary positions only). 0202 pctemp3 = nan(length(xidx), zidx, 1); 0203 for oidx = 1:length(xidx) 0204 pctemp3(oidx, :) = coarse.(fields{v}).data(xidx(oidx), yidx(oidx), 1:zidx, t); 0205 end 0206 % pcolor(lon(xidx(1):xidx(2), yidx(1):yidx(2)), lat(xidx(1):xidx(2), yidx(1):yidx(2)), pctemp3(:, :, 1)); shading flat; axis('equal', 'tight'); colorbar 0207 % hold on 0208 % plot(fvlon, fvlat, 'ro') 0209 % Remove NaNs since we've assumed later on that everything above 0210 % some value is NaN instead (currently 1.26e29). Make the value 0211 % that we replace the NaNs with larger since we check for larger 0212 % values when deleting invalid values. 0213 pctemp3(isnan(pctemp3)) = 1.27e29; 0214 0215 % Preallocate the intermediate results array. 0216 itempz = nan(nf, nz); 0217 0218 for j = 1:nz 0219 if j == zidx + 1 0220 if ftbverbose 0221 fprintf(' skipping remaining coarse data layers (%d) as their depths are below the FVCOM open boundary depths\n', nz - zidx) 0222 end 0223 continue 0224 elseif j > zidx 0225 continue 0226 else 0227 if ftbverbose 0228 fprintf(' coarse data layer %d of %d\n', j, nz) 0229 end 0230 end 0231 0232 % Now extract the relevant layer from the 3D subsets. 0233 tpctemp2 = pctemp3(:, j); 0234 0235 % Create and apply a mask to remove values outside the domain. 0236 % This inevitably flattens the arrays, but it shouldn't be a 0237 % problem since we'll be searching for the closest values in 0238 % such a manner as is appropriate for an unstructured grid 0239 % (i.e. we're assuming the coarse data is irregularly sampled 0240 % and interpolating with a triangulation). 0241 mask = tpctemp2 > 1.26e29; 0242 0243 % We need to do some more checks for the data which has been 0244 % saved via Python. This is a sort of bounds check to eliminate 0245 % unrealistic data. Warn if we actually delete anything. 0246 switch fields{v} 0247 case 'salinity' 0248 mask_alt = tpctemp2 < 0; 0249 if min(mask_alt(:)) == 1 && warned(1) 0250 warned(1) = false; 0251 warning('Removing negative salinities from the coarse data.') 0252 end 0253 case 'temperature' 0254 mask_alt = tpctemp2 < -20; 0255 if min(mask_alt(:)) == 1 && warned(2) 0256 warned(2) = false; 0257 warning('Removing temperature values below -20 celsius from the coarse data.') 0258 end 0259 case 'ssh' 0260 mask_alt = tpctemp2 < -20; 0261 if min(mask_alt(:)) == 1 && warned(3) 0262 warned(3) = false; 0263 warning('Removing sea surface height values below -20m from the coarse data.') 0264 end 0265 case {'u', 'v'} 0266 mask_alt = tpctemp2 > 100 | tpctemp2 < -100 ; 0267 if min(mask_alt(:)) == 1 && warned(4) 0268 warned(3) = false; 0269 warning('Removing non-tidal velocities above/below +/-100m/s from the coarse data.') 0270 end 0271 otherwise 0272 % Some other variable we won't mask. 0273 mask_alt = true(size(tpctemp2)); 0274 end 0275 mask = logical(~(~mask .* ~mask_alt)); 0276 clear mask_alt 0277 tpctemp2(mask) = []; 0278 0279 % Also apply the masks to the position arrays so we can't even 0280 % find positions outside the domain, effectively meaning if a 0281 % value is outside the domain, the nearest value to the current 0282 % boundary position will be used. 0283 % Extract the longitude and latitude data for the open boundary 0284 % region. 0285 tlon = lon(obc_coarse_idx); 0286 tlat = lat(obc_coarse_idx); 0287 tlon(mask) = []; 0288 tlat(mask) = []; 0289 0290 % If the coarse depths are now below the maximum depth in the 0291 % model domain, we'll have no valid data, so skip the 0292 % interpolation and just leave the NaNs in the itempz array. 0293 if isempty(tlon) && isempty(tlat) 0294 continue 0295 end 0296 0297 % Preallocate the intermediate results array. 0298 itempobc = nan(nf, 1); 0299 0300 % Speed up the tightest loop with a parallelized loop. 0301 parfor i = 1:nf 0302 fx = fvlon(i); 0303 fy = fvlat(i); 0304 0305 [dist, ii] = sort(sqrt((tlon - fx).^2 + (tlat - fy).^2)); 0306 % Get the n nearest positions from coarse data (more? 0307 % fewer?). 0308 np = 16; 0309 if length(ii) < np 0310 % Reset np to the number of points we actually have. 0311 np = length(ii); 0312 end 0313 ixy = ii(1:np); 0314 0315 % If the minimum distance away is greater than five times 0316 % the coarse grid resolution, skip this point at this 0317 % vertical level. If this is not done, we can be grabbing 0318 % data from a significant distance away 0319 if min(dist) > 5 * hdx 0320 continue 0321 end 0322 0323 % Get the variables into static variables for the 0324 % parallelisation. 0325 plon = double(tlon(ixy)); 0326 plat = double(tlat(ixy)); 0327 ptemp = tpctemp2(ixy); 0328 0329 % Use a triangulation to do the horizontal interpolation if 0330 % we have enough points, otherwise take the mean of the two 0331 % values. If the triangulation fails (when the points are 0332 % in a line), do an inverse distance weighed interpolation 0333 % instead 0334 if length(plon) < 3 0335 itempobc(i) = mean(ptemp); 0336 else 0337 try 0338 tritemp = scatteredInterpolant(plon, plat, ptemp, 'natural', 'none'); 0339 itempobc(i) = tritemp(fx, fy); 0340 catch err 0341 if strcmp(err.identifier, 'MATLAB:subsassignnumelmismatch') 0342 warning(['Scatter points failed the', ... 0343 ' triangular interpolation. Falling', ... 0344 ' back to inverse distance weighted', ... 0345 ' interpolation.']) 0346 % Use the inverse distance weighted mean of the 0347 % values for the interpolated value (the values 0348 % are in a straight line and can't be 0349 % interpolated with a triangular 0350 % interpolation). 0351 w = 1 ./ sqrt((plon - fx).^2 + (plat - fy).^2); 0352 w = w / max(w); 0353 itempobc(i) = sum(ptemp .* w) ./ sum(w); 0354 else 0355 error(err.message) 0356 end 0357 end 0358 end 0359 0360 % Check we're interpolating at the surface correctly. 0361 % figure(100) 0362 % if i == 1 0363 % clf 0364 % pcolor(coarse.lon - (hdx/2), coarse.lat - (hdx/2), coarse.(fields{v}).data(:, :, j, t)); colorbar; hold on 0365 % shading flat 0366 % hold on 0367 % end 0368 % plot(tlon(ixy), tlat(ixy), 'ko') 0369 % plot(fx, fy, 'rx') 0370 % plot(tlon(ixy(1)), tlat(ixy(1)), 'gs') 0371 % scatter(tlon(ixy(1)), tlat(ixy(1)), 20, ptemp(1), 'filled', 'markeredgecolor', 'r') 0372 % scatter(fx, fy, 20, ptemp(1), 'filled', 'markeredgecolor', 'g') 0373 % disp(ptemp(1)) 0374 % caxis([-0.1, 0.1]) 0375 % axis('equal', 'tight') 0376 % axis([fx - 0.5, fx + 0.5, fy - 0.5, fy + 0.5]) 0377 % legend('Data', 'Candidate valid', 'FVCOM node', 'Nearest valid', 'Interpolated', 'Location', 'NorthOutside', 'Orientation', 'Horizontal') 0378 % legend('BoxOff') 0379 % pause(0.1) 0380 0381 if isnan(itempobc(i)) 0382 % Use the surface layer as the canonical land mask and 0383 % check that the issue here is not just that the 0384 % current open boundary position is shallower than this 0385 % layer's depth. In the case where we're at the 0386 % surface, we always want a value, so use the closest 0387 % value, otherwise we can skip this data and leave it 0388 % as NaN. The vertical interpolation will strip out the 0389 % NaN depths so we shouldn't have any problems from 0390 % that perspective. 0391 0392 % I used to check we were on land, but really, 0393 % itempobc(i) will only equal NaN if we're on land for 0394 % this layer. This is only a problem when we're at the 0395 % surface as we always need at least one value to do 0396 % the vertical interpolation. So, check if we're at the 0397 % surface and if so, grab the nearest point. Otherwise, 0398 % leave the NaN in place. 0399 %[~, jj] = min(sqrt((lon(:) - fx).^2 + (lat(:) - fy).^2)); 0400 %[ir, ic] = ind2sub(size(lon), jj); 0401 %if landmask(ir, ic) == 1 && j == 1 0402 if j == 1 0403 %fprintf('Sea surface or on land (j = %i, lon: %.5f, lat: %.5f)\n', j, lon(ir, ic), lat(ir, ic)) 0404 itempobc(i) = tpctemp2(ii(1)); 0405 0406 %clf 0407 %pcolor(coarse.lon, coarse.lat, landmask * 1); colorbar; hold on 0408 %plot(lon(ir, ic), lat(ir, ic), 'ko') 0409 %plot(fx, fy, 'rx') 0410 %plot(tlon(ii(1)), tlat(ii(1)), 'gs') 0411 %axis('square') 0412 %axis([fx - 1.5, fx + 1.5, fy - 1.5, fy + 1.5]) 0413 %legend('Land mask', 'Mask test', 'FVCOM node', 'Nearest valid', 'Location', 'NorthOutside', 'Orientation', 'Horizontal') 0414 %legend('BoxOff') 0415 end 0416 end 0417 end 0418 0419 % Put the results in the intermediate array. 0420 itempz(:, j) = itempobc; 0421 0422 end 0423 0424 % If you want to check the interpolation has worked properly: 0425 % xx = repmat(fvlon, [1, nz]); 0426 % yy = repmat(fvlat, [1, nz]); 0427 % zz = repmat(-coarse.Depth.data, [1, nf])'; 0428 % dd = nanmax(hdepth, [], 3); 0429 % cc = itempz; 0430 % mm = ~isnan(cc); 0431 % ffx = repmat(fvlon, [1, fz]); 0432 % ffy = repmat(fvlat, [1, fz]); 0433 % ff = Mobj.siglayz(oNodes, :); 0434 % figure(10) 0435 % clf 0436 % scatter3(xx(mm), yy(mm), zz(mm), 40, cc(mm), 'filled') 0437 % hold on 0438 % scatter3(coarse.lon(:), coarse.lat(:), -dd(:), 40, 'c.') 0439 % scatter3(ffx(:), ffy(:), ff(:), 'k.') 0440 % colorbar 0441 % zlim([-300, 0]) 0442 0443 % Now we've interpolated in space, we can interpolate the z-values 0444 % to the sigma depths. 0445 0446 fvtempz = nan(nf, fz); 0447 for pp = 1:nf 0448 % Get the FVCOM depths at this position 0449 tfz = sigma(grid_ids(pp), :); % grid_ids is oNodes or oElems 0450 0451 % The coarse data is unusual in that the depths are fixed and 0452 % data are only saved at the depths shallower than the 0453 % bathymetry. As such, we get NaN values below the water depth 0454 % and we need to filter those out here. 0455 0456 % Find the coarse depths which cover the modelled depth range. 0457 tpz = -coarse.Depth.data; 0458 % Mask the coarse depths with the data array at this position. 0459 mm = isnan(itempz(pp, :)); 0460 tpz(mm) = []; 0461 0462 % If coarse has a single value, just repeat it across all depth 0463 % values. 0464 if length(tpz) == 1 0465 fvtempz(pp, :) = repmat(itempz(pp, ~mm), [1, length(tfz)]); 0466 else 0467 % To ensure we get the full vertical expression of the 0468 % vertical profiles, we need to normalise the coarse and 0469 % FVCOM depths to the same range. This is because in 0470 % instances where FVCOM depths are shallower (e.g. in 0471 % coastal regions), if we don't normalise the depths, we 0472 % end up truncating the vertical profile. This approach 0473 % ensures we always use the full vertical profile, but 0474 % we're potentially squeezing it into a smaller depth. 0475 A = max(tpz); 0476 B = min(tpz); 0477 C = max(tfz); 0478 D = min(tfz); 0479 norm_tpz = (((D - C) * (tpz - A)) / (B - A)) + C; 0480 0481 % Get the temperature and salinity values for this position 0482 % and interpolate down the water column (from coarse to 0483 % FVCOM). 0484 if any(~isnan(norm_tpz)) 0485 fvtempz(pp, :) = interp1(norm_tpz, itempz(pp, ~mm), tfz, 'pchip', 'extrap'); 0486 0487 %figure(800); 0488 %clf 0489 %plot(itempz(pp, ~mm), tpz, 'r-o') 0490 %hold on 0491 %plot(fvtempz(pp, :), tfz, 'k-x') 0492 %legend('coarse', 'FVCOM') 0493 else 0494 warning('Should never see this... ') % because we test for NaNs when fetching the values. 0495 warning('FVCOM boundary position at %f, %f is outside the coarse domain. Skipping.', fvlon(pp), fvlat(pp)) 0496 continue 0497 end 0498 end 0499 end 0500 0501 % Find and remove NaNs. 0502 parfor pp = 1:fz 0503 test = fvtempz(:, pp); 0504 if any(isnan(test)) 0505 igood = ~isnan(test); 0506 ftri = scatteredInterpolant(fvlon(igood), fvlat(igood), test(igood), 'nearest', 'nearest'); 0507 fvtempz(:, pp) = ftri(fvlon,fvlat); 0508 end 0509 end 0510 0511 % The horizontally- and vertically-interpolated values in the final 0512 % FVCOM results array. 0513 fvtemp(:, :, t) = fvtempz; 0514 0515 if ftbverbose 0516 fprintf('done.\n') 0517 end 0518 end 0519 % Dump the data into a temporary structure. 0520 fvcom.(fields{v}) = fvtemp; 0521 end % loop through variables 0522 if ftbverbose 0523 toc 0524 end 0525 0526 fvfields = fieldnames(fvcom); 0527 for s = 1:length(fvfields) 0528 switch fvfields{s} 0529 case 'temperature' 0530 Mobj.temperature = fvcom.temperature; 0531 case 'salinity' 0532 Mobj.salt = fvcom.salinity; 0533 case 'u' 0534 Mobj.u = fvcom.u; 0535 case 'v' 0536 Mobj.v = fvcom.v; 0537 case 'density' 0538 Mobj.rho1 = fvcom.density; 0539 otherwise 0540 Mobj.(fvfields{s}) = fvcom.(fvfields{s}); 0541 warning('Unrecognised variable %s. Leaving name as originally given', fvfields{s}) 0542 end 0543 end 0544 0545 if isfield(coarse, 'time') 0546 Mobj.ts_times = coarse.time; 0547 end 0548 0549 if ftbverbose 0550 fprintf('end : %s\n', subname) 0551 end 0552 0553 0554 %% 0555 % Plot a vertical profile for a boundary node (for my Irish Sea case, this 0556 % is one of the ones along the Celtic Sea boundary). Also plot the 0557 % distribution of interpolated values over the coarse data. Add the location 0558 % of the vertical profile (both FVCOM and coarse) to the plot. 0559 0560 % nn = 128; % open boundary index 0561 % tt = 1; % time index 0562 % fvz = 1; % fvcom depth index (currently 1-20) 0563 % hyz = 1; % coarse depth index (1-33) 0564 % 0565 % % Find the coarse seabed indices 0566 % % [~, hyz] = nanmax(hdepth, [], 3); 0567 % 0568 % % Get the corresponding indices for the coarse data 0569 % [~, idx] = min(sqrt((lon(:) - fvlon(nn)).^2 + (lat(:) - fvlat(nn)).^2)); 0570 % [xidx, yidx] = ind2sub(size(lon), idx); 0571 % 0572 % zidx = isfinite(hdepth(xidx, yidx, :)); 0573 % hz = 1:nz; 0574 % 0575 % % close all 0576 % 0577 % figure(100) 0578 % clf 0579 % subplot(2,2,1) 0580 % plot(Mobj.temperature(nn, :, tt), Mobj.siglayz(oNodes(nn), :), 'x-') 0581 % xlabel('Temperature (^{\circ}C)') 0582 % ylabel('Depth (m)') 0583 % title('FVCOM') 0584 % 0585 % subplot(2,2,2) 0586 % plot(squeeze(coarse.temperature.data(xidx, yidx, zidx, tt)), squeeze(-hdepth(xidx, yidx, zidx)), 'rx-') 0587 % xlabel('Temperature (^{\circ}C)') 0588 % ylabel('Depth (m)') 0589 % title('coarse') 0590 % 0591 % subplot(2,2,3) 0592 % plot(Mobj.temperature(nn, :, tt), 1:fz, 'x-') 0593 % xlabel('Temperature (^{\circ}C)') 0594 % ylabel('Array index') 0595 % title('FVCOM') 0596 % 0597 % subplot(2,2,4) 0598 % plot(squeeze(coarse.temperature.data(xidx, yidx, zidx, tt)), hz(zidx), 'rx-') 0599 % xlabel('Temperature (^{\circ}C)') 0600 % ylabel('Array index') 0601 % title('coarse') 0602 % 0603 % % Figure to check everything's as we'd expect. Plot first time step with 0604 % % the coarse surface temperature as a background with the interpolated 0605 % % boundary node surface values on top. 0606 % 0607 % figure(200) 0608 % clf 0609 % % Plot coarse surface data (last sigma layer) 0610 % dx = mean(diff(coarse.lon(:))); 0611 % dy = mean(diff(coarse.lat(:))); 0612 % temp = coarse.temperature.data(:, :, :, tt); 0613 % pcolor(coarse.lon - (dx / 2), coarse.lat - (dy / 2), ... 0614 % squeeze(temp(:, :, hyz))) 0615 % shading flat 0616 % axis('equal', 'tight') 0617 % daspect([1.5, 1, 1]) 0618 % hold on 0619 % % Add the interpolated surface data (first sigma layer) 0620 % scatter(Mobj.lon(oNodes), Mobj.lat(oNodes), 40, Mobj.temperature(:, fvz, tt), 'filled', 'MarkerEdgeColor', 'k') 0621 % axis([min(Mobj.lon(oNodes)), max(Mobj.lon(oNodes)), min(Mobj.lat(oNodes)), max(Mobj.lat(oNodes))]) 0622 % caxis([3, 13]) 0623 % plot(lon(xidx, yidx), lat(xidx, yidx), 'rs') % polcoms is all backwards 0624 % plot(Mobj.lon(oNodes(nn)), Mobj.lat(oNodes(nn)), 'wo') 0625 % colorbar