UFO
ufo_scatwind_neutralmetoffice_mod Module Reference

Fortran module for Met Office scatwind neutral wind forward operator. More...

Data Types

type  ufo_scatwind_neutralmetoffice
 Fortran derived type for neutral wind. More...
 

Functions/Subroutines

subroutine ufo_scatwind_neutralmetoffice_setup (self, f_conf)
 
subroutine ufo_scatwind_neutralmetoffice_simobs (self, geovals, obss, nvars, nlocs, hofx)
 Neutral wind forward operator for the Met Office system. More...
 
subroutine ops_scatwind_forwardmodel (za, u, v, ustr, oblen, seaice, orog, ycalc, cdr10)
 Scatterometer forward model. More...
 
subroutine ops_scatwind_phi_m_sea (recip_l_mo, z_uv, z0m, phi_m)
 Calculate the integrated froms of the Monin-Obukhov stability functions for surface exchanges. More...
 

Variables

character(len=maxvarlen), dimension(7), parameter geovars_default = (/ var_u, var_v, var_zi, var_sfc_ifrac, var_sfc_geomz, var_sea_fric_vel, var_obk_length /)
 

Detailed Description

Fortran module for Met Office scatwind neutral wind forward operator.

This code introduces the Met Office observation operator for scatterometer wind data. We assimilate the data as a "neutral" 10m wind, i.e. where the effects of atmospheric stability are neglected. For each observation we calculate the momentum roughness length using the Charnock relation. We then calculate the Monin-Obukhov stability function for momentum, integrated to the model's lowest wind level. The calculations are dependant upon on whether we have stable or unstable conditions according to the Obukhov Length. The neutral 10m wind components are then calculated from the lowest model level winds.

Author
J.Cotton (Met Office)
Date
22/12/2020: Created

Function/Subroutine Documentation

◆ ops_scatwind_forwardmodel()

subroutine ufo_scatwind_neutralmetoffice_mod::ops_scatwind_forwardmodel ( real(kind_real), dimension(:), intent(in)  za,
real(kind_real), dimension(:), intent(in)  u,
real(kind_real), dimension(:), intent(in)  v,
real(kind_real), intent(in)  ustr,
real(kind_real), intent(in)  oblen,
real(kind_real), intent(in)  seaice,
real(kind_real), intent(in)  orog,
real(kind_real), dimension(:), intent(inout)  ycalc,
real(kind_real), intent(inout)  cdr10 
)
private

Scatterometer forward model.

The main steps are as follows:

  • For each observation we calculate the momentum roughness length using the Charnock relation (Charnock, 1955) modified to include low-wind conditions (Smith, 1988):

    \[ z_{0m(sea)} = \frac{0.11\nu}{1.0 \times 10^{-5} + u_*} + \alpha_{ch} \frac{u_*^2}{g} \]

    where \(\nu = 14 \times 10^{-6}\) ms-1 is the dynamic viscosity of air, \(u_{*}\) is the friction velocity, \(\alpha_{ch} =0.018\) is the charnock parameter, and \(g\) is the acceleration due to gravity.
  • Call a subroutine to calculate the Monin-Obukhov stability function for momentum integrated to the lowest model level, \(\Phi_{m}\). The calculations are dependant upon on whether we have stable or unstable conditions according to the Obukhov Length, \(L\).
  • The neutral 10m wind components, \(v_{10n}\) are then calculated from the lowest model level winds, \(v_1\), as

    \[ \textbf{v}_{10n}=\frac{\ln\left(\left(10+z_{0m}\right)/z_{0m}\right)} {\Phi_m\left(L,z_1+z_{0m},z_{0m}\right)} \textbf{v}_1 \]

References:

Charnock, H. (1955). Wind stress on a water surface. Quart. J. Royal Meteorol. Soc., 81, 639-640. Smith, R. N. B. (1988). Coefficient for sea surface wind stress, heat flux and wind profiles as a function of wind speed and temperature. J. Geophys. Res., 93, 15467-15472.

Author
Met Office
Date
22/12/2020: Created
Parameters
[in]zaheights of rho levs
[in]uModel eastward wind profile
[in]vModel northward wind profile
[in]ustrModel friction velocity
[in]oblenModel obukhov length
[in]seaiceModel sea ice fraction
[in]orogModel orography
[in,out]ycalcModel equivalent of the obs
[in,out]cdr1010m interpolation coefficients

Definition at line 201 of file ufo_scatwind_neutralmetoffice_mod.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ops_scatwind_phi_m_sea()

subroutine ufo_scatwind_neutralmetoffice_mod::ops_scatwind_phi_m_sea ( real, intent(in)  recip_l_mo,
real, intent(in)  z_uv,
real, intent(in)  z0m,
real, intent(out)  phi_m 
)

Calculate the integrated froms of the Monin-Obukhov stability functions for surface exchanges.

The main steps are as follows:

  • In neutral conditions we have the logarithmic profile:

    \[ \Phi_{mn} = \ln\left(\frac{z_1 + z_{0m}}{z_{0m}}\right) \]

  • In stable conditions the stability functions of Beljaars and Holtslag (1991) are used:

    \[ \Phi_{m} = \Phi_{mn} + a\left(\zeta_1 - \zeta_{0m}\right) + b\left( \left(\zeta_1 - \frac{c}{d}\right) \exp\left(-d\zeta_1\right) - \left(\zeta_{0m} - \frac{c}{d}\right) \exp\left(-d\zeta_{0m}\right) \right) \]

    with \(\zeta_1=(z_1+z_{0m})/L\), \(\zeta_{0m}=z_{0m}/L\), and a=1, b=2/3, c=5, d=0.35.
  • In unstable conditions the Dyer and Hicks forms (Dyer, 1974) are used:

    \[ \Phi_{m} = \Phi_{mn} - 2\ln\left(\frac{1+X_1}{1+X_0}\right) - \ln\left(\frac{1+X_1^2}{1+X_0^2}\right) + 2\left(\tan^{-1}(X_1) - \tan^{-1}(X_0)\right) \]

    with \(X_1=(1-16\zeta_1)^{1/4}\) and \(X_0=(1-16\zeta_{0m})^{1/4}\) and

References:

Beljaars, A. C. M. and Holtslag, A. A. M. (1991). Flux parameterisation over land surfaces for atmospheric models. J. Appl. Meteor., 30, 327-341.

Dyer, A. J., 1974: A review of flux-profile relationships. Bound. Layer Meteor., 7, 363-372.

Author
Met Office
Date
22/12/2020: Created
Parameters
[in]recip_l_moReciprocal of Monin-Obukhov length (m^-1).
[in]z_uvHeight of wind level above roughness height(m).
[in]z0mRoughness length for momentum (m).
[out]phi_mStability function for momentum.

Definition at line 402 of file ufo_scatwind_neutralmetoffice_mod.F90.

Here is the caller graph for this function:

◆ ufo_scatwind_neutralmetoffice_setup()

subroutine ufo_scatwind_neutralmetoffice_mod::ufo_scatwind_neutralmetoffice_setup ( class(ufo_scatwind_neutralmetoffice), intent(inout)  self,
type(fckit_configuration), intent(in)  f_conf 
)
private

Definition at line 61 of file ufo_scatwind_neutralmetoffice_mod.F90.

◆ ufo_scatwind_neutralmetoffice_simobs()

subroutine ufo_scatwind_neutralmetoffice_mod::ufo_scatwind_neutralmetoffice_simobs ( class(ufo_scatwind_neutralmetoffice), 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 
)

Neutral wind forward operator for the Met Office system.

Author
Met Office
Date
22/12/2020: Created
Parameters
[in]selfThe object in which this operator is contained
[in]nlocsThe number of variables and locations
[in]geovalsThe model values, interpolated to the obsevation locations
[in,out]hofxThe output model equivalent of the observations
[in]obssThe observations, and meta-data for those observations

Definition at line 79 of file ufo_scatwind_neutralmetoffice_mod.F90.

Here is the call graph for this function:

Variable Documentation

◆ geovars_default

character(len=maxvarlen), dimension(7), parameter ufo_scatwind_neutralmetoffice_mod::geovars_default = (/ var_u, var_v, var_zi, var_sfc_ifrac, var_sfc_geomz, var_sea_fric_vel, var_obk_length /)

Definition at line 48 of file ufo_scatwind_neutralmetoffice_mod.F90.