Home > fvcom_prepro > get_POLCOMS_meanflow.m

get_POLCOMS_meanflow

PURPOSE ^

Read mean flow from the PML POLCOMS-ERSEM NetCDF AMM model output files

SYNOPSIS ^

function Mobj = get_POLCOMS_meanflow(Mobj, files)

DESCRIPTION ^

 Read mean flow from the PML POLCOMS-ERSEM NetCDF AMM model output files
 and interpolate onto the open boundaries in Mobj.

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

 DESCRIPTION:
    Interpolate u and v flow vectors to calculate daily mean flow at the
    FVCOM open boundaries at all sigma levels.

 INPUT:
   Mobj    = MATLAB mesh structure which must contain:
               - Mobj.lon, Mobj.lat - node coordinates (lat/long).
               - Mobj.obc_elements - list of open boundary element
               inidices.
               - Mobj.nObcElements - number of elements in each open boundary.
   files   = Cell array of PML POLCOMS-ERSEM NetCDF file(s) in which 4D
             variables of u and v velocity components (called 'ucurD' and
             'vcurD') exist. Their shape should be (y, x, sigma, time).

 NOTES:

   - If you supply multiple files in files, 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 etc.). 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 six new fields (mf_time, meanflow_u,
           meanflow_v, meanflow_ubar, meanflow_ubar and velocity).
           meanflow_ubar, meanflow_vbar and velocity have sizes of
           (sum(Mobj.nObcNodes), time); meanflow_u and meanflow_v have
           sizes (sum(Mobj.nObcNodes), siglay, time). The time dimension
           is determined based on the input NetCDF file. The mf_times
           variable is just the input file times in Modified Julian Day.

 EXAMPLE USAGE
    Mobj = get_POLCOMS_meanflow(Mobj, files)

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

 Revision history
    2013-02-20 First version.
    2013-02-26 Add interpolation of the u and v vectors separately and
    then calculate the interpolated velocity at the end.
    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-02-28 Change the interpolation to occur on the open boundary
    elements rather than on the open boundary nodes. This requires a
    couple of extra steps in the pre-processing (notably running
    find_boundary_elements) as well as transferring the sigma depths to
    the element centres.

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

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function Mobj = get_POLCOMS_meanflow(Mobj, files)
0002 % Read mean flow from the PML POLCOMS-ERSEM NetCDF AMM model output files
0003 % and interpolate onto the open boundaries in Mobj.
0004 %
0005 % function Mobj = get_POLCOMS_meanflow(Mobj, ts, polcoms_bathy, varlist)
0006 %
0007 % DESCRIPTION:
0008 %    Interpolate u and v flow vectors to calculate daily mean flow at the
0009 %    FVCOM open boundaries at all sigma levels.
0010 %
0011 % INPUT:
0012 %   Mobj    = MATLAB mesh structure which must contain:
0013 %               - Mobj.lon, Mobj.lat - node coordinates (lat/long).
0014 %               - Mobj.obc_elements - list of open boundary element
0015 %               inidices.
0016 %               - Mobj.nObcElements - number of elements in each open boundary.
0017 %   files   = Cell array of PML POLCOMS-ERSEM NetCDF file(s) in which 4D
0018 %             variables of u and v velocity components (called 'ucurD' and
0019 %             'vcurD') exist. Their shape should be (y, x, sigma, time).
0020 %
0021 % NOTES:
0022 %
0023 %   - If you supply multiple files in files, 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 etc.). The fourth dimension is
0028 %       time.
0029 %       - The order of the files given should be chronological.
0030 %
0031 %   - The NetCDF files used here are those from the PML POLCOMS-ERSEM model
0032 %   output.
0033 %
0034 % OUTPUT:
0035 %    Mobj = MATLAB structure in which six new fields (mf_time, meanflow_u,
0036 %           meanflow_v, meanflow_ubar, meanflow_ubar and velocity).
0037 %           meanflow_ubar, meanflow_vbar and velocity have sizes of
0038 %           (sum(Mobj.nObcNodes), time); meanflow_u and meanflow_v have
0039 %           sizes (sum(Mobj.nObcNodes), siglay, time). The time dimension
0040 %           is determined based on the input NetCDF file. The mf_times
0041 %           variable is just the input file times in Modified Julian Day.
0042 %
0043 % EXAMPLE USAGE
0044 %    Mobj = get_POLCOMS_meanflow(Mobj, files)
0045 %
0046 % Author(s):
0047 %    Pierre Cazenave (Plymouth Marine Laboratory)
0048 %
0049 % Revision history
0050 %    2013-02-20 First version.
0051 %    2013-02-26 Add interpolation of the u and v vectors separately and
0052 %    then calculate the interpolated velocity at the end.
0053 %    2013-02-27 Change the vertical interpolation to be scaled within the
0054 %    POLCOMS-ERSEM depth range for the current node. The net result is that
0055 %    the vertical profiles are squashed or stretched to fit within the
0056 %    FVCOM depths. This means the full profile structure is maintained in
0057 %    the resulting FVCOM boundary input despite the differing depths at the
0058 %    FVCOM boundary node.
0059 %    2013-02-28 Change the interpolation to occur on the open boundary
0060 %    elements rather than on the open boundary nodes. This requires a
0061 %    couple of extra steps in the pre-processing (notably running
0062 %    find_boundary_elements) as well as transferring the sigma depths to
0063 %    the element centres.
0064 %
0065 %==========================================================================
0066 
0067 subname = 'get_POLCOMS_meanflow';
0068 
0069 global ftbverbose;
0070 if ftbverbose
0071     fprintf('\n')
0072     fprintf(['begin : ' subname '\n'])
0073 end
0074 
0075 varlist = {'lon', 'lat', 'ucurD', 'vcurD', 'time', 'depth', 'pdepthD'};
0076 
0077 pc = get_POLCOMS_netCDF(files, varlist);
0078 
0079 [~, ~, nz, nt] = size(pc.ucurD.data);
0080 
0081 [lon, lat] = meshgrid(pc.lon.data, pc.lat.data);
0082 
0083 % We need positions of the open boundary element centroids (i.e. the MCE
0084 % centroid) at which to calculate the mean flow velocity. As such, we need
0085 % to first calculate a number of parameters at the element centres from the
0086 % nodal values we typically collect (e.g. positions, sigma layers etc.). Do
0087 % all that here.
0088 lonc = nodes2elems(Mobj.lon, Mobj);
0089 latc = nodes2elems(Mobj.lat, Mobj);
0090 % For the sigma levels, we need to wrap the conversion in a loop.
0091 siglayzc = nan([Mobj.nElems, size(Mobj.siglayz, 2)]);
0092 for zz = 1:size(Mobj.siglayz, 2)
0093     siglayzc(:, zz) = nodes2elems(Mobj.siglayz(:, zz), Mobj);
0094 end
0095 
0096 oElements = [Mobj.read_obc_elements{:}];
0097 
0098 % Find the FVCOM positions for all the open boundary elements.
0099 fvlon = lonc(oElements);
0100 fvlat = latc(oElements);
0101 
0102 % Number of open boundary elements
0103 ne = sum(Mobj.nObcElements);
0104 % Number of sigma layers.
0105 fz = size(siglayzc, 2);
0106 
0107 fvmfu = nan(ne, fz, nt); % FVCOM interpolated flow vector components
0108 fvmfv = nan(ne, fz, nt); % FVCOM interpolated flow vector components
0109 
0110 if ftbverbose
0111     tic
0112 end
0113 
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     pcu3 = pc.ucurD.data(:, :, :, t);
0120     pcv3 = pc.vcurD.data(:, :, :, t);
0121     pcz3 = pc.depth.data(:, :, :, t);
0122 
0123     % Preallocate the intermediate results arrays.
0124     iuz = nan(ne, nz);
0125     ivz = nan(ne, nz);
0126     izz = nan(ne, nz);
0127 
0128     % We need to create a mask to eliminate land values and apply it to the
0129     % depth averaged values.
0130     mask = squeeze(pc.depth.data(:, :, end, t))' >= 0;
0131 
0132     % Mask the longitude and latitude here (since it's depth independent,
0133     % there's no point doing it inside the next loop).
0134     tlon = lon;
0135     tlat = lat;
0136     tlon(mask) = [];
0137     tlat(mask) = [];
0138 
0139     for z = 1:nz
0140         % Turns out we do need vertical velocity for the mean flow, so
0141         % iterate through each vertical layer before interpolating
0142         % horizontally. The vertical interpolation happens after the
0143         % horizontal has been completed. Transpose the arrays to be (x, y)
0144         % rather than (y, x).
0145         pcu2 = pcu3(:, :, z)';
0146         pcv2 = pcv3(:, :, z)';
0147         pcz2 = pcz3(:, :, z)';
0148 
0149         % Flatten the arrays through the masking.
0150         pcu2(mask) = [];
0151         pcv2(mask) = [];
0152         pcz2(mask) = [];
0153 
0154         % Preallocate the intermediate results arrays.
0155         iuobc = nan(ne, 1);
0156         ivobc = nan(ne, 1);
0157         izobc = nan(ne, 1);
0158 
0159         % Speed up the tightest loops with a parallelized loop.
0160         parfor i = 1:ne
0161             % Now we can do each position within the current depth layer.
0162 
0163             fx = fvlon(i);
0164             fy = fvlat(i);
0165 
0166             [~, ii] = sort(sqrt((tlon - fx).^2 + (tlat - fy).^2));
0167             % Get the n nearest elements from PML POLCOMS-ERSEM data (more?
0168             % fewer?).
0169             ixy = ii(1:16);
0170 
0171             % Get the variables into static variables for the
0172             % parallelisation.
0173             plon = tlon(ixy);
0174             plat = tlat(ixy);
0175             pu = pcu2(ixy);
0176             pv = pcv2(ixy);
0177             pz = pcz2(ixy);
0178 
0179             % Use a triangulation to do the horizontal interpolation.
0180             triu = TriScatteredInterp(plon', plat', pu', 'natural');
0181             triv = TriScatteredInterp(plon', plat', pv', 'natural');
0182             triz = TriScatteredInterp(plon', plat', pz', 'natural');
0183             iuobc(i) = triu(fx, fy);
0184             ivobc(i) = triv(fx, fy);
0185             izobc(i) = triz(fx, fy);
0186 
0187             % Check if we have NaNs (mostly if the position is outside the
0188             % model domain).
0189             if isnan(iuobc(i)) || isnan(ivobc(i))
0190                 warning('FVCOM boundary element at %f, %f is outside the PML POLCOMS-ERSEM domain. Setting to the closest PML POLCOMS-ERSEM value.', fx, fy)
0191                 iuobc(i) = pcu2(ii(1));
0192                 ivobc(i) = pcv2(ii(1));
0193                 izobc(i) = pcz2(ii(1));
0194             end
0195         end
0196 
0197         % Put the results in this intermediate array.
0198         iuz(:, z) = iuobc;
0199         ivz(:, z) = ivobc;
0200         izz(:, z) = izobc;
0201     end
0202 
0203     % Now we've interpolated in space, we can interpolate the z-values
0204     % to the sigma depths.
0205 
0206     % Preallocate the output arrays
0207     fvuz = nan(ne, fz);
0208     fvvz = nan(ne, fz);
0209 
0210     for pp = 1:ne
0211         % Get the FVCOM depths at this element
0212         tfz = siglayzc(oElements(pp), :);
0213         % Now get the interpolated PML POLCOMS-ERSEM depth at this element
0214         tpz = izz(pp, :);
0215 
0216         % To ensure we get the full vertical expression of the vertical
0217         % profiles, we need to normalise the POLCOMS-ERSEM depths to the
0218         % FVCOM range for the current element. This is because in instances
0219         % where FVCOM depths are shallower (e.g. in coastal regions), if we
0220         % don't normalise the depths, we end up truncating the vertical
0221         % profile. This approach ensures we always use the full vertical
0222         % profile, but we're potentially squeezing it into a smaller depth.
0223         A = max(tpz);
0224         B = min(tpz);
0225         C = max(tfz);
0226         D = min(tfz);
0227         norm_tpz = (((D - C) * (tpz - A)) / (B - A)) + C;
0228 
0229         % Get the u and v velocity values for this elment and interpolate
0230         % down the water column (from PML POLCOMS-ERSEM to FVCOM). I
0231         % had originally planned to use csaps for the vertical
0232         % interplation/subsampling at each location. However, the demo
0233         % of csaps in the MATLAB documentation makes the interpolation
0234         % look horrible (shaving off extremes). interp1 provides a
0235         % range of interpolation schemes, of which pchip seems to do a
0236         % decent job of the interpolation (at least qualitatively).
0237         if ~isnan(tpz)
0238             fvuz(pp, :) = interp1(norm_tpz, iuz(pp, :), tfz, 'linear', 'extrap');
0239             fvvz(pp, :) = interp1(norm_tpz, ivz(pp, :), tfz, 'linear', 'extrap');
0240         else
0241             warning('Should never see this... ') % because we test for NaNs when fetching the values.
0242             warning('FVCOM boundary element at %f, %f is outside the PML POLCOMS-ERSEM domain. Skipping.', fvlon(pp), fvlat(pp))
0243             continue
0244         end
0245     end
0246 
0247     % The horizontally- and vertically-interpolated values in the final
0248     % FVCOM results array.
0249     fvmfu(:, :, t) = fvuz;
0250     fvmfv(:, :, t) = fvvz;
0251 
0252     if ftbverbose
0253         fprintf('done.\n')
0254     end
0255 end
0256 if ftbverbose
0257     toc
0258 end
0259 
0260 % Stick the values in the mesh structure.
0261 Mobj.meanflow_u = fvmfu;
0262 Mobj.meanflow_v = fvmfv;
0263 
0264 % Now we have the 3D arrays, create depth averaged velocities too
0265 Mobj.meanflow_ubar = squeeze(mean(Mobj.meanflow_u, 2));
0266 Mobj.meanflow_vbar = squeeze(mean(Mobj.meanflow_v, 2));
0267 
0268 % Depth averaged velocity
0269 Mobj.velocity = squeeze(mean(sqrt(Mobj.meanflow_u.^2 + Mobj.meanflow_v.^2), 2));
0270 
0271 % Convert the current times to Modified Julian Day (this is a bit ugly).
0272 pc.time.all = strtrim(regexp(pc.time.units, 'since', 'split'));
0273 pc.time.datetime = strtrim(regexp(pc.time.all{end}, ' ', 'split'));
0274 pc.time.ymd = str2double(strtrim(regexp(pc.time.datetime{1}, '-', 'split')));
0275 pc.time.hms = str2double(strtrim(regexp(pc.time.datetime{2}, ':', 'split')));
0276 
0277 Mobj.mf_times = greg2mjulian(...
0278     pc.time.ymd(1), ...
0279     pc.time.ymd(2), ...
0280     pc.time.ymd(3), ...
0281     pc.time.hms(1), ...
0282     pc.time.hms(2), ...
0283     pc.time.hms(3)) + (pc.time.data / 3600 / 24);
0284 
0285 if ftbverbose
0286     fprintf(['end   : ' subname '\n'])
0287 end
0288 
0289 % Check the interpolation along the boundary at the last time step. N.B.
0290 % Since the depths are scaled from the POLCOMS range to the FVCOM range,
0291 % the profiles won't match perfectly in these plots. This is because the
0292 % interpolated FVCOM profiles use the full water column from the POLCOMS
0293 % data rather than truncating it at the FVCOM depth.
0294 
0295 % % Map the open boundaries with the depth averaged velocity as a background
0296 % figure(1)
0297 % clf
0298 % imagesc(pc.lon.data, pc.lat.data, mean(sqrt(pcu3.^2 + pcv3.^2), 3)')
0299 % shading flat
0300 % axis('equal', 'tight')
0301 % set(gca,'YDir','normal')
0302 % colorbar
0303 % caxis([0, 0.05])
0304 % hold on
0305 % plot(fvlon, fvlat, 'wo')
0306 % axis([min(fvlon) - 0.1, max(fvlon) + 0.1, min(fvlat) - 0.1, max(fvlat) + 0.1])
0307 %
0308 % for i = 1:ne
0309 %
0310 %     % Add the current element's position and value as a scatter point to
0311 %     % the map plot.
0312 %     scatter(fvlon(i), fvlat(i), 50, Mobj.velocity(i, end), 'filled')
0313 %
0314 %     % Do vertical profiles of u, v and velocity.
0315 %     figure(2)
0316 %     clf
0317 %
0318 %     subplot(3,1,1)
0319 %     % FVCOM vertical profile. Has to be (i, :, end) because the
0320 %     % corresponding POLCOMS data isn't stored as a function of time (i.e.
0321 %     % iuz, ivz and izz are all for the last time step only).
0322 %     fvuz = Mobj.meanflow_u(i, :, end);
0323 %     fvz = siglayzc(oElements(i), :);
0324 %     plot(fvuz, fvz)
0325 %     hold on
0326 %     % The interpolated POLCOMS vertical profile (last time step only)
0327 %     plot(iuz(i, :), izz(i, :), 'g')
0328 %     % The depth-averaged velocity (again, last time step only)
0329 %     fvubar = Mobj.meanflow_ubar(i, end);
0330 %     plot([fvubar, fvubar], [min(fvz), max(fvz)], 'k')
0331 %     xlim([-0.1, 0.2])
0332 %     ylim([-100, 0])
0333 %     xlabel('Mean flow u-velocity (m^{3}s^{-1})')
0334 %     ylabel('Depth (m)')
0335 %     title('u-component')
0336 %     legend('FVCOM', 'POLCOMS', 'Mean FVCOM', 'Location', 'SouthEast')
0337 %     legend('boxoff')
0338 %
0339 %     subplot(3,1,2)
0340 %     % FVCOM vertical profile. Has to be (i, :, end) because the
0341 %     % corresponding POLCOMS data isn't stored as a function of time (i.e.
0342 %     % iuz, ivz and izz are all for the last time step only).
0343 %     fvvz = Mobj.meanflow_v(i, :, end);
0344 %     fvz = siglayzc(oElements(i), :);
0345 %     plot(fvvz, fvz)
0346 %     hold on
0347 %     % The interpolated POLCOMS vertical profile (last time step only)
0348 %     plot(ivz(i, :), izz(i, :), 'g')
0349 %     % The depth-averaged velocity (again, last time step only)
0350 %     fvvbar = Mobj.meanflow_vbar(i, end);
0351 %     plot([fvvbar, fvvbar], [min(fvz), max(fvz)], 'k')
0352 %     xlim([-0.1, 0.2])
0353 %     ylim([-100, 0])
0354 %     xlabel('Mean flow v-velocity (m^{3}s^{-1})')
0355 %     ylabel('Depth (m)')
0356 %     title('v-component')
0357 %     legend('FVCOM', 'POLCOMS', 'Mean FVCOM', 'Location', 'SouthEast')
0358 %     legend('boxoff')
0359 %
0360 %     subplot(3,1,3)
0361 %     % FVCOM vertical profile. Has to be (i, :, end) because the
0362 %     % corresponding POLCOMS data isn't stored as a function of time (i.e.
0363 %     % iuz, ivz and izz are all for the last time step only).
0364 %     fvvelz = sqrt(Mobj.meanflow_u(i, :, end).^2 + Mobj.meanflow_v(i, :, end).^2);
0365 %     fvz = siglayzc(oElements(i), :);
0366 %     plot(fvvelz, fvz)
0367 %     hold on
0368 %     % The interpolated POLCOMS vertical profile (last time step only)
0369 %     plot(sqrt(iuz(i, :).^2 + ivz(i, :).^2), izz(i, :), 'g')
0370 %     % The depth-averaged velocity (again, last time step only)
0371 %     fvvelbar = Mobj.velocity(i, end);
0372 %     plot([fvvelbar, fvvelbar], [min(fvz), max(fvz)], 'k')
0373 %     xlim([-0.1, 0.2])
0374 %     ylim([-100, 0])
0375 %     xlabel('Mean flow velocity (m^{3}s^{-1})')
0376 %     ylabel('Depth (m)')
0377 %     title('velocity')
0378 %     legend('FVCOM', 'POLCOMS', 'Mean FVCOM', 'Location', 'SouthEast')
0379 %     legend('boxoff')
0380 %     pause(0.1)
0381 % end
0382

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