global_cycle
1.4.0
|
This is a limited point version of surface program. More...
Go to the source code of this file.
Data Types | |
module | sfccyc_module |
This program runs in two different modes: More... | |
Functions/Subroutines | |
subroutine | albocn (albclm, slmask, albomx, len) |
Set the albedo at open water points. More... | |
subroutine | analy (lugb, iy, im, id, ih, fh, len, lsoil, slmask, fntsfa, fnweta, fnsnoa, fnzora, fnalba, fnaisa, fntg3a, fnscva, fnsmca, fnstca, fnacna, fnvega, fnveta, fnsota, fnvmna, fnvmxa, fnslpa, fnabsa, tsfanl, wetanl, snoanl, zoranl, albanl, aisanl, tg3anl, cvanl, cvbanl, cvtanl, smcanl, stcanl, slianl, scvanl, acnanl, veganl, vetanl, sotanl, alfanl, tsfan0, vmnanl, vmxanl, slpanl, absanl, kpdtsf, kpdwet, kpdsno, kpdsnd, kpdzor, kpdalb, kpdais, kpdtg3, kpdscv, kpdacn, kpdsmc, kpdstc, kpdveg, kprvet, kpdsot, kpdalf, kpdvmn, kpdvmx, kpdslp, kpdabs, irttsf, irtwet, irtsno, irtzor, irtalb, irtais, irttg3, irtscv, irtacn, irtsmc, irtstc, irtveg, irtvet, irtsot, irtalf, irtvmn, irtvmx, irtslp, irtabs, imsk, jmsk, slmskh, outlat, outlon, gaus, blno, blto, me, lanom) |
Read analysis fields. More... | |
subroutine | anomint (tsfan0, tsfclm, tsfcl0, tsfanl, len) |
Add initial SST anomaly to date interpolated climatology. More... | |
subroutine | clima (lugb, iy, im, id, ih, fh, len, lsoil, slmask, fntsfc, fnwetc, fnsnoc, fnzorc, fnalbc, fnaisc, fntg3c, fnscvc, fnsmcc, fnstcc, fnacnc, fnvegc, fnvetc, fnsotc, fnvmnc, fnvmxc, fnslpc, fnabsc, tsfclm, tsfcl2, wetclm, snoclm, zorclm, albclm, aisclm, tg3clm, cvclm, cvbclm, cvtclm, cnpclm, smcclm, stcclm, sliclm, scvclm, acnclm, vegclm, vetclm, sotclm, alfclm, vmnclm, vmxclm, slpclm, absclm, kpdtsf, kpdwet, kpdsno, kpdzor, kpdalb, kpdais, kpdtg3, kpdscv, kpdacn, kpdsmc, kpdstc, kpdveg, kpdvet, kpdsot, kpdalf, tsfcl0, kpdvmn, kpdvmx, kpdslp, kpdabs, deltsfc, lanom, imsk, jmsk, slmskh, outlat, outlon, gaus, blno, blto, me, lprnt, iprnt, fnalbc2, ialb, tile_num_ch, i_index, j_index) |
Driver routine that reads in climatological data for a given time, and, if necessary, interpolates it to the model grid. More... | |
subroutine | count (slimsk, sno, ijmax) |
Counts the number of model points that are snow covered land, snow-free land, open water, sea ice, and snow covered sea ice. More... | |
subroutine | dayoyr (iyr, imo, idy, ldy) |
Compute day of the year based on month and day. More... | |
subroutine | filanl (tsfanl, tsfan2, wetanl, snoanl, zoranl, albanl, aisanl, tg3anl, cvanl, cvbanl, cvtanl, cnpanl, smcanl, stcanl, slianl, scvanl, veganl, vetanl, sotanl, alfanl, sihanl, sicanl, vmnanl, vmxanl, slpanl, absanl, tsfclm, tsfcl2, wetclm, snoclm, zorclm, albclm, aisclm, tg3clm, cvclm, cvbclm, cvtclm, cnpclm, smcclm, stcclm, sliclm, scvclm, vegclm, vetclm, sotclm, alfclm, sihclm, sicclm, vmnclm, vmxclm, slpclm, absclm, len, lsoil) |
Fill in analysis arrays with climatology before reading analysis data. More... | |
subroutine | filfcs (tsffcs, wetfcs, snofcs, zorfcs, albfcs, tg3fcs, cvfcs, cvbfcs, cvtfcs, cnpfcs, smcfcs, stcfcs, slifcs, aisfcs, vegfcs, vetfcs, sotfcs, alffcs, sihfcs, sicfcs, vmnfcs, vmxfcs, slpfcs, absfcs, tsfanl, wetanl, snoanl, zoranl, albanl, tg3anl, cvanl, cvbanl, cvtanl, cnpanl, smcanl, stcanl, slianl, aisanl, veganl, vetanl, sotanl, alfanl, sihanl, sicanl, vmnanl, vmxanl, slpanl, absanl, len, lsoil) |
Fill in model grid guess arrays with analysis values if this is a dead start. More... | |
subroutine | fixrda (lugb, fngrib, kpds5, slmask, iy, im, id, ih, fh, gdata, len, iret, imsk, jmsk, slmskh, gaus, blno, blto, outlat, outlon, me) |
Read in grib1 analysis data for the requested date and horizontally interpolate to the model grid. More... | |
subroutine | fixrdc (lugb, fngrib, kpds5, kpds7, mon, slmask, gdata, len, iret, imsk, jmsk, slmskh, gaus, blno, blto, outlat, outlon, me) |
Read in grib1 climatology data for a specific month and and horizontally interpolate to the model grid. More... | |
subroutine | fixrdc_tile (filename_raw, tile_num_ch, i_index, j_index, kpds, var, mon, npts, me) |
Reads in climatological data on the model grid tile for a given month. More... | |
subroutine | fixrdg (lugb, idim, jdim, fngrib, kpds5, gdata, gaus, blno, blto, me) |
Read a GRIB1 file. More... | |
subroutine | ga2la (gauin, imxin, jmxin, regout, imxout, jmxout, wlon, rnlat, rlnout, rltout, gaus, blno, blto) |
Interpolation from lat/lon or gaussian grid to a lat/lon grid. More... | |
subroutine | gaulat (gaul, k) |
Calculate gaussian latitudes. More... | |
subroutine | getarea (kgds, dlat, dlon, rslat, rnlat, wlon, elon, ijordr, me) |
For a given GRIB1 grid description section array, determine some grid specifications. More... | |
subroutine | getscv (snofld, scvfld, len) |
Set snow cover flag based on snow depth. More... | |
subroutine | getsmc (wetfld, len, lsoil, smcfld, me) |
Set soil moisture from soil wetness. More... | |
subroutine | getstc (tsffld, tg3fld, slifld, len, lsoil, stcfld, tsfimx) |
Set soil temperature and sea ice column temperature. More... | |
subroutine | hmskrd (lugb, imsk, jmsk, fnmskh, kpds5, slmskh, gausm, blnmsk, bltmsk, me) |
Read a high-resolution land mask. More... | |
subroutine | la2ga (regin, imxin, jmxin, rinlon, rinlat, rlon, rlat, inttyp, gauout, len, lmask, rslmsk, slmask, outlat, outlon, me) |
Interpolate data from a lat/lon grid to the model grid. More... | |
subroutine | landtyp (vegtype, soiltype, slptype, slmask, len) |
Set vegetation, soil and slope type at undefined model points. More... | |
subroutine | maxmin (f, imax, kmax) |
Compute the maxmimum and minimum of a field. More... | |
subroutine | merge (len, lsoil, iy, im, id, ih, fh, deltsfc, sihfcs, sicfcs, vmnfcs, vmxfcs, slpfcs, absfcs, tsffcs, wetfcs, snofcs, zorfcs, albfcs, aisfcs, cvfcs, cvbfcs, cvtfcs, cnpfcs, smcfcs, stcfcs, slifcs, vegfcs, vetfcs, sotfcs, alffcs, sihanl, sicanl, vmnanl, vmxanl, slpanl, absanl, tsfanl, tsfan2, wetanl, snoanl, zoranl, albanl, aisanl, cvanl, cvbanl, cvtanl, cnpanl, smcanl, stcanl, slianl, veganl, vetanl, sotanl, alfanl, ctsfl, calbl, caisl, csnol, csmcl, czorl, cstcl, cvegl, ctsfs, calbs, caiss, csnos, csmcs, czors, cstcs, cvegs, ccv, ccvb, ccvt, ccnp, cvetl, cvets, csotl, csots, calfl, calfs, csihl, csihs, csicl, csics, cvmnl, cvmns, cvmxl, cvmxs, cslpl, cslps, cabsl, cabss, irttsf, irtwet, irtsno, irtzor, irtalb, irtais, irttg3, irtscv, irtacn, irtsmc, irtstc, irtveg, irtvmn, irtvmx, irtslp, irtabs, irtvet, irtsot, irtalf, landice, me) |
Blend the model forecast (or first guess) fields with the analysis/climatology. More... | |
subroutine | monitr (lfld, fld, slimsk, sno, ijmax) |
Determine the maximum and minimum values of a surface field at snow-free and snow covered land, open water, and snow-free and snow covered sea ice. More... | |
subroutine | netcdf_err (error) |
Print the error message for a given netCDF return code. More... | |
subroutine | newice (slianl, slifcs, tsfanl, tsffcs, len, lsoil, sihnew, sicnew, sihanl, sicanl, albanl, snoanl, zoranl, smcanl, stcanl, albsea, snosea, zorsea, smcsea, smcice, tsfmin, tsfice, albice, zorice, tgice, rla, rlo, me) |
Adjust surface fields when ice melts or forms. More... | |
subroutine | qcbyfc (tsffcs, snofcs, qctsfs, qcsnos, qctsfi, len, lsoil, snoanl, aisanl, slianl, tsfanl, albanl, zoranl, smcanl, smcclm, tsfsmx, albomx, zoromx, me) |
Quality control analysis fields using the first guess. More... | |
subroutine | qcmxice (glacir, amxice, len, me) |
Quality control maximum ice extent. More... | |
subroutine | qcmxmn (ttl, fld, slimsk, sno, iceflg, fldlmx, fldlmn, fldomx, fldomn, fldimx, fldimn, fldjmx, fldjmn, fldsmx, fldsmn, epsfld, rla, rlo, len, mode, percrit, lgchek, me) |
Range check a field. More... | |
subroutine | qcsice (ais, glacir, amxice, aicice, aicsea, sllnd, slmask, rla, rlo, len, me) |
Check the sea ice cover mask against the land-sea mask. More... | |
subroutine | qcsli (slianl, slifcs, len, me) |
Check consistency between the forecast and analysis land-sea-ice mask. More... | |
subroutine | qcsnow (snoanl, slmask, aisanl, glacir, len, snoval, landice, me) |
Quality control snow at the model points. More... | |
subroutine | rof01 (aisfld, len, op, crit) |
Round a field up to one or down to zero. More... | |
subroutine | scale (fld, len, scl) |
Multiply a field by a scaling factor. More... | |
subroutine | setlsi (slmask, aisfld, len, aicice, slifld) |
Set land-sea-ice mask at sea ice. More... | |
subroutine | setrmsk (kpds5, slmask, igaul, jgaul, wlon, rnlat, data, imax, jmax, rlnout, rltout, lmask, rslmsk, gaus, blno, blto, kgds1, kpds4, lbms) |
Set the mask for the input data. More... | |
subroutine | setzro (fld, eps, len) |
Set a field to zero if it is less than a threshold. More... | |
subroutine | sfccycle (lugb, len, lsoil, sig1t, deltsfc, iy, im, id, ih, fh, rla, rlo, slmask, orog, orog_uf, use_ufo, nst_anl, sihfcs, sicfcs, sitfcs, swdfcs, slcfcs, vmnfcs, vmxfcs, slpfcs, absfcs, tsffcs, snofcs, zorfcs, albfcs, tg3fcs, cnpfcs, smcfcs, stcfcs, slifcs, aisfcs, vegfcs, vetfcs, sotfcs, alffcs, cvfcs, cvbfcs, cvtfcs, me, nlunit, sz_nml, input_nml_file, ialb, isot, ivegsrc, tile_num_ch, i_index, j_index) |
Surface cycling driver routine. More... | |
subroutine | snodpth (scvanl, slianl, tsfanl, snoclm, glacir, snwmax, snwmin, landice, len, snoanl, me) |
Estimate snow depth at glacial, land and sea ice points. More... | |
subroutine | snodpth2 (glacir, snwmax, snoanl, len, me) |
Ensure deep snow pack at permanent glacial points. More... | |
subroutine | snosfc (snoanl, tsfanl, tsfsmx, len, me) |
Check skin temperature at points with snow. More... | |
subroutine | subst (data, imax, jmax, dlon, dlat, ijordr) |
Take an array of data on a lat/lon based grid and rearrange it so the corner point is in the 'lower left' and adjacent points in the 'i' direction are consecutive. More... | |
subroutine | tsfcor (tsfc, orog, slmask, umask, len, rlapse) |
Adjust skin temperature or SST for terrain. More... | |
subroutine | usesgt (sig1t, slianl, tg3anl, len, lsoil, tsfanl, stcanl, tsfimx) |
Set soil temperature and sea ice column temperature for a dead start. More... | |
This is a limited point version of surface program.
Definition in file sfcsub.F.
module sfccyc_module |
This program runs in two different modes:
if the date of the analysis does not match given iy,im,id,ih, (and fh), the program searches an old analysis by going back 6 hours, then 12 hours, then one day upto nrepmx days (parameter statement in the subrotine fixrd. now defined as 8). this allows the user to provide non-daily analysis to be used. if matching field is not found, the forecast guess will be used.
use of a combined earlier surface analyses and current analysis is not allowed (as was done in the old version for snow analysis in which old snow analysis is used in combination with initial guess), except for sea surface temperature. for sst anolmaly interpolation, you need to set lanom=.true. and must provide sst analysis at initial time.
if you want to do complex merging of past and present surface field analysis, you need to create a separate file that contains daily surface field.
for a dead start, do not supply fnbgsi or set fnbgsi=' '
subroutine albocn | ( | real (kind=kind_io8), dimension(len,4) | albclm, |
real (kind=kind_io8), dimension(len) | slmask, | ||
real (kind=kind_io8) | albomx, | ||
integer | len | ||
) |
Set the albedo at open water points.
[in,out] | albclm | Albedo. |
[in] | slmask | The land-sea-ice mask. |
[in] | albomx | The albedo at open water. |
[in] | len | Number of model points to process. |
Definition at line 6625 of file sfcsub.F.
Referenced by sfccycle().
subroutine analy | ( | lugb, | |
iy, | |||
im, | |||
id, | |||
ih, | |||
real (kind=kind_io8) | fh, | ||
len, | |||
lsoil, | |||
real (kind=kind_io8), dimension(len) | slmask, | ||
character*500 | fntsfa, | ||
character*500 | fnweta, | ||
character*500 | fnsnoa, | ||
character*500 | fnzora, | ||
character*500 | fnalba, | ||
character*500 | fnaisa, | ||
character*500 | fntg3a, | ||
character*500 | fnscva, | ||
character*500 | fnsmca, | ||
character*500 | fnstca, | ||
character*500 | fnacna, | ||
character*500 | fnvega, | ||
character*500 | fnveta, | ||
character*500 | fnsota, | ||
fnvmna, | |||
fnvmxa, | |||
fnslpa, | |||
fnabsa, | |||
real (kind=kind_io8), dimension(len) | tsfanl, | ||
real (kind=kind_io8), dimension(len) | wetanl, | ||
real (kind=kind_io8), dimension(len) | snoanl, | ||
real (kind=kind_io8), dimension(len) | zoranl, | ||
real (kind=kind_io8), dimension(len,4) | albanl, | ||
real (kind=kind_io8), dimension(len) | aisanl, | ||
real (kind=kind_io8), dimension(len) | tg3anl, | ||
real (kind=kind_io8), dimension (len) | cvanl, | ||
real (kind=kind_io8), dimension(len) | cvbanl, | ||
real (kind=kind_io8), dimension(len) | cvtanl, | ||
real (kind=kind_io8), dimension(len,lsoil) | smcanl, | ||
real (kind=kind_io8), dimension(len,lsoil) | stcanl, | ||
real (kind=kind_io8), dimension(len) | slianl, | ||
real (kind=kind_io8), dimension(len) | scvanl, | ||
real (kind=kind_io8), dimension(len) | acnanl, | ||
real (kind=kind_io8), dimension(len) | veganl, | ||
real (kind=kind_io8), dimension(len) | vetanl, | ||
real (kind=kind_io8), dimension(len) | sotanl, | ||
real (kind=kind_io8), dimension(len,2) | alfanl, | ||
real (kind=kind_io8), dimension(len) | tsfan0, | ||
vmnanl, | |||
vmxanl, | |||
slpanl, | |||
absanl, | |||
kpdtsf, | |||
kpdwet, | |||
kpdsno, | |||
kpdsnd, | |||
kpdzor, | |||
integer, dimension(4) | kpdalb, | ||
kpdais, | |||
kpdtg3, | |||
kpdscv, | |||
kpdacn, | |||
kpdsmc, | |||
kpdstc, | |||
kpdveg, | |||
kprvet, | |||
kpdsot, | |||
integer, dimension(2) | kpdalf, | ||
kpdvmn, | |||
kpdvmx, | |||
kpdslp, | |||
kpdabs, | |||
irttsf, | |||
irtwet, | |||
integer | irtsno, | ||
integer | irtzor, | ||
integer | irtalb, | ||
integer | irtais, | ||
integer | irttg3, | ||
integer | irtscv, | ||
integer | irtacn, | ||
integer | irtsmc, | ||
integer | irtstc, | ||
integer | irtveg, | ||
integer | irtvet, | ||
integer | irtsot, | ||
integer | irtalf, | ||
irtvmn, | |||
irtvmx, | |||
irtslp, | |||
irtabs, | |||
imsk, | |||
jmsk, | |||
real (kind=kind_io8), dimension(imsk,jmsk) | slmskh, | ||
real (kind=kind_io8), dimension(len) | outlat, | ||
real (kind=kind_io8), dimension(len) | outlon, | ||
logical | gaus, | ||
real (kind=kind_io8) | blno, | ||
real (kind=kind_io8) | blto, | ||
integer | me, | ||
logical | lanom | ||
) |
Read analysis fields.
[in] | lugb | Fortran unit number for analysis files. |
[in] | iy | Cycle year. |
[in] | im | Cycle month. |
[in] | id | Cycle day. |
[in] | ih | Cycle hour. |
[in] | fh | Forecast hour. |
[in] | len | Number of model points to process. |
[in] | lsoil | Number of soil layers. |
[in] | slmask | Model land-sea mask. |
[in] | fntsfa | SST analysis file. |
[in] | fnweta | Soil wetness analysis file. |
[in] | fnsnoa | Snow analysis file |
[in] | fnzora | Roughness length analysis file. |
[in] | fnalba | Snow-free albedo analysis file. |
[in] | fnaisa | Sea ice mask analysis file. |
[in] | fntg3a | Soil substrate analysis file. |
[in] | fnscva | Snow cover analysis file. |
[in] | fnsmca | Soil moisture analysis file. |
[in] | fnstca | Soil temperature analysis file. |
[in] | fnacna | Sea ice concentration analysis file. |
[in] | fnvega | Vegetation greenness analysis file. |
[in] | fnveta | Vegetation type analysis file. |
[in] | fnsota | Soil type analysis file. |
[in] | fnvmna | Minimum vegetation greenness analysis file. |
[in] | fnvmxa | Maximum vegetation greenness analysis file. |
[in] | fnslpa | Slope type analysis file. |
[in] | fnabsa | Maximum snow albedo analysis file. |
[out] | tsfanl | Skin temperature/SST analysis on model grid. |
[out] | wetanl | Soil wetness analysis on model grid. |
[out] | snoanl | Snow analysis on model grid. |
[out] | zoranl | Roughness length analysis on model grid. |
[out] | albanl | Snow-free albedo analysis on model grid. |
[out] | aisanl | Sea ice mask analysis on model grid. |
[out] | tg3anl | Soil substrate analysis on model grid. |
[out] | cvanl | Convective cloud cover analysis on model grid. |
[out] | cvbanl | Convective cloud base analysis on model grid. |
[out] | cvtanl | Convective cloud top analysis on model grid. |
[out] | smcanl | Soil moisture analysis on model grid. |
[out] | stcanl | Soil temperature analysis on model grid. |
[out] | slianl | Not used. |
[out] | scvanl | Snow cover analysis on model grid. |
[out] | acnanl | Sea ice concentration analysis on model grid. |
[out] | veganl | Vegetation greenness analysis on model grid. |
[out] | vetanl | Vegetation type analysis on model grid. |
[out] | sotanl | Soil type analysis on model grid. |
[out] | alfanl | Analysis of fraction for strongly and weakly zenith angle dependent albedo on model grid. |
[out] | tsfan0 | SST analysis at forecast hour 0 on model grid. |
[out] | vmnanl | Minimum vegetation greenness analysis on model grid. |
[out] | vmxanl | Maximum vegetation greenness analysis on model grid. |
[out] | slpanl | Slope type analysis on model grid. |
[out] | absanl | Maximum snow albedo analysis on model grid. |
[in] | kpdtsf | Grib parameter number of skin temperature/SST. |
[in] | kpdwet | Grib parameter number of soil wetness. |
[in] | kpdsno | Grib parameter number of liquid equivalent snow depth. |
[in] | kpdsnd | Grib parameter number of physical snow depth. |
[in] | kpdzor | Grib parameter number of roughness length. |
[in] | kpdalb | Grib parameter number of snow-free albedo. |
[in] | kpdais | Grib parameter number of sea ice mask. |
[in] | kpdtg3 | Grib parameter number of soil substrate temperature. |
[in] | kpdscv | Grib parameter number of snow cover. |
[in] | kpdacn | Grib parameter number of sea ice concentration. |
[in] | kpdsmc | Grib parameter number of soil moisture. |
[in] | kpdstc | Grib parameter number of soil temperature. |
[in] | kpdveg | Grib parameter number of vegetation greenness. |
[in] | kprvet | Grib parameter number of vegetation type. |
[in] | kpdsot | Grib parameter number of soil type. |
[in] | kpdalf | Grib parameter number for fraction for strongly and weakly zenith angle dependent albedo. |
[in] | kpdvmn | Grib parameter number of minimum vegetation greenness. |
[in] | kpdvmx | Grib parameter number of maximum vegetation greenness. |
[in] | kpdslp | Grib parameter number of slope type. |
[in] | kpdabs | Grib parameter number of maximum snow albedo. |
[out] | irttsf | Return code from read of skin temperature/SST analysis file. |
[out] | irtwet | Return code from read of soil wetness analysis file. |
[out] | irtsno | Return code from read of snow analysis file. |
[out] | irtzor | Return code from read of roughness length file. |
[out] | irtalb | Return code from read of snow-free albedo analysis file. |
[out] | irtais | Return code from read of ice mask analysis file. |
[out] | irttg3 | Return code from read of soil substrate temperature analysis file. |
[out] | irtscv | Return code from read of snow cover analysis file. |
[out] | irtacn | Return code from read of sea ice concentration analysis file. |
[out] | irtsmc | Return code from read of soil moisture analysis file. |
[out] | irtstc | Return code from read of soil temperature analysis file. |
[out] | irtveg | Return code from read of vegetation greenness analysis file. |
[out] | irtvet | Return code from read of vegetation type analysis file. |
[out] | irtsot | Return code from read of soil type analysis file. |
[out] | irtalf | Return code from read of file containing fraction for strongly and weakly zenith angle dependent albedo. |
[out] | irtvmn | Return code from read of minimum vegetation greenness analysis file. |
[out] | irtvmx | Return code from read of maximum vegetation greenness analysis file. |
[out] | irtslp | Return code from read of slope type analysis file. |
[out] | irtabs | Return code from read of maximum snow albedo analysis file. |
[in] | imsk | 'i' dimension of the high-res mask used for analysis data without a bitmap. |
[in] | jmsk | 'j' dimension of the high-res mask used for analysis data without a bitmap. |
[in] | slmskh | The high-resolution mask used for analysis data without a bitmap. |
[in] | outlat | Model latitudes |
[in] | outlon | Model longitudes |
[in] | gaus | When true, the high-res mask is on a gaussian grid. |
[in] | blno | Corner point longitude of the high-res mask. |
[in] | blto | Corner point latitude of the high-res mask. |
[in] | me | MPI task number. |
[in] | lanom | When true, do sst anomaly interpolation. |
Definition at line 4072 of file sfcsub.F.
References fixrda().
Referenced by sfccycle().
subroutine anomint | ( | real (kind=kind_io8), dimension(len) | tsfan0, |
real (kind=kind_io8), dimension(len) | tsfclm, | ||
real (kind=kind_io8), dimension(len) | tsfcl0, | ||
real (kind=kind_io8), dimension(len) | tsfanl, | ||
integer | len | ||
) |
Add initial SST anomaly to date interpolated climatology.
[in] | tsfan0 | Skin temperature/SST analysis at initial time. |
[in] | tsfclm | Skin temperature/SST climatology. |
[in] | tsfcl0 | Skin temperature/SST climatology at initial time. |
[out] | tsfanl | Updated skin temperature/SST. |
[in] | len | Number of model points to process. |
Definition at line 7699 of file sfcsub.F.
Referenced by sfccycle().
subroutine clima | ( | integer | lugb, |
integer | iy, | ||
integer | im, | ||
integer | id, | ||
integer | ih, | ||
real (kind=kind_io8) | fh, | ||
integer | len, | ||
integer | lsoil, | ||
real (kind=kind_io8), dimension(len) | slmask, | ||
character*500 | fntsfc, | ||
character*500 | fnwetc, | ||
character*500 | fnsnoc, | ||
character*500 | fnzorc, | ||
character*500 | fnalbc, | ||
character*500 | fnaisc, | ||
character*500 | fntg3c, | ||
character*500 | fnscvc, | ||
character*500 | fnsmcc, | ||
character*500 | fnstcc, | ||
character*500 | fnacnc, | ||
character*500 | fnvegc, | ||
character*500 | fnvetc, | ||
character*500 | fnsotc, | ||
character*500 | fnvmnc, | ||
character*500 | fnvmxc, | ||
character*500 | fnslpc, | ||
character*500 | fnabsc, | ||
real (kind=kind_io8), dimension(len) | tsfclm, | ||
real (kind=kind_io8), dimension(len) | tsfcl2, | ||
real (kind=kind_io8), dimension(len) | wetclm, | ||
real (kind=kind_io8), dimension(len) | snoclm, | ||
real (kind=kind_io8), dimension(len) | zorclm, | ||
real (kind=kind_io8), dimension(len,4) | albclm, | ||
real (kind=kind_io8), dimension(len) | aisclm, | ||
real (kind=kind_io8), dimension(len) | tg3clm, | ||
real (kind=kind_io8), dimension (len) | cvclm, | ||
real (kind=kind_io8), dimension(len) | cvbclm, | ||
real (kind=kind_io8), dimension(len) | cvtclm, | ||
real (kind=kind_io8), dimension(len) | cnpclm, | ||
real (kind=kind_io8), dimension(len,lsoil) | smcclm, | ||
real (kind=kind_io8), dimension(len,lsoil) | stcclm, | ||
real (kind=kind_io8), dimension(len) | sliclm, | ||
real (kind=kind_io8), dimension(len) | scvclm, | ||
real (kind=kind_io8), dimension(len) | acnclm, | ||
real (kind=kind_io8), dimension(len) | vegclm, | ||
real (kind=kind_io8), dimension(len) | vetclm, | ||
real (kind=kind_io8), dimension(len) | sotclm, | ||
real (kind=kind_io8), dimension(len,2) | alfclm, | ||
real (kind=kind_io8), dimension(len) | vmnclm, | ||
real (kind=kind_io8), dimension(len) | vmxclm, | ||
real (kind=kind_io8), dimension(len) | slpclm, | ||
real (kind=kind_io8), dimension(len) | absclm, | ||
integer | kpdtsf, | ||
integer | kpdwet, | ||
integer | kpdsno, | ||
integer | kpdzor, | ||
integer, dimension(4) | kpdalb, | ||
integer | kpdais, | ||
integer | kpdtg3, | ||
integer | kpdscv, | ||
integer | kpdacn, | ||
integer | kpdsmc, | ||
integer | kpdstc, | ||
integer | kpdveg, | ||
integer | kpdvet, | ||
integer | kpdsot, | ||
integer, dimension(2) | kpdalf, | ||
real (kind=kind_io8), dimension(len) | tsfcl0, | ||
integer | kpdvmn, | ||
integer | kpdvmx, | ||
integer | kpdslp, | ||
integer | kpdabs, | ||
real (kind=kind_io8) | deltsfc, | ||
logical | lanom, | ||
integer | imsk, | ||
integer | jmsk, | ||
real (kind=kind_io8), dimension(imsk,jmsk) | slmskh, | ||
real (kind=kind_io8), dimension(len) | outlat, | ||
real (kind=kind_io8), dimension(len) | outlon, | ||
logical | gaus, | ||
real (kind=kind_io8) | blno, | ||
real (kind=kind_io8) | blto, | ||
integer | me, | ||
logical | lprnt, | ||
integer | iprnt, | ||
character*500 | fnalbc2, | ||
integer | ialb, | ||
character(len=*), intent(in) | tile_num_ch, | ||
integer, dimension(len), intent(in) | i_index, | ||
integer, dimension(len), intent(in) | j_index | ||
) |
Driver routine that reads in climatological data for a given time, and, if necessary, interpolates it to the model grid.
[in] | lugb | Fortran unit number of GRIB1 climatological data. |
[in] | iy | Cycle year. |
[in] | im | Cycle month. |
[in] | id | Cycle day. |
[in] | ih | Cycle hour. |
[in] | fh | Forecast hour. |
[in] | len | Number of model points to process. |
[in] | lsoil | Number of soil layers. |
[in] | slmask | Model land-sea-ice mask. |
[in] | fntsfc | Climatological SST file. |
[in] | fnwetc | Climatological soil wetness file. |
[in] | fnsnoc | Climatological snow depth file. |
[in] | fnzorc | Climatological roughness length file. Or 'igbp' to use igbp vegetation type lookup table. Or 'sib' to use sib vegeation type lookup table. |
[in] | fnalbc | Climatological snow-free albedo file. |
[in] | fnaisc | Climatological sea ice mask file. |
[in] | fntg3c | Climatological soil substrate temperature file. |
[in] | fnscvc | Climatological snow cover file. |
[in] | fnsmcc | Climatological soil moisture file. |
[in] | fnstcc | Climatological soil temperature file. |
[in] | fnacnc | Climatological sea ice concentration file. |
[in] | fnvegc | Climatological vegetation greenness file. |
[in] | fnvetc | Climatological vegetation type file. |
[in] | fnsotc | Climatological soil type file. |
[in] | fnvmnc | Climatological minimum vegetation greenness file. |
[in] | fnvmxc | Climatological maximum vegetation greenness file. |
[in] | fnslpc | Climatological slope type file. |
[in] | fnabsc | Climatological maximum snow albedo file. |
[out] | tsfclm | Climatological skin temperature/SST on model grid. |
[out] | tsfcl2 | Climatological skin temperature/SST on model grid at time minus deltsfc. |
[out] | wetclm | Climatological soil wetness on model grid. |
[out] | snoclm | Climatological liquid equivalent snow depth on model grid. |
[out] | zorclm | Climatological roughness length on model grid. |
[out] | albclm | Climatological snow-free albedo on model grid. |
[out] | aisclm | Climatological sea ice mask on model grid. |
[out] | tg3clm | Climatological soil substrate temperature on model grid. |
[out] | cvclm | Climatological convective cloud cover on model grid. |
[out] | cvbclm | Climatological convective cloud base on model grid. |
[out] | cvtclm | Climatological convective cloud top on model grid. |
[out] | cnpclm | Climatological canopy water content on model grid. |
[out] | smcclm | Climatological soil moisture on model grid. |
[out] | stcclm | Climatologcial soil temperature on model grid. |
[out] | sliclm | Climatological model land-sea-ice mask. |
[out] | scvclm | Climatological snow cover on model grid. |
[out] | acnclm | Climatological sea ice concentration on model grid. |
[out] | vegclm | Climatological vegetation greenness on model grid. |
[out] | vetclm | Climatological vegetation type on model grid. |
[out] | sotclm | Climatological soil type on model grid. |
[out] | alfclm | Climatological fraction for strongly and weakly zenith angle dependent albedo on model grid. |
[out] | vmnclm | Climatological minimum vegetation greenness on model grid. |
[out] | vmxclm | Climatological maximum vegetation greenness on model grid. |
[out] | slpclm | Climatological slope type on model grid. |
[out] | absclm | Climatological maximum snow albedo on model grid. |
[in] | kpdtsf | Grib parameter number of skin temperature/SST. |
[in] | kpdwet | Grib parameter number of soil wetness. |
[in] | kpdsno | Grib parameter number of liquid equivalent snow depth. |
[in] | kpdzor | Grib parameter number of roughness length. |
[in] | kpdalb | Grib parameter number of snow-free albedo. |
[in] | kpdais | Grib parameter number of sea ice mask. |
[in] | kpdtg3 | Grib parameter number of soil substrate temperature. |
[in] | kpdscv | Grib parameter number of snow cover. |
[in] | kpdacn | Grib parameter number of sea ice concentration. |
[in] | kpdsmc | Grib parameter number of soil moisture. |
[in] | kpdstc | Grib parameter number of soil temperature. |
[in] | kpdveg | Grib parameter number of vegetation greenness. |
[in] | kpdvet | Grib parameter number of vegetation type. |
[in] | kpdsot | Grib parameter number of soil type. |
[in] | kpdalf | Grib parameter number for fraction for strongly and weakly zenith angle dependent albedo. |
[in] | tsfcl0 | Climatological SST at forecast hour 0. |
[in] | kpdvmn | Grib parameter number of minimum vegetation greenness. |
[in] | kpdvmx | Grib parameter number of maximum vegetation greenness. |
[in] | kpdslp | Grib parameter number of slope type. |
[in] | kpdabs | Grib parameter number of maximum snow albedo. |
[in] | deltsfc | Cycling frequency in hours. |
[in] | lanom | When true, do sst anomaly interpolation. |
[in] | imsk | 'i' dimension of the high-res mask used for climatological data without a bitmap. |
[in] | jmsk | 'j' dimension of the high-res mask used for climatological data without a bitmap. |
[in] | slmskh | The high-resolution mask used for climatological data without a bitmap. |
[in] | outlat | Model latitudes |
[in] | outlon | Model longitudes |
[in] | gaus | When true, the high-res mask is on a gaussian grid. |
[in] | blno | Corner point longitude of the high-res mask. |
[in] | blto | Corner point latitude of the high-res mask. |
[in] | me | MPI task number. |
[in] | lprnt | Turn of diagnostic print. |
[in] | iprnt | Index of diagnotic print point. |
[in] | fnalbc2 | File containing climatological fraction for strongly and weakly zenith angle dependent albedo. |
[in] | ialb | Use modis albedo when '1'. Use brigleb when '0'. |
[in] | tile_num_ch | Model tile number to process. |
[in] | i_index | The 'i' indices of the model grid to process. |
[in] | j_index | The 'j' indices of the model grid to process. |
Definition at line 7828 of file sfcsub.F.
References fixrdc(), and fixrdc_tile().
Referenced by sfccycle().
subroutine count | ( | real (kind=kind_io8), dimension(1) | slimsk, |
real (kind=kind_io8), dimension(1) | sno, | ||
integer | ijmax | ||
) |
Counts the number of model points that are snow covered land, snow-free land, open water, sea ice, and snow covered sea ice.
[in] | slimsk | The land-sea-ice mask. |
[in] | sno | Snow |
[in] | ijmax | Number of model points to process. |
Definition at line 2626 of file sfcsub.F.
Referenced by fixrdc_tile(), and read_write_data::read_salclm_gfs_nc().
subroutine dayoyr | ( | integer | iyr, |
integer | imo, | ||
integer | idy, | ||
integer | ldy | ||
) |
subroutine filanl | ( | real (kind=kind_io8), dimension(len) | tsfanl, |
real (kind=kind_io8), dimension(len) | tsfan2, | ||
real (kind=kind_io8), dimension(len) | wetanl, | ||
real (kind=kind_io8), dimension(len) | snoanl, | ||
real (kind=kind_io8), dimension(len) | zoranl, | ||
real (kind=kind_io8), dimension(len,4) | albanl, | ||
real (kind=kind_io8), dimension(len) | aisanl, | ||
real (kind=kind_io8), dimension(len) | tg3anl, | ||
real (kind=kind_io8), dimension (len) | cvanl, | ||
real (kind=kind_io8), dimension(len) | cvbanl, | ||
real (kind=kind_io8), dimension(len) | cvtanl, | ||
real (kind=kind_io8), dimension(len) | cnpanl, | ||
real (kind=kind_io8), dimension(len,lsoil) | smcanl, | ||
real (kind=kind_io8), dimension(len,lsoil) | stcanl, | ||
real (kind=kind_io8), dimension(len) | slianl, | ||
real (kind=kind_io8), dimension(len) | scvanl, | ||
real (kind=kind_io8), dimension(len) | veganl, | ||
real (kind=kind_io8), dimension(len) | vetanl, | ||
real (kind=kind_io8), dimension(len) | sotanl, | ||
real (kind=kind_io8), dimension(len,2) | alfanl, | ||
sihanl, | |||
sicanl, | |||
vmnanl, | |||
vmxanl, | |||
slpanl, | |||
absanl, | |||
real (kind=kind_io8), dimension(len) | tsfclm, | ||
real (kind=kind_io8), dimension(len) | tsfcl2, | ||
real (kind=kind_io8), dimension(len) | wetclm, | ||
real (kind=kind_io8), dimension(len) | snoclm, | ||
real (kind=kind_io8), dimension(len) | zorclm, | ||
real (kind=kind_io8), dimension(len,4) | albclm, | ||
real (kind=kind_io8), dimension(len) | aisclm, | ||
real (kind=kind_io8), dimension(len) | tg3clm, | ||
real (kind=kind_io8), dimension (len) | cvclm, | ||
real (kind=kind_io8), dimension(len) | cvbclm, | ||
real (kind=kind_io8), dimension(len) | cvtclm, | ||
real (kind=kind_io8), dimension(len) | cnpclm, | ||
real (kind=kind_io8), dimension(len,lsoil) | smcclm, | ||
real (kind=kind_io8), dimension(len,lsoil) | stcclm, | ||
real (kind=kind_io8), dimension(len) | sliclm, | ||
real (kind=kind_io8), dimension(len) | scvclm, | ||
real (kind=kind_io8), dimension(len) | vegclm, | ||
real (kind=kind_io8), dimension(len) | vetclm, | ||
real (kind=kind_io8), dimension(len) | sotclm, | ||
real (kind=kind_io8), dimension(len,2) | alfclm, | ||
sihclm, | |||
sicclm, | |||
vmnclm, | |||
vmxclm, | |||
slpclm, | |||
absclm, | |||
integer | len, | ||
integer | lsoil | ||
) |
Fill in analysis arrays with climatology before reading analysis data.
[out] | tsfanl | Skin temperature/SST analysis on model grid. |
[out] | tsfan2 | Skin temperature/SST analysis on model grid at time minus deltsfc. |
[out] | wetanl | Soil wetness analysis on model grid. |
[out] | snoanl | Liquid equivalent snow depth analysis on model grid. |
[out] | zoranl | Roughness length analysis on model grid. |
[out] | albanl | Snow-free albedo analysis on model grid. |
[out] | aisanl | Sea ice mask analysis on model grid. |
[out] | tg3anl | Soil substrate temperature analysis on model grid. |
[out] | cvanl | Convective cloud cover analysis on model grid. |
[out] | cvbanl | Convective cloud base analysis on model grid. |
[out] | cvtanl | Convective cloud top analysis on model grid. |
[out] | cnpanl | Canopy water content analysis on model grid. |
[out] | smcanl | Soil moisture analysis on model grid. |
[out] | stcanl | Soil temperature analysis on model grid. |
[out] | slianl | Land-sea-ice mask analysis on model grid. |
[out] | scvanl | Snow cover analysis on model grid. |
[out] | veganl | Vegetation greenness analysis on model grid. |
[out] | vetanl | Vegetation type analysis on model grid. |
[out] | sotanl | Soil type analysis on model grid. |
[out] | alfanl | Fraction for strongly and weakly zenith angle dependent albedo analysis on model grid. |
[out] | sihanl | Sea ice depth analysis on model grid. |
[out] | sicanl | Sea ice concentration analysis on model grid. |
[out] | vmnanl | Minimum vegetation greenness analysis on model grid. |
[out] | vmxanl | Maximum vegetation greenness analysis on model grid. |
[out] | slpanl | Slope type analysis on model grid. |
[out] | absanl | Maximum snow albedo analysis on model grid. |
[in] | tsfclm | Climatological skin temperature/SST on model grid. |
[in] | tsfcl2 | Climatological skin temperature/SST on model grid at time minus deltsfc. |
[in] | wetclm | Climatological soil wetness on model grid. |
[in] | snoclm | Climatological liquid equivalent snow depth on model grid. |
[in] | zorclm | Climatological roughness length on model grid. |
[in] | albclm | Climatological snow-free albedo on model grid. |
[in] | aisclm | Climatological sea ice mask on model grid. |
[in] | tg3clm | Climatological soil substrate temperature on model grid. |
[in] | cvclm | Climatological convective cloud cover on model grid. |
[in] | cvbclm | Climatological convective cloud base on model grid. |
[in] | cvtclm | Climatological convective cloud top on model grid. |
[in] | cnpclm | Climatological canopy water content on model grid. |
[in] | smcclm | Climatological soil moisture on model grid. |
[in] | stcclm | Climatologcial soil temperature on model grid. |
[in] | sliclm | Climatological model land-sea-ice mask. |
[in] | scvclm | Climatological snow cover on model grid. |
[in] | vegclm | Climatological vegetation greenness on model grid. |
[in] | vetclm | Climatological vegetation type on model grid. |
[in] | sotclm | Climatological soil type on model grid. |
[in] | alfclm | Climatological fraction for strongly and weakly zenith angle dependent albedo on model grid. |
[in] | sihclm | Climatological sea ice depth on the model grid. |
[in] | sicclm | Climatological sea ice concentration on the model grid. |
[in] | vmnclm | Climatological minimum vegetation greenness on model grid. |
[in] | vmxclm | Climatological maximum vegetation greenness on model grid. |
[in] | slpclm | Climatological slope type on model grid. |
[in] | absclm | Climatological maximum snow albedo on model grid. |
[in] | len | Number of model points to process. |
[in] | lsoil | Number of soil layers. |
Definition at line 3852 of file sfcsub.F.
Referenced by sfccycle().
subroutine filfcs | ( | real (kind=kind_io8), dimension(len) | tsffcs, |
real (kind=kind_io8), dimension(len) | wetfcs, | ||
real (kind=kind_io8), dimension(len) | snofcs, | ||
real (kind=kind_io8), dimension(len) | zorfcs, | ||
real (kind=kind_io8), dimension(len,4) | albfcs, | ||
real (kind=kind_io8), dimension(len) | tg3fcs, | ||
real (kind=kind_io8), dimension (len) | cvfcs, | ||
real (kind=kind_io8), dimension(len) | cvbfcs, | ||
real (kind=kind_io8), dimension(len) | cvtfcs, | ||
real (kind=kind_io8), dimension(len) | cnpfcs, | ||
real (kind=kind_io8), dimension(len,lsoil) | smcfcs, | ||
real (kind=kind_io8), dimension(len,lsoil) | stcfcs, | ||
real (kind=kind_io8), dimension(len) | slifcs, | ||
real (kind=kind_io8), dimension(len) | aisfcs, | ||
real (kind=kind_io8), dimension(len) | vegfcs, | ||
real (kind=kind_io8), dimension(len) | vetfcs, | ||
real (kind=kind_io8), dimension(len) | sotfcs, | ||
real (kind=kind_io8), dimension(len,2) | alffcs, | ||
sihfcs, | |||
sicfcs, | |||
vmnfcs, | |||
vmxfcs, | |||
slpfcs, | |||
absfcs, | |||
real (kind=kind_io8), dimension(len) | tsfanl, | ||
real (kind=kind_io8), dimension(len) | wetanl, | ||
real (kind=kind_io8), dimension(len) | snoanl, | ||
real (kind=kind_io8), dimension(len) | zoranl, | ||
real (kind=kind_io8), dimension(len,4) | albanl, | ||
real (kind=kind_io8), dimension(len) | tg3anl, | ||
real (kind=kind_io8), dimension (len) | cvanl, | ||
real (kind=kind_io8), dimension(len) | cvbanl, | ||
real (kind=kind_io8), dimension(len) | cvtanl, | ||
real (kind=kind_io8), dimension(len) | cnpanl, | ||
real (kind=kind_io8), dimension(len,lsoil) | smcanl, | ||
real (kind=kind_io8), dimension(len,lsoil) | stcanl, | ||
real (kind=kind_io8), dimension(len) | slianl, | ||
real (kind=kind_io8), dimension(len) | aisanl, | ||
real (kind=kind_io8), dimension(len) | veganl, | ||
real (kind=kind_io8), dimension(len) | vetanl, | ||
real (kind=kind_io8), dimension(len) | sotanl, | ||
real (kind=kind_io8), dimension(len,2) | alfanl, | ||
sihanl, | |||
sicanl, | |||
vmnanl, | |||
vmxanl, | |||
slpanl, | |||
absanl, | |||
integer | len, | ||
integer | lsoil | ||
) |
Fill in model grid guess arrays with analysis values if this is a dead start.
All fields are on the model grid.
[out] | tsffcs | First guess skin temperature/SST analysis. |
[out] | wetfcs | First guess soil wetness. |
[out] | snofcs | First guess liquid equivalent snow depth. |
[out] | zorfcs | First guess roughness length. |
[out] | albfcs | First guess snow-free albedo. |
[out] | tg3fcs | First guess soil substrate temperature. |
[out] | cvfcs | First guess cloud cover. |
[out] | cvbfcs | First guess convective cloud bottom. |
[out] | cvtfcs | First guess convective cloud top. |
[out] | cnpfcs | First guess canopy water content. |
[out] | smcfcs | First guess soil moisture. |
[out] | stcfcs | First guess soil temperature. |
[out] | slifcs | First guess land-sea-ice mask. |
[out] | aisfcs | First guess sea ice mask. |
[out] | vegfcs | First guess vegetation greenness. |
[out] | vetfcs | First guess vegetation type. |
[out] | sotfcs | First guess soil type. |
[out] | alffcs | First guess of fraction for strongly and weakly zenith angle dependent albedo. |
[out] | sihfcs | First guess sea ice depth. |
[out] | sicfcs | First guess sea ice concentration. |
[out] | vmnfcs | First guess minimum greenness fraction. |
[out] | vmxfcs | First guess maximum greenness fraction. |
[out] | slpfcs | First guess slope type. |
[out] | absfcs | First guess maximum snow albedo analysis. |
[in] | tsfanl | Skin temperature/SST analysis. |
[in] | wetanl | Soil wetness analysis. |
[in] | snoanl | Liquid equivalent snow depth analysis. |
[in] | zoranl | Roughness length analysis. |
[in] | albanl | Snow-free albedo analysis. |
[in] | tg3anl | Soil substrate temperature analysis. |
[in] | cvanl | Convective cloud cover analysis. |
[in] | cvbanl | Convective cloud base analysis. |
[in] | cvtanl | Convective cloud top analysis. |
[in] | cnpanl | Canopy water content analysis. |
[in] | smcanl | Soil moisture analysis. |
[in] | stcanl | Soil temperature analysis. |
[in] | slianl | Land-sea-ice mask analysis. |
[in] | aisanl | Sea ice mask analysis. |
[in] | veganl | Vegetation greenness analysis. |
[in] | vetanl | Vegetation type analysis. |
[in] | sotanl | Soil type analysis. |
[in] | alfanl | Analysis of fraction for strongly and weakly zenith angle dependent albedo. |
[in] | sihanl | Sea ice depth analysis. |
[in] | sicanl | Sea ice concentration analysis. |
[in] | vmnanl | Minimum greenness fraction analysis. |
[in] | vmxanl | Maximum greenness fraction analysis. |
[in] | slpanl | Slope type analysis. |
[in] | absanl | Maximum snow albedo analysis. |
[in] | len | Number of model points to process. |
[in] | lsoil | Number of soil layers. |
Definition at line 4793 of file sfcsub.F.
Referenced by sfccycle().
subroutine fixrda | ( | integer | lugb, |
character*500 | fngrib, | ||
integer | kpds5, | ||
real (kind=kind_io8), dimension(len) | slmask, | ||
integer | iy, | ||
integer | im, | ||
integer | id, | ||
integer | ih, | ||
real (kind=kind_io8) | fh, | ||
real (kind=kind_io8), dimension(len) | gdata, | ||
integer | len, | ||
integer | iret, | ||
integer | imsk, | ||
integer | jmsk, | ||
real (kind=kind_io8), dimension(imsk,jmsk) | slmskh, | ||
logical | gaus, | ||
real (kind=kind_io8) | blno, | ||
real (kind=kind_io8) | blto, | ||
real (kind=kind_io8), dimension(len) | outlat, | ||
real (kind=kind_io8), dimension(len) | outlon, | ||
integer | me | ||
) |
Read in grib1 analysis data for the requested date and horizontally interpolate to the model grid.
If data not found for the requested date, a backwards search is performed.
[in] | lugb | Fortran unit number of the grib1 analysis data file. |
[in] | fngrib | The name of the grib1 analysis data file. |
[in] | kpds5 | The grib1 parameter number of the requested data. |
[in] | slmask | Model land-sea-ice mask. |
[in] | iy | Year. |
[in] | im | Month. |
[in] | id | Day. |
[in] | ih | Hour. |
[in] | fh | Forecast hour. |
[out] | gdata | The analysis data interpolated to the model grid. |
[in] | len | The number of model points to process. |
[in] | iret | Return code. '0' is success. |
[in] | imsk | 'i' dimension of the high-res mask used for input analysis data without a bitmap. |
[in] | jmsk | 'j' dimension of the high-res mask used for input analysis data without a bitmap. |
[in] | slmskh | The high-resolution mask used for input analysis data witouth a bitmap. |
[in] | gaus | When true, the high-res mask is on a gaussian grid. |
[in] | blno | Corner point longitude of the high-res mask. |
[in] | blto | Corner point latitude of the high-res mask. |
[in] | outlat | Model latitudes. |
[in] | outlon | Model longitudes. |
[in] | me | MPI task number |
Definition at line 9406 of file sfcsub.F.
References getarea(), la2ga(), setrmsk(), and subst().
Referenced by analy().
subroutine fixrdc | ( | integer | lugb, |
character*500 | fngrib, | ||
integer | kpds5, | ||
integer, intent(in) | kpds7, | ||
integer | mon, | ||
real (kind=kind_io8), dimension(len) | slmask, | ||
real (kind=kind_io8), dimension(len) | gdata, | ||
integer | len, | ||
integer | iret, | ||
integer | imsk, | ||
integer | jmsk, | ||
real (kind=kind_io8), dimension(imsk,jmsk) | slmskh, | ||
logical | gaus, | ||
real (kind=kind_io8) | blno, | ||
real (kind=kind_io8) | blto, | ||
real (kind=kind_io8), dimension(len) | outlat, | ||
real (kind=kind_io8), dimension(len) | outlon, | ||
integer | me | ||
) |
Read in grib1 climatology data for a specific month and and horizontally interpolate to the model grid.
[in] | lugb | Fortran unit number of the grib1 climatology data file. |
[in] | fngrib | The name of the grib1 climatology data file. |
[in] | kpds5 | The grib1 parameter number of the requested data. |
[in] | kpds7 | The grib1 level indicator of the requested data. |
[in] | mon | The requested month. |
[in] | slmask | Model land-sea-ice mask. |
[out] | gdata | The climatology data interpolated to the model grid. |
[in] | len | The number of model points to process. |
[in] | iret | Return code. '0' is success. |
[in] | imsk | 'i' dimension of the high-res mask used for input data without a bitmap. |
[in] | jmsk | 'j' dimension of the high-res mask used for input data without a bitmap. |
[in] | slmskh | The high-resolution mask used for input data witouth a bitmap. |
[in] | gaus | When true, the high-res mask is on a gaussian grid. |
[in] | blno | Corner point longitude of the high-res mask. |
[in] | blto | Corner point latitude of the high-res mask. |
[in] | outlat | Model latitudes. |
[in] | outlon | Model longitudes. |
[in] | me | MPI task number |
Definition at line 9205 of file sfcsub.F.
References getarea(), la2ga(), setrmsk(), and subst().
Referenced by clima(), and sfccycle().
subroutine fixrdc_tile | ( | character(len=*), intent(in) | filename_raw, |
character(len=*), intent(in) | tile_num_ch, | ||
integer, dimension(npts), intent(in) | i_index, | ||
integer, dimension(npts), intent(in) | j_index, | ||
integer, intent(in) | kpds, | ||
real(kind_io8), dimension(npts), intent(out) | var, | ||
integer, intent(in) | mon, | ||
integer, intent(in) | npts, | ||
integer, intent(in) | me | ||
) |
Reads in climatological data on the model grid tile for a given month.
Unlike fixrdc, this routine assumes the input data is already interpolated to the model tile. The climatological data must be NetCDF.
[in] | filename_raw | The name of the climatological file. |
[in] | tile_num_ch | The tile number to be processed. |
[in] | i_index | The 'i' indices of the model grid to process. |
[in] | j_index | The 'j' indices of the model grid to process. |
[in] | kpds | Parameter code for the data. |
[out] | var | The climatological data on the model grid. |
[in] | mon | Which month of data to read. |
[in] | npts | Number of model points to process. |
[in] | me | MPI task number. |
Definition at line 9024 of file sfcsub.F.
References count(), and netcdf_err().
Referenced by clima().
subroutine fixrdg | ( | integer | lugb, |
integer | idim, | ||
integer | jdim, | ||
character*(*) | fngrib, | ||
integer | kpds5, | ||
real (kind=kind_io8), dimension(idim*jdim) | gdata, | ||
logical | gaus, | ||
real (kind=kind_io8) | blno, | ||
real (kind=kind_io8) | blto, | ||
integer | me | ||
) |
Read a GRIB1 file.
Return the requested data record and some grid specifications.
[in] | lugb | Fortran unit number of the grib file. |
[in,out] | idim | "i" dimension of the data. |
[in,out] | jdim | "j" dimension of the data. |
[in] | fngrib | Name of the grib file. |
[in] | kpds5 | The grib1 parameter number for the requested field. |
[out] | gdata | The requested data. |
[out] | gaus | When true, grid is gaussian. |
[out] | blno | Corner point longitude of grid. |
[out] | blto | Corner point latitude of grid. |
[in] | me | MPI task number. |
Definition at line 2840 of file sfcsub.F.
Referenced by hmskrd().
subroutine ga2la | ( | real (kind=kind_io8), dimension (imxin,jmxin) | gauin, |
integer | imxin, | ||
integer | jmxin, | ||
real (kind=kind_io8), dimension(imxout,jmxout) | regout, | ||
integer | imxout, | ||
integer | jmxout, | ||
real (kind=kind_io8) | wlon, | ||
real (kind=kind_io8) | rnlat, | ||
real (kind=kind_io8), dimension(imxout) | rlnout, | ||
real (kind=kind_io8), dimension(jmxout) | rltout, | ||
logical | gaus, | ||
real (kind=kind_io8) | blno, | ||
real (kind=kind_io8) | blto | ||
) |
Interpolation from lat/lon or gaussian grid to a lat/lon grid.
[in] | gauin | Input data. |
[in] | imxin | 'i' dimension of input data. |
[in] | jmxin | 'j' dimension of input data. |
[out] | regout | Output data. |
[in] | imxout | 'i' dimension of output data. |
[in] | jmxout | 'j' dimension of output data. |
[in] | wlon | Longitude of west boundary of input data. |
[in] | rnlat | Latitude of north row of input data. |
[in] | rlnout | Longitudes on output data grid. |
[in] | rltout | Latitudes on output data grid. |
[in] | gaus | Is input data on gaussian grid? |
[in] | blno | Corner point longitude of input data. |
[in] | blto | Corner point latitude of input data. |
Definition at line 7388 of file sfcsub.F.
References num_parthds().
Referenced by setrmsk().
subroutine gaulat | ( | real (kind=kind_io8), dimension(k) | gaul, |
integer | k | ||
) |
subroutine getarea | ( | integer, dimension(22) | kgds, |
real (kind=kind_io8) | dlat, | ||
real (kind=kind_io8) | dlon, | ||
real (kind=kind_io8) | rslat, | ||
real (kind=kind_io8) | rnlat, | ||
real (kind=kind_io8) | wlon, | ||
real (kind=kind_io8) | elon, | ||
logical | ijordr, | ||
integer | me | ||
) |
For a given GRIB1 grid description section array, determine some grid specifications.
[in] | kgds | GRIB1 grid description section array. |
[out] | dlat | Grid resolution in the 'j' direction. |
[out] | dlon | Grid resolution in the 'i' direction. |
[out] | rslat | Latitude of the southern row of data. |
[out] | rnlat | Latitude of the northern row of data. |
[out] | wlon | Longitude of the western column of data. |
[out] | elon | Longitude of the eastern column of data. |
[out] | ijordr | When false, adjacent points in the 'i' direction are consecutive. Otherwise, 'j' points are consecutive. |
[in] | me | MPI rank number. |
subroutine getscv | ( | real (kind=kind_io8), dimension(len) | snofld, |
real (kind=kind_io8), dimension(len) | scvfld, | ||
integer | len | ||
) |
Set snow cover flag based on snow depth.
[in] | snofld | Snow depth. |
[out] | scvfld | Snow cover. |
[in] | len | Number of model points to process. |
Definition at line 6463 of file sfcsub.F.
Referenced by sfccycle().
subroutine getsmc | ( | real (kind=kind_io8), dimension(len) | wetfld, |
integer | len, | ||
integer | lsoil, | ||
real (kind=kind_io8), dimension(len,lsoil) | smcfld, | ||
integer | me | ||
) |
Set soil moisture from soil wetness.
[in] | wetfld | Soil wetness |
[in] | len | Number of model points to process. |
[in] | lsoil | Number of soil layers. |
[out] | smcfld | Soil moisture. |
[in] | me | MPI task number. |
Definition at line 6529 of file sfcsub.F.
Referenced by sfccycle().
subroutine getstc | ( | real (kind=kind_io8), dimension(len) | tsffld, |
real (kind=kind_io8), dimension(len) | tg3fld, | ||
real (kind=kind_io8), dimension(len) | slifld, | ||
integer | len, | ||
integer | lsoil, | ||
real (kind=kind_io8), dimension(len,lsoil) | stcfld, | ||
real (kind=kind_io8) | tsfimx | ||
) |
Set soil temperature and sea ice column temperature.
[in] | tsffld | Skin temperature/SST. |
[in] | tg3fld | Soil substrate temperature. |
[in] | slifld | Land-sea-ice mask. |
[in] | len | Number of model points to process. |
[in] | lsoil | Number of soil layers. |
[out] | stcfld | Soil/sea ice column temperature. |
[in] | tsfimx | Freezing point of sea water. |
Definition at line 6487 of file sfcsub.F.
Referenced by sfccycle(), and usesgt().
subroutine hmskrd | ( | integer | lugb, |
integer | imsk, | ||
integer | jmsk, | ||
character*500 | fnmskh, | ||
integer | kpds5, | ||
real (kind=kind_io8), dimension(mdata) | slmskh, | ||
logical | gausm, | ||
real (kind=kind_io8) | blnmsk, | ||
real (kind=kind_io8) | bltmsk, | ||
integer | me | ||
) |
Read a high-resolution land mask.
It will be used as a surrogate for GRIB input data without a bitmap. This is NOT the model land mask. This mask file is GRIB1.
[in] | lugb | Fortran unit number of the mask file. |
[out] | imsk | 'i' dimension of the mask. |
[out] | jmsk | 'j' dimension of the mask. |
[in] | fnmskh | Name of the mask file. |
[in] | kpds5 | GRIB1 parameter number for mask. |
[out] | slmskh | The high-resolution mask. |
[out] | gausm | When true, mask is on a gaussian grid. |
[out] | blnmsk | Corner point longitude of the mask grid. |
[out] | bltmsk | Corner point latitude of the mask grid. |
[in] | me | MPI task number. |
Definition at line 2792 of file sfcsub.F.
References fixrdg().
Referenced by sfccycle().
subroutine la2ga | ( | real (kind=kind_io8), dimension (imxin,jmxin) | regin, |
integer | imxin, | ||
integer | jmxin, | ||
real (kind=kind_io8), dimension(imxin) | rinlon, | ||
real (kind=kind_io8), dimension(jmxin) | rinlat, | ||
real (kind=kind_io8) | rlon, | ||
real (kind=kind_io8) | rlat, | ||
integer | inttyp, | ||
real (kind=kind_io8), dimension(len) | gauout, | ||
integer | len, | ||
logical | lmask, | ||
real (kind=kind_io8), dimension(imxin,jmxin) | rslmsk, | ||
real (kind=kind_io8), dimension(len) | slmask, | ||
real (kind=kind_io8), dimension(len) | outlat, | ||
real (kind=kind_io8), dimension(len) | outlon, | ||
integer | me | ||
) |
Interpolate data from a lat/lon grid to the model grid.
[in] | regin | Input data. |
[in] | imxin | 'i' dimension of input data. |
[in] | jmxin | 'j' dimension of input data. |
[in] | rinlon | |
[in] | rinlat | |
[in] | rlon | |
[in] | rlat | |
[in] | inttyp | |
[out] | gauout | The interpolated data on the model grid. |
[in] | len | Number of model points to process. |
[in] | lmask | |
[in] | rslmsk | |
[in] | slmask | Model land-sea-ice mask. |
[in] | outlat | Latitudes on the model grid. |
[in] | outlon | Longitudes on the model grid. |
[in] | me | MPI task number. |
Definition at line 3210 of file sfcsub.F.
References num_parthds().
subroutine landtyp | ( | real (kind=kind_io8), dimension(len) | vegtype, |
real (kind=kind_io8), dimension(len) | soiltype, | ||
real (kind=kind_io8), dimension(len) | slptype, | ||
real (kind=kind_io8), dimension(len) | slmask, | ||
integer | len | ||
) |
Set vegetation, soil and slope type at undefined model points.
[in,out] | vegtype | Vegetation type |
[in,out] | soiltype | Soil type |
[in,out] | slptype | Slope type |
[in] | slmask | Land-sea-ice mask. |
[in] | len | Number of model points to process. |
Definition at line 7643 of file sfcsub.F.
Referenced by sfccycle().
subroutine maxmin | ( | real (kind=kind_io8), dimension(imax,kmax) | f, |
integer | imax, | ||
integer | kmax | ||
) |
subroutine merge | ( | integer | len, |
integer | lsoil, | ||
integer | iy, | ||
integer | im, | ||
integer | id, | ||
integer | ih, | ||
real (kind=kind_io8) | fh, | ||
real (kind=kind_io8) | deltsfc, | ||
real (kind=kind_io8), dimension(len) | sihfcs, | ||
real (kind=kind_io8), dimension(len) | sicfcs, | ||
real (kind=kind_io8), dimension(len) | vmnfcs, | ||
real (kind=kind_io8), dimension(len) | vmxfcs, | ||
real (kind=kind_io8), dimension(len) | slpfcs, | ||
real (kind=kind_io8), dimension(len) | absfcs, | ||
real (kind=kind_io8), dimension(len) | tsffcs, | ||
real (kind=kind_io8), dimension(len) | wetfcs, | ||
real (kind=kind_io8), dimension(len) | snofcs, | ||
real (kind=kind_io8), dimension(len) | zorfcs, | ||
real (kind=kind_io8), dimension(len,4) | albfcs, | ||
real (kind=kind_io8), dimension(len) | aisfcs, | ||
real (kind=kind_io8), dimension (len) | cvfcs, | ||
real (kind=kind_io8), dimension(len) | cvbfcs, | ||
real (kind=kind_io8), dimension(len) | cvtfcs, | ||
real (kind=kind_io8), dimension(len) | cnpfcs, | ||
real (kind=kind_io8), dimension(len,lsoil) | smcfcs, | ||
real (kind=kind_io8), dimension(len,lsoil) | stcfcs, | ||
real (kind=kind_io8), dimension(len) | slifcs, | ||
real (kind=kind_io8), dimension(len) | vegfcs, | ||
real (kind=kind_io8), dimension(len) | vetfcs, | ||
real (kind=kind_io8), dimension(len) | sotfcs, | ||
real (kind=kind_io8), dimension(len,2) | alffcs, | ||
real (kind=kind_io8), dimension(len) | sihanl, | ||
real (kind=kind_io8), dimension(len) | sicanl, | ||
real (kind=kind_io8), dimension(len) | vmnanl, | ||
real (kind=kind_io8), dimension(len) | vmxanl, | ||
real (kind=kind_io8), dimension(len) | slpanl, | ||
real (kind=kind_io8), dimension(len) | absanl, | ||
real (kind=kind_io8), dimension(len) | tsfanl, | ||
real (kind=kind_io8), dimension(len) | tsfan2, | ||
real (kind=kind_io8), dimension(len) | wetanl, | ||
real (kind=kind_io8), dimension(len) | snoanl, | ||
real (kind=kind_io8), dimension(len) | zoranl, | ||
real (kind=kind_io8), dimension(len,4) | albanl, | ||
real (kind=kind_io8), dimension(len) | aisanl, | ||
real (kind=kind_io8), dimension (len) | cvanl, | ||
real (kind=kind_io8), dimension(len) | cvbanl, | ||
real (kind=kind_io8), dimension(len) | cvtanl, | ||
real (kind=kind_io8), dimension(len) | cnpanl, | ||
real (kind=kind_io8), dimension(len,lsoil) | smcanl, | ||
real (kind=kind_io8), dimension(len,lsoil) | stcanl, | ||
real (kind=kind_io8), dimension(len) | slianl, | ||
real (kind=kind_io8), dimension(len) | veganl, | ||
real (kind=kind_io8), dimension(len) | vetanl, | ||
real (kind=kind_io8), dimension(len) | sotanl, | ||
real (kind=kind_io8), dimension(len,2) | alfanl, | ||
real (kind=kind_io8) | ctsfl, | ||
real (kind=kind_io8) | calbl, | ||
real (kind=kind_io8) | caisl, | ||
real (kind=kind_io8) | csnol, | ||
real (kind=kind_io8), dimension(lsoil) | csmcl, | ||
real (kind=kind_io8) | czorl, | ||
real (kind=kind_io8), dimension(lsoil) | cstcl, | ||
real (kind=kind_io8) | cvegl, | ||
real (kind=kind_io8) | ctsfs, | ||
real (kind=kind_io8) | calbs, | ||
real (kind=kind_io8) | caiss, | ||
real (kind=kind_io8) | csnos, | ||
real (kind=kind_io8), dimension(lsoil) | csmcs, | ||
real (kind=kind_io8) | czors, | ||
real (kind=kind_io8), dimension(lsoil) | cstcs, | ||
real (kind=kind_io8) | cvegs, | ||
real (kind=kind_io8) | ccv, | ||
real (kind=kind_io8) | ccvb, | ||
real (kind=kind_io8) | ccvt, | ||
real (kind=kind_io8) | ccnp, | ||
real (kind=kind_io8) | cvetl, | ||
real (kind=kind_io8) | cvets, | ||
real (kind=kind_io8) | csotl, | ||
real (kind=kind_io8) | csots, | ||
real (kind=kind_io8) | calfl, | ||
real (kind=kind_io8) | calfs, | ||
real (kind=kind_io8) | csihl, | ||
real (kind=kind_io8) | csihs, | ||
real (kind=kind_io8) | csicl, | ||
real (kind=kind_io8) | csics, | ||
real (kind=kind_io8) | cvmnl, | ||
real (kind=kind_io8) | cvmns, | ||
real (kind=kind_io8) | cvmxl, | ||
real (kind=kind_io8) | cvmxs, | ||
real (kind=kind_io8) | cslpl, | ||
real (kind=kind_io8) | cslps, | ||
real (kind=kind_io8) | cabsl, | ||
real (kind=kind_io8) | cabss, | ||
integer | irttsf, | ||
integer | irtwet, | ||
integer | irtsno, | ||
integer | irtzor, | ||
integer | irtalb, | ||
integer | irtais, | ||
integer | irttg3, | ||
integer | irtscv, | ||
integer | irtacn, | ||
integer | irtsmc, | ||
integer | irtstc, | ||
integer | irtveg, | ||
integer | irtvmn, | ||
integer | irtvmx, | ||
integer | irtslp, | ||
integer | irtabs, | ||
integer | irtvet, | ||
integer | irtsot, | ||
integer | irtalf, | ||
logical, intent(in) | landice, | ||
integer | me | ||
) |
Blend the model forecast (or first guess) fields with the analysis/climatology.
The forecast/first guess variables are named with "fcs". The analysis/climatology variables are named with "anl". On output, the blended fields are stored in the "anl" variables.
[in] | len | Number of model points to process. |
[in] | lsoil | Number of soil layers. |
[in] | iy | Cycle year. |
[in] | im | Cycle month. |
[in] | id | Cycle day. |
[in] | ih | Cycle hour. |
[in] | fh | Forecast hour. |
[in] | deltsfc | Cycling frequency in hours. |
[in] | sihfcs | First guess sea ice depth. |
[in] | sicfcs | First guess sea ice concentration. |
[in] | vmnfcs | First guess minimum vegetation greenness. |
[in] | vmxfcs | First guess maximum vegetation greenness. |
[in] | slpfcs | First guess slope type. |
[in] | absfcs | First guess maximum snow albedo. |
[in] | tsffcs | First guess skin temperature/SST. |
[in] | wetfcs | First guess soil wetness. |
[in] | snofcs | First guess liquid equivalent snow depth. |
[in] | zorfcs | First guess roughness length. |
[in] | albfcs | First guess snow free albedo. |
[in] | aisfcs | First guess ice. |
[in] | cvfcs | First guess convective cloud cover. |
[in] | cvbfcs | First guess convective cloud bottom. |
[in] | cvtfcs | First guess convective cloud top. |
[in] | cnpfcs | First guess canopy water content. |
[in] | smcfcs | First guess soil moisture. |
[in] | stcfcs | First guess soil/ice temperature. |
[in] | slifcs | First guess land-sea-ice mask. |
[in] | vegfcs | First guess vegetation greenness. |
[in] | vetfcs | First guess vegetation type. |
[in] | sotfcs | First guess soil type. |
[in] | alffcs | First guess strong/weak zenith angle dependent albedo. |
[in,out] | sihanl | Blended sea ice depth. |
[in,out] | sicanl | Blended sea ice concentration. |
[in,out] | vmnanl | Blended minimum vegetation greenness. |
[in,out] | vmxanl | Blended maximum vegetation greenness. |
[in,out] | slpanl | Blended slope type. |
[in,out] | absanl | Blended maximum snow albedo. |
[in,out] | tsfanl | Blended skin temperature/SST. |
[in,out] | tsfan2 | Not used. |
[in,out] | wetanl | Blended soil wetness. |
[in,out] | snoanl | Blended liquid equivalent snow depth. |
[in,out] | zoranl | Blended roughness length. |
[in,out] | albanl | Blended snow-free albedo. |
[in,out] | aisanl | Blended ice. |
[in,out] | cvanl | Blended convective cloud cover. |
[in,out] | cvbanl | Blended convective cloud base. |
[in,out] | cvtanl | Blended convective cloud top. |
[in,out] | cnpanl | Blended canopy moisture content. |
[in,out] | smcanl | Blended soil moisture. |
[in,out] | stcanl | Blended soil/ice temperature. |
[in,out] | slianl | Blended land-sea-ice mask. |
[in,out] | veganl | Blended vegetation greenness. |
[in,out] | vetanl | Blended vegetation type. |
[in,out] | sotanl | Blended soil type. |
[in,out] | alfanl | Blended strong/weak zenith angle dependent albedo. |
[in] | ctsfl | Merging coefficient for skin temperature. |
[in] | calbl | Merging coefficient for snow-free albedo at land points. |
[in] | caisl | Merging coefficient for sea ice at land points. |
[in] | csnol | Merging coefficient for snow at land points. |
[in] | csmcl | Merging coefficient for soil moisture at land points. |
[in] | czorl | Merging coefficient for roughness length at land points. |
[in] | cstcl | Merging coefficient for soil temperature at land points. |
[in] | cvegl | Merging coefficient for vegetation greenness at land points. |
[in] | ctsfs | Merging coefficient for SST. |
[in] | calbs | Merging coefficient for snow-free albedo at water points. |
[in] | caiss | Merging coefficient for sea ice at water points. |
[in] | csnos | Merging coefficient for snow at water points. |
[in] | csmcs | Merging coefficient for soil moisture at water points. |
[in] | czors | Merging coefficient for roughness length at water points. |
[in] | cstcs | Merging coefficient for sea ice temperature at water points. |
[in] | cvegs | Merging coefficient for vegetation greenness at water points. |
[in] | ccv | Merging coefficient for convective cloud cover. |
[in] | ccvb | Merging coefficient for covective cloud bottom. |
[in] | ccvt | Merging coefficient for covective cloud top. |
[in] | ccnp | Merging coefficient for canopy moisture. |
[in] | cvetl | Merging coefficient for vegetation type at land points. |
[in] | cvets | Merging coefficient for vegetation type at water points. |
[in] | csotl | Merging coefficient for soil type at land points. |
[in] | csots | Merging coefficient for soil type at water points. |
[in] | calfl | Merging coefficient for strong/weak zenith angle dependent albedo at land points. |
[in] | calfs | Merging coefficient for strong/weak zenith angle dependent albedo at water points. |
[in] | csihl | Merging coefficient for sea ice depth at land points. |
[in] | csihs | Merging coefficient for sea ice depth at water points. |
[in] | csicl | Merging coefficient for sea ice concentration at land points. |
[in] | csics | Merging coefficient for sea ice concentration at water points. |
[in] | cvmnl | Merging coefficient for minimum vegetation greenness at land points. |
[in] | cvmns | Merging coefficient for minimum vegetation greenness at water points. |
[in] | cvmxl | Merging coefficient for maximum vegetation greenness at land points. |
[in] | cvmxs | Merging coefficient for maximum vegetation greenness at water points. |
[in] | cslpl | Merging coefficient for slope type at land points. |
[in] | cslps | Merging coefficient for slope type at water points. |
[in] | cabsl | Merging coefficient for maximum snow albedo at land points. |
[in] | cabss | Merging coefficient for maximum snow albedo at water points. |
[in] | irttsf | Return code from read of skin temperature/SST analysis file. |
[in] | irtwet | Return code from read of soil wetness analysis file. |
[in] | irtsno | Return code from read of snow analysis file. |
[in] | irtzor | Return code from read of roughness length file. |
[in] | irtalb | Return code from read of snow-free albedo analysis file. |
[in] | irtais | Return code from read of ice mask analysis file. |
[in] | irttg3 | Return code from read of soil substrate temperature analysis file. |
[in] | irtscv | Return code from read of snow cover analysis file. |
[in] | irtacn | Return code from read of sea ice concentration analysis file. |
[in] | irtsmc | Return code from read of soil moisture analysis file. |
[in] | irtstc | Return code from read of soil temperature analysis file. |
[in] | irtveg | Return code from read of vegetation greenness analysis file. |
[in] | irtvmn | Return code from read of minimum vegetation greenness analysis file. |
[in] | irtvmx | Return code from read of maximum vegetation greenness analysis file. |
[in] | irtslp | Return code from read of slope type analysis file. |
[in] | irtabs | Return code from read of maximum snow albedo analysis file. |
[in] | irtvet | Return code from read of vegetation type analysis file. |
[in] | irtsot | Return code from read of soil type analysis file. |
[in] | irtalf | Return code from read of file containing fraction for strongly and weakly zenith angle dependent albedo. |
[in] | landice | Permanent land ice flag. |
[in] | me | MPI task number. |
Definition at line 5192 of file sfcsub.F.
References num_parthds().
Referenced by sfccycle().
subroutine monitr | ( | character(len=*) | lfld, |
real (kind=kind_io8), dimension(ijmax) | fld, | ||
real (kind=kind_io8), dimension(ijmax) | slimsk, | ||
real (kind=kind_io8), dimension(ijmax) | sno, | ||
integer | ijmax | ||
) |
Determine the maximum and minimum values of a surface field at snow-free and snow covered land, open water, and snow-free and snow covered sea ice.
[in] | lfld | The name of the surface field to monitory. |
[in] | fld | The surface field to monitor. |
[in] | slimsk | Land-sea-ice mask. |
[in] | sno | Snow. |
[in] | ijmax | Number of model points to process. |
Definition at line 2689 of file sfcsub.F.
Referenced by sfccycle().
subroutine netcdf_err | ( | integer, intent(in) | error | ) |
Print the error message for a given netCDF return code.
[in] | error |
Definition at line 9166 of file sfcsub.F.
Referenced by fixrdc_tile(), read_write_data::read_data(), read_write_data::read_gsi_data(), read_write_data::read_lat_lon_orog(), and read_write_data::write_data().
subroutine newice | ( | real (kind=kind_io8), dimension(len) | slianl, |
real (kind=kind_io8), dimension(len) | slifcs, | ||
real (kind=kind_io8), dimension(len) | tsfanl, | ||
real (kind=kind_io8), dimension(len) | tsffcs, | ||
integer | len, | ||
integer | lsoil, | ||
sihnew, | |||
sicnew, | |||
real (kind=kind_io8), dimension(len) | sihanl, | ||
real (kind=kind_io8), dimension(len) | sicanl, | ||
real (kind=kind_io8), dimension(len,4) | albanl, | ||
real (kind=kind_io8), dimension(len) | snoanl, | ||
real (kind=kind_io8), dimension(len) | zoranl, | ||
real (kind=kind_io8), dimension(len,lsoil) | smcanl, | ||
real (kind=kind_io8), dimension(len,lsoil) | stcanl, | ||
real (kind=kind_io8) | albsea, | ||
real (kind=kind_io8) | snosea, | ||
real (kind=kind_io8) | zorsea, | ||
real (kind=kind_io8) | smcsea, | ||
real (kind=kind_io8) | smcice, | ||
real (kind=kind_io8) | tsfmin, | ||
real (kind=kind_io8) | tsfice, | ||
real (kind=kind_io8) | albice, | ||
real (kind=kind_io8) | zorice, | ||
real (kind=kind_io8) | tgice, | ||
real (kind=kind_io8), dimension(len) | rla, | ||
real (kind=kind_io8), dimension(len) | rlo, | ||
integer | me | ||
) |
Adjust surface fields when ice melts or forms.
[in] | slianl | Land-sea-ice mask. |
[in] | slifcs | First guess land-sea-ice mask. |
[out] | tsfanl | Skin temperature/SST. |
[in] | tsffcs | First guess skin temperature/SST. |
[in] | len | Number of model points to process. |
[in] | lsoil | Number of soil layers. |
[in] | sihnew | Sea ice depth for new ice. |
[in] | sicnew | Sea ice concentration for new ice. |
[out] | sihanl | Sea ice depth. |
[out] | sicanl | Sea ice concentration. |
[out] | albanl | Snow-free albedo. |
[out] | snoanl | Liquid equivalent snow depth. |
[out] | zoranl | Roughness length. |
[out] | smcanl | Soil moisture |
[out] | stcanl | Soil temperature |
[in] | albsea | Albedo at open water. |
[in] | snosea | Snow at open water. |
[in] | zorsea | Roughness length at open water. |
[in] | smcsea | Soil moisture at open water. |
[in] | smcice | Soil moisture at ice. |
[in] | tsfmin | SST at open water. |
[in] | tsfice | Skin temperature at ice. |
[in] | albice | Ice albedo. |
[in] | zorice | Roughness length of ice. |
[in] | tgice | Freezing point of salt water. |
[in] | rla | Model latitudes. |
[in] | rlo | Model longitudes. |
[in] | me | MPI task number. |
Definition at line 5676 of file sfcsub.F.
Referenced by sfccycle().
subroutine qcbyfc | ( | real (kind=kind_io8), dimension(len) | tsffcs, |
real (kind=kind_io8), dimension(len) | snofcs, | ||
real (kind=kind_io8) | qctsfs, | ||
real (kind=kind_io8) | qcsnos, | ||
real (kind=kind_io8) | qctsfi, | ||
integer | len, | ||
integer | lsoil, | ||
real (kind=kind_io8), dimension(len) | snoanl, | ||
real (kind=kind_io8), dimension(len) | aisanl, | ||
real (kind=kind_io8), dimension(len) | slianl, | ||
real (kind=kind_io8), dimension(len) | tsfanl, | ||
real (kind=kind_io8), dimension(len,4) | albanl, | ||
real (kind=kind_io8), dimension(len) | zoranl, | ||
real (kind=kind_io8), dimension(len,lsoil) | smcanl, | ||
real (kind=kind_io8), dimension(len,lsoil) | smcclm, | ||
real (kind=kind_io8) | tsfsmx, | ||
real (kind=kind_io8) | albomx, | ||
real (kind=kind_io8) | zoromx, | ||
integer | me | ||
) |
Quality control analysis fields using the first guess.
[in] | tsffcs | First guess skin temperature/SST. |
[in] | snofcs | First guess snow. |
[in] | qctsfs | Surface temperature above which no snow allowed. |
[in] | qcsnos | Snow depth above which snow must exits. |
[in] | qctsfi | SST above which sea ice is not allowed. |
[in] | len | Number of model points to process. |
[in] | lsoil | Number of soil layers. |
[out] | snoanl | Snow analysis. |
[in] | aisanl | Ice mask. |
[in] | slianl | Land-sea-ice mask. |
[out] | tsfanl | Skin temperature/SST analysis. |
[in] | albanl | Snow-free albedo analysis. |
[in] | zoranl | Roughness length analysis. |
[in] | smcanl | Soil moisture analysis. |
[in] | smcclm | Soil moisture climatology. |
[in] | tsfsmx | Not used. |
[in] | albomx | Snow-free albedo at open water. |
[in] | zoromx | Roughness length at open water. |
[in] | me | MPI task number. |
Definition at line 6766 of file sfcsub.F.
Referenced by sfccycle().
subroutine qcmxice | ( | real (kind=kind_io8), dimension(len) | glacir, |
real (kind=kind_io8), dimension(len) | amxice, | ||
integer | len, | ||
integer | me | ||
) |
Quality control maximum ice extent.
[in] | glacir | Glacial flag |
[out] | amxice | Maximum ice extent. |
[in] | len | Number of model points to process. |
[in] | me | MPI task number. |
Definition at line 6649 of file sfcsub.F.
Referenced by sfccycle().
subroutine qcmxmn | ( | character*8 | ttl, |
real (kind=kind_io8), dimension(len) | fld, | ||
real (kind=kind_io8), dimension(len) | slimsk, | ||
real (kind=kind_io8), dimension(len) | sno, | ||
logical, dimension(len) | iceflg, | ||
real (kind=kind_io8) | fldlmx, | ||
real (kind=kind_io8) | fldlmn, | ||
real (kind=kind_io8) | fldomx, | ||
real (kind=kind_io8) | fldomn, | ||
real (kind=kind_io8) | fldimx, | ||
real (kind=kind_io8) | fldimn, | ||
real (kind=kind_io8) | fldjmx, | ||
real (kind=kind_io8) | fldjmn, | ||
real (kind=kind_io8) | fldsmx, | ||
real (kind=kind_io8) | fldsmn, | ||
real (kind=kind_io8) | epsfld, | ||
real (kind=kind_io8), dimension(len) | rla, | ||
real (kind=kind_io8), dimension(len) | rlo, | ||
integer | len, | ||
integer | mode, | ||
real (kind=kind_io8) | percrit, | ||
logical | lgchek, | ||
integer | me | ||
) |
Range check a field.
The check is a function of whether the point is land, water or ice.
[in] | ttl | Name of field. |
[in,out] | fld | The field array to be ranged checked. |
[in] | slimsk | The land-sea-ice mask. |
[in] | sno | Snow. |
[in] | iceflg | When true, check ice points. |
[in] | fldlmx | Maximum allowable value at snow-free land. |
[in] | fldlmn | Minimum allowable value at snow-free land. |
[in] | fldomx | Maximum allowable value at open water. |
[in] | fldomn | Minimum allowable value at open water. |
[in] | fldimx | Maximum allowable value at snow-free ice. |
[in] | fldimn | Minimum allowable value at snow-free ice. |
[in] | fldjmx | Maximum allowable value at snow covered ice. |
[in] | fldjmn | Minimum allowable value at snow covered ice. |
[in] | fldsmx | Maximum allowable value at snow covered land. |
[in] | fldsmn | Minimum allowable value at snow covered land. |
[in] | epsfld | Difference from the max/min allowable value at which the field is updated. |
[in] | rla | Latitude of the points to process. |
[in] | rlo | Longitude of the points to process. |
[in] | len | Number of points to process. |
[in] | mode | When '1', update the field. When not '1', run routine in diagnostic mode. |
[in] | percrit | Critical percentage of 'bad' points required for abort. |
[in] | lgchek | When true, abort when a critical percentage of points are outside acceptable range. |
[in] | me | MPI task number. |
Definition at line 6021 of file sfcsub.F.
References num_parthds().
Referenced by sfccycle().
subroutine qcsice | ( | real (kind=kind_io8), dimension(len) | ais, |
real (kind=kind_io8), dimension(len) | glacir, | ||
real (kind=kind_io8), dimension(len) | amxice, | ||
real (kind=kind_io8) | aicice, | ||
real (kind=kind_io8) | aicsea, | ||
real (kind=kind_io8) | sllnd, | ||
real (kind=kind_io8), dimension(len) | slmask, | ||
real (kind=kind_io8), dimension(len) | rla, | ||
real (kind=kind_io8), dimension(len) | rlo, | ||
integer | len, | ||
integer | me | ||
) |
Check the sea ice cover mask against the land-sea mask.
[out] | ais | Sea ice cover mask. |
[in] | glacir | Glacial flag |
[in] | amxice | Maximum ice extent. |
[in] | aicice | Ice indicator. |
[in] | aicsea | Open water indicator. |
[in] | sllnd | Land indicator. |
[in] | slmask | Land-sea-ice mask. |
[in] | rla | Model latitudes |
[in] | rlo | Model longitudes |
[in] | len | Number of model points to process. |
[in] | me | MPI task number. |
Definition at line 5844 of file sfcsub.F.
Referenced by sfccycle().
subroutine qcsli | ( | real (kind=kind_io8), dimension(len) | slianl, |
real (kind=kind_io8), dimension(len) | slifcs, | ||
integer | len, | ||
integer | me | ||
) |
Check consistency between the forecast and analysis land-sea-ice mask.
[in] | slianl | Analysis mask. |
[in,out] | slifcs | Forecast mask. |
[in] | len | Number of model points to process. |
[in] | me | MPI task number. |
Definition at line 6678 of file sfcsub.F.
Referenced by sfccycle().
subroutine qcsnow | ( | real (kind=kind_io8), dimension(len) | snoanl, |
real (kind=kind_io8), dimension(len) | slmask, | ||
real (kind=kind_io8), dimension(len) | aisanl, | ||
real (kind=kind_io8), dimension(len) | glacir, | ||
integer | len, | ||
real (kind=kind_io8) | snoval, | ||
logical, intent(in) | landice, | ||
integer | me | ||
) |
Quality control snow at the model points.
[out] | snoanl | Model snow to be qc'd. |
[in] | slmask | Land-sea-ice mask. |
[in] | aisanl | Ice mask. |
[in] | glacir | Permanent glacial point when '1'. |
[in] | len | Number of model points to process. |
[in] | snoval | Minimum snow depth at glacial points. |
[in] | landice | Is the point permanent land ice? |
[in] | me | MPI task number. |
Definition at line 5783 of file sfcsub.F.
Referenced by sfccycle().
subroutine rof01 | ( | real (kind=kind_io8), dimension(len) | aisfld, |
integer | len, | ||
character*2 | op, | ||
real (kind=kind_io8) | crit | ||
) |
Round a field up to one or down to zero.
[in,out] | aisfld | The field to adjust. |
[in] | len | Number of model points to process. |
[in] | op | Operation on field. |
[in] | crit | Critical value above/below the field is rounded. |
Definition at line 4890 of file sfcsub.F.
Referenced by setrmsk(), and sfccycle().
subroutine scale | ( | real (kind=kind_io8), dimension(len) | fld, |
integer | len, | ||
real (kind=kind_io8) | scl | ||
) |
Multiply a field by a scaling factor.
[in,out] | fld | The field to scale. |
[in] | len | Number of model points to process. |
[in] | scl | Scaling factor. |
Definition at line 5978 of file sfcsub.F.
Referenced by sfccycle().
subroutine setlsi | ( | real (kind=kind_io8), dimension(len) | slmask, |
real (kind=kind_io8), dimension(len) | aisfld, | ||
integer | len, | ||
real (kind=kind_io8) | aicice, | ||
real (kind=kind_io8), dimension(len) | slifld | ||
) |
Set land-sea-ice mask at sea ice.
[in] | slmask | The model mask on input. |
[in] | aisfld | The ice mask on the model grid. |
[in] | len | Number of model points to process. |
[in] | aicice | Ice concentration theshold for an ice point. |
[out] | slifld | The model mask with sea ice added. |
Definition at line 5953 of file sfcsub.F.
Referenced by sfccycle().
subroutine setrmsk | ( | integer | kpds5, |
real (kind=kind_io8), dimension(igaul,jgaul) | slmask, | ||
integer | igaul, | ||
integer | jgaul, | ||
real (kind=kind_io8) | wlon, | ||
real (kind=kind_io8) | rnlat, | ||
real (kind=kind_io8), dimension(imax,jmax) | data, | ||
integer | imax, | ||
integer | jmax, | ||
real (kind=kind_io8), dimension(imax) | rlnout, | ||
real (kind=kind_io8), dimension(jmax) | rltout, | ||
logical | lmask, | ||
real (kind=kind_io8), dimension(imax,jmax) | rslmsk, | ||
logical | gaus, | ||
real (kind=kind_io8) | blno, | ||
real (kind=kind_io8) | blto, | ||
integer | kgds1, | ||
integer, intent(in) | kpds4, | ||
logical*1, dimension(imax,jmax), intent(in) | lbms | ||
) |
Set the mask for the input data.
(Not the model mask). The mask is determined from the input data bitmap or by interpolating a high-resolution mask to the input data grid. Note: not all data has a mask.
[in] | kpds5 | Grib parameter number. |
[in] | slmask | High-resolution mask. |
[in] | igaul | "i" dimension of the high-res mask. |
[in] | jgaul | "j" dimension of the high-res mask. |
[in] | wlon | Longitude of 'west' boundary of input data. |
[in] | rnlat | Latitude of north row of input data. |
[in] | data | The input data. |
[in] | imax | "i" dimension of input grid. |
[in] | jmax | "j" dimension of input grid. |
[out] | rlnout | Longitudes on input data grid. |
[out] | rltout | Latitudes on input data grid. |
[out] | lmask | True, when input data has a mask. |
[out] | rslmsk | The mask of the input data grid. |
[in] | gaus | Is high-resolution mask on a gaussian grid. |
[in] | blno | Corner point longitude of the high-res mask. |
[in] | blto | Corner point longitude of the high-res mask. |
[in] | kgds1 | Grib indicator for grid type. |
[in] | kpds4 | Grib indicator for bitmap. |
[in] | lbms | Bitmap of input data. |
Definition at line 6888 of file sfcsub.F.
subroutine setzro | ( | real (kind=kind_io8), dimension(len) | fld, |
real (kind=kind_io8) | eps, | ||
integer | len | ||
) |
Set a field to zero if it is less than a threshold.
[in,out] | fld | Field to set. |
[in] | eps | Threshold. |
[in] | len | Number of model points to process. |
Definition at line 6445 of file sfcsub.F.
Referenced by sfccycle().
subroutine sfccycle | ( | integer | lugb, |
integer | len, | ||
integer | lsoil, | ||
real (kind=kind_io8), dimension(len) | sig1t, | ||
real (kind=kind_io8) | deltsfc, | ||
integer | iy, | ||
integer | im, | ||
integer | id, | ||
integer | ih, | ||
real (kind=kind_io8) | fh, | ||
real (kind=kind_io8), dimension(len) | rla, | ||
real (kind=kind_io8), dimension(len) | rlo, | ||
real (kind=kind_io8), dimension(len) | slmask, | ||
real (kind=kind_io8), dimension(len) | orog, | ||
real (kind=kind_io8), dimension(len) | orog_uf, | ||
logical | use_ufo, | ||
logical | nst_anl, | ||
real (kind=kind_io8), dimension(len) | sihfcs, | ||
real (kind=kind_io8), dimension(len) | sicfcs, | ||
real (kind=kind_io8), dimension(len) | sitfcs, | ||
real (kind=kind_io8), dimension(len) | swdfcs, | ||
real (kind=kind_io8), dimension(len,lsoil) | slcfcs, | ||
real (kind=kind_io8), dimension(len) | vmnfcs, | ||
real (kind=kind_io8), dimension(len) | vmxfcs, | ||
real (kind=kind_io8), dimension(len) | slpfcs, | ||
real (kind=kind_io8), dimension(len) | absfcs, | ||
real (kind=kind_io8), dimension(len) | tsffcs, | ||
real (kind=kind_io8), dimension(len) | snofcs, | ||
real (kind=kind_io8), dimension(len) | zorfcs, | ||
real (kind=kind_io8), dimension(len,4) | albfcs, | ||
real (kind=kind_io8), dimension(len) | tg3fcs, | ||
real (kind=kind_io8), dimension(len) | cnpfcs, | ||
real (kind=kind_io8), dimension(len,lsoil) | smcfcs, | ||
real (kind=kind_io8), dimension(len,lsoil) | stcfcs, | ||
real (kind=kind_io8), dimension(len) | slifcs, | ||
real (kind=kind_io8), dimension(len) | aisfcs, | ||
real (kind=kind_io8), dimension(len) | vegfcs, | ||
real (kind=kind_io8), dimension(len) | vetfcs, | ||
real (kind=kind_io8), dimension(len) | sotfcs, | ||
real (kind=kind_io8), dimension(len,2) | alffcs, | ||
real (kind=kind_io8), dimension (len) | cvfcs, | ||
real (kind=kind_io8), dimension(len) | cvbfcs, | ||
real (kind=kind_io8), dimension(len) | cvtfcs, | ||
integer | me, | ||
integer | nlunit, | ||
integer | sz_nml, | ||
character(len=*), dimension(sz_nml), intent(in) | input_nml_file, | ||
integer | ialb, | ||
integer | isot, | ||
integer | ivegsrc, | ||
character(len=*), intent(in) | tile_num_ch, | ||
integer, dimension(len), intent(in) | i_index, | ||
integer, dimension(len), intent(in) | j_index | ||
) |
Surface cycling driver routine.
Update 'first guess' surface fields with either analysis or climatological data.
[in] | lugb | Fortran unit number used to read data files. |
[in] | len | Number of model points to process. |
[in] | lsoil | Number of soil layers. |
[in] | sig1t | Sigma level 1 temperature for dead start. |
[in] | deltsfc | Cycling frequency in hours. |
[in] | iy | Year of initial state. |
[in] | im | Month of initial state. |
[in] | id | Day of initial state. |
[in] | ih | Hour of initial state. |
[in] | fh | Forecast hour. |
[in] | rla | Latitude of model points. |
[in] | rlo | Longitude of model points. |
[in] | slmask | Model land-sea mask without ice flag. |
[in] | orog | Filtered orography on model grid. |
[in] | orog_uf | Unfiltered orography on model grid. |
[in] | use_ufo | When true, adjust sst and soil substrate temperature for differences between filtered and unfiltered terrain. |
[in] | nst_anl | Set to true when using nst model. |
[in,out] | sihfcs | Sea ice thickness on model grid. |
[in,out] | sicfcs | Sea ice concentration on model grid. |
[in,out] | sitfcs | Sea ice skin temperature on model grid. |
[in,out] | swdfcs | Physical snow depth on model grid. |
[in,out] | slcfcs | Liquid portion of volumentric soil moisture on model grid. |
[in,out] | vmnfcs | Minimum greenness fraction of model grid. |
[in,out] | vmxfcs | Maximum greenness fraction of model grid. |
[in,out] | slpfcs | Slope type on model grid. |
[in,out] | absfcs | Maximum snow albedo on model grid. |
[in,out] | tsffcs | Skin temperature/SST on model grid. |
[in,out] | snofcs | Liquid equivalent snow depth on model grid. |
[in,out] | zorfcs | Roughness length on model grid. |
[in,out] | albfcs | Snow-free albedo on model grid. |
[in,out] | tg3fcs | Soil substrate temperature on model grid. |
[in,out] | cnpfcs | Canopy moisture content on model grid. |
[in,out] | smcfcs | Total volumetric soil moisture on model grid. |
[in,out] | stcfcs | Soil/ice depth temperature on model grid. |
[in,out] | slifcs | Model land-sea mask including ice flag. |
[in,out] | aisfcs | Model ice mask. |
[in,out] | vegfcs | Vegetation greenness on model grid. |
[in,out] | vetfcs | Vegetation type on model grid. |
[in,out] | sotfcs | Soil type on model grid. |
[in,out] | alffcs | Fraction for strongly and weakly zenith angle dependent albedo on model grid. |
[in,out] | cvfcs | Convective cloud cover on model grid. |
[in,out] | cvbfcs | Convective cloud base on model grid. |
[in,out] | cvtfcs | Convective cloud top on model grid. |
[in] | me | MPI task number. |
[in] | nlunit | Program namelist unit number. |
[in] | sz_nml | Dimension of input_nml_file. |
[in] | input_nml_file | Name of program namelist file. |
[in] | ialb | Use modis albedo when '1'. Use brigleb when '0'. |
[in] | isot | When '1', use statsgo soil type. When '0' use zobler soil type. |
[in] | ivegsrc | When '1', use igbp vegetation type. When '2' use sib vegetation type. |
[in] | tile_num_ch | Model tile number to process. |
[in] | i_index | The 'i' indices of the model grid to process. |
[in] | j_index | The 'j' indices of the model grid to process. |
Definition at line 175 of file sfcsub.F.
References albocn(), analy(), anomint(), clima(), filanl(), filfcs(), fixrdc(), getscv(), getsmc(), getstc(), hmskrd(), landtyp(), merge(), monitr(), newice(), qcbyfc(), qcmxice(), qcmxmn(), qcsice(), qcsli(), qcsnow(), rof01(), scale(), setlsi(), setzro(), snodpth(), snodpth2(), snosfc(), tsfcor(), and usesgt().
Referenced by sfcdrv().
subroutine snodpth | ( | real (kind=kind_io8), dimension(len) | scvanl, |
real (kind=kind_io8), dimension(len) | slianl, | ||
real (kind=kind_io8), dimension(len) | tsfanl, | ||
real (kind=kind_io8), dimension(len) | snoclm, | ||
real (kind=kind_io8), dimension(len) | glacir, | ||
real (kind=kind_io8) | snwmax, | ||
real (kind=kind_io8) | snwmin, | ||
logical, intent(in) | landice, | ||
integer | len, | ||
real (kind=kind_io8), dimension(len) | snoanl, | ||
integer | me | ||
) |
Estimate snow depth at glacial, land and sea ice points.
[in] | scvanl | Snow cover. |
[in] | slianl | Land-sea-ice mask. |
[in] | tsfanl | Skin temperature. |
[in] | snoclm | Climatological snow depth. |
[in] | glacir | Permanent glacial point when '1'. |
[in] | snwmax | Maximum snow depth. |
[in] | snwmin | Minimum snow depth. |
[in] | landice | When true, point is permanent land ice. |
[in] | len | Number of model points to process. |
[out] | snoanl | Snow depth. |
[in] | me | MPI task number. |
Definition at line 4978 of file sfcsub.F.
Referenced by sfccycle().
subroutine snodpth2 | ( | real (kind=kind_io8), dimension(len) | glacir, |
real (kind=kind_io8) | snwmax, | ||
real (kind=kind_io8), dimension(len) | snoanl, | ||
integer | len, | ||
integer | me | ||
) |
Ensure deep snow pack at permanent glacial points.
[in] | glacir | Glacial flag |
[in] | snwmax | Deep snow pack depth |
[in,out] | snoanl | Model snow |
[in] | len | Number of model points to process. |
[in] | me | MPI task number. |
Definition at line 9730 of file sfcsub.F.
Referenced by sfccycle().
subroutine snosfc | ( | real (kind=kind_io8), dimension(len) | snoanl, |
real (kind=kind_io8), dimension(len) | tsfanl, | ||
real (kind=kind_io8) | tsfsmx, | ||
integer | len, | ||
integer | me | ||
) |
Check skin temperature at points with snow.
If it is above a threshold, reset it to that value.
[in] | snoanl | Snow. |
[in,out] | tsfanl | Skin temperature at the model points. |
[in] | tsfsmx | Maximum allowable skin temperature at snow points. |
[in] | len | Number of model points to process. |
[in] | me | MPI rank |
Definition at line 6593 of file sfcsub.F.
Referenced by sfccycle().
subroutine subst | ( | real (kind=kind_io8), dimension(imax,jmax) | data, |
integer | imax, | ||
integer | jmax, | ||
real (kind=kind_io8) | dlon, | ||
real (kind=kind_io8) | dlat, | ||
logical | ijordr | ||
) |
Take an array of data on a lat/lon based grid and rearrange it so the corner point is in the 'lower left' and adjacent points in the 'i' direction are consecutive.
[in,out] | data | The data to be adjusted. |
[in] | imax | 'i' dimension of data. |
[in] | jmax | 'j' dimension of data. |
[in] | dlon | Delta longitude of the data. |
[in] | dlat | Delta latitude of the data. |
[in] | ijordr | When false, adjacent points in the 'i' direction are consecutive. Otherwise, 'j' points are consecutive. |
subroutine tsfcor | ( | real (kind=kind_io8), dimension(len) | tsfc, |
real (kind=kind_io8), dimension(len) | orog, | ||
real (kind=kind_io8), dimension(len) | slmask, | ||
real (kind=kind_io8) | umask, | ||
integer | len, | ||
real (kind=kind_io8) | rlapse | ||
) |
Adjust skin temperature or SST for terrain.
[in,out] | tsfc | Skin temperature/SST |
[in] | orog | Orography height. |
[in] | slmask | Model land-sea mask. |
[in] | umask | When '0' adjust SST, when '1' adjust skin temperature. |
[in] | len | Number of model points to process. |
[in] | rlapse | Standard atmospheric lapse rate. |
Definition at line 4948 of file sfcsub.F.
Referenced by sfccycle().
subroutine usesgt | ( | real (kind=kind_io8), dimension(len) | sig1t, |
real (kind=kind_io8), dimension(len) | slianl, | ||
real (kind=kind_io8), dimension(len) | tg3anl, | ||
integer | len, | ||
integer | lsoil, | ||
real (kind=kind_io8), dimension(len) | tsfanl, | ||
real (kind=kind_io8), dimension(len,lsoil) | stcanl, | ||
real (kind=kind_io8) | tsfimx | ||
) |
Set soil temperature and sea ice column temperature for a dead start.
[in] | sig1t | Sigma level 1 temperature for dead start. |
[in] | slianl | Land-sea-ice mask. |
[in] | tg3anl | Soil substrate temperature. |
[in] | len | Number of model points to process. |
[in] | lsoil | Number of soil layers. |
[out] | tsfanl | Skin temperature. |
[out] | stcanl | Soil/sea ice column temperature. |
[in] | tsfimx | Freezing point of sea water. |
Definition at line 6560 of file sfcsub.F.
References getstc().
Referenced by sfccycle().