Home > fvcom_prepro > get_ERA_forcing.m

get_ERA_forcing

PURPOSE ^

Reads in ERA-20C netCDF files and outputs a struct containing the

SYNOPSIS ^

function era = get_ERA_forcing(Mobj, files, varargin)

DESCRIPTION ^

 Reads in ERA-20C netCDF files and outputs a struct containing the
 requested variables for a given spatial domain and time period.

 ERA = get_ERA_forcing(Mobj, files, 'varlist', {'lon', 'lat', 'ssr'})

 DESCRIPTION:
   For the given netCDF(s), load the data into a struct for the variables
   in varlist. If varlist is omitted, all variables are loaded.

 INPUT:
   Mobj - Mesh object containing the following fields:
       * nativeCoords - 'spherical' or 'cartesian'
       * lon, lat or x, y - node coordinates (depending on nativeCoords)
   files - array of paths to the ERA-20C netCDF files
   Optional keyword/value pairs:
   'varlist' - give an array of variables to extract from the netCDF
       files. Missing variables will cause an error.
   'clip' - set to true for clipping the data to the model domain, false
       to load the full data in the netCDFs.

 OUTPUT:
   era - struct containing the variables specified in varlist.

 NOTES:
   This script has been written against the netCDFs of the ECMWF-ERA20C
   data created by the python script ecmwf-era20c.py.

 Author(s)
   Pierre Cazenave (Plymouth Marine Laboratory)

 Revision history:
   2015-08-13 First version based loosely on read_ERA_wind.m in the
   fvcom-toolbox.

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

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function era = get_ERA_forcing(Mobj, files, varargin)
0002 % Reads in ERA-20C netCDF files and outputs a struct containing the
0003 % requested variables for a given spatial domain and time period.
0004 %
0005 % ERA = get_ERA_forcing(Mobj, files, 'varlist', {'lon', 'lat', 'ssr'})
0006 %
0007 % DESCRIPTION:
0008 %   For the given netCDF(s), load the data into a struct for the variables
0009 %   in varlist. If varlist is omitted, all variables are loaded.
0010 %
0011 % INPUT:
0012 %   Mobj - Mesh object containing the following fields:
0013 %       * nativeCoords - 'spherical' or 'cartesian'
0014 %       * lon, lat or x, y - node coordinates (depending on nativeCoords)
0015 %   files - array of paths to the ERA-20C netCDF files
0016 %   Optional keyword/value pairs:
0017 %   'varlist' - give an array of variables to extract from the netCDF
0018 %       files. Missing variables will cause an error.
0019 %   'clip' - set to true for clipping the data to the model domain, false
0020 %       to load the full data in the netCDFs.
0021 %
0022 % OUTPUT:
0023 %   era - struct containing the variables specified in varlist.
0024 %
0025 % NOTES:
0026 %   This script has been written against the netCDFs of the ECMWF-ERA20C
0027 %   data created by the python script ecmwf-era20c.py.
0028 %
0029 % Author(s)
0030 %   Pierre Cazenave (Plymouth Marine Laboratory)
0031 %
0032 % Revision history:
0033 %   2015-08-13 First version based loosely on read_ERA_wind.m in the
0034 %   fvcom-toolbox.
0035 %
0036 %==========================================================================
0037 
0038 subname = 'get_ERA_forcing';
0039 
0040 global ftbverbose
0041 if ftbverbose
0042     fprintf('\nbegin : %s \n', subname)
0043 end
0044 
0045 assert(nargin >= 2, 'Incorrect number of input arguments (2 minimum).')
0046 
0047 varlist = [];
0048 clip = false;
0049 for var = 1:2:length(varargin)
0050     switch varargin{var}
0051         case 'varlist'
0052             varlist = varargin{var + 1};
0053         case 'clip'
0054             clip = varargin{var + 1};
0055     end
0056 end
0057 
0058 nf = length(files);
0059 
0060 for f = 1:nf
0061     
0062     assert(exist(files{f}, 'file') == 2, 'file: %s does not exist', files{f})
0063     
0064     if ftbverbose
0065         fprintf('Loading file %s (%i of %i)\n', files{f}, f, nf)
0066     end
0067 
0068     if isempty(varlist)
0069         ncdetails = ncinfo(files{f});
0070         varlist = {ncdetails.Variables.Name};
0071     end
0072 
0073     % Only do the spatial data on the first file (the data are global).
0074     if f == 1
0075         lon = double(ncread(files{f}, 'lon'));
0076         lat = double(ncread(files{f}, 'lat'));
0077         lonvector = unique(lon);
0078         latvector = unique(lat);
0079         dx = unique(diff(lonvector));
0080         dy = mean(unique(diff(latvector)));
0081 
0082         if clip
0083             % To save on memory, cull the data which doesn't come from the
0084             % domain we've specified (with a 2 element buffer).
0085             if min(Mobj.lon) < 0 && max(Mobj.lon) < 0
0086                 % Easy, just shift by 360.
0087                 idx_lon = find(lonvector > (min(Mobj.lon) - (2 * dx)) + 360 ...
0088                     & lonvector < (max(Mobj.lon) + (2 * dx)) + 360);
0089             elseif min(Mobj.lon) < 0 && max(Mobj.lon) > 0
0090                 % This is the tricky one. We'll do two passes to extract
0091                 % the western chunk first (min(Mobj.lon) + 360 to 360),
0092                 % then the eastern chunk (0 - max(Mobj.lon))
0093                 idx_lon{1} = find(lonvector >= (min(Mobj.lon) - (2 * dx)) + 360)';
0094                 idx_lon{2} = find(lonvector < (max(Mobj.lon) + (2 * dx)))';
0095                 idx_lon = [idx_lon{:}];
0096             else
0097                 % Dead easy, we're in the eastern hemisphere, so nothing
0098                 % too strenuous here.
0099                 idx_lon = find(lonvector > (min(Mobj.lon) - (2 * dx)) ...
0100                     & lonvector < (max(Mobj.lon) + (2 * dx)));
0101             end
0102             % Latitudes are easy because there's only one system.
0103             idx_lat = find(latvector > (min(Mobj.lat) - (2 * dy)) & latvector < (max(Mobj.lat) + (2 * dy)));
0104 
0105         else
0106             idx_lon = 1:length(lonvector);
0107             idx_lat = 1:length(latvector);
0108         end
0109         % Make a grid of the domain
0110         [era.lon, era.lat] = deal(lonvector(idx_lon), latvector(idx_lat));
0111         era.lon(era.lon > 180) = era.lon(era.lon > 180) - 360;
0112     end
0113 
0114     if ftbverbose
0115         fprintf('Variables:\n')
0116     end
0117     for v = 1:length(varlist)
0118         if ftbverbose
0119             fprintf('\t%s... ', varlist{v})
0120         end
0121         vardump = double(ncread(files{f}, varlist{v}));
0122         era.(varlist{v}).data = vardump(idx_lon, idx_lat, :);
0123         era.(varlist{v}).lat = latvector(idx_lat);
0124         era.(varlist{v}).lon = lonvector(idx_lon);
0125         era.(varlist{v}).time = double(ncread(files{f}, 'time'));
0126         clear vardump
0127         if ftbverbose
0128             fprintf('done.\n')
0129         end
0130     end
0131 end
0132 
0133 if ftbverbose
0134     fprintf('end   : %s \n', subname)
0135 end

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