


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 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 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_tsobc(Mobj, ts, depth)
Author(s):
Pierre Cazenave (Plymouth Marine Laboratory)
Revision history
2013-02-07 First version based on get_POLCOMS_tsobc.m.
==========================================================================

0001 function [Mobj,pc] = get_POLCOMS_tsobc_gcoms(Mobj, ts,inputConf) 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 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 three new fields (called temperature, 0036 % salinity and ts_time). temperature and salinity have sizes 0037 % (sum(Mobj.nObcNodes), sigma, time). The time dimension is 0038 % determined based on the input NetCDF file. The ts_time variable 0039 % is just the input file times in Modified Julian Day. 0040 % 0041 % EXAMPLE USAGE 0042 % Mobj = get_POLCOMS_tsobc(Mobj, ts, depth) 0043 % 0044 % Author(s): 0045 % Pierre Cazenave (Plymouth Marine Laboratory) 0046 % 0047 % Revision history 0048 % 2013-02-07 First version based on get_POLCOMS_tsobc.m. 0049 % 0050 %========================================================================== 0051 wasOpened = false; 0052 if license('test', 'Distrib_Computing_Toolbox') 0053 % We have the Parallel Computing Toolbox, so launch a bunch of workers. 0054 if matlabpool('size') == 0 0055 % Force pool to be local in case we have remote pools available. 0056 matlabpool open local 0057 wasOpened = true; 0058 end 0059 end 0060 0061 subname = 'get_POLCOMS_tsobc_gcoms'; 0062 0063 global ftbverbose; 0064 if ftbverbose 0065 fprintf('\n') 0066 fprintf(['begin : ' subname '\n']) 0067 fprintf(['Expecting two files : ' ts{1} '\n' ts{2} '\n']) 0068 end 0069 0070 varlist = {'lon', 'lat', 'ETW', 'x1X', 'time'}; 0071 0072 varlistdmeaneco3d = {'depth', 'pdepth','N3n','N1p','N5s','N4n'}; 0073 % varlistdmeaneco3d = {'depth', 'pdepth'}; 0074 % 4D variables 0075 varnames_time={'ETW', 'x1X','N3n','N1p','N5s','N4n' } 0076 varnames_fvcom={'temperature', 'salt','N3n','N1p','N5s','N4n'} 0077 0078 % Data format: 0079 % 0080 % pc.ETWD.data and pc.x1XD.data are y, x, sigma, time 0081 % 0082 dd = 1; 0083 catstr = 'pc=catstruct('; 0084 for ff = 1:2:length(ts)-1 0085 pc = get_POLCOMS_netCDF(ts{ff}, varlist); 0086 % Convert the current times to Modified Julian Day (this is a bit ugly). 0087 pc.time.all = strtrim(regexp(pc.time.units, 'since', 'split')); 0088 pc.time.datetime = strtrim(regexp(pc.time.all{end}, ' ', 'split')); 0089 pc.time.ymd = str2double(strtrim(regexp(pc.time.datetime{1}, '-', 'split'))); 0090 pc.time.hms = str2double(strtrim(regexp(datestr(datenum(pc.time.ymd),'HH:MM:SS'), ':', 'split'))); 0091 pc.time.data = datenum(... 0092 pc.time.ymd(1), ... 0093 pc.time.ymd(2), ... 0094 pc.time.ymd(3), ... 0095 pc.time.hms(1), ... 0096 pc.time.hms(2), ... 0097 pc.time.hms(3)) + (pc.time.data / 3600 / 24); 0098 0099 pcdmean = get_POLCOMS_netCDF(ts{ff+1}, varlistdmeaneco3d); 0100 eval(['dump',num2str(dd),'=catstruct(pc,pcdmean);']) 0101 0102 catstr = [catstr,'dump',num2str(dd),',']; 0103 dd = dd+1; 0104 end 0105 fieldN=fieldnames(dump1); 0106 for fdn=1:length(fieldN) 0107 catstr = ''; 0108 for fn=1:dd-1 0109 % -- time is always the last field!! 0110 catstr = [catstr,'dump',num2str(fn),'.',(fieldN{fdn}),'.data,']; 0111 end 0112 eval(['[nx,ny,nz,nt]=size(dump1.',(fieldN{fdn}),'.data);']) 0113 ND=find([nx,ny,nz,nt]==length(dump1.time.data)) 0114 if isempty(ND) 0115 % no time dimension... do not concatenate, use one of the dumps 0116 eval(['pc.',(fieldN{fdn}),'.data=dump1.',(fieldN{fdn}),'.data;']) 0117 else 0118 disp( ['pc.',(fieldN{fdn}),'.data=cat(ND,',catstr(1:end-1),');']) 0119 eval(['pc.',(fieldN{fdn}),'.data=cat(ND,',catstr(1:end-1),');']) 0120 end 0121 end 0122 clear pcdmean dump* 0123 % find time limit for range of interest 0124 d0=datenum(inputConf.startDate); 0125 d1=datenum(inputConf.endDate); 0126 % Include one day before and one day after the required time range 0127 d0=d0-1;d1=d1+1; 0128 igood=find(pc.time.data >= d0 & pc.time.data <= d1); 0129 0130 [yyyy,mm,dd,HH,MM,SS]=datevec(pc.time.data(igood)) 0131 Mobj.ts_times = greg2mjulian(yyyy,mm,dd,HH,MM,SS); 0132 0133 0134 for fdn=1:length(varnames_time) 0135 pc.(varnames_time{fdn}).data=pc.(varnames_time{fdn}).data(:,:,:,igood); 0136 end 0137 % progress to interpolation 0138 0139 [~, ~, nz, nt] = size(pc.ETW.data); 0140 0141 % Make rectangular arrays for the nearest point lookup. 0142 [lon, lat] = meshgrid(pc.lon.data, pc.lat.data); 0143 0144 fvlon = Mobj.lon(Mobj.obc_nodes(Mobj.obc_nodes ~= 0)); 0145 fvlat = Mobj.lat(Mobj.obc_nodes(Mobj.obc_nodes ~= 0)); 0146 0147 % Number of boundary nodes 0148 nf = sum(Mobj.nObcNodes); 0149 % Number of sigma layers. 0150 fz = size(Mobj.siglayz, 2); 0151 0152 fvtemp = nan(nf, fz, nt); % FVCOM interpolated temperatures 0153 0154 if ftbverbose 0155 tic 0156 end 0157 for fdn=1:length(varnames_time) 0158 0159 for t = 1:nt 0160 if ftbverbose 0161 fprintf('%s : %i of %i timesteps... ', subname, t, nt) 0162 end 0163 % Get the current 3D array of PML POLCOMS-ERSEM results. 0164 pctemp3 = pc.(varnames_time{fdn}).data(:, :, :, t); 0165 % pcsalt3 = pc.x1X.data(:, :, :, t); 0166 0167 % Preallocate the intermediate results arrays. 0168 itempz = nan(nf, nz); 0169 idepthz = nan(nf, nz); 0170 0171 for j = 1:nz 0172 % Now extract the relevant layer from the 3D subsets. Transpose the 0173 % data to be (x, y) rather than (y, x). 0174 pctemp2 = pctemp3(:, :, j)'; 0175 pcdepth2 = squeeze(pc.depth.data(:, :, j, t))'; 0176 0177 % Create new arrays which will be flattened when masking (below). 0178 tpctemp2 = pctemp2; 0179 tpcdepth2 = pcdepth2; 0180 tlon = lon; 0181 tlat = lat; 0182 0183 % Create and apply a mask to remove values outside the domain. This 0184 % inevitably flattens the arrays, but it shouldn't be a problem 0185 % since we'll be searching for the closest values in such a manner 0186 % as is appropriate for an unstructured grid (i.e. we're assuming 0187 % the PML POLCOMS-ERSEM data is irregularly spaced). 0188 mask = tpcdepth2 < -10e+10; 0189 tpctemp2(mask) = []; 0190 tpcdepth2(mask) = []; 0191 % Also apply the masks to the position arrays so we can't even find 0192 % positions outside the domain, effectively meaning if a value is 0193 % outside the domain, the nearest value to the boundary node will 0194 % be used. 0195 tlon(mask) = []; 0196 tlat(mask) = []; 0197 0198 % Preallocate the intermediate results arrays. 0199 itempobc = nan(nf, 1); 0200 idepthobc = nan(nf, 1); 0201 0202 % Speed up the tightest loop with a parallelized loop. 0203 parfor i = 1:nf 0204 % for i = 1:nf 0205 % Now we can do each position within the 2D layer. 0206 0207 fx = fvlon(i); 0208 fy = fvlat(i); 0209 0210 [~, ii] = sort(sqrt((tlon - fx).^2 + (tlat - fy).^2)); 0211 % Get the n nearest nodes from PML POLCOMS-ERSEM data (more? 0212 % fewer?). 0213 ixy = ii(1:16); 0214 0215 % Get the variables into static variables for the 0216 % parallelisation. 0217 plon = tlon(ixy); 0218 plat = tlat(ixy); 0219 ptemp = tpctemp2(ixy); 0220 pdepth = tpcdepth2(ixy); 0221 0222 % Use a triangulation to do the horizontal interpolation. 0223 tritemp = TriScatteredInterp(plon', plat', ptemp', 'linear'); 0224 triz = TriScatteredInterp(plon', plat', pdepth', 'linear'); 0225 itempobc(i) = tritemp(fx, fy); 0226 idepthobc(i) = triz(fx, fy); 0227 0228 % Check all three, though if one is NaN, they all will be. 0229 if isnan(itempobc(i)) || isnan(idepthobc(i)) 0230 warning('FVCOM boundary node at %f, %f is outside the PML POLCOMS-ERSEM domain. Setting to the closest PML POLCOMS-ERSEM value.', fx, fy) 0231 itempobc(i) = tpctemp2(ii(1)); 0232 idepthobc(i) = tpcdepth2(ii(1)); 0233 end 0234 end 0235 0236 % Put the results in this intermediate array. 0237 itempz(:, j) = itempobc; 0238 idepthz(:, j) = idepthobc; 0239 end 0240 0241 % Now we've interpolated in space, we can interpolate the z-values 0242 % to the sigma depths. 0243 oNodes = Mobj.obc_nodes(Mobj.obc_nodes ~= 0); 0244 0245 % Preallocate the output arrays 0246 fvtempz = nan(nf, fz); 0247 0248 for pp = 1:nf 0249 % Get the FVCOM depths at this node 0250 tfz = Mobj.siglayz(oNodes(pp), :); 0251 % Now get the interpolated PML POLCOMS-ERSEM depth at this node 0252 tpz = idepthz(pp, :); 0253 0254 % To ensure we get the full vertical expression of the vertical 0255 % profiles, we need to normalise the POLCOMS-ERSEM and FVCOM 0256 % depths to the same range. This is because in instances where 0257 % FVCOM depths are shallower (e.g. in coastal regions), if we 0258 % don't normalise the depths, we end up truncating the vertical 0259 % profile. This approach ensures we always use the full 0260 % vertical profile, but we're potentially squeezing it into a 0261 % smaller depth. 0262 A = max(tpz); 0263 B = min(tpz); 0264 C = max(tfz); 0265 D = min(tfz); 0266 norm_tpz = (((D - C) * (tpz - A)) / (B - A)) + C; 0267 0268 % Get the temperature and salinity values for this node and 0269 % interpolate down the water column (from PML POLCOMS-ERSEM to 0270 % FVCOM). I had originally planned to use csaps for the 0271 % vertical interplation/subsampling at each location. However, 0272 % the demo of csaps in the MATLAB documentation makes the 0273 % interpolation look horrible (shaving off extremes). interp1 0274 % provides a range of interpolation schemes, of which pchip 0275 % seems to do a decent job of the interpolation (at least 0276 % qualitatively). 0277 if ~isnan(tpz) 0278 fvtempz(pp, :) = interp1(norm_tpz, itempz(pp, :), tfz, 'linear', 'extrap'); 0279 else 0280 warning('Should never see this... ') % because we test for NaNs when fetching the values. 0281 warning('FVCOM boundary node at %f, %f is outside the PML POLCOMS-ERSEM domain. Skipping.', fvlon(pp), fvlat(pp)) 0282 continue 0283 end 0284 end 0285 0286 % The horizontally- and vertically-interpolated values in the final 0287 % FVCOM results array. 0288 fvtemp(:, :, t) = fvtempz; 0289 0290 if ftbverbose 0291 fprintf('done.\n') 0292 end 0293 end 0294 Mobj.(varnames_fvcom{fdn}) = fvtemp; 0295 end 0296 0297 if ftbverbose 0298 toc 0299 end 0300 0301 0302 0303 if ftbverbose 0304 fprintf(['end : ' subname '\n']) 0305 end 0306 if wasOpened 0307 matlabpool close 0308 end 0309 return