Home > fvcom_prepro > wrf2fvcom.m

wrf2fvcom

PURPOSE ^

Performs minor fixes to the output of wrf_to_fvcom in the FVCOM source

SYNOPSIS ^

function wrf2fvcom(wrffile, fvfile)

DESCRIPTION ^

 Performs minor fixes to the output of wrf_to_fvcom in the FVCOM source
 directory to netCDF.

 wrf2fvcom(ncfile, wrffile)

 DESCRIPTION:
   Take FVCOM routine wrf_to_fvcom netCDF output and apply minor fixes to
   the variables.

   The fixes are:
       - convert pressure from millibars to Pascals

 INPUT:
   wrffile - path to the wrf_to_fvcom netCDF file.
   fvfile - path to the FVCOM forcing file to create.

 OUTPUT:
   FVCOM forcing netCDF file (fvfile).

 EXAMPLE USAGE:
   wrf2fvcom('wrf_to_fvcom_output.nc', 'casename_v01_wnd.nc')

 NOTES:
   This currently only supports regularly gridded output (i.e. WRF native
   grid). Support for writing out unstructured data may come in the
   future.

 Author(s):
   Dmitry Aleynik (Scottish Association for Marine Science)
   Pierre Cazenave (Plymouth Marine Laboratory)

 Revision history:
   2015-10-23 - First version developed with reference to Dmitry Aleynik's
   dump_meteo_fvgrid.m function.

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

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function wrf2fvcom(wrffile, fvfile)
0002 % Performs minor fixes to the output of wrf_to_fvcom in the FVCOM source
0003 % directory to netCDF.
0004 %
0005 % wrf2fvcom(ncfile, wrffile)
0006 %
0007 % DESCRIPTION:
0008 %   Take FVCOM routine wrf_to_fvcom netCDF output and apply minor fixes to
0009 %   the variables.
0010 %
0011 %   The fixes are:
0012 %       - convert pressure from millibars to Pascals
0013 %
0014 % INPUT:
0015 %   wrffile - path to the wrf_to_fvcom netCDF file.
0016 %   fvfile - path to the FVCOM forcing file to create.
0017 %
0018 % OUTPUT:
0019 %   FVCOM forcing netCDF file (fvfile).
0020 %
0021 % EXAMPLE USAGE:
0022 %   wrf2fvcom('wrf_to_fvcom_output.nc', 'casename_v01_wnd.nc')
0023 %
0024 % NOTES:
0025 %   This currently only supports regularly gridded output (i.e. WRF native
0026 %   grid). Support for writing out unstructured data may come in the
0027 %   future.
0028 %
0029 % Author(s):
0030 %   Dmitry Aleynik (Scottish Association for Marine Science)
0031 %   Pierre Cazenave (Plymouth Marine Laboratory)
0032 %
0033 % Revision history:
0034 %   2015-10-23 - First version developed with reference to Dmitry Aleynik's
0035 %   dump_meteo_fvgrid.m function.
0036 %
0037 %==========================================================================
0038 
0039 global ftbverbose;
0040 subname = 'wrf2fvcom';
0041 if ftbverbose
0042     fprintf('\nbegin : %s \n', subname)
0043     fprintf('Read WRF output from %s... ', wrffile)
0044 end
0045 
0046 % Read the WRF output and build a struct with the relevant fields in it.
0047 wrf_fields = {'Shortwave' ,'Net_Heat', ...
0048     'U10', 'V10', ...
0049     'Precipitation', 'Evaporation', 'SLP', ...
0050     'XLAT', 'XLONG', 'Times'};
0051 
0052 for ff = 1:length(wrf_fields)
0053     wrf.(wrf_fields{ff}) = ncread(wrffile, wrf_fields{ff});
0054 end
0055 
0056 % Get the required dimensions.
0057 finfo = ncinfo(wrffile);
0058 dimnames = {finfo.Dimensions.Name};
0059 for dd = 1:length(dimnames);
0060     switch dimnames{dd}
0061         case 'west_east'
0062             wrf.west_east = finfo.Dimensions(dd).Length;
0063         case 'south_north'
0064             wrf.south_north = finfo.Dimensions(dd).Length;
0065         case 'DateStrLen'
0066             wrf.DateStrLen = finfo.Dimensions(dd).Length;
0067     end
0068 end
0069 
0070 % Sort out times.
0071 wrf.ntimes = length(wrf.Times(1, :));
0072 wrf.mtime = datenum(wrf.Times');
0073 
0074 if ftbverbose
0075     fprintf('done.\n')
0076 end
0077 
0078 grid_name = 'wrf_grid';
0079 grid_source = 'wrf grid (structured) surface forcing';
0080 CoordVar = 'lat lon';
0081 
0082 % -------------------------------------------------------------------------
0083 % Define the netCDF parameters
0084 %--------------------------------------------------------------------------
0085 
0086 if ftbverbose
0087     fprintf('Export WRF to FVCOM input file %s... ', fvfile)
0088 end
0089 
0090 % Define global attributes.
0091 nc.type = 'FVCOM METEO FORCING FILE';
0092 nc.title = 'WRF model forcing';
0093 nc.history = sprintf('File created with %s from the MATLAB fvcom-toolbox', subname);
0094 nc.source =  grid_source;
0095 nc.START_DATE = datestr(wrf.mtime(1, 1), 31);
0096 nc.END_DATE = datestr(wrf.mtime(end, 1), 31);
0097 nc.file = fvfile;
0098 nc.Conventions = 'CF-1.0';
0099 
0100 % Populate the global attributes.
0101 ncid = netcdf.create(fvfile, 'clobber');
0102 netcdf.putAtt(ncid, netcdf.getConstant('NC_GLOBAL'), 'type', nc.type);
0103 netcdf.putAtt(ncid, netcdf.getConstant('NC_GLOBAL'), 'title', nc.title);
0104 netcdf.putAtt(ncid, netcdf.getConstant('NC_GLOBAL'), 'history', nc.history);
0105 netcdf.putAtt(ncid, netcdf.getConstant('NC_GLOBAL'), 'source', nc.source);
0106 netcdf.putAtt(ncid, netcdf.getConstant('NC_GLOBAL'), 'START_DATE', nc.START_DATE);
0107 netcdf.putAtt(ncid, netcdf.getConstant('NC_GLOBAL'), 'END_DATE', nc.END_DATE);
0108 netcdf.putAtt(ncid, netcdf.getConstant('NC_GLOBAL'), 'file', nc.file);
0109 netcdf.putAtt(ncid, netcdf.getConstant('NC_GLOBAL'), 'Conventions', nc.Conventions);
0110 
0111 sn_idim = netcdf.defDim(ncid, 'south_north', wrf.south_north);
0112 we_idim = netcdf.defDim(ncid, 'west_east', wrf.west_east);
0113 DateStrLen_idim = netcdf.defDim(ncid, 'DateStrLen', wrf.DateStrLen);
0114 time_idim = netcdf.defDim(ncid, 'time', 0);
0115 
0116 % Variables
0117 varid_Times = netcdf.defVar(ncid, 'Times', 'NC_CHAR', [DateStrLen_idim, time_idim]);
0118 netcdf.putAtt(ncid, varid_Times, 'description', 'UTC');
0119 
0120 varid_XLAT = netcdf.defVar(ncid, 'XLAT', 'NC_FLOAT', [we_idim, sn_idim]);
0121 varid_XLONG = netcdf.defVar(ncid, 'XLONG', 'NC_FLOAT', [we_idim, sn_idim]);
0122 varid_Shortwave = netcdf.defVar(ncid, 'Shortwave', 'NC_FLOAT', [we_idim, sn_idim, time_idim]);
0123 varid_Net_Heat = netcdf.defVar(ncid, 'Net_Heat', 'NC_FLOAT', [we_idim, sn_idim, time_idim]);
0124 varid_U10 = netcdf.defVar(ncid, 'U10', 'NC_FLOAT', [we_idim, sn_idim, time_idim]);
0125 varid_V10 = netcdf.defVar(ncid, 'V10', 'NC_FLOAT', [we_idim, sn_idim, time_idim]);
0126 varid_Precipitation = netcdf.defVar(ncid, 'Precipitation', 'NC_FLOAT', [we_idim, sn_idim, time_idim]);
0127 varid_Evaporation = netcdf.defVar(ncid, 'Evaporation', 'NC_FLOAT', [we_idim, sn_idim, time_idim]);
0128 varid_SLP = netcdf.defVar(ncid, 'air_pressure', 'NC_FLOAT',[we_idim, sn_idim, time_idim]);
0129 
0130 netcdf.putAtt(ncid, varid_XLAT, 'long_name', 'latitude');
0131 netcdf.putAtt(ncid, varid_XLAT, 'description', 'LATITUDE, SOUTH IS NEGATIVE');
0132 netcdf.putAtt(ncid, varid_XLAT, 'units', 'degrees_north');
0133 netcdf.putAtt(ncid, varid_XLAT, 'type', 'data');
0134 netcdf.putAtt(ncid, varid_XLONG, 'long_name', 'longitude');
0135 netcdf.putAtt(ncid, varid_XLONG, 'description', 'LONGITUDE, WEST IS NEGATIVE');
0136 netcdf.putAtt(ncid, varid_XLONG, 'units', 'degrees_east');
0137 netcdf.putAtt(ncid, varid_XLONG, 'type', 'data');
0138 
0139 netcdf.putAtt(ncid, varid_Shortwave, 'long_name', 'Shortwave, upward is negative');
0140 netcdf.putAtt(ncid, varid_Shortwave, 'units', 'W m-2');
0141 netcdf.putAtt(ncid, varid_Shortwave, 'grid', grid_name);
0142 netcdf.putAtt(ncid, varid_Shortwave, 'coordinates', CoordVar);
0143 netcdf.putAtt(ncid, varid_Shortwave, 'type', 'data');
0144 
0145 netcdf.putAtt(ncid, varid_Net_Heat, 'long_name', 'Sum of shortwave, longwave, sensible and latent heat fluxes, ocean lose heat is negative');
0146 netcdf.putAtt(ncid, varid_Net_Heat, 'units', 'W m-2');
0147 netcdf.putAtt(ncid, varid_Net_Heat, 'grid', grid_name);
0148 netcdf.putAtt(ncid, varid_Net_Heat, 'coordinates', CoordVar);
0149 netcdf.putAtt(ncid, varid_Net_Heat, 'type', 'data');
0150 
0151 netcdf.putAtt(ncid, varid_U10, 'long_name', 'Eastward Wind Velocity');
0152 netcdf.putAtt(ncid, varid_U10, 'description', 'U at 10 M');
0153 netcdf.putAtt(ncid, varid_U10, 'units', 'm s-1'  );
0154 netcdf.putAtt(ncid, varid_U10, 'grid', grid_name);
0155 netcdf.putAtt(ncid, varid_U10, 'type', 'data');
0156 
0157 netcdf.putAtt(ncid, varid_V10, 'long_name', 'Northward Wind Velocity');
0158 netcdf.putAtt(ncid, varid_V10, 'description', 'V at 10 M');
0159 netcdf.putAtt(ncid, varid_V10, 'units', 'm s-1'  );
0160 netcdf.putAtt(ncid, varid_V10, 'grid', grid_name);
0161 netcdf.putAtt(ncid, varid_V10, 'type', 'data');
0162 
0163 netcdf.putAtt(ncid, varid_Precipitation, 'long_name', 'Precipitation');
0164 netcdf.putAtt(ncid, varid_Precipitation, 'description', 'Precipitation, ocean lose water is negative');
0165 netcdf.putAtt(ncid, varid_Precipitation, 'units', 'm s-1'  );
0166 netcdf.putAtt(ncid, varid_Precipitation, 'grid', grid_name);
0167 netcdf.putAtt(ncid, varid_Precipitation, 'coordinates', CoordVar);
0168 netcdf.putAtt(ncid, varid_Precipitation, 'type', 'data');
0169 
0170 netcdf.putAtt(ncid, varid_Evaporation, 'long_name', 'Evaporation');
0171 netcdf.putAtt(ncid, varid_Evaporation, 'description', 'Evaporation, ocean lose water is negative');
0172 netcdf.putAtt(ncid, varid_Evaporation, 'units', 'm s-1'  );
0173 netcdf.putAtt(ncid, varid_Evaporation, 'grid', grid_name);
0174 netcdf.putAtt(ncid, varid_Evaporation, 'coordinates', CoordVar);
0175 netcdf.putAtt(ncid, varid_Evaporation, 'type', 'data');
0176 
0177 netcdf.putAtt(ncid, varid_SLP, 'long_name', 'Air Pressure');
0178 netcdf.putAtt(ncid, varid_SLP, 'description', 'Sea surface airpressure');
0179 netcdf.putAtt(ncid, varid_SLP, 'units', 'Pa'  );
0180 netcdf.putAtt(ncid, varid_SLP, 'grid', grid_name);
0181 netcdf.putAtt(ncid, varid_SLP, 'coordinates', CoordVar);
0182 netcdf.putAtt(ncid, varid_SLP, 'type', 'data');
0183 
0184 netcdf.endDef(ncid);
0185 
0186 % -------------------------------------------------------------------------
0187 % Dump the data
0188 % -------------------------------------------------------------------------
0189 
0190 for i = 1:wrf.ntimes
0191     cc = datestr(wrf.mtime(i), 31);
0192     netcdf.putVar(ncid, varid_Times, [0, i-1], [19, 1], cc);
0193 end
0194 
0195 netcdf.putVar(ncid, varid_XLONG, [0, 0], size(wrf.XLONG), wrf.XLONG)
0196 netcdf.putVar(ncid, varid_XLAT, [0, 0], size(wrf.XLAT), wrf.XLAT)
0197 
0198 netcdf.putVar(ncid, varid_Shortwave, [0, 0, 0], size(wrf.Shortwave), wrf.Shortwave); % shortwave W/m^2
0199 netcdf.putVar(ncid, varid_Net_Heat, [0, 0, 0], size(wrf.Net_Heat), wrf.Net_Heat); % net W/m^2
0200 netcdf.putVar(ncid, varid_U10, [0, 0, 0], size(wrf.U10), wrf.U10); % m/s
0201 netcdf.putVar(ncid, varid_V10, [0, 0, 0], size(wrf.V10), wrf.V10); % m/s
0202 netcdf.putVar(ncid, varid_Precipitation, [0, 0, 0], size(wrf.Precipitation), wrf.Precipitation); % m/s
0203 netcdf.putVar(ncid, varid_Evaporation, [0, 0, 0], size(wrf.Evaporation), wrf.Evaporation); % m/s
0204 netcdf.putVar(ncid, varid_SLP, [0, 0, 0], size(wrf.SLP), wrf.SLP * 100); % mbar -> Pa
0205 
0206 netcdf.close(ncid);
0207 
0208 if ftbverbose
0209     fprintf('done.\n')
0210     fprintf('end   : %s \n', subname)
0211 end

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