


Extract river data from the supplied river positions for the FVCOM
grid in Mobj from the NEMO data
get_NEMO_rivers(Mobj, dist_thresh)
DESCRIPTION:
For the positions in Mobj.rivers.positions, find the nearest
unstructured grid node and extract the river discharge from
Mobj.rivers.river_flux. To correct the flux values from the
NEMO data, we need the grid size data (in
Mobj.rivers.river_coordinates). The river positions must fall
within the specified distance (dist_thresh).
INPUT:
Mobj - MATLAB mesh object containing:
* have_lonlat - boolean to check for spherical coordinates.
* lon, lat - positions for the unstructured grid.
* tri - triangulation table for the unstructured grid.
* nVerts - number of nodes in the grid.
* read_obc_nodes - open boundary node IDs.
* rivers - river data struct with the following fields:
- positions - river positions in lon, lat.
- names - list of river names
- river_flux - path to the NEMO netCDF data file. Must contain
the area of each grid element for conversion of the fluxes.
dist_thresh - maximum distance away from a river node beyond
which the search for an FVCOM node is abandoned. Units in degrees.
The following keyword-argument pairs are also valid:
'model_year' - [optional] when giving climatology, a year must be
specified so that the time series can be anchored in time. The
returned time series will be 3 years long centred on the specified
year. Discharges will be repeated for the two additional years.
'dump_positions' - [optional] dump the NEMO river positions to the
specified text file.
'alternate_positions' - [optional] read in a CSV file with alternate
positions for the NEMO rivers. Supply a comma separated file with
the new values and then the old values as (xnew,ynew,xold,yold).
This is useful if you've manually moved the NEMO rivers onto more
realistic locations.
'remove_baltic' - [optional] remove the Baltic river inputs. Set to
true to remove; defaults to false.
OUTPUT:
Mobj.river_flux - volume flux at the nodes within the model domain.
Mobj.river_nh4 - ammonia
Mobj.river_no3 - notrate
Mobj.river_o - oxygen
Mobj.river_p - phosphate
Mobj.river_sio3 - silicate
Mobj.river_dic - dissolved inorganic carbon
Mobj.river_totalk - total alkalinity
Mobj.river_bioalk - bio-alkalinity
Mobj.river_nodes - node IDs for the rivers.
Mobj.river_names - river names which fall within the model domain. For
rivers where the discharge has been summed, the name is compoud,
with each contributing name separated by a hyphen (-).
Mobj.river_time - Modified Julian Day time series for the river
discharge data.
Mobj.river_nemo_location - river locations (NEMO positions).
EXAMPLE USAGE:
Mobj = get_NEMO_rivers(Mobj, 0.15)
To extract the NEMO river locations to file:
Mobj = get_NEMO_rivers(Mobj, 0.15, 'dump_positions', '/tmp/nemo.txt');
Author(s):
Pierre Cazenave (Plymouth Marine Laboratory)
Revision history:
2016-03-02 - First version based on get_nemo_rivers.m.
2016-05-03 - Add option to dump NEMO river locations to text file.
2016-06-06 - Remove temperature and salinity from the description as
they're really the Baltic Sea inputs only. Also read in the grid area
from the new format ERSEM file (variable dA). Add total alkalinity
to the outputs.
2016-08-10 - Add new option to use a file specifying alternative river
positions.
2017-10-04 - Add new option to remove the Baltic Sea inputs from the
river data.
==========================================================================

0001 function Mobj = get_NEMO_rivers(Mobj, dist_thresh, varargin) 0002 % Extract river data from the supplied river positions for the FVCOM 0003 % grid in Mobj from the NEMO data 0004 % 0005 % get_NEMO_rivers(Mobj, dist_thresh) 0006 % 0007 % DESCRIPTION: 0008 % For the positions in Mobj.rivers.positions, find the nearest 0009 % unstructured grid node and extract the river discharge from 0010 % Mobj.rivers.river_flux. To correct the flux values from the 0011 % NEMO data, we need the grid size data (in 0012 % Mobj.rivers.river_coordinates). The river positions must fall 0013 % within the specified distance (dist_thresh). 0014 % 0015 % INPUT: 0016 % Mobj - MATLAB mesh object containing: 0017 % * have_lonlat - boolean to check for spherical coordinates. 0018 % * lon, lat - positions for the unstructured grid. 0019 % * tri - triangulation table for the unstructured grid. 0020 % * nVerts - number of nodes in the grid. 0021 % * read_obc_nodes - open boundary node IDs. 0022 % * rivers - river data struct with the following fields: 0023 % - positions - river positions in lon, lat. 0024 % - names - list of river names 0025 % - river_flux - path to the NEMO netCDF data file. Must contain 0026 % the area of each grid element for conversion of the fluxes. 0027 % dist_thresh - maximum distance away from a river node beyond 0028 % which the search for an FVCOM node is abandoned. Units in degrees. 0029 % The following keyword-argument pairs are also valid: 0030 % 'model_year' - [optional] when giving climatology, a year must be 0031 % specified so that the time series can be anchored in time. The 0032 % returned time series will be 3 years long centred on the specified 0033 % year. Discharges will be repeated for the two additional years. 0034 % 'dump_positions' - [optional] dump the NEMO river positions to the 0035 % specified text file. 0036 % 'alternate_positions' - [optional] read in a CSV file with alternate 0037 % positions for the NEMO rivers. Supply a comma separated file with 0038 % the new values and then the old values as (xnew,ynew,xold,yold). 0039 % This is useful if you've manually moved the NEMO rivers onto more 0040 % realistic locations. 0041 % 'remove_baltic' - [optional] remove the Baltic river inputs. Set to 0042 % true to remove; defaults to false. 0043 % 0044 % OUTPUT: 0045 % Mobj.river_flux - volume flux at the nodes within the model domain. 0046 % Mobj.river_nh4 - ammonia 0047 % Mobj.river_no3 - notrate 0048 % Mobj.river_o - oxygen 0049 % Mobj.river_p - phosphate 0050 % Mobj.river_sio3 - silicate 0051 % Mobj.river_dic - dissolved inorganic carbon 0052 % Mobj.river_totalk - total alkalinity 0053 % Mobj.river_bioalk - bio-alkalinity 0054 % Mobj.river_nodes - node IDs for the rivers. 0055 % Mobj.river_names - river names which fall within the model domain. For 0056 % rivers where the discharge has been summed, the name is compoud, 0057 % with each contributing name separated by a hyphen (-). 0058 % Mobj.river_time - Modified Julian Day time series for the river 0059 % discharge data. 0060 % Mobj.river_nemo_location - river locations (NEMO positions). 0061 % 0062 % EXAMPLE USAGE: 0063 % Mobj = get_NEMO_rivers(Mobj, 0.15) 0064 % 0065 % To extract the NEMO river locations to file: 0066 % 0067 % Mobj = get_NEMO_rivers(Mobj, 0.15, 'dump_positions', '/tmp/nemo.txt'); 0068 % 0069 % Author(s): 0070 % Pierre Cazenave (Plymouth Marine Laboratory) 0071 % 0072 % Revision history: 0073 % 2016-03-02 - First version based on get_nemo_rivers.m. 0074 % 2016-05-03 - Add option to dump NEMO river locations to text file. 0075 % 2016-06-06 - Remove temperature and salinity from the description as 0076 % they're really the Baltic Sea inputs only. Also read in the grid area 0077 % from the new format ERSEM file (variable dA). Add total alkalinity 0078 % to the outputs. 0079 % 2016-08-10 - Add new option to use a file specifying alternative river 0080 % positions. 0081 % 2017-10-04 - Add new option to remove the Baltic Sea inputs from the 0082 % river data. 0083 % 0084 %========================================================================== 0085 0086 [~, subname] = fileparts(mfilename('fullpath')); 0087 0088 global ftbverbose 0089 if ftbverbose 0090 fprintf('\nbegin : %s\n', subname) 0091 end 0092 0093 % Check inputs 0094 if ~Mobj.have_lonlat 0095 error(['Require unstructured grid positions in lon/lat format ', ... 0096 'to compare against supplied river positions.']) 0097 end 0098 0099 yr = []; 0100 0101 % If we have only three arguments, we have to assume we've been given a 0102 % year for the climatology. Otherwise, we need to read the arguments based 0103 % on keyword-value pairs. Really, we want keyword-value pairs all the time, 0104 % so silently work when given three arguments and don't mention it in the 0105 % help. This is going to bite me at some point in the future, I'm sure. 0106 dump_positions = false; 0107 alt_positions = false; 0108 drop_baltic = false; 0109 if nargin == 3 0110 yr = varargin{1}; 0111 elseif nargin > 3 0112 for aa = 1:2:length(varargin) 0113 switch varargin{aa} 0114 case 'model_year' 0115 yr = varargin{aa + 1}; 0116 case 'dump_positions' 0117 dump_positions = true; 0118 position_file = varargin{aa + 1}; 0119 case 'alternate_positions' 0120 alt_positions = true; 0121 alternate_file = varargin{aa + 1}; 0122 case 'remove_baltic' 0123 drop_baltic = true; 0124 end 0125 end 0126 end 0127 0128 if ~isempty(yr) && ~isnumeric(yr) 0129 error('Trying to do climatology, but don''t have an anchor year. Supply one via the ''model_year'' keyword-value pair.') 0130 end 0131 0132 % Load the NEMO data. 0133 nemo.lon = ncread(Mobj.rivers.river_flux, 'x'); 0134 nemo.lat = ncread(Mobj.rivers.river_flux, 'y'); 0135 [nemo.LON, nemo.LAT] = meshgrid(nemo.lon, nemo.lat); 0136 nemo.time = ncread(Mobj.rivers.river_flux, 'time_counter'); 0137 nemo.flux = ncread(Mobj.rivers.river_flux, 'rorunoff'); 0138 nemo.nh4 = ncread(Mobj.rivers.river_flux, 'ronh4'); 0139 nemo.no3 = ncread(Mobj.rivers.river_flux, 'rono3'); 0140 nemo.o = ncread(Mobj.rivers.river_flux, 'roo'); 0141 nemo.p = ncread(Mobj.rivers.river_flux, 'rop'); 0142 nemo.sio3 = ncread(Mobj.rivers.river_flux, 'rosio2'); 0143 nemo.dic = ncread(Mobj.rivers.river_flux, 'rodic'); 0144 nemo.alt = ncread(Mobj.rivers.river_flux, 'rototalk'); 0145 nemo.bioalk = ncread(Mobj.rivers.river_flux, 'robioalk'); 0146 % nemo.temp = ncread(Mobj.rivers.river_flux, 'rotemper'); 0147 % nemo.salt = ncread(Mobj.rivers.river_flux, 'rosaline'); 0148 0149 if drop_baltic 0150 % Set the data at the Baltic indices to zero before we get too carried 0151 % away. Find the indices based on positions rather than hard-coding 0152 % them so we can still do this for newer NEMO river inputs (assuming 0153 % they still use the same approach for the Baltic). 0154 baltic.lon = [10.7777, 12.5555]; 0155 baltic.lat = [55.5998, 56.1331]; 0156 names = fieldnames(nemo); 0157 for pp = 1:length(baltic.lon) 0158 [~, lon_idx] = min(abs(nemo.lon - baltic.lon(pp))); 0159 [~, lat_idx] = min(abs(nemo.lat - baltic.lat(pp))); 0160 for n = 1:length(names) 0161 switch names{n} 0162 case {'lon', 'lat', 'LON', 'LAT', 'time'} 0163 continue 0164 end 0165 0166 % Drop the Baltic data (replace with zeros to match the other 0167 % non-river data in the netCDF). 0168 nemo.(names{n})(lon_idx, lat_idx, :) = 0; 0169 end 0170 end 0171 clearvars lon_idx lat_idx baltic names n 0172 end 0173 0174 % Now get the NEMO grid data. 0175 % nemo.e1t = ncread(Mobj.rivers.river_coordinates, 'e1t'); 0176 % nemo.e2t = ncread(Mobj.rivers.river_coordinates, 'e2t'); 0177 % Calculate the area for all elements in the grid. 0178 % nemo.area = nemo.e1t .* nemo.e2t; 0179 % Just use the new area variable instead of needing two files. 0180 nemo.area = ncread(Mobj.rivers.river_flux, 'dA'); 0181 0182 % NEMO does conversions to ERSEM units internally. Whilst this is easy in 0183 % some ways, it's not particularly transparent. So, instead, we'll do all 0184 % the conversions up front and then the data that get loaded into ERSEM are 0185 % already in the correct units. To summarise those conversions, we have: 0186 % 0187 % nutrients (nh4, no3, o, p, sio3) from grams/l to mmol/m^3 0188 % flux from kg/m^{2}/s to m^{3}/s (divide by freshwater density) 0189 % dic - no change whatsoever. 0190 0191 [~, ~, nt] = size(nemo.flux); 0192 0193 % Flux in NEMO is specified in kg/m^{2}/s. FVCOM wants m^{3}/s. Divide by 0194 % freshwater density to get m/s and then multiply by the area of each 0195 % element to get flux. 0196 nemo.flux = nemo.flux / 1000; 0197 % Now multiply by the relevant area to (finally!) get to m^{3}/s. 0198 nemo.flux = nemo.flux .* repmat(nemo.area, 1, 1, nt); 0199 % Set zero values to a very small number instead. 0200 tmp = nemo.flux; 0201 tmp(tmp==0) = 1E-8; 0202 0203 % Convert units from grams to millimoles where appropriate. 0204 nemo.nh4 = (nemo.nh4 / 14) *1000 ./ tmp; %g/s to mmol/m3 0205 nemo.no3 = (nemo.no3 / 14) *1000 ./ tmp;%g/s to mmol/m3 0206 nemo.o = (nemo.o / 16) *1000 ./ tmp; % Nemo oxygen concentrations are for O rather than O2 0207 nemo.p = (nemo.p / 35.5)*1000 ./ tmp;%g/s to mmol/m3 0208 nemo.sio3 = (nemo.sio3 / 28) *1000./ tmp;%g/2 to mmol/m3 0209 nemo.bioalk = nemo.bioalk./ tmp / 1000; % bioalk is in umol/s need umol/kg 0210 nemo.dic = nemo.dic./12./ tmp *1000; % dic is in gC/s need mmol/m3 0211 % total alkalinity is already in umol/Kg as expected by ERSEM. 0212 clear tmp 0213 0214 % Now we've got the data, use the flux data to find the indices of the 0215 % rivers in the arrays and extract those as time series in a format 0216 % suitable for writing out with write_FVCOM_river. 0217 mask = sum(nemo.flux, 3) ~= 0; 0218 % Get the indices so we can extract the time series. 0219 [mirow, micol] = ind2sub(size(mask), find(mask == true)); 0220 nr = length(mirow); 0221 0222 % Now do all the data. 0223 names = fieldnames(nemo); 0224 for n = 1:length(names) 0225 switch names{n} 0226 case {'lon', 'lat', 'LON', 'LAT', 'time'} 0227 continue 0228 end 0229 0230 % Preallocate and then fill with the time series for each valid river 0231 % location. 0232 nemo.rivers.(names{n}) = nan(nt, nr); 0233 for r = 1:nr 0234 nemo.rivers.(names{n})(:, r) = squeeze(nemo.(names{n})(mirow(r), micol(r), :)); 0235 end 0236 end 0237 0238 % Do the positions in the same way (with a loop) to avoid having to figure 0239 % out if using the mask collapses the arrays in teh same way as the find 0240 % did. 0241 nemo.rivers.positions = nan(nr, 2); 0242 nemo.rivers.names = {}; 0243 for r = 1:nr 0244 % For some reason, the output of meshgrid is the wrong way around, so 0245 % use the row and column indices the wrong way around too. 0246 nemo.rivers.positions(r, :) = [nemo.LON(micol(r), mirow(r)), ... 0247 nemo.LAT(micol(r), mirow(r))]; 0248 % We don't have sensible names for the NEMO rivers, so use the position 0249 % instead. Perhaps there exists a list of the names somewhere which I 0250 % could use. Yuri has asked Sarah Wakelin to see if she's got that 0251 % list. 0252 nemo.rivers.names{r, 1} = sprintf('river_%.6f-%.6f', ... 0253 nemo.rivers.positions(r, :)); 0254 end 0255 0256 if alt_positions 0257 % We've been given a file with alternate river positions in. Use that 0258 % for the river locations instead. 0259 f = fopen(alternate_file, 'r'); 0260 assert(f > 1, 'Error opening river locations file (%s)', alternate_file) 0261 new_positions = textscan(f, '%f%f%f%f%[^\n\r]', ... 0262 'Delimiter', ',', ... 0263 'HeaderLines', 1, ... 0264 'ReturnOnError', false); 0265 new_x = new_positions{1}; 0266 new_y = new_positions{2}; 0267 original_x = new_positions{3}; 0268 original_y = new_positions{4}; 0269 % Although in principle the "old" positions should be identical, odd 0270 % little precision issues can creep in and cause all sorts of problems. 0271 % So, we'll instead search for the nearest location and use that 0272 % instead. 0273 for ri = 1:length(original_x) 0274 xdiffs = original_x(ri) - nemo.rivers.positions(:, 1); 0275 ydiffs = original_y(ri) - nemo.rivers.positions(:, 2); 0276 [~, index] = min(sqrt(xdiffs.^2 - ydiffs.^2)); 0277 % Now we know which river we're updating, just replace the NEMO 0278 % positions read in from netCDF with our CSV versions. 0279 nemo.rivers.positions(index, :) = [new_x(ri), new_y(ri)]; 0280 end 0281 clear xdiffs ydiffs ri index 0282 end 0283 0284 % Separate the inputs into separate arrays. 0285 nemo_name = nemo.rivers.names; 0286 nemo_xy = nemo.rivers.positions; 0287 0288 if dump_positions 0289 % Add a header for GIS 0290 dlmwrite(position_file, 'lonDD,latDD', 'delimiter', ''); 0291 dlmwrite(position_file, nemo_xy, 'precision', '%0.6f', '-append'); 0292 end 0293 0294 fv_nr = length(nemo_name); 0295 0296 % Check each location in the NEMO positions against the grid in Mobj and 0297 % for the indices within the dist_thresh, load and extract the relevant 0298 % time series data. 0299 0300 vc = 0; % valid FVCOM boundary node counter 0301 0302 % We need to find the unstructured grid boundary nodes and exclude the open 0303 % boundary nodes from them. This will be our list of potential candidates 0304 % for the river nodes (i.e. the land coastline). 0305 [~, ~, ~, bnd] = connectivity([Mobj.lon, Mobj.lat], Mobj.tri); 0306 boundary_nodes = 1:Mobj.nVerts; 0307 boundary_nodes = boundary_nodes(bnd); 0308 coast_nodes = boundary_nodes(~ismember(boundary_nodes, [Mobj.read_obc_nodes{:}])); 0309 tlon = Mobj.lon(coast_nodes); 0310 tlat = Mobj.lat(coast_nodes); 0311 0312 fv_obc = nan; 0313 fv_names = cell(0); 0314 fv_location = [nan, nan]; 0315 % Initialise the flow array with a 366 day long time series of nans. This 0316 % array will be appended to (unless all rivers are outside the domain). 0317 % Only do this if we're doing climatology (signified by a non-empty year). 0318 skipped = 0; 0319 0320 % Preallocate all the arrays. 0321 for n = 1:length(names) 0322 switch names{n} 0323 case {'lon', 'lat', 'LON', 'LAT', 'time'} 0324 continue 0325 end 0326 fv.(names{n}) = nan(nt, nr); 0327 end 0328 0329 for ff = 1:fv_nr 0330 % Find the coastline node closest to this river. Don't bother with sqrt 0331 % for the distance threshold, instead just square the distance when 0332 % doing the comparison. This should increase performance, although 0333 % probably only marginally. 0334 fv_dist = (nemo_xy(ff, 1) - tlon).^2 + ... 0335 (nemo_xy(ff, 2) - tlat).^2; 0336 [c, idx] = min(fv_dist); 0337 if c > dist_thresh^2 0338 skipped = skipped + 1; 0339 continue 0340 else 0341 if ftbverbose 0342 fprintf('candidate river %s found (%f, %f)... ', nemo_name{ff}, nemo_xy(ff, 1), nemo_xy(ff, 2)) 0343 end 0344 end 0345 0346 vc = vc + 1; 0347 0348 % We need to make sure the element in which this node occurs does not 0349 % have two land boundaries (otherwise the model just fills up that 0350 % element because that element will always have a zero velocity). 0351 0352 % Find the other nodes which are joined to the node we've just found. 0353 % We don't need the column to get the other nodes in the element, only 0354 % the row is required. 0355 if ~isempty(coast_nodes(idx)) 0356 [row, ~] = find(Mobj.tri == coast_nodes(idx)); 0357 0358 if length(row) == 1 0359 % This is a bad node because it is a part of only one element. The 0360 % rivers need two adjacent elements to work reliably (?). So, we 0361 % need to repeat the process above until we find a node that's 0362 % connected to two elements. We'll try the other nodes in the 0363 % current element before searching the rest of the coastline (which 0364 % is computationally expensive). 0365 0366 % Remove the current node index from the list of candidates (i.e. 0367 % leave only the two other nodes in the element). 0368 mask = Mobj.tri(row, :) ~= coast_nodes(idx); 0369 n_tri = Mobj.tri(row, mask); 0370 0371 % Remove values which aren't coastline values (we don't want to set 0372 % the river node to an open water node). 0373 n_tri = intersect(n_tri, coast_nodes); 0374 0375 % Of the remaining nodes in the element, find the closest one to 0376 % the original river location (in fvcom_xy). 0377 [~, n_idx] = sort(sqrt( ... 0378 (nemo_xy(ff, 1) - Mobj.lon(n_tri)).^2 ... 0379 + (nemo_xy(ff, 2) - Mobj.lat(n_tri)).^2)); 0380 0381 [row_2, ~] = find(Mobj.tri == n_tri(n_idx(1))); 0382 if length(n_idx) > 1 0383 [row_3, ~] = find(Mobj.tri == n_tri(n_idx(2))); 0384 end 0385 % Closest first 0386 if length(row_2) > 1 0387 idx = find(coast_nodes == n_tri(n_idx(1))); 0388 % The other one (only if we have more than one node to consider). 0389 elseif length(n_idx) > 1 && length(row_3) > 1 0390 idx = find(coast_nodes == n_tri(n_idx(2))); 0391 % OK, we need to search across all the other coastline nodes. 0392 else 0393 % TODO: Implement a search of all the other coastline nodes. 0394 % My testing indicates that we never get here (at least for the 0395 % grids I've tested). I'd be interested to see the mesh which 0396 % does get here... 0397 continue 0398 end 0399 if ftbverbose 0400 fprintf('alternate node ') 0401 end 0402 end 0403 0404 % Add it to the list of valid rivers 0405 fv_obc(vc) = coast_nodes(idx); 0406 fv_names{vc} = nemo.rivers.names{ff}; 0407 fv_location(vc, :) = nemo.rivers.positions(ff, :); 0408 0409 % Add the current river data to the relevant arrays. 0410 for n = 1:length(names) 0411 switch names{n} 0412 case {'lon', 'lat', 'LON', 'LAT', 'time'} 0413 continue 0414 end 0415 fv.(names{n})(:, vc) = nemo.rivers.(names{n})(:, ff); 0416 end 0417 0418 if ftbverbose 0419 fprintf('added (%f, %f)\n', Mobj.lon(fv_obc(vc)), Mobj.lat(fv_obc(vc))) 0420 end 0421 end 0422 end 0423 0424 % Trim the data arrays for the rivers we've extracted since we preallocated 0425 % the array assuming we'd use all the rivers. 0426 for n = 1:length(names) 0427 switch names{n} 0428 case {'lon', 'lat', 'LON', 'LAT', 'time'} 0429 continue 0430 end 0431 fv.(names{n}) = fv.(names{n})(:, 1:vc); 0432 end 0433 0434 % Save a list of the field names in the FVCOM river data. 0435 fnames = fieldnames(fv); 0436 0437 % Assign the relevant arrays to the Mobj. Flux is added in the section 0438 % dealing with either climatology or time series data. 0439 Mobj.river_nodes = fv_obc; 0440 Mobj.river_names = fv_names'; 0441 Mobj.river_nemo_location = fv_location; 0442 Mobj.have_rivers = true; 0443 Mobj.nRivers = length(fv_obc); 0444 0445 % Dump the data into the relevant arrays in Mobj. 0446 for n = 1:length(fnames) 0447 new = sprintf('river_%s', fnames{n}); 0448 Mobj.(new) = fv.(fnames{n}); 0449 end 0450 0451 % Create a Modified Julian Day time series of the NEMO river data. Assume 0452 % it's a climatology, so generate the right time series based on that 0453 % assumption. 0454 if isempty(yr) && ~isnumeric(yr) 0455 error('For climatology, a year must be specified for the time series to be generated.') 0456 end 0457 0458 % Make three years of data starting from the year before the current one. 0459 % Do so accounting for leap years. 0460 daysinyr = [sum(eomday(yr - 1, 1:12)), ... 0461 sum(eomday(yr, 1:12)), ... 0462 sum(eomday(yr + 1, 1:12))]; 0463 % Offset the checkdate by one to add zero for the first day. Alternative 0464 % would be to specify the day as the end of the previous month (or a day 0465 % value of zero?). 0466 offsetdays = (1:sum(daysinyr)) - 1; 0467 mtime = datevec(datenum(yr - 1, 1, 1, 0, 0, 0) + offsetdays); 0468 Mobj.river_time = greg2mjulian(mtime(:, 1), mtime(:, 2), ... 0469 mtime(:, 3), mtime(:, 4), mtime(:, 5), mtime(:, 6)); 0470 0471 % Repeat the river data for the climatology before adding to the Mobj. 0472 for n = 1:length(fnames) 0473 new = sprintf('river_%s', fnames{n}); 0474 0475 % NEMO is insane. Old versions of the rivers data had enough days for 0476 % leap and non-leap years; new versions don't. So, let's pad the data 0477 % by a day in either case and that should work for data with both 365 0478 % and 366 times. In the case where we're using data with 366 values 0479 % already, this shouldn't make any difference as we're explicitly 0480 % indexing to the days in the year using daysinyr anyway. 0481 fv.(fnames{n}) = [fv.(fnames{n}); fv.(fnames{n})(end, :)]; 0482 0483 Mobj.(new) = [fv.(fnames{n})(1:daysinyr(1), :); ... 0484 fv.(fnames{n})(1:daysinyr(2), :); ... 0485 fv.(fnames{n})(1:daysinyr(3), :)]; 0486 end 0487 Mobj.rivers_orig = nemo.rivers; 0488 if ftbverbose 0489 fprintf('included %d of %d rivers (skipped %d)\n', ... 0490 fv_nr - skipped, fv_nr, skipped) 0491 fprintf('end : %s \n', subname) 0492 end