Home > fvcom_prepro > get_POLCOMS_tsrestart_NOCL.m

get_POLCOMS_tsrestart_NOCL

PURPOSE ^

Extract temperature and salinity forcing information from NOC Operational

SYNOPSIS ^

function restart = get_POLCOMS_tsrestart_NOCL(Mobj, inputConf)

DESCRIPTION ^

 Extract temperature and salinity forcing information from NOC Operational
 Tide Surge model output and interpolate onto the FVCOM domain for use in 
 a reproduced FVCOM restart file.

 function restart = get_POLCOMS_tsrestart_NOCL(Mobj, inputConf)

 DESCRIPTION:
    Interpolate temperature and salinity values onto the FVCOM nodes
    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 and/or Mobj.x, Mobj.y - node
                   coordinates.
   inputConf   = MATLAB structure which must contain: 
                   - inputConf.polcoms_ts - location of NOC Operational
                   Model output containing 4D variables of temperature
                   (tem) and salinity (sal). They should have dimensions
                   (x, y, sigma, time).
                   - inputConf.polcoms_z - location of NOC Operational
                   Model output containing 4D variables of bathymetry
                   (XXX) and sigma layer thickness (XXX).
                   - inputConf.startDate - start date and time for FVCOM
                   model run
 
 OUTPUT:
    restart = MATLAB structure which contains three fields (called
              temperature, salinity and ts_time). temperature and salinity
              have sizes (Mobj.nVerts, 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
    restart = get_POLCOMS_tsrestart_NOCL(Mobj, inputConf)

 Author(s):
    Pierre Cazenave (Plymouth Marine Laboratory)
    Karen Amoudry (National Oceanography Centre, Liverpool)

 PWC Revision history
    2013-01-09 First version based on the FVCOM shelf model
    get_POLCOMS_forcing.m script (i.e. not a function but a plain script).

 KJA Revision history:
    2014-01-15 First version, adapted from KJA's
    'get_POLCOMS_tsobc_NOCL.m', which in turn was based on PWC's
    'get_POLCOMS_tsobc.m'.

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

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function restart = get_POLCOMS_tsrestart_NOCL(Mobj, inputConf)
0002 % Extract temperature and salinity forcing information from NOC Operational
0003 % Tide Surge model output and interpolate onto the FVCOM domain for use in
0004 % a reproduced FVCOM restart file.
0005 %
0006 % function restart = get_POLCOMS_tsrestart_NOCL(Mobj, inputConf)
0007 %
0008 % DESCRIPTION:
0009 %    Interpolate temperature and salinity values onto the FVCOM nodes
0010 %    at all sigma levels.
0011 %
0012 % INPUT:
0013 %   Mobj        = MATLAB mesh structure which must contain:
0014 %                   - Mobj.siglayz - sigma layer depths for all model
0015 %                   nodes.
0016 %                   - Mobj.lon, Mobj.lat and/or Mobj.x, Mobj.y - node
0017 %                   coordinates.
0018 %   inputConf   = MATLAB structure which must contain:
0019 %                   - inputConf.polcoms_ts - location of NOC Operational
0020 %                   Model output containing 4D variables of temperature
0021 %                   (tem) and salinity (sal). They should have dimensions
0022 %                   (x, y, sigma, time).
0023 %                   - inputConf.polcoms_z - location of NOC Operational
0024 %                   Model output containing 4D variables of bathymetry
0025 %                   (XXX) and sigma layer thickness (XXX).
0026 %                   - inputConf.startDate - start date and time for FVCOM
0027 %                   model run
0028 %
0029 % OUTPUT:
0030 %    restart = MATLAB structure which contains three fields (called
0031 %              temperature, salinity and ts_time). temperature and salinity
0032 %              have sizes (Mobj.nVerts, sigma, time). The time dimension is
0033 %              determined based on the input NetCDF file. The ts_time
0034 %              variable is just the input file times in Modified Julian Day.
0035 %
0036 % EXAMPLE USAGE
0037 %    restart = get_POLCOMS_tsrestart_NOCL(Mobj, inputConf)
0038 %
0039 % Author(s):
0040 %    Pierre Cazenave (Plymouth Marine Laboratory)
0041 %    Karen Amoudry (National Oceanography Centre, Liverpool)
0042 %
0043 % PWC Revision history
0044 %    2013-01-09 First version based on the FVCOM shelf model
0045 %    get_POLCOMS_forcing.m script (i.e. not a function but a plain script).
0046 %
0047 % KJA Revision history:
0048 %    2014-01-15 First version, adapted from KJA's
0049 %    'get_POLCOMS_tsobc_NOCL.m', which in turn was based on PWC's
0050 %    'get_POLCOMS_tsobc.m'.
0051 %
0052 %==========================================================================
0053 
0054 subname = 'get_POLCOMS_tsrestart_NOCL';
0055 
0056 global ftbverbose;
0057 if ftbverbose
0058     fprintf('\n')
0059     fprintf(['begin : ' subname '\n'])
0060 end
0061 %%
0062 % Which variables do we want from the POLCOMS file?
0063 varlist = {'lon', 'lat', 'tem', 'sal', 'time'};
0064 
0065 % Create the POLCOMS filename based on the year and month of interest
0066 polcoms_ts = [inputConf.polcoms_ts,num2str(inputConf.startDate(1)),...
0067     '-',num2str(inputConf.startDate(2),'%02i'),'.nc'];
0068 
0069 % Get the results from the POLCOMS file
0070 nc = netcdf.open(polcoms_ts, 'NOWRITE');
0071 
0072 for var=1:numel(varlist)
0073     
0074     getVar = varlist{var};
0075     varid_pc = netcdf.inqVarID(nc, getVar);
0076     
0077     data = netcdf.getVar(nc, varid_pc, 'single');
0078     
0079     pc.(getVar).data = double(data);
0080     
0081     % Try to get some units (important for the calculation of MJD).
0082     try
0083         units = netcdf.getAtt(nc,varid_pc,'units');
0084     catch
0085         units = [];
0086     end
0087     pc.(getVar).units = units;
0088 end
0089 
0090 netcdf.close(nc)
0091     
0092 % Get rid of data outside the time we're interested in
0093 % Convert FVCOM timestep into AMM/S12 output (seconds since 20071101:000000)
0094 AMM_time = etime(inputConf.startDate,[2007,11,01,0,0,0]);
0095 keep_time = (pc.time.data == AMM_time+(12*60*60)); % Allowance for POLCOMS data being at noon, not midnight
0096 pc.tem.data = pc.tem.data(:,:,:,keep_time);
0097 pc.sal.data = pc.sal.data(:,:,:,keep_time);
0098 
0099 % Import the POLCOMS sigma info
0100 pc = get_POLCOMS_sigma(pc,inputConf);
0101 
0102 % Data format:
0103 %   pc.tem.data and pc.sal.data are x, y, sigma
0104 [~, ~, nz] = size(pc.tem.data);
0105 
0106 % Make rectangular arrays for the nearest point lookup.
0107 [plon, plat] = meshgrid(pc.lon.data, pc.lat.data);
0108 
0109 % Number of sigma layers.
0110 fz = size(Mobj.siglayz, 2);
0111 
0112 if ftbverbose
0113     tic
0114 end
0115 
0116 %  Flip the vertical layer dimension to make the POLCOMS data go from
0117 % surface to seabed to match its depth data and to match how FVCOM
0118 % works.
0119 pctemp3 = flipdim(pc.tem.data, 3);
0120 pcsalt3 = flipdim(pc.sal.data, 3);
0121 
0122 % Preallocate the FVCOM results arrays
0123 fvtemp = nan(Mobj.nVerts, fz); % FVCOM interpolated temperatures
0124 fvsal = nan(Mobj.nVerts, fz); % FVCOM interpolated salinities
0125 
0126 % Preallocate the intermediate results arrays.
0127 itempz = nan(Mobj.nVerts, nz);
0128 isalz = nan(Mobj.nVerts, nz);
0129 idepthz = nan(Mobj.nVerts,nz);
0130 
0131 for j = 1:nz
0132     % Transpose the POLCOMS data to be (x,y) rather than (y,x)
0133     pctemp2 = pctemp3(:, :, j)';
0134     pcsalt2 = pcsalt3(:, :, j)';
0135     pcdepth2 = squeeze(pc.depth.data(:, :, j))';
0136     
0137     % Reshape the arrays to allow the sort to work properly later
0138     tlon=reshape(plon,(size(plon,1)*size(plon,2)),1);
0139     tlat=reshape(plat,(size(plat,1)*size(plat,2)),1);
0140     pctemp2 = reshape(pctemp2,(size(pctemp2,1)*size(pctemp2,2)),1);
0141     pcsalt2 = reshape(pcsalt2,(size(pcsalt2,1)*size(pcsalt2,2)),1);
0142     pcdepth2 = reshape(pcdepth2,(size(pcdepth2,1)*size(pcdepth2,2)),1);
0143     
0144     % Find the points which aren't NaNs
0145     keeptemp = find(~isnan(pctemp2));
0146     keepsalt = find(~isnan(pcsalt2));
0147     keepdepth = find(~isnan(pcdepth2));
0148     
0149     keep = intersect(keeptemp,keepsalt);
0150     keep = intersect(keep, keepdepth);
0151     
0152     % Interpolate the POLCOMS results by (1) turning the non-NaN POLCOMS
0153     % results into a Tri object and (2) extracting the values of that Tri
0154     % Object at the FVCOM node locations.
0155     tritemp = TriScatteredInterp(tlon(keep),tlat(keep),pctemp2(keep),'natural');
0156     itemp = tritemp(Mobj.lon,Mobj.lat);
0157     
0158     trisalt = TriScatteredInterp(tlon(keep),tlat(keep),pcsalt2(keep),'natural');
0159     isalt = trisalt(Mobj.lon,Mobj.lat);
0160     
0161     tridepth = TriScatteredInterp(tlon(keep),tlat(keep),pcdepth2(keep),'natural');
0162     idepth = tridepth(Mobj.lon,Mobj.lat);
0163     
0164     % Put the results in this intermediate array.
0165     itempz(:, j) = itemp;
0166     isalz(:, j) = isalt;
0167     idepthz(:,j) = idepth;
0168 end
0169 
0170 % Now we've interpolated in space, we can interpolate the z-values
0171 % to the sigma depths.
0172 % for zi = 1:fz
0173     for pp = 1:Mobj.nVerts
0174         % Get the FVCOM depths at this node
0175         tfz = Mobj.siglayz(pp, :);
0176         % Now get the interpolated POLCOMS depth at this node
0177         tpz = idepthz(pp, :);
0178         
0179         % To ensure we get the full vertical expression of the vertical
0180         % profiles, we need to normalise the POLCOMS and FVCOM
0181         % depths to the same range. This is because in instances where
0182         % FVCOM depths are shallower (e.g. in coastal regions), if we
0183         % don't normalise the depths, we end up truncating the vertical
0184         % profile. This approach ensures we always use the full
0185         % vertical profile, but we're potentially squeezing it into a
0186         % smaller depth.
0187         A = max(tpz);
0188         B = min(tpz);
0189         C = max(tfz);
0190         D = min(tfz);
0191         norm_tpz = (((D - C) * (tpz - A)) / (B - A)) + C;
0192         
0193         % Get the temperature and salinity values for this node and
0194         % interpolate down the water column (from POLCOMS to FVCOM).
0195         % Change to 'pchip' to match PWC parent code.
0196         if ~isnan(tpz)
0197             fvtemp(pp, :) = interp1(norm_tpz, itempz(pp, :), tfz, 'pchip', 'extrap');
0198             fvsal(pp, :) = interp1(norm_tpz, isalz(pp, :), tfz, 'pchip', 'extrap');
0199         else
0200             warning('Should never see this... ') % because we test for NaNs when fetching the values.
0201             warning('FVCOM boundary node at %f, %f is outside the POLCOMS domain. Skipping.', Mobj.lon(pp), Mobj.lat(pp))
0202             continue
0203         end
0204     end
0205 % end
0206 
0207 if ftbverbose
0208     toc
0209 end
0210 
0211 % Convert NOC model output temperatures from Kelvin to Celsius
0212 fvtemp = fvtemp - 273.15;
0213 
0214 % Timeshift to match the expected FVCOM input times. The temperature and
0215 % salinity values are a daily average (midnight to midnight). They are
0216 % given a timestamp of the middle time in this period (noon). Each day's
0217 % value applies to the whole day (e.g. from midnight to midnight).
0218 % Therefore, we can shift the time to midnight at the start of the relevant
0219 % day and apply the value to the whole day. This makes FVCOM (and me)
0220 % happy.
0221 restart.ts_times = datenum(2007, 11, 0, 0, 0, 0) + ...    % Convert NOC model reference time to Matlab datenum
0222     ((pc.time.data(keep_time)+(12*3600)) / 3600 / 24);  % Add NOC model output time and 12 hour offset
0223 restart.ts_times = datevec(restart.ts_times);   % Convert times to datevec format
0224 
0225 % Final output variable names must match FVCOM restart file variable names
0226 restart.temp = fvtemp;
0227 restart.salinity = fvsal;
0228 
0229 if ftbverbose
0230     fprintf(['end   : ' subname '\n'])
0231 end

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