get_CFS_forcing.m 24.4 KB
Newer Older
1
function data = get_CFS_forcing(Mobj, modelTime, varargin)
Pierre Cazenave's avatar
Pierre Cazenave committed
2 3
% Get the required parameters from CFSv2 reanalysis products to force
% FVCOM.
4 5 6 7
%
% data = get_CFS_forcing(Mobj, modelTime)
%
% DESCRIPTION:
Pierre Cazenave's avatar
Pierre Cazenave committed
8 9
%   Using the NOAA OPeNDAP server, extract the necessary parameters to
%   create an FVCOM forcing file.
10 11 12 13 14 15 16
%
% 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
Pierre Cazenave's avatar
Pierre Cazenave committed
17
%   varargin - optional parameter/value pairs:
18
%       - list of variables to extract:
Pierre Cazenave's avatar
Pierre Cazenave committed
19
%           'varlist', {'tmp2m', 'uwnd', 'vwnd'}
20 21 22 23
%
% OUTPUT:
%   data - struct of the data necessary to force FVCOM. These can be
%   interpolated onto an unstructured grid in Mobj using grid2fvcom.m.
Pierre Cazenave's avatar
Pierre Cazenave committed
24
%   Contains vectors of the longitude and latitude data (lon, lat).
25 26
%
% The parameters which can be obtained from the NCEP data are:
Pierre Cazenave's avatar
Pierre Cazenave committed
27 28 29 30 31 32 33 34 35 36 37
%     - Net shortwave radiation surface (nswrs = uswrf - dswrf) [surface]
%     - Downward longwave radiation surface (dlwrf) [surface]
%     - Surface pressure (pressfc) [surface]
%     - u wind component (uwnd) [10m]
%     - v wind component (vwnd) [10m]
%     - Air temperature (tmp2m) [2m]
%     - Precipitation rate (prate) [surface]
%     - Specific humidity (q2m) [2m]
%     - Relative humidity (rhum) [2m] - calculated from q2m.
%     - Latent heat flux (lhtfl) [surface]
%     - Evaporation rate (Et) [surface]
38
%
Pierre Cazenave's avatar
Pierre Cazenave committed
39 40 41 42 43
% Relative humidity is calculated from specific humidity with the QAIR2RH
% function (see fvcom-toolbox/utilities). 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).
44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59
%
% 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:
Pierre Cazenave's avatar
Pierre Cazenave committed
60
%   2015-05-19 First version loosely based on get_NCEP_forcing.m.
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
%
%==========================================================================

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

98
% Create year and month arrays for the period we've been given.
99 100 101 102
[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);
103 104 105 106
[months, idx] = unique(months, 'stable');
years = years(idx);
nt = length(months);

107 108 109 110 111 112
for t = 1:nt
    month = months(t);
    year = years(t);
    if ftbverbose
        fprintf('Downloading for %04d/%02d\n', year, month)
    end
113

114 115 116
    % Set up a struct of the remote locations in which we're
    % interested.
    url = 'http://nomads.ncdc.noaa.gov/thredds/dodsC/cfsr1hr/';
Pierre Cazenave's avatar
Pierre Cazenave committed
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
    ncep.dlwsfc  = [url, ...
        sprintf('%04d%02d/dlwsfc.gdas.%04d%02d.grb2', year, month, ...
        year, month)];
    ncep.dswsfc  = [url, ...
        sprintf('%04d%02d/dswsfc.gdas.%04d%02d.grb2', year, month, ...
        year, month)];
    ncep.lhtfl   = [url, ...
        sprintf('%04d%02d/lhtfl.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.q2m     = [url, ...
        sprintf('%04d%02d/q2m.gdas.%04d%02d.grb2', year, month, ...
        year, month)];
    ncep.tmp2m   = [url, ...
        sprintf('%04d%02d/tmp2m.gdas.%04d%02d.grb2', year, month, ...
        year, month)];
    ncep.uswsfc  = [url, ...
        sprintf('%04d%02d/uswsfc.gdas.%04d%02d.grb2', year, month, ...
        year, month)];
    ncep.uwnd    = [url, ...
        sprintf('%04d%02d/wnd10m.gdas.%04d%02d.grb2', year, month, ...
        year, month)];
    ncep.vwnd    = [url, ...
        sprintf('%04d%02d/wnd10m.gdas.%04d%02d.grb2', year, month, ...
        year, month)];
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

    % We need variable names too since we can't store them as the keys in
    % ncep due to characters which MATLAB won't allow in fields (mainly -).
    names.dlwsfc = 'Downward_Long-Wave_Rad_Flux';
    names.dswsfc = 'Downward_Short-Wave_Rad_Flux';
    names.lhtfl = 'Latent_heat_net_flux';
    names.prate = 'Precipitation_rate';
    names.pressfc = 'Pressure';
    names.q2m = 'Specific_humidity';
    names.tmp2m = 'Temperature';
    names.uswsfc = 'Upward_Short-Wave_Rad_Flux';
    names.uwnd = 'U-component_of_wind';
    names.vwnd = '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
173

174 175 176 177 178 179 180 181 182 183 184
        % These are needed when catting the arrays together.
        if t == 1
            data.(fields{aa}).data = [];
            data.(fields{aa}).time = [];
            data.(fields{aa}).lat = [];
            data.(fields{aa}).lon = [];
        end
        scratch.(fields{aa}).data = [];
        scratch.(fields{aa}).time = [];
        scratch.(fields{aa}).lat = [];
        scratch.(fields{aa}).lon = [];
185

186
        % ncid_info = ncinfo(ncep.(fields{aa}));
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
        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 = netcdf.getVar(ncid, varid, 'double');
        if max(data_time) == 6
            % Precipitation data has times as 0-6 repeated for n days.
            % We need a sensible set of hours since the start of the
            % month for subsequent subsampling in time.
            data_time = 0:length(data_time) - 1;
        end

        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));
219
            end
220 221 222 223 224 225
        catch
            true;
        end
        if isempty(data_level_idx) % default to the first
            data_level_idx = 1;
        end
226

227 228 229 230
        % Time is in hours relative to the start of the month for CFSv2.
        timevec = datevec((data_time / 24) + datenum(year, month, 1, 0, 0, 0));

        % Get the data time and convert to Modified Julian Day.
231
        scratch.time = greg2mjulian(...
232 233 234 235 236 237
            timevec(:, 1), ...
            timevec(:, 2), ...
            timevec(:, 3), ...
            timevec(:, 4), ...
            timevec(:, 5), ...
            timevec(:, 6));
238 239 240 241 242 243 244 245 246 247 248 249
        % Clip the time to the given range. Because of some oddness with
        % some variables giving data beyond the end of the month whilst
        % others don't, set the limits in time for each month to be the
        % first/last day of the month or the modelTime start/end, whichever
        % is larger/smaller.
        startTime = max([modelTime(1), ...
            greg2mjulian(year, month, 1, 0, 0, 0)]);
        % Offset end by one day to capture the right number of days
        % (midnight falls at the beginning of the specified day).
        endTime = min([modelTime(end), ...
            greg2mjulian(year, month, eomday(year, month), 0, 0, 0) + 1]);
        data_time_mask = scratch.time >= startTime & scratch.time < endTime;
250
        data_time_idx = 1:size(scratch.time, 1);
251 252
        data_time_idx = data_time_idx(data_time_mask);
        if ~isempty(data_time_idx)
253
            scratch.time = scratch.time(data_time_mask);
254 255 256
        else
            % Reset the index to its original size. This is for data
            % with only a single time stamp which falls outside the
257
            % model time.
258 259
            if length(scratch.time) == 1
                data_time_idx = 1:size(scratch.time, 1);
260
            end
261 262
        end

263 264 265 266 267 268 269 270 271
        % Check the times.
        % [y, m, d, hh, mm, ss] = mjulian2greg(scratch.time);
        % fprintf('(%s - %s) ', ...
        %     datestr([y(1),m(1),d(1),hh(1),mm(1),ss(1)], ...
        %         'yyyy-mm-dd HH:MM:SS'), ...
        %     datestr([y(end),m(end),d(end),hh(end),mm(end),ss(end)], ...
        %         'yyyy-mm-dd HH:MM:SS'))
        % clearvars y m d hh mm ss oftv

272 273 274 275 276 277 278 279
        % 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;
Pierre Cazenave's avatar
Pierre Cazenave committed
280 281
            index_lon = find(data_lon.lon > extents(1) & ...
                data_lon.lon < extents(2));
282 283 284 285 286 287 288 289 290
        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.
Pierre Cazenave's avatar
Pierre Cazenave committed
291 292
            index_lon = find(data_lon.lon > extents(1) & ...
                data_lon.lon < extents(2));
293 294 295 296
        end

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

        % Get the data
        if iscell(index_lon)
301
            scratch.(fields{aa}).lon = data_lon.lon(cat(1, index_lon{:}));
302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325

            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)];
326 327
            end

Pierre Cazenave's avatar
Pierre Cazenave committed
328 329
            data_west.(fields{aa}).(fields{aa}) = ...
                netcdf.getVar(ncid, varid, start, count, 'double');
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

            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
            data_east.(fields{aa}).(fields{aa}) = ...
                netcdf.getVar(ncid, varid, start, count, 'double');

            scratch.(fields{aa}).(fields{aa}).(fields{aa}) = ...
                cat(1, ...
                data_west.(fields{aa}).(fields{aa}), ...
                data_east.(fields{aa}).(fields{aa}));

            % Merge the two sets of data together
            structfields = fieldnames(data_west.(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).
                        scratch.(fields{aa}).(structfields{ii}) = ...
                            [data_west.(fields{aa}).(structfields{ii}); ...
                            data_east.(fields{aa}).(structfields{ii})];
                    case fields{aa}
                        % This is the actual data.
                        scratch.(fields{aa}).(structfields{ii}) = ...
                            [rot90(data_west.(fields{aa}).(structfields{ii})), ...
                            rot90(data_east.(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 = ...
                                data_west.(fields{aa}).(structfields{ii}) - ...
                                data_east.(fields{aa}).(structfields{ii});
                            if range(tdata(:)) == 0
                                % They're the same data
                                scratch.(fields{aa}).(structfields{ii}) = ...
                                    data_west.(fields{aa}).(structfields{ii});
                            else
395
                                warning(['Unexpected data field and the', ...
396
                                    ' west and east halves don''t match.', ...
397
                                    ' Skipping.'])
398 399
                            end
                        catch
400
                            warning(['Unexpected data field and the', ...
401
                                ' west and east halves don''t match.', ...
402
                                ' Skipping.'])
403 404
                        end
                        clearvars tdata
405 406
                end
            end
407 408 409
            clearvars data_west data_east
        else
            % We have a straightforward data extraction
410
            scratch.(fields{aa}).lon = data_lon.lon(index_lon);
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

            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), ...
                min(index_lat), ...
                min(data_time_idx)] - 1;
            count = [...
                length(index_lon), ...
                length(index_lat), ...
                length(data_time_idx)];
            % The air data (NCEP version of this script) 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 the data_level_idx one instead (should be the first
            % level).
            try
                scratch.(fields{aa}).(fields{aa}).(fields{aa}) = ...
                    netcdf.getVar(ncid,varid,start,count,'double');
            catch
                start = [...
                    min(index_lon), ...
                    min(index_lat), ...
                    data_level_idx, ...
                    min(data_time_idx)] - 1;
                count = [...
                    length(index_lon), ...
                    length(index_lat), ...
                    1, ...
                    length(data_time_idx)];
Pierre Cazenave's avatar
Pierre Cazenave committed
445 446
                scratch.(fields{aa}).(fields{aa}) = ...
                    netcdf.getVar(ncid, varid, start, count, 'double');
447 448
            end

449 450
        end
        clearvars data_time* data_level_idx
451

452 453
        scratch.(fields{aa}).lon(scratch.(fields{aa}).lon > 180) = ...
            scratch.(fields{aa}).lon(scratch.(fields{aa}).lon > 180) - 360;
454

455
        datatmp = squeeze(scratch.(fields{aa}).(fields{aa}));
456

457 458 459
        % data.(fields{aa}).data = datatmp;
        data.(fields{aa}).data = cat(3, data.(fields{aa}).data, datatmp);
        % data.(fields{aa}).time = data.time;
460 461 462 463 464 465
        data.(fields{aa}).time = cat(1, data.(fields{aa}).time, scratch.time);
        % data.(fields{aa}).time = cat(1, data.(fields{aa}).time, ...
        %     squeeze(scratch.(fields{aa}).(fields{aa}).time));
        data.(fields{aa}).lat = scratch.(fields{aa}).lat;
        data.(fields{aa}).lon = scratch.(fields{aa}).lon;

466 467 468
        if ftbverbose
            if isfield(data, fields{aa})
                fprintf('done.\n')
469
            else
470
                fprintf('error!\n')
471 472
            end
        end
473
    end
474 475 476 477 478 479 480 481 482 483 484
end

% Now we have the data, we need to fix the averaging to be hourly instead
% of n-hourly, where n varies from 0 to 6. See
% http://rda.ucar.edu/datasets/ds094.1/#docs/FAQs_hrly_timeseries.html with
% the question "How can the individual one-hour averages be computed?".
fields = fieldnames(data);
for f = 1:length(fields)
    if isfield(data.(fields{f}), 'data')
        [~, ~, nt] = size(data.(fields{f}).data);
        fixed = data.(fields{f}).data;
485 486 487 488
        if ftbverbose
            fprintf('De-averaging the n-hourly %s data to hourly... ', ....
                fields{f})
        end
489 490 491 492 493 494 495 496 497 498

        for t = 1:6:nt
            % Fix the next 5 hours of data. Assume 0th hour is just the
            % original data - since the formula multiplies by the n-1 hour,
            % if we want the first hour's worth of data, then the second
            % term in the formula with multiply by zero, so the formula is
            % essentially only using the first term, which is just the data
            % at n (i.e. 0).
            for n = 1:5
                if t + n <= nt
499 500 501
                    fixed(:, :, t + n) = ...
                        (n * data.(fields{f}).data(:, :, t + n)) - ...
                        ((n - 1) * data.(fields{f}).data(:, :, t + n - 1));
502 503 504
                end
            end
        end
505 506
        data.(fields{f}).data = fixed;
        clearvars fixed
507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611
        if ftbverbose; fprintf('done.\n'); end
    end
end

% Calculate the net long and shortwave radiation fluxes.
if isfield(data, 'uswsfc') && isfield(data, 'dswsfc')
    data.nswsfc.data = data.uswsfc.data - data.dswsfc.data;
    data.nswsfc.time = data.uswsfc.time;
    data.nswsfc.lon = data.uswsfc.lon;
    data.nswsfc.lat = data.uswsfc.lat;
end

% 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') && isfield(data, 'lhtfl')
    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 (four including pres and press) alternatively
        % gridded data.
        case {'topo', 'rhum', '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;

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

% Convert specific humidity to relative humidity.
if isfield(data, 'q2m') && isfield(data, 'tmp2m') && isfield(data, 'pressfc')
    % Convert pressure from Pascals to millibars. Save relative humidity as
    % percentage. Convert specific humidity to percent too.
    data.rhum.data = 100 * qair2rh(data.q2m.data, data.tmp2m.data, data.pressfc.data / 100);
end
if isfield(data, 'q2m')
    data.q2m.data = 100 * data.q2m.data;
end

% Make sure all the data we have downloaded are the same shape as the
% longitude and latitude arrays.
for aa = 1:length(fields)
    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 on a ', ...
                        'different horizontal grid?'])
                end
            end
        else
            warning('Variable %s requested but not downloaded.', fields{aa})
        end
612 613 614 615 616
    end
end

% Have a look at some data.
% [X, Y] = meshgrid(data.lon, data.lat);
617
% for i = 1:size(data.uwnd.data, 3)
618 619 620 621 622 623 624 625 626 627 628 629 630
%     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