


Extract temperature and salinity boundary forcing information from NOC
Operation Tide Surge model output.
function Mobj = get_POLCOMS_tsobc_NOCL(Mobj, inputConf)
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 and/or Mobj.x, Mobj.y - node
coordinates.
- Mobj.obc_nodes - list of open boundary node inidices.
- Mobj.nObcNodes - number of nodes in each open
boundary.
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
- inputConf.endDate - end date and time for FVCOM
model run
OUTPUT:
Mobj = MATLAB structure in which three new fields (called temperature,
salinity and ts_time). 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_forcing_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:
2013-02-05 Adapted from PWC's script to fit NOCL file formats.
2013-08-01 Fixed bug which transposed arrays and resulted in incorrect
outputs. Incorporated PWC's bugfixes from his 'get_POLCOMS_tsobc.m'
function. His notes: "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."
2013-10-14 Added support for timeseries which cross multiple months
(and thus multiple POLCOMS files).
==========================================================================

0001 function Mobj = get_POLCOMS_tsobc_NOCL(Mobj, inputConf) 0002 % Extract temperature and salinity boundary forcing information from NOC 0003 % Operation Tide Surge model output. 0004 % 0005 % function Mobj = get_POLCOMS_tsobc_NOCL(Mobj, inputConf) 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 0014 % nodes. 0015 % - Mobj.lon, Mobj.lat and/or Mobj.x, Mobj.y - node 0016 % coordinates. 0017 % - Mobj.obc_nodes - list of open boundary node inidices. 0018 % - Mobj.nObcNodes - number of nodes in each open 0019 % boundary. 0020 % inputConf = MATLAB structure which must contain: 0021 % - inputConf.polcoms_ts - location of NOC Operational 0022 % Model output containing 4D variables of temperature 0023 % (tem) and salinity (sal). They should have dimensions 0024 % (x, y, sigma, time). 0025 % - inputConf.polcoms_z - location of NOC Operational 0026 % Model output containing 4D variables of bathymetry 0027 % (XXX) and sigma layer thickness (XXX). 0028 % - inputConf.startDate - start date and time for FVCOM 0029 % model run 0030 % - inputConf.endDate - end date and time for FVCOM 0031 % model run 0032 % 0033 % OUTPUT: 0034 % Mobj = MATLAB structure in which three new fields (called temperature, 0035 % salinity and ts_time). 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_forcing_NOCL(Mobj, inputConf) 0042 % 0043 % Author(s): 0044 % Pierre Cazenave (Plymouth Marine Laboratory) 0045 % Karen Amoudry (National Oceanography Centre, Liverpool) 0046 % 0047 % PWC Revision history 0048 % 2013-01-09 First version based on the FVCOM shelf model 0049 % get_POLCOMS_forcing.m script (i.e. not a function but a plain script). 0050 % 0051 % KJA Revision history: 0052 % 2013-02-05 Adapted from PWC's script to fit NOCL file formats. 0053 % 2013-08-01 Fixed bug which transposed arrays and resulted in incorrect 0054 % outputs. Incorporated PWC's bugfixes from his 'get_POLCOMS_tsobc.m' 0055 % function. His notes: "2013-06-03 Fix some bugs in the way the open 0056 % boundary node values were stored (the order in which they were stored 0057 % did not match the order of the nodes in Casename_obc.dat). Also fix 0058 % the order of the vertically interpolated values so that FVCOM starts 0059 % at the surface instead of mirroring POLCOMS' approach (where the first 0060 % value is the seabed). The effect of these two fixes (nodes and 0061 % vertical) should match what FVCOM expects." 0062 % 2013-10-14 Added support for timeseries which cross multiple months 0063 % (and thus multiple POLCOMS files). 0064 % 0065 %========================================================================== 0066 0067 subname = 'get_POLCOMS_tsobc_NOCL'; 0068 0069 global ftbverbose; 0070 if ftbverbose 0071 fprintf('\n') 0072 fprintf(['begin : ' subname '\n']) 0073 end 0074 %% 0075 varlist = {'lon', 'lat', 'tem', 'sal', 'time'}; 0076 0077 % Create an array of hourly timesteps, ensuring the output time series is 0078 % at least as long as the FVCOM model run time. 0079 timesteps = datevec(datenum(inputConf.startDate):1/24:datenum(inputConf.endDate)+1); 0080 0081 % Find the number of months in the timeseries 0082 [months,ia]=unique(timesteps(:,1:2),'rows'); 0083 0084 for i=1:size(months,1) 0085 % Create the POLCOMS filename based on the year and month of interest 0086 polcoms_ts = [inputConf.polcoms_ts,num2str(timesteps(ia(i),1)),... 0087 '-',num2str(timesteps(ia(i),2),'%02i'),'.nc']; 0088 0089 % Get the results from the POLCOMS file 0090 nc = netcdf.open(polcoms_ts, 'NOWRITE'); 0091 0092 for var=1:numel(varlist) 0093 0094 getVar = varlist{var}; 0095 varid_pc = netcdf.inqVarID(nc, getVar); 0096 0097 data = netcdf.getVar(nc, varid_pc, 'single'); 0098 0099 if i==1 % If it's the first file, get the data 0100 pc.(getVar).data = double(data); 0101 elseif ~strcmp(getVar,'time') && ~strcmp(getVar,'lat') && ~strcmp(getVar,'lon') 0102 % if it's the second or greater file, concatecate the data 0103 pc.(getVar).data = cat(4,pc.(getVar).data,double(data)); 0104 elseif strcmp(getVar,'time') 0105 % if it's the second or greater file, concatecate the data 0106 % time is concatenated differently 0107 pc.(getVar).data = cat(1,pc.(getVar).data,double(data)); 0108 end 0109 0110 % Try to get some units (important for the calculation of MJD). 0111 try 0112 units = netcdf.getAtt(nc,varid_pc,'units'); 0113 catch 0114 units = []; 0115 end 0116 pc.(getVar).units = units; 0117 end 0118 0119 netcdf.close(nc) 0120 0121 end 0122 0123 % Get rid of data outside the time we're interested in 0124 % Convert FVCOM timestep into AMM/S12 output (seconds since 20071101:000000) 0125 AMM_time = etime(timesteps,repmat([2007,11,01,0,0,0],size(timesteps,1),1)); 0126 keep_time = (pc.time.data >= AMM_time(1) & pc.time.data <= AMM_time(end)); 0127 pc.tem.data = pc.tem.data(:,:,:,keep_time); 0128 pc.sal.data = pc.sal.data(:,:,:,keep_time); 0129 0130 % Import the POLCOMS sigma info 0131 pc = get_POLCOMS_sigma(pc,inputConf); 0132 0133 % Data format: 0134 % 0135 % pc.tem.data and pc.sal.data are x, y, sigma, time 0136 % 0137 [~, ~, nz, nt] = size(pc.tem.data); 0138 0139 % Make rectangular arrays for the nearest point lookup. 0140 [lon, lat] = meshgrid(pc.lon.data, pc.lat.data); 0141 0142 % Change the way the nodes are listed to match the order in the 0143 % Casename_obc.dat file. 0144 tmpObcNodes = Mobj.obc_nodes'; 0145 oNodes = tmpObcNodes(tmpObcNodes ~= 0)'; 0146 0147 % Find the nearest POLCOMS point to each point in the FVCOM open boundaries 0148 fvlon = Mobj.lon(oNodes); 0149 fvlat = Mobj.lat(oNodes); 0150 0151 % Number of boundary nodes 0152 nf = sum(Mobj.nObcNodes); 0153 % Number of sigma layers. 0154 fz = size(Mobj.siglayz, 2); 0155 0156 % itemp = nan(nf, nz, nt); % POLCOMS interpolated temperatures 0157 % isal = nan(nf, nz, nt); % POLCOMS interpolated salinities 0158 fvtemp = nan(nf, fz, nt); % FVCOM interpolated temperatures 0159 fvsal = nan(nf, fz, nt); % FVCOM interpolated salinities 0160 0161 if ftbverbose 0162 tic 0163 end 0164 %% 0165 for t = 1:nt 0166 % Get the current 3D array of POLCOMS results. 0167 pctemp3 = pc.tem.data(:, :, :, t); 0168 pcsalt3 = pc.sal.data(:, :, :, t); 0169 0170 % Flip the vertical layer dimension to make the POLCOMS data go from 0171 % surface to seabed to match its depth data and to match how FVCOM 0172 % works. 0173 pctemp3 = flipdim(pctemp3, 3); 0174 pcsalt3 = flipdim(pcsalt3, 3); 0175 0176 % Preallocate the intermediate results arrays. 0177 itempz = nan(nf, nz); 0178 isalz = nan(nf, nz); 0179 idepthz = nan(nf,nz); 0180 for j = 1:nz 0181 % Now extract the relevant layer from the 3D subsets. Transpose the 0182 % data to be (x, y) rather than (y, x). 0183 pctemp2 = pctemp3(:, :, j)'; 0184 pcsalt2 = pcsalt3(:, :, j)'; 0185 pcdepth2 = squeeze(pc.depth.data(:, :, j))'; 0186 % Create new arrays which will be flattened when masking (below). 0187 tpctemp2 = pctemp2; 0188 tpcsalt2 = pcsalt2; 0189 tpcdepth2 = pcdepth2; 0190 tlon = lon; 0191 tlat = lat; 0192 0193 % Create and apply a mask to remove values outside the domain. This 0194 % inevitably flattens the arrays, but it shouldn't be a problem 0195 % since we'll be searching for the closest values in such a manner 0196 % as is appropriate for an unstructured grid (i.e. we're assuming 0197 % the POLCOMS data is irregularly spaced). 0198 mask = tpcdepth2 < -20000; 0199 tpctemp2(mask) = []; 0200 tpcsalt2(mask) = []; 0201 tpcdepth2(mask) = []; 0202 % Also apply the masks to the position arrays so we can't even find 0203 % positions outside the domain, effectively meaning if a value is 0204 % outside the domain, the nearest value to the boundary node will 0205 % be used. 0206 tlon(mask) = []; 0207 tlat(mask) = []; 0208 % Reshape the arrays to allow the sort to work properly later 0209 tlon=reshape(tlon,(size(tlon,1)*size(tlon,2)),1); 0210 tlat=reshape(tlat,(size(tlat,1)*size(tlat,2)),1); 0211 0212 % Preallocate the intermediate results arrays. 0213 itempobc = nan(nf, 1); 0214 isalobc = nan(nf, 1); 0215 idepthobc = nan(nf, 1); 0216 0217 % Speed up the tightest loop with a parallelized loop. 0218 for i = 1:nf 0219 % Now we can do each position within the 2D layer. 0220 0221 fx = fvlon(i); 0222 fy = fvlat(i); 0223 0224 [~, ii] = sort(sqrt((tlon - fx).^2 + (tlat - fy).^2)); 0225 % Get the n nearest nodes from POLCOMS (more? fewer?). 0226 ixy = ii(1:16); 0227 0228 % Get the variables into static variables for the 0229 % parallelisation. 0230 plon = tlon(ixy); 0231 plat = tlat(ixy); 0232 ptemp = tpctemp2(ixy); 0233 psal = tpcsalt2(ixy); 0234 pdepth = tpcdepth2(ixy); 0235 0236 % Use a triangulation to do the horizontal interpolation. 0237 tritemp = TriScatteredInterp(plon, plat, ptemp, 'natural'); 0238 trisal = TriScatteredInterp(plon, plat, psal, 'natural'); 0239 triz = TriScatteredInterp(plon, plat, pdepth, 'natural'); 0240 itempobc(i) = tritemp(fx, fy); 0241 isalobc(i) = trisal(fx, fy); 0242 idepthobc(i) = triz(fx, fy); 0243 0244 % Check all three, though if one is NaN, they all will be. 0245 if isnan(itempobc(i)) || isnan(isalobc(i)) || isnan(idepthobc(i)) 0246 if ftbverbose 0247 warning('FVCOM boundary node at %f, %f is outside the POLCOMS domain. Setting to the closest POLCOMS value.', fx, fy) 0248 end 0249 p = 1; 0250 while isnan(tpctemp2(ii(p))) 0251 p = p+1; 0252 end 0253 itempobc(i) = tpctemp2(ii(p)); 0254 p = 1; 0255 while isnan(tpcsalt2(ii(p))) 0256 p = p+1; 0257 end 0258 isalobc(i) = tpcsalt2(ii(p)); 0259 p = 1; 0260 while isnan(tpcdepth2(ii(p))) 0261 p = p+1; 0262 end 0263 idepthobc(i) = tpcdepth2(ii(p)); 0264 end 0265 end 0266 0267 % Put the results in this intermediate array. 0268 itempz(:, j) = itempobc; 0269 isalz(:, j) = isalobc; 0270 idepthz(:,j) = idepthobc; 0271 end 0272 0273 % Now we've interpolated in space, we can interpolate the z-values 0274 % to the sigma depths. 0275 oNodes = Mobj.obc_nodes(Mobj.obc_nodes ~= 0); 0276 for zi = 1:fz 0277 0278 % Preallocate the output arrays 0279 fvtempz = nan(nf, fz); 0280 fvsalz = nan(nf, fz); 0281 0282 for pp = 1:nf 0283 % Get the FVCOM depths at this node 0284 tfz = Mobj.siglayz(oNodes(pp), :); 0285 % Now get the interpolated POLCOMS depth at this node 0286 tpz = idepthz(pp, :); 0287 0288 % To ensure we get the full vertical expression of the vertical 0289 % profiles, we need to normalise the POLCOMS and FVCOM 0290 % depths to the same range. This is because in instances where 0291 % FVCOM depths are shallower (e.g. in coastal regions), if we 0292 % don't normalise the depths, we end up truncating the vertical 0293 % profile. This approach ensures we always use the full 0294 % vertical profile, but we're potentially squeezing it into a 0295 % smaller depth. 0296 A = max(tpz); 0297 B = min(tpz); 0298 C = max(tfz); 0299 D = min(tfz); 0300 norm_tpz = (((D - C) * (tpz - A)) / (B - A)) + C; 0301 0302 % Get the temperature and salinity values for this node and 0303 % interpolate down the water column (from POLCOMS to FVCOM). 0304 % Change to 'pchip' to match PWC parent code. 0305 if ~isnan(tpz) 0306 fvtempz(pp, :) = interp1(norm_tpz, itempz(pp, :), tfz, 'pchip', 'extrap'); 0307 fvsalz(pp, :) = interp1(norm_tpz, isalz(pp, :), tfz, 'pchip', 'extrap'); 0308 else 0309 warning('Should never see this... ') % because we test for NaNs when fetching the values. 0310 warning('FVCOM boundary node at %f, %f is outside the POLCOMS domain. Skipping.', fvlon(pp), fvlat(pp)) 0311 continue 0312 end 0313 end 0314 end 0315 0316 % The horizontally-interpolated values in the final results array. 0317 % itemp(:, :, t) = itempz; 0318 % isal(:, :, t) = isalz; 0319 % The horizontally- and vertically-interpolated values in the final 0320 % FVCOM results array. 0321 fvtemp(:, :, t) = fvtempz; 0322 fvsal(:, :, t) = fvsalz; 0323 end 0324 0325 if ftbverbose 0326 toc 0327 end 0328 0329 %% 0330 % Convert NOC model output temperatures from Kelvin to Celsius 0331 fvtemp = fvtemp - 273.15; 0332 0333 % Timeshift to match the expected FVCOM input times. The temperature and 0334 % salinity values are a daily average (midnight to midnight). They are 0335 % given a timestamp of the middle time in this period (noon). Each day's 0336 % value applies to the whole day (e.g. from midnight to midnight). 0337 % Therefore, we can shift the time to midnight at the start of the relevant 0338 % day and apply the value to the whole day. This makes FVCOM (and me) 0339 % happy. 0340 Mobj.ts_times = greg2mjulian(2007, 11, 0, 0, 0, 0) + ... % Convert NOC model reference time to MJD 0341 ((pc.time.data(keep_time)+(12*3600)) / 3600 / 24); % Add NOC model output time and 12 hour offset 0342 0343 Mobj.temperature = fvtemp; 0344 Mobj.salt = fvsal; 0345 0346 if ftbverbose 0347 fprintf(['end : ' subname '\n']) 0348 end