get_CFS_forcing.m 22.3 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498
function data = get_CFS_forcing(Mobj, modelTime, varargin)
% Get the required parameters from CFSv2 products to force FVCOM.
%
% data = get_CFS_forcing(Mobj, modelTime)
%
% DESCRIPTION:
%   Using OPeNDAP, extract the necessary parameters to create an FVCOM
%   forcing file. Requires the air_sea toolbox.
%
% INPUT:
%   Mobj - MATLAB mesh object. Must contain fields:
%       lon, lat    - array of longitude and latitudes.
%       have_lonlat - boolean to signify whether coordinates are spherical
%                   or cartesian.
%   modelTime - Modified Julian Date start and end times
%   varargin - parameter/value pairs
%       - list of variables to extract:
%           'varlist', {'nshf', 'uwnd', 'vwnd'}
%
% OUTPUT:
%   data - struct of the data necessary to force FVCOM. These can be
%   interpolated onto an unstructured grid in Mobj using grid2fvcom.m.
%
% The parameters which can be obtained from the NCEP data are:
%     - u wind component (uwnd)
%     - v wind component (vwnd)
%     - Downward longwave radiation surface (dlwrf)
%     - Net shortwave radiation surface (nswrs = uswrf - dswrf)
%     - Air temperature (air)
%     - Relative humidity (rhum)
%     - Precipitation rate (prate)
%     - Surface pressure (pres or press)
%     - Latent heat flux (lhtfl)
%     - Potential evaporation rate (pevpr)
%
% In addition to these, the momentum flux (tau) is calculated from wind
% data. Precipitation is converted from kg/m^2/s to m/s. Evaporation (Et)
% is calculated from the mean daily latent heat net flux (lhtfl) at the
% surface. Precipitation-evaporation is also created (P_E).
%
% EXAMPLE USAGE:
%   To download the default set of data (see list above):
%
%       forcing = get_CFS_forcing(Mobj, [51345, 51376]);
%
%   To only download wind data:
%
%       forcing = get_CFS_forcing(Mobj, [51345, 51376], 'varlist', {'uwnd', 'vwnd'});
%
% Author(s)
%   Pierre Cazenave (Plymouth Marine Laboratory)
%   Ricardo Torres (Plymouth Marine Laboratory)
%   Rory O'Hara Murray (Marine Scotland Science)
%
% Revision history:
%   2015-05-19 First version based on get_NCEP_forcing.m.
%
%==========================================================================

subname = 'get_CFS_forcing';

global ftbverbose;
if ftbverbose
    fprintf('\nbegin : %s\n', subname)
end

% Parse the input arguments
varlist = [];
if nargin > 2
    for a = 1:2:nargin - 2
        switch varargin{a}
            case 'varlist'
                varlist = varargin{a + 1};
        end
    end
end
if ftbverbose
    fprintf('Extracting CFSv2 data.\n')
end

% Get the extent of the model domain (in spherical)
if ~Mobj.have_lonlat
    error('Need spherical coordinates to extract the forcing data')
else
    % Add a buffer of one grid cell in latitude and two in longitude to
    % make sure the model domain is fully covered by the extracted data.
    [dx, dy] = deal(0.5, 0.5); % approximate CFSv2 resolution in degrees
    extents = [min(Mobj.lon(:)) - (2 * dx), ...
        max(Mobj.lon(:)) + (2 * dx), ...
        min(Mobj.lat(:)) - dy, ...
        max(Mobj.lat(:)) + dy];
end

% if modelTime(end) - modelTime(1) > 365
%     error('Can''t (yet) process more than a year at a time.')
% end
[yyyy, mm, dd, HH, MM, SS] = mjulian2greg(modelTime);
dates = datenum([yyyy; mm; dd; HH; MM; SS]');
serial = dates(1):dates(2);
[years, months, ~, ~, ~, ~] = datevec(serial);
years = unique(years, 'stable');
months = unique(months, 'stable');
ny = length(years);
nm = length(months);

for y = 1:ny
    year = years(y);
    for m = 1:nm
        month = months(m);
        % Set up a struct of the remote locations in which we're
        % interested.

        url = 'http://nomads.ncdc.noaa.gov/thredds/dodsC/cfsr1hr/';
        % Get the forcing data.
        ncep.tmp2m   = [url, sprintf('%04d%02d/tmp2m.gdas.%04d%02d.grb2', year, month, year, month)];
        ncep.wnd     = [url, sprintf('%04d%02d/wnd10m.gdas.%04d%02d.grb2', year, month, year, month)];
        ncep.q2m     = [url, sprintf('%04d%02d/q2m.gdas.%04d%02d.grb2', year, month, year, month)];
        ncep.prate   = [url, sprintf('%04d%02d/prate.gdas.%04d%02d.grb2', year, month, year, month)];
        ncep.pressfc = [url, sprintf('%04d%02d/pressfc.gdas.%04d%02d.grb2', year, month, year, month)];
        ncep.lhtfl   = [url, sprintf('%04d%02d/lhtfl.gdas.%04d%02d.grb2', year, month, year, month)];
        ncep.dswsfc  = [url, sprintf('%04d%02d/dswsfc.gdas.%04d%02d.grb2', year, month, year, month)];
        ncep.uswsfc  = [url, sprintf('%04d%02d/uswsfc.gdas.%04d%02d.grb2', year, month, year, month)];
        ncep.dlwsfc  = [url, sprintf('%04d%02d/dlwsfc.gdas.%04d%02d.grb2', year, month, year, month)];

        names.tmp2m = 'Temperature';
        names.uwnd = 'U-component_of_wind';
        names.uwnd = 'V-component_of_wind';
        
        fields = fieldnames(ncep);

        for aa = 1:length(fields)        
            % We've been given a list of variables to do, so skip those that
            % aren't in the list.
            if ~isempty(varlist) && max(strcmp(fields{aa}, varlist)) ~= 1
                continue
            end

            if ftbverbose
                fprintf('getting ''%s'' data... ', fields{aa})
            end

            data.(fields{aa}).data = [];
            data.(fields{aa}).time = [];
            data.(fields{aa}).lat = [];
            data.(fields{aa}).lon = [];
            data_attributes.(fields{aa}) = [];

            %ncid_info = ncinfo(ncep.(fields{aa}));
            ncid = netcdf.open(ncep.(fields{aa}));

            % If you don't know what it contains, start by using the
            % 'netcdf.inq' operation:
            %[numdims,numvars,numglobalatts,unlimdimid] = netcdf.inq(ncid);
            % Time is in hours since the start of the month. We want
            % sensible times, so we'll have to offset at some point.
            varid = netcdf.inqVarID(ncid, 'time');
            data_time.time = netcdf.getVar(ncid, varid, 'double');

            varid = netcdf.inqVarID(ncid,'lon');
            data_lon.lon = netcdf.getVar(ncid,varid,'double');
            varid = netcdf.inqVarID(ncid,'lat');
            data_lat.lat = netcdf.getVar(ncid,varid,'double');
            % Some of the NCEP Reanalysis 2 data are 4D, but with a single
            % vertical level (e.g. uwnd, vwnd, air, rhum).
            data_level_idx = [];
            try % not all data have a 'level', so fail gracefully here.
                varid = netcdf.inqVarID(ncid, 'level');
                data_level.level = netcdf.getVar(ncid, varid, 'double');
                if length(data_level.level) > 1
                    % Assume we've got rhum and we want humidity at the sea
                    % surface (1013 millibars (or hPa)). As such, ZQQ must be
                    % 0.0 in the FVCOM model namelist. Find the closest level
                    % to pressure at 1 standard atmosphere.
                    [~, data_level_idx] = min(abs(data_level.level - 1013));
                end
            catch
                true;
            end
            if isempty(data_level_idx) % default to the first
                data_level_idx = 1;
            end

            % Time is in hours relative to the start of the month.
            timevec = datevec((data_time.time / 24) + datenum(year, month, 1, 0, 0, 0));

            % Get the data time and convert to Modified Julian Day.
            data.time = greg2mjulian(...
                timevec(:, 1), ...
                timevec(:, 2), ...
                timevec(:, 3), ...
                timevec(:, 4), ...
                timevec(:, 5), ...
                timevec(:, 6));
            % Clip the time to the given range
            data_time_mask = data.time >= modelTime(1) & data.time <= modelTime(end);
            data_time_idx = 1:size(data.time, 1);
            data_time_idx = data_time_idx(data_time_mask);

            % Check the times
            %[yyyy,mm,dd,hh,MM,ss] = mjulian2greg(data.time(1))
            %[yyyy,mm,dd,hh,MM,ss] = mjulian2greg(data.time(end))
            % Get the data in two goes, once for the end of the grid (west of
            % Greenwich), once for the beginning (east of Greenwich), and then
            % stick the two bits together.
            clear index_lon index_lat
            if extents(1) < 0 && extents(2) < 0
                % This is OK, we can just shunt the values by 360.
                extents(1) = extents(1) + 360;
                extents(2) = extents(2) + 360;
                index_lon = find(data_lon.lon > extents(1) & data_lon.lon < extents(2));
            elseif extents(1) < 0 && extents(2) > 0
                % This is the tricky one. We'll do two passes to extract the
                % western chunk first (extents(1)+360 to 360), then the eastern
                % chunk (0-extents(2)).
                index_lon{1} = find(data_lon.lon >= extents(1) + 360);
                index_lon{2} = find(data_lon.lon <= extents(2));
            else
                % Dead easy, we're in the eastern hemisphere, so nothing too
                % strenuous here.
                index_lon = find(data_lon.lon > extents(1) & data_lon.lon < extents(2));
            end

            % Latitude is much more straightforward
            index_lat = find(data_lat.lat > extents(3) & data_lat.lat < extents(4));
            data.(fields{aa}).lat = data_lat.lat(index_lat);

            % Get the data
            if iscell(index_lon)
                data.(fields{aa}).lon = data_lon.lon(cat(1,index_lon{:}));

                varid = netcdf.inqVarID(ncid, names.(fields{aa}));

                [~, ~, dimids, ~] = netcdf.inqVar(ncid,varid);
                if length(dimids) == 4
                    start = [min(index_lon{1}), min(index_lat), data_level_idx, min(data_time_idx)] - 1;
                    count = [length(index_lon{1}), length(index_lat), length(data_level_idx), length(data_time_idx)];
                elseif length(dimids) == 3
                    start = [min(index_lon{1}), min(index_lat), min(data_time_idx)] - 1;
                    count = [length(index_lon{1}), length(index_lat), length(data_time_idx)];
                end

                data1_west.(fields{aa}).(fields{aa}) = netcdf.getVar(ncid, varid, start, count, 'double');

                if length(dimids) == 4
                    start = [min(index_lon{2}), min(index_lat), data_level_idx, min(data_time_idx)] - 1;
                    count = [length(index_lon{2}), length(index_lat), length(data_level_idx), length(data_time_idx)];
                elseif length(dimids) == 3
                    start = [min(index_lon{2}), min(index_lat), min(data_time_idx)] - 1;
                    count = [length(index_lon{2}), length(index_lat), length(data_time_idx)];
                end
                data1_east.(fields{aa}).(fields{aa}) = netcdf.getVar(ncid, varid, start, count, 'double');

                data1.(fields{aa}).(fields{aa}).(fields{aa}) = ...
                    cat(1, data1_west.(fields{aa}).(fields{aa}), data1_east.(fields{aa}).(fields{aa}));

                % Merge the two sets of data together
                structfields = fieldnames(data1_west.(fields{aa}).(fields{aa}));
                for ii = 1:length(structfields)
                    switch structfields{ii}
                        case 'lon'
                            % Only the longitude and the actual data need
                            % sticking together, but each must be done
                            % along a different axis (lon is a vector, the
                            % data is an array).
                            data1.(fields{aa}).(fields{aa}).(structfields{ii}) = ...
                                [data1_west.(fields{aa}).(fields{aa}).(structfields{ii});data1_east.(fields{aa}).(fields{aa}).(structfields{ii})];
                        case fields{aa}
                            % This is the actual data
                            data1.(fields{aa}).(fields{aa}).(structfields{ii}) = ...
                                [data1_west.(fields{aa}).(fields{aa}).(structfields{ii}),data1_east.(fields{aa}).(fields{aa}).(structfields{ii})];
                        otherwise
                            % Assume the data are the same in both arrays.
                            % A simple check of the range of values in the
                            % difference between the two arrays should show
                            % whether they're the same or not. If they are,
                            % use the western values, otherwise, warn about
                            % the differences. It might be the data are
                            % relatively unimportant anyway (i.e. not used
                            % later on).
                            try
                                tdata = data1_west.(fields{aa}).(fields{aa}).(structfields{ii}) - data1_east.(fields{aa}).(fields{aa}).(structfields{ii});
                                if range(tdata(:)) == 0
                                    % They're the same data
                                    data1.(fields{aa}).(fields{aa}).(structfields{ii}) = ...
                                        data1_west.(fields{aa}).(fields{aa}).(structfields{ii});
                                else
                                    warning('Unexpected data field and the west and east halves don''t match. Skipping.')
                                end
                            catch
                                warning('Unexpected data field and the west and east halves don''t match. Skipping.')
                            end
                            clear tdata
                    end
                end
            else
                % We have a straightforward data extraction
                data.(fields{aa}).lon = data_lon.lon(index_lon);

                varid = netcdf.inqVarID(ncid,(fields{aa}));
                % [varname,xtype,dimids,natts] = netcdf.inqVar(ncid,varid);
                % [~,length1] = netcdf.inqDim(ncid,dimids(1))
                % [~,length2] = netcdf.inqDim(ncid,dimids(2))
                % [~,length3] = netcdf.inqDim(ncid,dimids(3))
                start=[min(index_lon)-1,min(index_lat)-1,min(data_time_idx)-1];
                count=[length(index_lon),length(index_lat),length(data_time_idx)];
                % The air data was failing with a three long start and count
                % array, so try first without (to retain original behaviour for
                % other potentially unaffected variables) but fall back to
                % getting only the first level (start = 0, count = 1).
                try
                    data1.(fields{aa}).(fields{aa}).(fields{aa}) = netcdf.getVar(ncid,varid,start,count,'double');
                catch
                    start=[min(index_lon)-1,min(index_lat)-1,0,min(data_time_idx)-1];
                    count=[length(index_lon),length(index_lat),1,length(data_time_idx)];
                    data1.(fields{aa}).(fields{aa}).(fields{aa}) = netcdf.getVar(ncid,varid,start,count,'double');
                end

            end

            datatmp = squeeze(data1.(fields{aa}).(fields{aa}).(fields{aa}));
            datatmp = (datatmp * data_attributes.(fields{aa}).(fields{aa}).scale_factor) + data_attributes.(fields{aa}).(fields{aa}).add_offset;

            % Fix the longitude ranges for all data.
            data.(fields{aa}).lon(data.(fields{aa}).lon > 180) = ...
                data.(fields{aa}).lon(data.(fields{aa}).lon > 180) - 360;

            data.(fields{aa}).data = datatmp;
            data.(fields{aa}).time = data.time;
            data.(fields{aa}).unpacked_valid_range = ...
                data_attributes.(fields{aa}).(fields{aa}).unpacked_valid_range;
            %     data.(fields{aa}).time = cat(1, data.(fields{aa}).time, squeeze(data1.(fields{aa}).(fields{aa}).time));
            %     data.(fields{aa}).lat = squeeze(data1.(fields{aa}).(fields{aa}).lat);
            %     data.(fields{aa}).lon = squeeze(data1.(fields{aa}).(fields{aa}).lon);

            % Replace values outside the specified actual range with NaNs. For the
            % majority of the variables, this shouldn't ever really generate values
            % of NaN since the coverage is global (land and sea). This did crop up
            % as a problem with the pevpr data (which is land only). In some ways,
            % if something fails later on (e.g. the interpolation) because there's
            % NaNs, that should be a wakeup call to check what's going on with the
            % data.
            if isfield(data_attributes.(fields{aa}).(fields{aa}), 'actual_range')
                actual_min = data_attributes.(fields{aa}).(fields{aa}).actual_range(1);
                actual_max = data_attributes.(fields{aa}).(fields{aa}).actual_range(2);
                mask = data.(fields{aa}).data < actual_min | data.(fields{aa}).data > actual_max;
                data.(fields{aa}).data(mask) = NaN;
            end

            if ftbverbose
                if isfield(data, fields{aa})
                    fprintf('done.\n')
                else
                    fprintf('error!\n')
                end
            end
        end

        % Calculate the net long and shortwave radiation fluxes.
        if isfield(data, 'ulwrf') && isfield(data, 'uswrf') && isfield(data, 'dlwrf') && isfield(data, 'dswrf')
            vars = {'nswrs', 'nlwrs'};
            up = {'uswrf', 'ulwrf'};
            down = {'dswrf', 'dlwrf'};
            for i = 1:length(vars)
                data.(vars{i}).data = data.(up{i}).data - data.(down{i}).data;
                data.(vars{i}).time = data.(up{i}).time;
                data.(vars{i}).lon = data.(up{i}).lon;
                data.(vars{i}).lat = data.(up{i}).lat;
            end
        end

        % Now we have some data, we need to create some additional parameters
        % required by FVCOM.

        % Convert precipitation from kg/m^2/s to m/s (required by FVCOM) by
        % dividing by freshwater density (kg/m^3).
        if isfield(data, 'prate')
            data.prate.data = data.prate.data/1000;
        end

        % Evaporation can be approximated by:
        %
        %   E(m/s) = lhtfl/Llv/rho
        %
        % where:
        %
        %   lhtfl   = "Mean daily latent heat net flux at the surface"
        %   Llv     = Latent heat of vaporization (approx to 2.5*10^6 J kg^-1)
        %   rho     = 1025 kg/m^3
        %
        if isfield(data, 'prate')
            Llv = 2.5*10^6;
            rho = 1025; % using a typical value for seawater.
            Et = data.lhtfl.data/Llv/rho;
            data.P_E.data = data.prate.data - Et;
            % Evaporation and precipitation need to have the same sign for FVCOM (ocean
            % losing water is negative in both instances). So, flip the evaporation
            % here.
            data.Et.data = -Et;
        end

        % Get the fields we need for the subsequent interpolation Find the position
        % of a sensibly sized array (i.e. not 'topo', 'rhum' or 'pres').
        for vv = 1:length(fields)
            if ~isempty(varlist) && max(strcmp(fields{vv}, varlist)) ~= 1
                continue
            end

            switch fields{vv}
                % Set ii in each instance in case we've been told to only use
                % one of the three alternatively gridded data.
                case 'topo'
                    ii = vv;
                    continue
                case 'rhum'
                    ii = vv;
                    continue
                case {'pres', 'press'}
                    ii = vv;
                    continue
                otherwise
                    % We've got one, so stop looking.
                    ii = vv;
                    break
            end
        end
        data.lon = data.(fields{ii}).lon;
        data.lon(data.lon > 180) = data.lon(data.lon > 180) - 360;
        data.lat = data.(fields{ii}).lat;

        % Create a land mask from the pevpr data (if it's been extracted).
        if isfield(data, 'pevpr')
            % Find any value less than or equal to the valid maximum across all
            % time steps.
            data.land_mask = max(data.pevpr.data <= data.pevpr.unpacked_valid_range(2), [], 3);
        end

        % Convert temperature to degrees Celsius (from Kelvin)
        if isfield(data, 'air')
            data.air.data = data.air.data - 273.15;
        end

        % Make sure all the data we have downloaded is the same shape as the
        % longitude and latitude arrays. This is complicated by the fact the NCEP
        % surface products (e.g. rhum, pres) are on a different grid from the rest
        % (e.g. uwnd).
        for aa = 1:length(fields)
            %     if strcmpi(fields{aa}, 'dswrf') || strcmpi(fields{aa}, 'dlwrf') || strcmpi(fields{aa}, 'uswrf') || strcmpi(fields{aa}, 'ulwrf')
            %         % But only if we haven't been given a list of variables to fetch.
            %         if nargin ~= 3
            %             continue
            %         end
            %     end

            if ~isempty(varlist) && max(strcmp(fields{aa}, varlist)) ~= 1
                % We've been given a list of variables to extract, so skip those
                % that aren't in that list
                continue
            else
                if isfield(data, fields{aa})
                    [px, py] = deal(length(data.(fields{aa}).lon), length(data.(fields{aa}).lat));
                    [ncx, ncy, ~] = size(data.(fields{aa}).data);
                    if ncx ~= px || ncy ~= py
                        data.(fields{aa}).data = permute(data.(fields{aa}).data, [2, 1, 3]);

                        % Check everything's OK now.
                        [ncx, ncy, ~] = size(data.(fields{aa}).data);
                        if ncx ~= px || ncy ~= py
                            error('Unable to resize data arrays to match position data orientation. Are these data NCEP surface data (i.e. on a different horizontal grid?)')
                        else
                            if ftbverbose
                                fprintf('Matching %s data dimensions to position arrays\n', fields{aa})
                            end
                        end
                    end
                else
                    warning('Variable %s requested but not downloaded?', fields{aa})
                end
            end
        end
    end
end

% Have a look at some data.
% [X, Y] = meshgrid(data.lon, data.lat);
% for i=1:size(data.uwnd.data, 3)
%     figure(1)
%     clf
%     uv = sqrt(data.uwnd.data(:, :, i).^2 + data.vwnd.data(:, :, i).^2);
%     pcolor(X, Y, uv')
%     shading flat
%     axis('equal','tight')
%     pause(0.1)
% end

if ftbverbose
    fprintf('end   : %s\n', subname)
end