chgres_cube  1.12.0
 All Data Structures Files Functions Variables
atmosphere Module Reference

Process atmospheric fields. More...

Public Member Functions

subroutine, public atmosphere_driver (localpet)
 Driver routine to process for atmospheric fields. More...
 
subroutine create_atm_b4adj_esmf_fields
 Create target grid field objects to hold data before vertical interpolation. More...
 
subroutine, public read_vcoord_info
 Reads model vertical coordinate definition file (as specified by namelist variable vcoord_file_target_grid). More...
 
subroutine terp3 (IM, IXZ1, IXQ1, IXZ2, IXQ2, NM, NXQ1, NXQ2, KM1, KXZ1, KXQ1, Z1, Q1, KM2, KXZ2, KXQ2, Z2, Q2)
 Cubically interpolate in one dimension. More...
 

Private Member Functions

subroutine cleanup_all_target_atm_data
 Cleanup target grid atmospheric field objects. More...
 
subroutine cleanup_target_atm_b4adj_data
 Cleanup atmospheric field (before adjustment) objects. More...
 
subroutine compute_zh
 Compute vertical level height. More...
 
subroutine convert_winds_to_uv
 Convert 3-d component winds to u and v. More...
 
subroutine create_atm_esmf_fields
 Create target grid field objects. More...
 
subroutine horiz_interp_thomp_mp_climo
 Horizontally interpolate thompson microphysics data to the target model grid. More...
 
subroutine newpr1 (localpet)
 Computes 3-D pressure given an adjusted surface pressure. More...
 
subroutine newps (localpet)
 Computes adjusted surface pressure given a new terrain height. More...
 
subroutine rsearch (IM, KM1, IXZ1, KXZ1, Z1, KM2, IXZ2, KXZ2, Z2, IXL2, KXL2, L2)
 Search for a surrounding real interval. More...
 
subroutine vintg
 Vertically interpolate upper-air fields. More...
 
subroutine vintg_thomp_mp_climo
 Vertically interpolate atmospheric fields to target FV3 grid. More...
 
subroutine vintg_wam (YEAR, MONTH, DAY, HOUR, PF)
 Vertically extend model top into thermosphere for whole atmosphere model. More...
 

Private Attributes

type(esmf_field) dzdt_b4adj_target_grid
 vertical vel before vert adj More...
 
type(esmf_field) pres_b4adj_target_grid
 3-d pres before terrain adj More...
 
type(esmf_field) pres_target_grid
 3-d pressure More...
 
type(esmf_field) ps_b4adj_target_grid
 sfc pres before terrain adj More...
 
type(esmf_field) qnifa_climo_b4adj_target_grid
 number concentration of ice friendly aerosols before vert adj More...
 
type(esmf_field) qnwfa_climo_b4adj_target_grid
 number concentration of water friendly aerosols before vert adj More...
 
type(esmf_field) temp_b4adj_target_grid
 temp before vert adj More...
 
type(esmf_field) terrain_interp_to_target_grid
 Input grid terrain interpolated to target grid. More...
 
type(esmf_field) thomp_pres_climo_b4adj_target_grid
 pressure of each level on target grid More...
 
type(esmf_field), dimension(:),
allocatable 
tracers_b4adj_target_grid
 tracers before vert adj More...
 
type(esmf_field) xwind_b4adj_target_grid
 x-component wind, before vert adj More...
 
type(esmf_field) xwind_s_target_grid
 x-component wind, 'south' edge More...
 
type(esmf_field) xwind_target_grid
 x-component wind, grid box center More...
 
type(esmf_field) xwind_w_target_grid
 x-component wind, 'west' edge More...
 
type(esmf_field) ywind_b4adj_target_grid
 y-component wind, before vert adj More...
 
type(esmf_field) ywind_s_target_grid
 y-component wind, 'south' edge More...
 
type(esmf_field) ywind_target_grid
 y-component wind, grid box center More...
 
type(esmf_field) ywind_w_target_grid
 y-component wind, 'west' edge More...
 
type(esmf_field) zwind_b4adj_target_grid
 z-component wind, before vert adj More...
 
type(esmf_field) zwind_s_target_grid
 z-component wind, 'south' edge More...
 
type(esmf_field) zwind_target_grid
 z-component wind, grid box center More...
 
type(esmf_field) zwind_w_target_grid
 z-component wind, 'west' edge More...
 

Detailed Description

Process atmospheric fields.

Horizontally interpolate from input to target FV3 grid using ESMF regridding. Adjust surface pressure according to terrain differences between input and target grid. Vertically interpolate to target FV3 grid vertical levels. Processing based on the spectral GFS version of CHGRES.

For variables "b4adj" indicates fields on the target grid before vertical adjustment. "target" indicates data on target grid. "input" indicates data on input grid. "_s" indicates fields on the 'south' edge of the grid box. "_w" indicate fields on the 'west' edge of the grid box. Otherwise, fields are at the center of the grid box.

Author
George Gayno NCEP/EMC

Definition at line 19 of file atmosphere.F90.

Member Function/Subroutine Documentation

subroutine atmosphere::cleanup_all_target_atm_data ( )
private

Cleanup target grid atmospheric field objects.

Author
George Gayno

Definition at line 2385 of file atmosphere.F90.

References atmosphere_target_data::cleanup_atmosphere_target_data().

Referenced by atmosphere_driver().

subroutine atmosphere::cleanup_target_atm_b4adj_data ( )
private

Cleanup atmospheric field (before adjustment) objects.

Author
George Gayno

Definition at line 2358 of file atmosphere.F90.

Referenced by atmosphere_driver().

subroutine atmosphere::compute_zh ( )
private

Compute vertical level height.

Author
George Gayno

Definition at line 2268 of file atmosphere.F90.

References utilities::error_handler().

Referenced by atmosphere_driver().

subroutine atmosphere::convert_winds_to_uv ( )
private

Convert 3-d component winds to u and v.

Author
George Gayno

Definition at line 779 of file atmosphere.F90.

References utilities::error_handler().

Referenced by atmosphere_driver().

subroutine atmosphere::create_atm_b4adj_esmf_fields ( )

Create target grid field objects to hold data before vertical interpolation.

These will be defined with the same number of vertical levels as the input grid.

Author
George Gayno

Definition at line 494 of file atmosphere.F90.

References utilities::error_handler().

Referenced by atmosphere_driver().

subroutine atmosphere::create_atm_esmf_fields ( )
private

Create target grid field objects.

Author
George Gayno

Definition at line 586 of file atmosphere.F90.

References utilities::error_handler().

Referenced by atmosphere_driver().

subroutine atmosphere::horiz_interp_thomp_mp_climo ( )
private

Horizontally interpolate thompson microphysics data to the target model grid.

Author
George Gayno

Definition at line 1291 of file atmosphere.F90.

References thompson_mp_climo_data::cleanup_thomp_mp_climo_input_data(), and utilities::error_handler().

Referenced by atmosphere_driver().

subroutine atmosphere::newpr1 ( integer, intent(in)  localpet)
private

Computes 3-D pressure given an adjusted surface pressure.

program history log: 2005-04-11 Hann-Ming Henry Juang hybrid sigma, sigma-p, and sigma-

  • PRGMMR: Henry Juang ORG: W/NMC23 DATE: 2005-04-11
  • PRGMMR: Fanglin Yang ORG: W/NMC23 DATE: 2006-11-28
  • PRGMMR: S. Moorthi ORG: NCEP/EMC DATE: 2006-12-12
  • PRGMMR: S. Moorthi ORG: NCEP/EMC DATE: 2007-01-02

    INPUT ARGUMENT LIST: IM INTEGER NUMBER OF POINTS TO COMPUTE KM INTEGER NUMBER OF LEVELS IDVC INTEGER VERTICAL COORDINATE ID (1 FOR SIGMA AND 2 FOR HYBRID) IDSL INTEGER TYPE OF SIGMA STRUCTURE (1 FOR PHILLIPS OR 2 FOR MEAN) NVCOORD INTEGER NUMBER OF VERTICAL COORDINATES VCOORD REAL (KM+1,NVCOORD) VERTICAL COORDINATE VALUES FOR IDVC=1, NVCOORD=1: SIGMA INTERFACE FOR IDVC=2, NVCOORD=2: HYBRID INTERFACE A AND B FOR IDVC=3, NVCOORD=3: JUANG GENERAL HYBRID INTERFACE AK REAL (KM+1) HYBRID INTERFACE A BK REAL (KM+1) HYBRID INTERFACE B PS REAL (IX) SURFACE PRESSURE (PA) OUTPUT ARGUMENT LIST: PM REAL (IX,KM) MID-LAYER PRESSURE (PA) DP REAL (IX,KM) LAYER DELTA PRESSURE (PA)

Parameters
[in]localpetESMF local persistent execution thread
Author
Hann Ming Henry Juang, Juang, Fanglin Yang, S. Moorthi

Definition at line 948 of file atmosphere.F90.

References utilities::error_handler().

Referenced by atmosphere_driver().

subroutine atmosphere::newps ( integer, intent(in)  localpet)
private

Computes adjusted surface pressure given a new terrain height.

Computes a new surface pressure given a new orography. The new pressure is computed assuming a hydrostatic balance and a constant temperature lapse rate. Below ground, the lapse rate is assumed to be -6.5 k/km.

program history log:

  • 91-10-31 mark iredell
  • 2018-apr adapt for fv3. george gayno
Parameters
[in]localpetESMF local persistent execution thread
Author
Mark Iredell, George Gayno
Date
92-10-31

Definition at line 1059 of file atmosphere.F90.

References utilities::error_handler().

Referenced by atmosphere_driver().

subroutine, public atmosphere::read_vcoord_info ( )

Reads model vertical coordinate definition file (as specified by namelist variable vcoord_file_target_grid).

Author
George Gayno

Definition at line 1256 of file atmosphere.F90.

References utilities::error_handler().

Referenced by atmosphere_driver().

subroutine atmosphere::rsearch ( integer, intent(in)  IM,
integer, intent(in)  KM1,
integer, intent(in)  IXZ1,
integer, intent(in)  KXZ1,
real(esmf_kind_r8), dimension(1+(im-1)*ixz1+(km1-1)*kxz1), intent(in)  Z1,
integer, intent(in)  KM2,
integer, intent(in)  IXZ2,
integer, intent(in)  KXZ2,
real(esmf_kind_r8), dimension(1+(im-1)*ixz2+(km2-1)*kxz2), intent(in)  Z2,
integer, intent(in)  IXL2,
integer, intent(in)  KXL2,
integer, dimension(1+(im-1)*ixl2+(km2-1)*kxl2), intent(out)  L2 
)
private

Search for a surrounding real interval.

This subprogram searches monotonic sequences of real numbers for intervals that surround a given search set of real numbers. The sequences may be monotonic in either direction; the real numbers may be single or double precision; the input sequences and sets and the output locations may be arbitrarily dimensioned.

If the array z1 is dimensioned (im,km1), then the skip numbers are ixz1=1 and kxz1=im; if it is dimensioned (km1,im), then the skip numbers are ixz1=km1 and kxz1=1; if it is dimensioned (im,jm,km1), then the skip numbers are ixz1=1 and kxz1=im*jm; etcetera. Similar examples apply to the skip numbers for z2 and l2.

Returned values of 0 or km1 indicate that the given search value is outside the range of the sequence.

If a search value is identical to one of the sequence values then the location returned points to the identical value. If the sequence is not strictly monotonic and a search value is identical to more than one of the sequence values, then the location returned may point to any of the identical values.

to be exact, for each i from 1 to im and for each k from 1 to km2, z=z2(1+(i-1)*ixz2+(k-1)*kxz2) is the search value and l=l2(1+(i-1)*ixl2+(k-1)*kxl2) is the location returned. if l=0, then z is less than the start point z1(1+(i-1)*ixz1) for ascending sequences (or greater than for descending sequences). if l=km1, then z is greater than or equal to the end point z1(1+(i-1)*ixz1+(km1-1)*kxz1) for ascending sequences (or less than or equal to for descending sequences). otherwise z is between the values z1(1+(i-1)*ixz1+(l-1)*kxz1) and z1(1+(i-1)*ixz1+(l-0)*kxz1) and may equal the former.

Parameters
[in]iminteger number of sequences to search
[in]km1integer number of points in each sequence
[in]ixz1integer sequence skip number for z1
[in]kxz1integer point skip number for z1
[in]z1real (1+(im-1)*ixz1+(km1-1)*kxz1) sequence values to search (z1 must be monotonic in either direction)
[in]km2integer number of points to search for in each respective sequence
[in]ixz2integer sequence skip number for z2
[in]kxz2integer point skip number for z2
[in]z2real (1+(im-1)*ixz2+(km2-1)*kxz2) set of values to search for (z2 need not be monotonic)
[in]ixl2integer sequence skip number for l2
[in]kxl2integer point skip number for l2
[out]l2integer (1+(im-1)*ixl2+(km2-1)*kxl2) interval locations having values from 0 to km1 (z2 will be between z1(l2) and z1(l2+1))
Author
Mark Iredell
Date
98-05-01

Definition at line 2220 of file atmosphere.F90.

Referenced by terp3().

subroutine atmosphere::terp3 ( integer  IM,
integer  IXZ1,
integer  IXQ1,
integer  IXZ2,
integer  IXQ2,
integer  NM,
integer  NXQ1,
integer  NXQ2,
integer  KM1,
integer  KXZ1,
integer  KXQ1,
real(esmf_kind_r8), dimension(1+(im-1)*ixz1+(km1-1)*kxz1)  Z1,
real(esmf_kind_r8), dimension(1+(im-1)*ixq1+(km1-1)*kxq1+(nm-1)*nxq1)  Q1,
integer  KM2,
integer  KXZ2,
integer  KXQ2,
real(esmf_kind_r8), dimension(1+(im-1)*ixz2+(km2-1)*kxz2)  Z2,
real(esmf_kind_r8), dimension(1+(im-1)*ixq2+(km2-1)*kxq2+(nm-1)*nxq2)  Q2 
)

Cubically interpolate in one dimension.

Interpolate field(s) in one dimension along the column(s). The interpolation is cubic lagrangian with a monotonic constraint in the center of the domain. In the outer intervals it is linear. Outside the domain, fields are held constant.

PROGRAM HISTORY LOG:

  • 98-05-01 MARK IREDELL
  • 1999-01-04 IREDELL USE ESSL SEARCH
Parameters
[in]iminteger number of columns
[in]ixz1integer column skip number for z1
[in]ixq1integer column skip number for q1
[in]ixz2integer column skip number for z2
[in]ixq2integer column skip number for q2
[in]nminteger number of fields per column
[in]nxq1integer field skip number for q1
[in]nxq2integer field skip number for q2
[in]km1integer number of input points
[in]kxz1integer point skip number for z1
[in]kxq1integer point skip number for q1
[in]z1real (1+(im-1)*ixz1+(km1-1)*kxz1) input coordinate values in which to interpolate (z1 must be strictly monotonic in either direction)
[in]q1real (1+(im-1)*ixq1+(km1-1)*kxq1+(nm-1)*nxq1) input fields to interpolate
[in]km2integer number of output points
[in]kxz2integer point skip number for z2
[in]kxq2integer point skip number for q2
[in]z2real (1+(im-1)*ixz2+(km2-1)*kxz2) output coordinate values to which to interpolate (z2 need not be monotonic)
[out]q2real (1+(im-1)*ixq2+(km2-1)*kxq2+(nm-1)*nxq2) output interpolated fields
Author
Mark Iredell
Date
98-05-01

Definition at line 2026 of file atmosphere.F90.

References rsearch().

Referenced by vintg(), and vintg_thomp_mp_climo().

subroutine atmosphere::vintg ( )
private

Vertically interpolate upper-air fields.

Vertically interpolate upper-air fields. Wind, temperature, humidity and other tracers are interpolated. The interpolation is cubic lagrangian in log pressure with a monotonic constraint in the center of the domain. In the outer intervals it is linear in log pressure. Outside the domain, fields are generally held constant, except for temperature and humidity below the input domain, where the temperature lapse rate is held fixed at -6.5 k/km and the relative humidity is held constant. This routine expects fields ordered from bottom to top of atmosphere.

Author
Mark Iredell
Date
92-10-31

Definition at line 1763 of file atmosphere.F90.

References utilities::error_handler(), and terp3().

Referenced by atmosphere_driver().

subroutine atmosphere::vintg_thomp_mp_climo ( )
private

Vertically interpolate atmospheric fields to target FV3 grid.

Vertically interpolate thompson microphysics climo tracers to the target model levels.

Author
George Gayno

Definition at line 1403 of file atmosphere.F90.

References utilities::error_handler(), and terp3().

Referenced by atmosphere_driver().

subroutine atmosphere::vintg_wam ( integer, intent(in)  YEAR,
integer, intent(in)  MONTH,
integer, intent(in)  DAY,
integer, intent(in)  HOUR,
character(*), intent(in)  PF 
)
private

Vertically extend model top into thermosphere for whole atmosphere model.

Use climatological data to extent model top into thermosphere for temperature and consoder primary compositions of neutral atmosphere in term of specific values of oxygen, single oxygen, and ozone.

Parameters
[in]yearinitial year
[in]monthinitial month
[in]dayinitial day
[in]hourinitial hour
[in]pfpath to MSIS2.1 parm file
Author
Hann-Ming Henry Juang NCEP/EMC

Definition at line 1527 of file atmosphere.F90.

References utilities::error_handler(), and gettemp().

Referenced by atmosphere_driver().

Field Documentation

type(esmf_field) atmosphere::dzdt_b4adj_target_grid
private

vertical vel before vert adj

Definition at line 80 of file atmosphere.F90.

type(esmf_field) atmosphere::pres_b4adj_target_grid
private

3-d pres before terrain adj

Definition at line 84 of file atmosphere.F90.

type(esmf_field) atmosphere::pres_target_grid
private

3-d pressure

Definition at line 83 of file atmosphere.F90.

type(esmf_field) atmosphere::ps_b4adj_target_grid
private

sfc pres before terrain adj

Definition at line 82 of file atmosphere.F90.

type(esmf_field) atmosphere::qnifa_climo_b4adj_target_grid
private

number concentration of ice friendly aerosols before vert adj

Definition at line 102 of file atmosphere.F90.

type(esmf_field) atmosphere::qnwfa_climo_b4adj_target_grid
private

number concentration of water friendly aerosols before vert adj

Definition at line 104 of file atmosphere.F90.

type(esmf_field) atmosphere::temp_b4adj_target_grid
private

temp before vert adj

Definition at line 85 of file atmosphere.F90.

type(esmf_field) atmosphere::terrain_interp_to_target_grid
private

Input grid terrain interpolated to target grid.

Definition at line 86 of file atmosphere.F90.

type(esmf_field) atmosphere::thomp_pres_climo_b4adj_target_grid
private

pressure of each level on target grid

Definition at line 106 of file atmosphere.F90.

type(esmf_field), dimension(:), allocatable atmosphere::tracers_b4adj_target_grid
private

tracers before vert adj

Definition at line 81 of file atmosphere.F90.

type(esmf_field) atmosphere::xwind_b4adj_target_grid
private

x-component wind, before vert adj

Definition at line 90 of file atmosphere.F90.

type(esmf_field) atmosphere::xwind_s_target_grid
private

x-component wind, 'south' edge

Definition at line 93 of file atmosphere.F90.

type(esmf_field) atmosphere::xwind_target_grid
private

x-component wind, grid box center

Definition at line 87 of file atmosphere.F90.

type(esmf_field) atmosphere::xwind_w_target_grid
private

x-component wind, 'west' edge

Definition at line 96 of file atmosphere.F90.

type(esmf_field) atmosphere::ywind_b4adj_target_grid
private

y-component wind, before vert adj

Definition at line 91 of file atmosphere.F90.

type(esmf_field) atmosphere::ywind_s_target_grid
private

y-component wind, 'south' edge

Definition at line 94 of file atmosphere.F90.

type(esmf_field) atmosphere::ywind_target_grid
private

y-component wind, grid box center

Definition at line 88 of file atmosphere.F90.

type(esmf_field) atmosphere::ywind_w_target_grid
private

y-component wind, 'west' edge

Definition at line 97 of file atmosphere.F90.

type(esmf_field) atmosphere::zwind_b4adj_target_grid
private

z-component wind, before vert adj

Definition at line 92 of file atmosphere.F90.

type(esmf_field) atmosphere::zwind_s_target_grid
private

z-component wind, 'south' edge

Definition at line 95 of file atmosphere.F90.

type(esmf_field) atmosphere::zwind_target_grid
private

z-component wind, grid box center

Definition at line 89 of file atmosphere.F90.

type(esmf_field) atmosphere::zwind_w_target_grid
private

z-component wind, 'west' edge

Definition at line 98 of file atmosphere.F90.


The documentation for this module was generated from the following file: