website fvcom.docx
TRANSCRIPT
https://github.com/pwcazenave/fvcom-toolbox
Add new example script to create inputs for FVCOM (including heating, tides and open boundaries). master 20150319commitcd5f46cee7ebf5d9315d73523532eb1f44ad3c3b1 parent7c5ca4cpwcazenaveauthoredon Mar 19UnifiedSplitShowing1 changed filewith566 additionsand0 deletions.View566examples/example_fvcom_inputs.m@@ -0,0 +1,566 @@
+% Create the necessary files for FVCOM from a given SMS input unstructured
+% grid.
+%
+% This requires:
+% - my fork of the fvcom-toolbox
+% > https://gitlab.em.pml.ac.uk/pica/fvcom-toolbox or
+% > https://github.com/pwcazenave/fvcom-toolbox
+% - the Tide Model Driver toolbox and data
+% > http://polaris.esr.org/ptm_index.html
+% - the air-sea toolbox
+% > http://woodshole.er.usgs.gov/operations/sea-mat/air_sea-html
+%
+% Author(s):
+% Pierre Cazenave (Plymouth Marine Laboratory)
+%
+% Revision history:
+%
+% 2015-03-19 Example script to generate FVCOM inputs from TPXO, NCEP and
+% HYCOM data sources.
+
+matlabrc
+close all
+clc
+
+global ftbverbose
+ftbverbose = 1; % be noisy
+
+addpath('/users/modellers/pica/Code/fvcom-toolbox/utilities')
+addpath('/users/modellers/pica/Code/fvcom-toolbox/fvcom_prepro/')
+addpath('/users/modellers/rito/matlab/air-sea')
+conf.obc_tides.dirTMD = '/users/modellers/pica/Code/MATLAB/toolboxes/TMD2.03/';
+addpath(conf.obc_tides.dirTMD)
+addpath(fullfile(conf.obc_tides.dirTMD, 'FUNCTIONS'))
+
+%%%------------------------------------------------------------------------
+%%% INPUT CONFIGURATION
+%%%
+%%% Define the model input parameters here e.g. number of tidal components,
+%%% type of boundary forcing, estimated current velocity, sponge layer
+%%% coefficient and radius, forcing source etc.
+%%%
+%%%------------------------------------------------------------------------
+
+conf.base = '/data/medusa/pica/models/FVCOM/pml-tamar/run/';
+
+% Which version of FVCOM are we using (for the forcing file formats)?
+conf.FVCOM_version = '3.2.1';
+
+%%%------------------------------------------------------------------------
+%%% Time stuff
+%%%------------------------------------------------------------------------
+
+% Model time ([Y, M, D, h, m, s])
+conf.modelYear = 2013;
+conf.startDate = [conf.modelYear, 10, 01, 00, 00, 00];
+conf.endDate = [conf.modelYear, 11, 01, 00, 00, 00];
+
+% Time sampling for the surface forcing interpolation and the open boundary
+% tidal forcing (in hours)
+conf.sampling.surface = 1;
+conf.sampling.tides = 5 / 60;
+
+%%%------------------------------------------------------------------------
+%%% Spatial stuff
+%%%------------------------------------------------------------------------
+
+% Case name for the model inputs and outputs
+conf.casename = 'grid_v01';
+
+conf.coordType = 'cartesian'; % 'cartesian' or 'spherical'
+% Input grid UTM Zone (if applicable)
+conf.utmZone = {'30 U'}; % syntax for utm2deg
+
+% Option to smooth the bathymetry data.
+conf.smoothBathy = 'no'; % 'yes' or 'no'.
+if strcmpi(conf.smoothBathy, 'yes')
+ % Set the smoothing factor and number of iterations (see smoothmesh).
+ conf.smoothFactors = [0.5, 4]; % [factor, iterations]
+end
+
+% Sigma layer definition file.
+conf.sigma_file = fullfile(conf.base, 'input/configs/', conf.casename, 'sigma_geom.dat');
+
+% Give some names to the boundaries. This must match the number of node
+% strings defined in SMS. Ideally, the order of the names should match the
+% order in which the boundaries were made in SMS.
+conf.boundaryNames = {'South'}; % Number is (relatively) important!
+
+%%%------------------------------------------------------------------------
+%%% Model constants and stuff
+%%%------------------------------------------------------------------------
+
+% Type of open boundary treatment (see Table 6.1 (?) in the manual).
+% 1 - Active (ASL): sea level is specified at the open boundary.
+% 2 - Clamped: zeta = 0 at the boundary (Beardsley and Haidvogel, 1981).
+% 3 - Implicit Gravity Wave Radiation.
+% 4 - Partial Clamped Gravity Wave Radiation (Blumberg and Kantha, 1985).
+% 5 - Explicit Orlanski Radiation (Orlanski, 1976; Chapman, 1985)
+conf.obc_type = 1;
+
+% Sponge layer parameters
+conf.sponge.radius = 10000; % in metres
+conf.sponge.coeff = 0.0001;
+
+% z0 value in metres (for uniform) or 'random'.
+conf.bedRoughness = 0.03; % or 0.015, 0.025 or 0.03 - Davies and Furnes (1980) shelf model
+
+% Estimated velocity (m/s) and tidal range (m) for time step estimate
+conf.estVel = 2;
+conf.estRange = 5;
+
+%%%------------------------------------------------------------------------
+%%% Tides and stuff
+%%%------------------------------------------------------------------------
+
+% Open boundary nodal forcing type
+% 'z' for predicted surface elevation
+% 'fvcom' for FVCOM modelled surface elevation
+% 'phase-amp' for amplitudes and phases
+% 'polpred' for POLPRED amplitudes and phases
+conf.obc_tides.forcing = 'z';
+
+% How many tidal constituents do we actually want to use at the model
+% boundaries? Case sensitive (M2 != m2).
+% conf.obc_tides.components = {'M2','S2','N2','K2','K1','O1','P1','Q1','Mf','Mm','Ssa','M4','MS4','MN4'};
+conf.obc_tides.components = {'M2','S2','N2','K2','K1','O1','P1','Q1','M4'};
+
+% Location of the TMD model description file. Fall back to the global model
+% if the regional model doesn't work (NaNs in the time series). Regional
+% model used here is one suggested by Dima Aleynik and requires
+% ftp://ftp.oce.orst.edu/dist/tides/regional/ES.tar.Z to be extracted into
+% the TMD toolbox directory.
+conf.obc_tides.globalModel = '/users/modellers/pica/Code/MATLAB/toolboxes/TMD2.03/DATA/Model_tpxo7.2';
+conf.obc_tides.model = '/users/modellers/pica/Code/MATLAB/toolboxes/TMD2.03/DATA/Model_ES2008';
+
+%%%------------------------------------------------------------------------
+%%% Forcing and stuff
+%%%------------------------------------------------------------------------
+
+% Open boundary temperatures (string for source or number for constant).
+conf.obc_temp = 'HYCOM';
+
+% Open boundary salinities (string for source or number for constant).
+conf.obc_salt = 'HYCOM';
+
+% Surface heat fluxes.
+conf.surface_heat = 'NCEP';
+
+%%%------------------------------------------------------------------------
+%%% END OF INPUT CONFIGURATION
+%%%------------------------------------------------------------------------
+
+%% Process the files needed for all months.
+
+% Read the input mesh and bathymetry. Also creates the data necessary for
+% the Coriolis correction in FVCOM.
+Mobj = read_sms_mesh(...
+ '2dm', fullfile(conf.base, 'raw_data', [conf.casename, '.2dm']),...
+ 'bath', fullfile(conf.base, 'raw_data', [conf.casename, '.dat']),...
+ 'coordinate', conf.coordType, 'addCoriolis', true);
+
+% Clean out some unneeded fields.
+Mobj = rmfield(Mobj, 'riv_nodes');
+
+% Add grid metrics.
+Mobj = setup_metrics(Mobj);
+
+% Smooth the bathymetry if desired.
+if strcmpi(conf.smoothBathy, 'yes')
+ Mobj = setup_metrics(Mobj);
+ Mobj.h = smoothfield(Mobj.h, Mobj, ...
+ conf.smoothFactors(1), conf.smoothFactors(2));
+ % smoothfield2 is really inappropriate for bathymetry data.
+ % Mobj.h = smoothfield2(Mobj.h,Mobj,inputconf.smoothFactors(2));
+end
+
+% Create a Coriolis file from the bathy which varies with latitude. Given
+% the size of the domain, this is probably necessary. First need to convert
+% the UTM coordinates to lat/long to be able to calculate it, if
+% appropriate.
+if Mobj.have_lonlat == 0
+ utmZones = cellfun(@(x) repmat(x, length(Mobj.x), 1), conf.utmZone, 'uni', false);
+ [Mobj.lat, Mobj.lon] = utm2deg(Mobj.x, Mobj.y, utmZones{1});
+ Mobj.have_lonlat = true;
+ clear utmZones
+end
+Mobj = add_coriolis(Mobj, 'uselatitude');
+
+% Parse the open boundary nodes and add accordingly.
+if Mobj.have_strings
+ for i = 1:size(Mobj.read_obc_nodes, 2)
+ nodeList = double(cell2mat(Mobj.read_obc_nodes(i)));
+ Mobj = add_obc_nodes_list(Mobj, nodeList, conf.boundaryNames{i}, conf.obc_type);
+ end
+ clear nodeList
+end
+
+% Create a sponge layer
+if Mobj.have_strings
+ for i = 1:size(Mobj.read_obc_nodes, 2)
+ % Get the list of nodes in this boundary
+ nodeList = double(cell2mat(Mobj.read_obc_nodes(i)));
+ Mobj = add_sponge_nodes_list(Mobj,nodeList,...
+ [conf.boundaryNames{i}, ' sponge'], conf.sponge.radius,...
+ conf.sponge.coeff);
+ clear nodeList
+ end
+end
+assert(length(Mobj.sponge_nodes(Mobj.sponge_nodes ~= 0)) == sum(Mobj.nObcNodes), 'Number of sponge nodes does not match number of open boundary nodes.')
+clear i
+
+% Get the sigma depths in order to interpolate from the POLCOMS depths
+Mobj = read_sigma(Mobj, conf.sigma_file);
+
+% Do the bed roughness (uniform or variable)
+if ~ischar(conf.bedRoughness)
+ % Create a constant roughness z0 file
+ Mobj.z0 = repmat(conf.bedRoughness, 1, Mobj.nElems);
+elseif ischar(conf.bedRoughness) && strcmpi(conf.bedRoughness, 'random')
+ % Create a random bed roughness distribution with a maximum value of
+ % ~0.1 mm (value will be converted to metres).
+ Mobj.z0 = ((abs(randn(1,Mobj.nElems))/5)*0.1)/1000; % convert to metres
+else
+ fprintf('Unrecognised bed roughness type.\nSpecify a size (in m) or ''random'' for random bed roughness.\n')
+end
+
+% Estimate model time step. Supply estimated velocity (m/s) and tidal
+% range (m) after the mesh object.
+Mobj = estimate_ts(Mobj, conf.estVel, conf.estRange);
+
+%% Prepare the data in month long sections and output to subdirectories
+
+% Create an array of the days in this year's months to use to get the
+% Modified Julian Day start and end date range.
+daysOfMonths = eomday(conf.modelYear, 1:12);
+% Get the MJD of the start of the year. We'll use this within the loop to
+% add the days of the months we're iterating through.
+s0 = greg2mjulian(conf.modelYear, 1, 1, 0, 0, 0);
+
+for mm = conf.startDate(2):conf.endDate(2)
+ % Get the Modified Julian Day for the start and end of the month we're
+ % on at the moment. We'll pad by few days either way to give us a bit
+ % of leeway. Only do this within the year (i.e. don't pad to before the
+ % start of the year because the get_NCEP_forcing script can't get data
+ % across year boundaries.
+ if mm == 1
+ dOffsets = [0, 4];
+ elseif mm == 12
+ dOffsets = [2, 0];
+ else
+ dOffsets = [2, 4];
+ end
+ conf.startDateMJD = s0 + sum(daysOfMonths(1:mm)) - daysOfMonths(mm) - dOffsets(1);
+ conf.endDateMJD = s0 + sum(daysOfMonths(1:mm)) + dOffsets(2);
+ conf.time.tides = ...
+ conf.startDateMJD:conf.sampling.tides/24:conf.endDateMJD;
+ [sYr, sMon, sDay, sHr, sMin, sSec] = mjulian2greg(conf.startDateMJD);
+ [eYr, eMon, eDay, eHr, eMin, eSec] = mjulian2greg(conf.endDateMJD);
+ conf.startDate = [sYr, sMon, sDay, sHr, sMin, sSec];
+ conf.endDate = [eYr, eMon, eDay, eHr, eMin, eSec];
+ clear sYr sMon sDay sHr sMin sSec eYr eMon eDay eHr eMin eSec
+
+ % Output directory will contain some duplicate files (e.g. model grid
+ % etc.) but this makes things easier to manage. One subdirectory per
+ % month of the year in question.
+ conf.outbase = fullfile(conf.base, ...
+ 'input/configs/', ...
+ conf.casename, ...
+ sprintf('%04d/%02d', conf.modelYear, mm));
+ % Make the output directory if it doesn't exist.
+ if exist(conf.outbase, 'dir') ~= 7
+ mkdir(conf.outbase)
+ end
+
+ % Generate a surface elevation time series for open boundary forcing.
+ if strcmpi(conf.obc_tides.forcing, 'z') && exist('TMD') == 2 %#ok
+ % Use tmd_tide_pred to predict surface elevations for a given time
+ % range.
+
+ % Change to the TMD directory so the tidal generation works. We'll
+ % change back at the end.
+ oldDir = pwd;
+ cd(conf.obc_tides.dirTMD) % for TPXO to work
+
+ % Add the tidal components to the Mobj.
+ Mobj.Components = conf.obc_tides.components;
+
+ % Create a time series in MATLAB datenum format with ten minute
+ % inputs. First need to go from MJD to gregorian, then from
+ % gregorian to MATLAB dates.
+ conf.time.tidesMJD = datenum(...
+ conf.startDateMJD):...
+ conf.sampling.tides/24:...
+ datenum(conf.endDateMJD);
+ [tmpYY, tmpMM, tmpDD, tmphh, tmpmm, tmpss] = ...
+ mjulian2greg(conf.time.tidesMJD);
+ conf.time.tidesTPXO = ...
+ datenum(tmpYY, tmpMM, tmpDD, tmphh, tmpmm, tmpss);
+ clear tmpYY tmpMM tmpDD tmphh tmpmm tmpss
+
+ % Get the indices to use the tidal constituents defined in
+ % conf.obc_tides.components for TPXO (which requires a
+ % numerical array of the constituents to be used). The order of the
+ % TPXO constituents is M2, S2, N2, K2, K1, O1, P1, Q1, MF, MM, M4,
+ % MS4, MN4.
+ tpxoConsts = {'M2', 'S2', 'N2', 'K2', 'K1', 'O1', 'P1', 'Q1', ...
+ 'MF', 'MM', 'M4', 'MS4', 'MN4'};
+ tIndUse = nan(length(Mobj.Components), 1);
+ tInd = 1:length(tpxoConsts);
+ for i=1:length(Mobj.Components)
+ tPos = tInd(strcmp(Mobj.Components{i}, tpxoConsts));
+ if ~isempty(tPos)
+ tIndUse(i) = tPos;
+ else
+ warning('Supplied constituent (%s) is not present in the TPXO data', Mobj.Components{i}) %#ok
+ end
+ end
+ % Tidy up a bit
+ clear c tpxoConsts tPos tInd
+ tIndUse = tIndUse(~isnan(tIndUse));
+
+ % We can't just use tmd_tide_pred to do all the surface elevations
+ % at once. Instead, the useful approaches are:
+ %
+ % 1. Time series at a single location
+ % 2. Map of a given time step at all locations
+ %
+ % Since I'm likely to have many more time steps than locations,
+ % it's probably best to do the time series at each location than
+ % all the locations and a single time step.
+ %
+ % The order of the surface elevations in Mobj.surfaceElevation
+ % should reflect the order of the open boundary node IDs as FVCOM
+ % assumes they just map directly. So, rather than iterate through
+ % each position, we need to get the position based on the list of
+ % node IDs (Mobj.obc_nodes, without the zeros and in order of each
+ % boundary).
+ tmpObcNodes = Mobj.obc_nodes';
+ ObcNodes = tmpObcNodes(tmpObcNodes~=0)';
+ clear tmpObcNodes
+ surfaceElevation = nan(size(ObcNodes,2), size(conf.time.tidesTPXO, 2));
+ parfor i = 1:size(ObcNodes, 2)
+ % Get the current location (from the node ID)
+ currLon = Mobj.lon(ObcNodes(i));
+ currLat = Mobj.lat(ObcNodes(i));
+ if ftbverbose
+ fprintf('Position %i of %i (%.3f %.3f)... \n', i, size(ObcNodes, 2), currLon, currLat)
+ end
+ [surfaceElevation(i, :), ~] = ...
+ tmd_tide_pred(conf.obc_tides.model, ...
+ conf.time.tidesTPXO, currLat, currLon, ...
+ 'z', tIndUse);
+ if isnan(surfaceElevation(i, :))
+ % Try the global model instead.
+ [surfaceElevation(i, :), ~] = ...
+ tmd_tide_pred(conf.obc_tides.globalModel, ...
+ conf.time.tidesTPXO, currLat, currLon, ...
+ 'z', tIndUse);
+ end
+ end
+ Mobj.surfaceElevation = surfaceElevation;
+ % Tidy up some more
+ clear tIndUse obc_lat obc_lon ObcNodes currLon currLat surfaceElevation
+
+ if ftbverbose; fprintf('done.\n'); end
+
+ % Change back to the directory we were in before we needed to do
+ % the TMD stuff.
+ cd(oldDir); clear oldDir
+
+ end
+
+ % Write out all the initial output files.
+ try
+ % Grid
+ write_FVCOM_grid(Mobj, fullfile(conf.outbase, ...
+ [conf.casename, '_grd.dat']));
+
+ % Bathymetry
+ write_FVCOM_bath(Mobj, fullfile(conf.outbase, ...
+ [conf.casename, '_dep.dat']));
+
+ % Coriolis
+ write_FVCOM_cor(Mobj, fullfile(conf.outbase, ...
+ [conf.casename, '_cor.dat']));
+
+ % Open boundaries
+ write_FVCOM_obc(Mobj, fullfile(conf.outbase, ...
+ [conf.casename, '_obc.dat']))
+
+ % Sponge file
+ write_FVCOM_sponge(Mobj, fullfile(conf.outbase, ...
+ [conf.casename, '_spg.dat']))
+
+ % Bed roughness (constant or variable (see above)).
+ write_FVCOM_z0(Mobj.z0, fullfile(conf.outbase, ...
+ [conf.casename, '_z0=', ...
+ num2str(conf.bedRoughness), '.nc']), 'bottom roughness');
+
+ % Time series wave stations
+ write_FVCOM_stations(Mobj, fullfile(conf.outbase, ...
+ [conf.casename,'_station.dat']));
+
+ % Sigma file
+ copyfile(conf.sigma_file, conf.outbase)
+
+ catch err
+ rethrow(err)
+ end
+
+ % Get the surface heating data.
+ if strcmpi(conf.surface_heat, 'NCEP')
+ % Use the OPeNDAP NCEP script (get_NCEP_forcing.m) to get the
+ % following parameters:
+ % - Downward longwave radiation surface (dlwrs) [W/m^2]
+ % - Downward shortwave radiation surface (dswrs) [W/m^2]
+ % - Air temperature (air) [celsius]
+ % - Relative humidity (rhum) [%]
+ % - Sea level pressure (pres) [Pa]
+ %
+ % The script converts the NCEP data from the OPeNDAP server from
+ % longitudes in the 0 to 360 range to the -180 to 180 range. It
+ % also subsets for the right region (defined by Mobj.lon and
+ % Mobj.lat).
+ heating = get_NCEP_forcing(Mobj, ...
+ [conf.startDateMJD, conf.endDateMJD], ...
+ 'varlist', {'dlwrf', 'dswrf', 'air', 'rhum', 'pres'});
+
+ heating.domain_cols = length(heating.lon);
+ heating.domain_rows = length(heating.lat);
+ if isfield(heating, 'rhum')
+ heating.domain_cols_alt = length(heating.rhum.lon);
+ heating.domain_rows_alt = length(heating.rhum.lat);
+ end
+
+ % Convert the small subdomain into cartesian coordinates. We need
+ % to do this twice because some of the NCEP data are on different
+ % grids (e.g. sea level pressure, relative humidity etc.).
+ tmpZone = regexpi(conf.utmZone,'\ ','split');
+ [tmpLon, tmpLat] = meshgrid(heating.lon, heating.lat);
+ [heating.x, heating.y] = wgs2utm(tmpLat(:), tmpLon(:), str2double(char(tmpZone{1}(1))), 'N');
+ if isfield(heating, 'rhum')
+ [tmpLon2, tmpLat2] = meshgrid(heating.rhum.lon, heating.rhum.lat);
+ [heating.xalt, heating.yalt] = wgs2utm(tmpLat2(:), tmpLon2(:), str2double(char(tmpZone{1}(1))), 'N');
+ end
+ clear tmpLon tmpLat tmpLon2 tmpLat2 tmpZone
+ % Create arrays of the x and y positions.
+ heating.x = reshape(heating.x, heating.domain_rows, heating.domain_cols);
+ heating.y = reshape(heating.y, heating.domain_rows, heating.domain_cols);
+ if isfield(heating, 'rhum')
+ heating.xalt = reshape(heating.xalt, heating.domain_rows_alt, heating.domain_cols_alt);
+ heating.yalt = reshape(heating.yalt, heating.domain_rows_alt, heating.domain_cols_alt);
+ end
+
+ [heating.lon, heating.lat] = meshgrid(heating.lon, heating.lat);
+
+ heating = rmfield(heating, {'domain_rows', 'domain_cols'});
+ if isfield(heating, 'rhum')
+ heating = rmfield(heating, {'domain_rows_alt', 'domain_cols_alt'});
+ end
+
+ % Get rid of the alternative arrays as we don't need those anymore.
+ if isfield(heating, 'rhum')
+ heating = rmfield(heating, {'xalt', 'yalt'});
+ end
+
+ % Interpolate the data onto the FVCOM unstructured grid.
+ tic
+ interpfields = {'dswrf', 'dlwrf', 'pres', 'air', 'rhum', ...
+ 'time', 'lon', 'lat', 'x', 'y'};
+ heating_interp = grid2fvcom(Mobj, interpfields, heating);
+ if ftbverbose
+ fprintf('Elapsed interpolation time: %.2f minutes\n', toc / 60)
+ end
+
+ % Write out the surface forcing data to netCDF.
+ if exist('heating_interp', 'var')
+ heatBase = fullfile(conf.outbase, conf.casename);
+ write_FVCOM_forcing(Mobj, heatBase, heating_interp, ...
+ 'NCEP atmospheric forcing data', ...
+ conf.FVCOM_version);
+ end
+
+ end
+
+ % Now we need some boundary temperature and salinity conditions.
+ if any(strcmpi('HYCOM', {conf.obc_temp, conf.obc_salt}))
+ % Use HYCOM data for the boundary forcing.
+
+ % Offset the times to give us a bit of wiggle room.
+ modelTime = [conf.startDateMJD - dOffsets(1), conf.endDateMJD + dOffsets(2)];
+ hycom = get_HYCOM_forcing(Mobj, modelTime);
+
+ % Interpolate the 4D HYCOM data on the FVCOM vertical grid at the
+ % open boundaries.
+ Mobj = get_HYCOM_tsobc(Mobj, hycom);
+
+ write_FVCOM_tsobc(fullfile(conf.outbase, conf.casename), ...
+ Mobj.ts_times, ...
+ size(Mobj.temperature, 2), ...
+ Mobj.temperature, ...
+ Mobj.salt,...
+ Mobj)
+ end
+
+ % Add the daily SSH data on top of the predicted tidal elevations if
+ % we're using both predicted elevations and HYCOM data.
+ if strcmpi(conf.obc_tides.forcing, 'z') && any(strcmpi('HYCOM', {conf.obc_temp, conf.obc_salt}))
+
+ if ~isfield(hycom, 'ssh')
+ warning('No sea surface height field in the HYCOM data.')
+ elseif isfield(hycom, 'ssh') && ~strcmpi(conf.obc_tides.forcing, 'fvcom')
+ % Add the SSH from HYCOM to the predicted surface elevation.
+ % Just do a linear interpolation between the daily values and
+ % find the nearest point in the HYCOM grid (don't worry about
+ % doing fancy interpolations).
+ if Mobj.have_strings
+ % Get the list of open boundary nodes.
+ tmpObcNodes = Mobj.obc_nodes';
+ oNodes = tmpObcNodes(tmpObcNodes ~= 0)';
+
+ fvlon = Mobj.lon(oNodes);
+ fvlat = Mobj.lat(oNodes);
+
+ % Loop through all the nodes and find the nearest SSH
+ % values.
+ assert(length(fvlon) == length(fvlat), 'Inconsistent number of coordinates for the open boundary.')
+
+ for p = 1:length(fvlon)
+ % Find the closest HYCOM position and interpolate the
+ % SSH to the same sampling as the predicted elevations.
+ fx = fvlon(p);
+ fy = fvlat(p);
+ [~, jj] = sort(sqrt((hycom.lon(:) - fx).^2 + (hycom.lat(:) - fy).^2));
+ % Try the indices in order of proximity until we come
+ % across the first one that isn't the no data value
+ % (anything in excess of 1.26e29). Limit the search to
+ % the first 100 nearest HYCOM locations.
+ [ir, ic] = ind2sub(size(hycom.lon), jj(1));
+ cc = 0;
+ while max(hycom.ssh.data(ir, ic, :)) > 1.26e29
+ cc = cc + 1;
+ [ir, ic] = ind2sub(size(hycom.lon), jj(cc));
+ if cc > 100
+ error('Couldn''t find sea surface height value within the 100 nearest elements in the HYCOM data.')
+ end
+ end
+
+ % Interpolate in time, use nearest value for space.
+ ssh = interp1(hycom.time, squeeze(hycom.ssh.data(ir, ic, :)), conf.time.tidesMJD);
+
+ % Add this node's sea surface height component to the
+ % predicted.
+ Mobj.surfaceElevation(p, :) = Mobj.surfaceElevation(p, :) + ssh;
+ end
+ clear ssh ir ic cc fx fy p fvlon fvlat oNodes tmpObcNodes
+ end
+ end
+ % Write out the TPXO predicted surface elevation.
+ ElevationFile = fullfile(conf.outbase, [conf.casename, '_elevtide.nc']);
+ write_FVCOM_elevtide(Mobj, conf.time.tidesMJD, ElevationFile, ...
+ 'Model surface elevation boundary input')
+ end
+
+end