Home > fvcom_prepro > get_POLCOMS_tsobc.m

get_POLCOMS_tsobc

PURPOSE ^

Read temperature and salinity from the PML POLCOMS-ERSEM NetCDF model

SYNOPSIS ^

function Mobj = get_POLCOMS_tsobc(Mobj, ts)

DESCRIPTION ^

 Read temperature and salinity from the PML POLCOMS-ERSEM NetCDF model
 output files and interpolate onto the open boundaries in Mobj.

 function Mobj = get_POLCOMS_tsobc(Mobj, ts, polcoms_bathy, varlist)

 DESCRIPTION:
    Interpolate temperature and salinity 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.lon, Mobj.lat - node coordinates (lat/long).
               - Mobj.obc_nodes - list of open boundary node inidices.
               - Mobj.nObcNodes - number of nodes in each open boundary.
   ts      = Cell array of PML POLCOMS-ERSEM NetCDF file(s) in which 4D
             variables of temperature and salinity (called 'ETWD' and
             'x1XD') exist. Their shape should be (y, x, sigma, time).

 NOTES:

   - If you supply multiple files in ts, there are a few assumptions:

       - Variables are only appended if there are 4 dimensions; fewer than
       that, and the values are assumed to be static across all the given
       files (e.g. longitude, latitude). The fourth dimension is time.
       - The order of the files given should be chronological.
 
   - The NetCDF files used here are those from the PML POLCOMS-ERSEM model
   output.

 OUTPUT:
    Mobj = MATLAB structure in which three new fields (called temperature,
           salinity and ts_times). temperature and salinity have sizes
           (sum(Mobj.nObcNodes), sigma, time). The time dimension is
           determined based on the input NetCDF file. The ts_time variable
           is just the input file times in Modified Julian Day.

 EXAMPLE USAGE
    Mobj = get_POLCOMS_tsobc(Mobj, ts)

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

 Revision history
    2013-02-07 First version.
    2013-02-27 Change the vertical interpolation to be scaled within the
    POLCOMS-ERSEM depth range for the current node. The net result is that
    the vertical profiles are squashed or stretched to fit within the
    FVCOM depths. This means the full profile structure is maintained in
    the resulting FVCOM boundary input despite the differing depths at the
    FVCOM boundary node.
    2013-06-03 Fix some bugs in the way the open boundary node values were
    stored (the order in which they were stored did not match the order of
    the nodes in Casename_obc.dat). Also fix the order of the vertically
    interpolated values so that FVCOM starts at the surface instead of
    mirroring POLCOMS' approach (where the first value is the seabed). The
    effect of these two fixes (nodes and vertical) should match what FVCOM
    expects. Also add a set of figures (commented out) at the end for
    diagnostic purposes.

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

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function Mobj = get_POLCOMS_tsobc(Mobj, ts)
0002 % Read temperature and salinity from the PML POLCOMS-ERSEM NetCDF model
0003 % output files and interpolate onto the open boundaries in Mobj.
0004 %
0005 % function Mobj = get_POLCOMS_tsobc(Mobj, ts, polcoms_bathy, varlist)
0006 %
0007 % DESCRIPTION:
0008 %    Interpolate temperature and salinity values onto the FVCOM open
0009 %    boundaries at all 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.lon, Mobj.lat - node coordinates (lat/long).
0015 %               - Mobj.obc_nodes - list of open boundary node inidices.
0016 %               - Mobj.nObcNodes - number of nodes in each open boundary.
0017 %   ts      = Cell array of PML POLCOMS-ERSEM NetCDF file(s) in which 4D
0018 %             variables of temperature and salinity (called 'ETWD' and
0019 %             'x1XD') exist. Their shape should be (y, x, sigma, time).
0020 %
0021 % NOTES:
0022 %
0023 %   - If you supply multiple files in ts, there are a few assumptions:
0024 %
0025 %       - Variables are only appended if there are 4 dimensions; fewer than
0026 %       that, and the values are assumed to be static across all the given
0027 %       files (e.g. longitude, latitude). The fourth dimension is time.
0028 %       - The order of the files given should be chronological.
0029 %
0030 %   - The NetCDF files used here are those from the PML POLCOMS-ERSEM model
0031 %   output.
0032 %
0033 % OUTPUT:
0034 %    Mobj = MATLAB structure in which three new fields (called temperature,
0035 %           salinity and ts_times). temperature and salinity have sizes
0036 %           (sum(Mobj.nObcNodes), sigma, time). The time dimension is
0037 %           determined based on the input NetCDF file. The ts_time variable
0038 %           is just the input file times in Modified Julian Day.
0039 %
0040 % EXAMPLE USAGE
0041 %    Mobj = get_POLCOMS_tsobc(Mobj, ts)
0042 %
0043 % Author(s):
0044 %    Pierre Cazenave (Plymouth Marine Laboratory)
0045 %
0046 % Revision history
0047 %    2013-02-07 First version.
0048 %    2013-02-27 Change the vertical interpolation to be scaled within the
0049 %    POLCOMS-ERSEM depth range for the current node. The net result is that
0050 %    the vertical profiles are squashed or stretched to fit within the
0051 %    FVCOM depths. This means the full profile structure is maintained in
0052 %    the resulting FVCOM boundary input despite the differing depths at the
0053 %    FVCOM boundary node.
0054 %    2013-06-03 Fix some bugs in the way the open boundary node values were
0055 %    stored (the order in which they were stored did not match the order of
0056 %    the nodes in Casename_obc.dat). Also fix the order of the vertically
0057 %    interpolated values so that FVCOM starts at the surface instead of
0058 %    mirroring POLCOMS' approach (where the first value is the seabed). The
0059 %    effect of these two fixes (nodes and vertical) should match what FVCOM
0060 %    expects. Also add a set of figures (commented out) at the end for
0061 %    diagnostic purposes.
0062 %
0063 %==========================================================================
0064 
0065 subname = 'get_POLCOMS_tsobc';
0066 
0067 global ftbverbose;
0068 if ftbverbose
0069     fprintf('\n')
0070     fprintf(['begin : ' subname '\n'])
0071 end
0072 
0073 if license('test', 'Distrib_Computing_Toolbox')
0074     % We have the Parallel Computing Toolbox, so launch a bunch of workers.
0075     if isempty(gcp('nocreate'))
0076         % Force pool to be local in case we have remote pools available.
0077         parpool('local');
0078     end
0079 end
0080 
0081 varlist = {'lon', 'lat', 'ETWD', 'x1XD', 'time', 'depth', 'pdepthD'};
0082 
0083 % Data format:
0084 %
0085 %   pc.ETWD.data and pc.x1XD.data are y, x, sigma, time
0086 %
0087 pc = get_POLCOMS_netCDF(ts, varlist);
0088 
0089 [~, ~, nz, nt] = size(pc.ETWD.data);
0090 
0091 % Make rectangular arrays for the nearest point lookup.
0092 [lon, lat] = meshgrid(pc.lon.data, pc.lat.data);
0093 
0094 %oNodes = Mobj.obc_nodes(Mobj.obc_nodes ~= 0);
0095 % Change the way the nodes are listed to match the order in the
0096 % Casename_obc.dat file.
0097 tmpObcNodes = Mobj.obc_nodes';
0098 oNodes = tmpObcNodes(tmpObcNodes ~= 0)';
0099 
0100 fvlon = Mobj.lon(oNodes);
0101 fvlat = Mobj.lat(oNodes);
0102 
0103 % Number of boundary nodes
0104 nf = sum(Mobj.nObcNodes);
0105 % Number of sigma layers.
0106 fz = size(Mobj.siglayz, 2);
0107 
0108 fvtemp = nan(nf, fz, nt); % FVCOM interpolated temperatures
0109 fvsal = nan(nf, fz, nt); % FVCOM interpolated salinities
0110 
0111 if ftbverbose
0112     tic
0113 end
0114 for t = 1:nt
0115     if ftbverbose
0116         fprintf('%s : %i of %i timesteps... ', subname, t, nt)
0117     end
0118     % Get the current 3D array of PML POLCOMS-ERSEM results.
0119     pctemp3 = pc.ETWD.data(:, :, :, t);
0120     pcsalt3 = pc.x1XD.data(:, :, :, t);
0121 
0122     % Preallocate the intermediate results arrays.
0123     itempz = nan(nf, nz);
0124     isalz = nan(nf, nz);
0125     idepthz = nan(nf, nz);
0126 
0127     for j = 1:nz
0128         % Now extract the relevant layer from the 3D subsets. Transpose the
0129         % data to be (x, y) rather than (y, x).
0130         pctemp2 = pctemp3(:, :, j)';
0131         pcsalt2 = pcsalt3(:, :, j)';
0132         pcdepth2 = squeeze(pc.depth.data(:, :, j, t))';
0133 
0134         % Create new arrays which will be flattened when masking (below).
0135         tpctemp2 = pctemp2;
0136         tpcsalt2 = pcsalt2;
0137         tpcdepth2 = pcdepth2;
0138         tlon = lon;
0139         tlat = lat;
0140 
0141         % Create and apply a mask to remove values outside the domain. This
0142         % inevitably flattens the arrays, but it shouldn't be a problem
0143         % since we'll be searching for the closest values in such a manner
0144         % as is appropriate for an unstructured grid (i.e. we're assuming
0145         % the PML POLCOMS-ERSEM data is irregularly spaced).
0146         mask = tpcdepth2 > 20000;
0147         tpctemp2(mask) = [];
0148         tpcsalt2(mask) = [];
0149         tpcdepth2(mask) = [];
0150         % Also apply the masks to the position arrays so we can't even find
0151         % positions outside the domain, effectively meaning if a value is
0152         % outside the domain, the nearest value to the boundary node will
0153         % be used.
0154         tlon(mask) = [];
0155         tlat(mask) = [];
0156 
0157         % Preallocate the intermediate results arrays.
0158         itempobc = nan(nf, 1);
0159         isalobc = nan(nf, 1);
0160         idepthobc = nan(nf, 1);
0161 
0162         % Speed up the tightest loop with a parallelized loop.
0163         parfor i = 1:nf
0164             % Now we can do each position within the 2D layer.
0165 
0166             fx = fvlon(i);
0167             fy = fvlat(i);
0168 
0169             [~, ii] = sort(sqrt((tlon - fx).^2 + (tlat - fy).^2));
0170             % Get the n nearest nodes from PML POLCOMS-ERSEM data (more?
0171             % fewer?).
0172             ixy = ii(1:16);
0173 
0174             % Get the variables into static variables for the
0175             % parallelisation.
0176             plon = tlon(ixy);
0177             plat = tlat(ixy);
0178             ptemp = tpctemp2(ixy);
0179             psal = tpcsalt2(ixy);
0180             pdepth = tpcdepth2(ixy);
0181 
0182             % Use a triangulation to do the horizontal interpolation.
0183             tritemp = TriScatteredInterp(plon', plat', ptemp', 'natural');
0184             trisal = TriScatteredInterp(plon', plat', psal', 'natural');
0185             triz = TriScatteredInterp(plon', plat', pdepth', 'natural');
0186             itempobc(i) = tritemp(fx, fy);
0187             isalobc(i) = trisal(fx, fy);
0188             idepthobc(i) = triz(fx, fy);
0189 
0190             % Check all three, though if one is NaN, they all will be.
0191             if isnan(itempobc(i)) || isnan(isalobc(i)) || isnan(idepthobc(i))
0192                 warning('FVCOM boundary node at %f, %f is outside the PML POLCOMS-ERSEM domain. Setting to the closest PML POLCOMS-ERSEM value.', fx, fy)
0193                 itempobc(i) = tpctemp2(ii(1));
0194                 isalobc(i) = tpcsalt2(ii(1));
0195                 idepthobc(i) = tpcdepth2(ii(1));
0196             end
0197         end
0198 
0199         % Put the results in the intermediate array.
0200         itempz(:, j) = itempobc;
0201         isalz(:, j) = isalobc;
0202         idepthz(:, j) = idepthobc;
0203 
0204     end
0205 
0206     % Now we've interpolated in space, we can interpolate the z-values
0207     % to the sigma depths.
0208 
0209     % Preallocate the output arrays
0210     fvtempz = nan(nf, fz);
0211     fvsalz = nan(nf, fz);
0212 
0213     for pp = 1:nf
0214         % Get the FVCOM depths at this node
0215         tfz = Mobj.siglayz(oNodes(pp), :);
0216         % Now get the interpolated PML POLCOMS-ERSEM depth at this node
0217         tpz = idepthz(pp, :);
0218 
0219         % To ensure we get the full vertical expression of the vertical
0220         % profiles, we need to normalise the POLCOMS-ERSEM and FVCOM
0221         % depths to the same range. This is because in instances where
0222         % FVCOM depths are shallower (e.g. in coastal regions), if we
0223         % don't normalise the depths, we end up truncating the vertical
0224         % profile. This approach ensures we always use the full
0225         % vertical profile, but we're potentially squeezing it into a
0226         % smaller depth.
0227         A = max(tpz);
0228         B = min(tpz);
0229         C = max(tfz);
0230         D = min(tfz);
0231         norm_tpz = (((D - C) * (tpz - A)) / (B - A)) + C;
0232 
0233         % Get the temperature and salinity values for this node and
0234         % interpolate down the water column (from PML POLCOMS-ERSEM to
0235         % FVCOM). I had originally planned to use csaps for the
0236         % vertical interplation/subsampling at each location. However,
0237         % the demo of csaps in the MATLAB documentation makes the
0238         % interpolation look horrible (shaving off extremes). interp1
0239         % provides a range of interpolation schemes, of which pchip
0240         % seems to do a decent job of the interpolation (at least
0241         % qualitatively).
0242         if ~isnan(tpz)
0243             fvtempz(pp, :) = interp1(norm_tpz, itempz(pp, :), tfz, 'pchip', 'extrap');
0244             fvsalz(pp, :) = interp1(norm_tpz, isalz(pp, :), tfz, 'pchip', 'extrap');
0245         else
0246             warning('Should never see this... ') % because we test for NaNs when fetching the values.
0247             warning('FVCOM boundary node at %f, %f is outside the PML POLCOMS-ERSEM domain. Skipping.', fvlon(pp), fvlat(pp))
0248             continue
0249         end
0250     end
0251 
0252     % The horizontally- and vertically-interpolated values in the final
0253     % FVCOM results array.
0254     fvtemp(:, :, t) = fvtempz;
0255     fvsal(:, :, t) = fvsalz;
0256 
0257     if ftbverbose
0258         fprintf('done.\n')
0259     end
0260 end
0261 if ftbverbose
0262     toc
0263 end
0264 
0265 Mobj.temperature = fvtemp;
0266 Mobj.salt = fvsal;
0267 
0268 % Convert the current times to Modified Julian Day (this is a bit ugly).
0269 pc.time.all = strtrim(regexp(pc.time.units, 'since', 'split'));
0270 pc.time.datetime = strtrim(regexp(pc.time.all{end}, ' ', 'split'));
0271 pc.time.ymd = str2double(strtrim(regexp(pc.time.datetime{1}, '-', 'split')));
0272 pc.time.hms = str2double(strtrim(regexp(pc.time.datetime{2}, ':', 'split')));
0273 
0274 Mobj.ts_times = greg2mjulian(...
0275     pc.time.ymd(1), ...
0276     pc.time.ymd(2), ...
0277     pc.time.ymd(3), ...
0278     pc.time.hms(1), ...
0279     pc.time.hms(2), ...
0280     pc.time.hms(3)) + (pc.time.data / 3600 / 24);
0281 
0282 if ftbverbose
0283     fprintf(['end   : ' subname '\n'])
0284 end
0285 
0286 
0287 %%
0288 % Plot a vertical profile for a boundary node (for my Irish Sea case, this
0289 % is one of the ones along the Celtic Sea boundary). Also plot the
0290 % distribution of interpolated values over the POLCOMS data. Add the
0291 % location of the vertical profile (both FVCOM and POLCOMS) to the plot.
0292 % nn = 55;   % open boundary index
0293 % tt = 1;    % time index
0294 %
0295 % % Get the corresponding indices for the POLCOMS data
0296 % [~, xidx] = min(abs(lon(1, :) - fvlon(nn)));
0297 % [~, yidx] = min(abs(lat(:, 1) - fvlat(nn)));
0298 %
0299 % % close all
0300 %
0301 % figure
0302 % clf
0303 % subplot(2,2,1)
0304 % plot(Mobj.temperature(nn, :, 1), Mobj.siglayz(oNodes(nn), :), 'x-')
0305 % xlabel('Temperature (^{\circ}C)')
0306 % ylabel('Depth (m)')
0307 % title('FVCOM')
0308 %
0309 % subplot(2,2,2)
0310 % % Although POLCOMS stores its temperature values from seabed to surface,
0311 % % the depths are stored surface to seabed. Nice.
0312 % plot(squeeze(pc.ETWD.data(xidx, yidx, :, 1)), squeeze(pc.depth.data(xidx, yidx, :, 1)), 'rx-')
0313 % xlabel('Temperature (^{\circ}C)')
0314 % ylabel('Depth (m)')
0315 % title('POLCOMS')
0316 %
0317 % subplot(2,2,3)
0318 % plot(Mobj.temperature(nn, :, tt), 1:fz, 'x-')
0319 % xlabel('Temperature (^{\circ}C)')
0320 % ylabel('Array index')
0321 % title('FVCOM')
0322 %
0323 % subplot(2,2,4)
0324 % plot(squeeze(pc.ETWD.data(xidx, yidx, :, tt)), 1:nz, 'rx-')
0325 % xlabel('Temperature (^{\circ}C)')
0326 % ylabel('Array index')
0327 % title('POLCOMS')
0328 %
0329 % % Figure to check everything's as we'd expect. Plot first time step with
0330 % % the POLCOMS surface temperature as a background with the interpolated
0331 % % boundary node surface values on top.
0332 %
0333 % figure
0334 % clf
0335 % % Plot POLCOMS surface data (last sigma layer)
0336 % dx = mean(diff(pc.lon.data));
0337 % dy = mean(diff(pc.lat.data));
0338 % pcolor(pc.lon.data - (dx / 2), pc.lat.data - (dy / 2), ...
0339 %     squeeze(pc.ETWD.data(:, :, 1, tt))')
0340 % shading flat
0341 % axis('equal', 'tight')
0342 % daspect([1.5, 1, 1])
0343 % hold on
0344 % % Add the interpolated surface data (first sigma layer)
0345 % scatter(Mobj.lon(oNodes), Mobj.lat(oNodes), 40, Mobj.temperature(:, 1, tt), 'filled', 'MarkerEdgeColor', 'k')
0346 % axis([min(Mobj.lon(oNodes)), max(Mobj.lon(oNodes)), min(Mobj.lat(oNodes)), max(Mobj.lat(oNodes))])
0347 % caxis([6, 20])
0348 % plot(lon(yidx, xidx), lat(yidx, xidx), 'rs') % polcoms is all backwards
0349 % plot(Mobj.lon(oNodes(nn)), Mobj.lat(oNodes(nn)), 'wo')
0350 % colorbar

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