UFO
ufo_sfcpcorrected_mod Module Reference

Fortran module for sfcpcorrected observation operator. More...

Data Types

type  ufo_sfcpcorrected
 Fortran derived type for the observation type. More...
 

Functions/Subroutines

subroutine ufo_sfcpcorrected_setup (self, f_conf)
 
subroutine ufo_sfcpcorrected_simobs (self, geovals, obss, nvars, nlocs, hofx)
 
subroutine da_intpsfc_prs (nobs, missing, P_o2m, H_o, P_o, H_m, TV_m, T_o, Q_o)
 \Conduct terrain height correction for surface pressure More...
 
subroutine da_intpsfc_prs_ukmo (nobs, missing, P_o2m, H_o, P_o, H_m, P_m, TV_2000, P_2000)
 \Conduct terrain height correction for surface pressure More...
 

Variables

integer, parameter max_string = 800
 
character(len=maxvarlen), dimension(5) geovars_list = (/ var_ps, var_geomz, var_sfc_geomz, var_tv, var_prs /)
 

Detailed Description

Fortran module for sfcpcorrected observation operator.

Function/Subroutine Documentation

◆ da_intpsfc_prs()

subroutine ufo_sfcpcorrected_mod::da_intpsfc_prs ( integer, intent(in)  nobs,
real(c_double), intent(in)  missing,
real(kind_real), dimension(nobs), intent(out)  P_o2m,
real(kind_real), dimension(nobs), intent(in)  H_o,
real(kind_real), dimension(nobs), intent(in)  P_o,
real(kind_real), dimension(nobs), intent(in)  H_m,
real(kind_real), dimension(nobs), intent(in)  TV_m,
real(kind_real), dimension(nobs), intent(in), optional  T_o,
real(kind_real), dimension(nobs), intent(in), optional  Q_o 
)

\Conduct terrain height correction for surface pressure

\This subroutine is based on a subroutine from WRFDA da_intpsfc_prs.inc file corresponding to sfc_assi_options = 1 in WRFDA's namelist

\Date: June 2019: Created \Method: hydrosatic equation

P_o2m = P_o * exp [-grav/rd * (H_m-H_o) / (TV_m + TV_o)/2)

Where: H_m = model surface height H_o = station height TV_m = virtual temperature at model surface height TV_o = virtual temperature at station height P_o2m = pressure interpolated from station height to model surface height P_o = pressure at station height grav = gravitational acceleration rd = gas constant per mole

Parameters
[in]nobstotal observation number
[out]p_o2mobserved PS at model sfc height
[in]p_oobserved Height and PS
[in]tv_mmodel sfc height and TV
[in]q_oobserbed T and Q

Definition at line 313 of file ufo_sfcpcorrected_mod.F90.

Here is the caller graph for this function:

◆ da_intpsfc_prs_ukmo()

subroutine ufo_sfcpcorrected_mod::da_intpsfc_prs_ukmo ( integer, intent(in)  nobs,
real(c_double), intent(in)  missing,
real(kind_real), dimension(nobs), intent(out)  P_o2m,
real(kind_real), dimension(nobs), intent(in)  H_o,
real(kind_real), dimension(nobs), intent(in)  P_o,
real(kind_real), dimension(nobs), intent(in)  H_m,
real(kind_real), dimension(nobs), intent(in)  P_m,
real(kind_real), dimension(nobs), intent(in)  TV_2000,
real(kind_real), dimension(nobs), intent(in)  P_2000 
)
private

\Conduct terrain height correction for surface pressure

\Reference: Ingleby,2013. UKMO Technical Report No: 582. Appendix 1.

\Method: integrate the hydrosatic equation dp/dz=-rho*g/RT to get P_m2o first, equation:

(P_m2o/P_m)=(T_m2o/T_m)** (grav/rd*L)

Where: P_m2o = model surface pressure at station height P_m = model surface pressure T_m = temperature at model surface height; derived from TV_2000 T_m2o = model surface temperature at station height grav = gravitational acceleration rd = gas constant per mole Lclr = constant lapse rate (0.0065 K/m)

To avoid dirunal/local variations, use TV_2000 (2000 m above the model surface height) instead of direct T_m

T_m = TV_2000 * (P_o / P_2000) ** (rd*L/grav)

Where: P_2000 = background pressure at 2000 m TV_2000 = background virtual temperature at 2000 m P_o = pressure at station height

Finally, in practice, adjust P_o to the model surface height using

P_o2m = P_o * (P_m / P_m2o)

Parameters
[in]nobstotal observation number
[out]p_o2mobserved PS at model sfc height
[in]p_oobserved Height and PS
[in]p_2000model Height, PS, TV at 2000 m, and P at 2000 m

Definition at line 379 of file ufo_sfcpcorrected_mod.F90.

Here is the caller graph for this function:

◆ ufo_sfcpcorrected_setup()

subroutine ufo_sfcpcorrected_mod::ufo_sfcpcorrected_setup ( class(ufo_sfcpcorrected), intent(inout)  self,
type(fckit_configuration), intent(in)  f_conf 
)
private

In the case where a user wants to specify the geoVaLs variable name of model height of vertical levels and/or the surface height. Example: MPAS is height but FV-3 uses geopotential_height.

Definition at line 40 of file ufo_sfcpcorrected_mod.F90.

◆ ufo_sfcpcorrected_simobs()

subroutine ufo_sfcpcorrected_mod::ufo_sfcpcorrected_simobs ( class(ufo_sfcpcorrected), intent(in)  self,
type(ufo_geovals), intent(in)  geovals,
type(c_ptr), intent(in), value  obss,
integer, intent(in)  nvars,
integer, intent(in)  nlocs,
real(c_double), dimension(nvars, nlocs), intent(inout)  hofx 
)

Definition at line 77 of file ufo_sfcpcorrected_mod.F90.

Here is the call graph for this function:

Variable Documentation

◆ geovars_list

character(len=maxvarlen), dimension(5) ufo_sfcpcorrected_mod::geovars_list = (/ var_ps, var_geomz, var_sfc_geomz, var_tv, var_prs /)

Definition at line 35 of file ufo_sfcpcorrected_mod.F90.

◆ max_string

integer, parameter ufo_sfcpcorrected_mod::max_string = 800
private

Definition at line 20 of file ufo_sfcpcorrected_mod.F90.