FV3 Bundle
qg_fields Module Reference

Handle fields for the QG model. More...

Data Types

type  qg_field
 Fortran derived type to hold QG fields. More...
 

Functions/Subroutines

subroutine, public create (self, geom, vars)
 Linked list implementation. More...
 
subroutine, public delete (self)
 
subroutine, public zeros (self)
 
subroutine, public dirac (self, c_conf)
 
subroutine, public random (self)
 
subroutine, public copy (self, rhs)
 
subroutine, public self_add (self, rhs)
 
subroutine, public self_schur (self, rhs)
 
subroutine, public self_sub (self, rhs)
 
subroutine, public self_mul (self, zz)
 
subroutine, public axpy (self, zz, rhs)
 
subroutine, public dot_prod (fld1, fld2, zprod)
 
subroutine, public add_incr (self, rhs)
 
subroutine, public diff_incr (lhs, x1, x2)
 
subroutine, public change_resol (fld, rhs)
 
subroutine, public read_file (fld, c_conf, vdate)
 
subroutine, public analytic_init (fld, geom, config, vdate)
 Analytic Initialization for the QG model. More...
 
subroutine calc_streamfunction (kx, ky, x, u, v, dx, dy)
 Calculate Streamfunction from u and v. More...
 
subroutine deriv_1d (df, f, n, delta, periodic)
 Finite Difference Derivative. More...
 
subroutine, public write_file (fld, c_conf, vdate)
 
subroutine, public gpnorm (fld, nf, pstat)
 
subroutine, public fldrms (fld, prms)
 
subroutine, public interp_tl (fld, locs, vars, gom)
 
subroutine, public interp_ad (fld, locs, vars, gom)
 
subroutine lin_weights (kk, delta1, delta2, k1, k2, w1, w2)
 
subroutine ug_size (self, ug)
 
subroutine, public ug_coord (self, ug, colocated)
 
subroutine, public field_to_ug (self, ug, colocated)
 
subroutine, public field_from_ug (self, ug)
 
character(len=length) function genfilename (c_conf, length, vdate)
 
integer function common_vars (x1, x2)
 
subroutine check_resolution (x1, x2)
 
subroutine check (self)
 
subroutine ncerr (info)
 

Variables

logical netcdfio = .true.
 
type(registry_t), public qg_field_registry
 Linked list interface - defines registry_t type. More...
 

Detailed Description

Handle fields for the QG model.

Function/Subroutine Documentation

◆ add_incr()

subroutine, public qg_fields::add_incr ( type(qg_field), intent(inout)  self,
type(qg_field), intent(in)  rhs 
)

Definition at line 392 of file qg_fields.F90.

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

◆ analytic_init()

subroutine, public qg_fields::analytic_init ( type(qg_field), intent(inout)  fld,
type(qg_geom), intent(inout)  geom,
type(c_ptr), intent(in)  config,
type(datetime), intent(inout)  vdate 
)

Analytic Initialization for the QG model.

analytic_init() initializes the qg field using one of several alternative idealized analytic models. This is intended to faciltitate testing by eliminating the need to read in the initial state from a file and by providing exact expressions to test interpolations. This function is activated by setting the initial.analytic_init field in the configuration file.

Initialization options that begin with "dcmip" refer to tests defined by the multi-institutional 2012 Dynamical Core Intercomparison Project and the associated Summer School, sponsored by NOAA, NSF, DOE, NCAR, and the University of Michigan.

Currently implemented options for analytic_init include:

  • baroclinic-instability: A two-layer flow that is baroclinically unstable in the presence of non-zero orography
  • dcmip-test-1-1: 3D deformational flow
  • dcmip-test-1-2: 3D Hadley-like meridional circulation

Relevant parameters specified in the State.StateGenerate section of the config file include:

  • analytic_init the type of initialization as described in analytic_init())
  • variables currently must be set to ["x","q","u","v"]
  • date a reference date and time for the state
  • top_layer_depth, bottom_layer_depth are needed in order to convert the normalized height coordinate into meters so it can be used in the analytic formulae
Author
M. Miesch (JCSDA, adapted from a pre-existing call to invent_state)
Date
March 7, 2018: Created
Warning
If you use the dcmip routines, be sure to include u and v in the variables list.
For dcmip initial conditions, this routine computes the streamfunction x and the potential vorticity q from the horizontal winds u and v. Since q (or alternatively x) is the state variable for the evolution of the QG model, only the non-divergent (vortical) component of the horizontal flow is captured. In other words, u and v are projected onto a flow field with zero horizontal divergence. Furthermore, this projection neglects geometric curvature factors. This projection is used to calculate the streamfunction X and the potential vorticity q. However, on exiting this function, the u and v components of the qg_field structure are still set to their original (non-projected) values. This is for use with the State interpolation test. The winds u and v are re-calculated from X before doing a forecast by the routine c_qg_prepare_integration() in c_qg_model.F90.
For all but the baroclinic-instability option, the PV here is computed assuming zero orography. For particular model configurations, the orography is defined in c_qg_setup() [c_qg_model.F90] and the PV is recalculated in c_qg_prepare_integration() [c_qg_model.F90] before the integration of the model commences.
Parameters
[in,out]fldFields
[in,out]geomGrid information
[in]configConfiguration
[in,out]vdateDateTime

Definition at line 707 of file qg_fields.F90.

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

◆ axpy()

subroutine, public qg_fields::axpy ( type(qg_field), intent(inout)  self,
real(kind=kind_real), intent(in)  zz,
type(qg_field), intent(in)  rhs 
)

Definition at line 333 of file qg_fields.F90.

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

◆ calc_streamfunction()

subroutine qg_fields::calc_streamfunction ( integer, intent(in)  kx,
integer, intent(in)  ky,
real(kind=kind_real), dimension(kx,ky,2), intent(out)  x,
real(kind=kind_real), dimension(kx,ky,2), intent(in)  u,
real(kind=kind_real), dimension(kx,ky,2), intent(in)  v,
real(kind=kind_real), intent(in)  dx,
real(kind=kind_real), intent(in)  dy 
)

Calculate Streamfunction from u and v.

calc_streamfunction() computes the streamfunction x from the velocities u and v. It was originally developed to support the dcmip initial conditions in analytic_init() but is potentially of broader applicability. The algorithm first computes the vertical vorticity and then inverts the resulting Laplacian for x.

Author
M. Miesch
Date
March 7, 2018: Created
Warning
This routine currently does not take into account latitudinal boundary conditions on x.






Parameters
[in]kxZonal grid dimension
[in]kyMeridional grid dimension
[out]xStreamfunction
[in]uzonal velocity (non-dimensional)
[in]vmeridional velocity (non-dimensional)
[in]dxZonal grid spacing (non-dimensional)
[in]dyMeridional grid spacing (non-dimensional)

Definition at line 895 of file qg_fields.F90.

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

◆ change_resol()

subroutine, public qg_fields::change_resol ( type(qg_field), intent(inout)  fld,
type(qg_field), intent(in)  rhs 
)

Definition at line 463 of file qg_fields.F90.

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

◆ check()

subroutine qg_fields::check ( type(qg_field), intent(in)  self)
private

Definition at line 1671 of file qg_fields.F90.

◆ check_resolution()

subroutine qg_fields::check_resolution ( type(qg_field), intent(in)  x1,
type(qg_field), intent(in)  x2 
)
private

Definition at line 1656 of file qg_fields.F90.

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

◆ common_vars()

integer function qg_fields::common_vars ( type(qg_field), intent(in)  x1,
type(qg_field), intent(in)  x2 
)

Definition at line 1633 of file qg_fields.F90.

Here is the caller graph for this function:

◆ copy()

subroutine, public qg_fields::copy ( type(qg_field), intent(inout)  self,
type(qg_field), intent(in)  rhs 
)

Definition at line 225 of file qg_fields.F90.

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

◆ create()

subroutine, public qg_fields::create ( type(qg_field), intent(inout)  self,
type(qg_geom), intent(in), target  geom,
type(qg_vars), intent(in)  vars 
)

Linked list implementation.

Definition at line 76 of file qg_fields.F90.

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

◆ delete()

subroutine, public qg_fields::delete ( type(qg_field), intent(inout)  self)

Definition at line 136 of file qg_fields.F90.

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

◆ deriv_1d()

subroutine qg_fields::deriv_1d ( real(kind=kind_real), dimension(:), intent(inout)  df,
real(kind=kind_real), dimension(:), intent(in)  f,
integer, intent(in)  n,
real(kind=kind_real), intent(in)  delta,
logical, intent(in), optional  periodic 
)
private

Finite Difference Derivative.

deriv() computes the first derivative of a function using a simple fourth-order, central-difference scheme

Author
M. Miesch
Date
March 7, 2018: Created
Parameters
[in]f1D function
[in]deltagrid spacing (assumed uniform)
[in]nNumber of grid points
[in,out]dfresult
[in]periodicSet to true if the domain is periodic

Definition at line 945 of file qg_fields.F90.

Here is the caller graph for this function:

◆ diff_incr()

subroutine, public qg_fields::diff_incr ( type(qg_field), intent(inout)  lhs,
type(qg_field), intent(in)  x1,
type(qg_field), intent(in)  x2 
)

Definition at line 424 of file qg_fields.F90.

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

◆ dirac()

subroutine, public qg_fields::dirac ( type(qg_field), intent(inout)  self,
type(c_ptr), intent(in)  c_conf 
)
Parameters
[in]c_confConfiguration

Definition at line 167 of file qg_fields.F90.

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

◆ dot_prod()

subroutine, public qg_fields::dot_prod ( type(qg_field), intent(in)  fld1,
type(qg_field), intent(in)  fld2,
real(kind=kind_real), intent(inout)  zprod 
)

Definition at line 355 of file qg_fields.F90.

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

◆ field_from_ug()

subroutine, public qg_fields::field_from_ug ( type(qg_field), intent(inout)  self,
type(unstructured_grid), intent(in)  ug 
)

Definition at line 1534 of file qg_fields.F90.

Here is the caller graph for this function:

◆ field_to_ug()

subroutine, public qg_fields::field_to_ug ( type(qg_field), intent(in)  self,
type(unstructured_grid), intent(inout)  ug,
integer, intent(in)  colocated 
)

Definition at line 1479 of file qg_fields.F90.

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

◆ fldrms()

subroutine, public qg_fields::fldrms ( type(qg_field), intent(in)  fld,
real(kind=kind_real), intent(out)  prms 
)

Definition at line 1126 of file qg_fields.F90.

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

◆ genfilename()

character(len=length) function qg_fields::genfilename ( type(c_ptr), intent(in)  c_conf,
integer, intent(in)  length,
type(datetime), intent(in)  vdate 
)
Parameters
[in]c_confConfiguration

Definition at line 1579 of file qg_fields.F90.

Here is the caller graph for this function:

◆ gpnorm()

subroutine, public qg_fields::gpnorm ( type(qg_field), intent(in)  fld,
integer, intent(in)  nf,
real(kind=kind_real), dimension(3, nf), intent(inout)  pstat 
)

Definition at line 1089 of file qg_fields.F90.

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

◆ interp_ad()

subroutine, public qg_fields::interp_ad ( type(qg_field), intent(inout)  fld,
type(qg_locs), intent(in)  locs,
type(qg_vars), intent(in)  vars,
type(qg_goms), intent(inout)  gom 
)

Definition at line 1253 of file qg_fields.F90.

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

◆ interp_tl()

subroutine, public qg_fields::interp_tl ( type(qg_field), intent(in)  fld,
type(qg_locs), intent(in)  locs,
type(qg_vars), intent(in)  vars,
type(qg_goms), intent(inout)  gom 
)

Definition at line 1166 of file qg_fields.F90.

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

◆ lin_weights()

subroutine qg_fields::lin_weights ( integer, intent(in)  kk,
real(kind=kind_real), intent(in)  delta1,
real(kind=kind_real), intent(in)  delta2,
integer, intent(out)  k1,
integer, intent(out)  k2,
real(kind=kind_real), intent(out)  w1,
real(kind=kind_real), intent(out)  w2 
)
private

Definition at line 1345 of file qg_fields.F90.

Here is the caller graph for this function:

◆ ncerr()

subroutine qg_fields::ncerr ( integer, intent(in)  info)
private
Parameters
[in]infoInfo index

Definition at line 1720 of file qg_fields.F90.

Here is the caller graph for this function:

◆ random()

subroutine, public qg_fields::random ( type(qg_field), intent(inout)  self)

Definition at line 208 of file qg_fields.F90.

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

◆ read_file()

subroutine, public qg_fields::read_file ( type(qg_field), intent(inout)  fld,
type(c_ptr), intent(in)  c_conf,
type(datetime), intent(inout)  vdate 
)
Parameters
[in,out]fldFields
[in]c_confConfiguration
[in,out]vdateDateTime

Definition at line 524 of file qg_fields.F90.

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

◆ self_add()

subroutine, public qg_fields::self_add ( type(qg_field), intent(inout)  self,
type(qg_field), intent(in)  rhs 
)

Definition at line 252 of file qg_fields.F90.

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

◆ self_mul()

subroutine, public qg_fields::self_mul ( type(qg_field), intent(inout)  self,
real(kind=kind_real), intent(in)  zz 
)

Definition at line 315 of file qg_fields.F90.

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

◆ self_schur()

subroutine, public qg_fields::self_schur ( type(qg_field), intent(inout)  self,
type(qg_field), intent(in)  rhs 
)

Definition at line 273 of file qg_fields.F90.

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

◆ self_sub()

subroutine, public qg_fields::self_sub ( type(qg_field), intent(inout)  self,
type(qg_field), intent(in)  rhs 
)

Definition at line 294 of file qg_fields.F90.

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

◆ ug_coord()

subroutine, public qg_fields::ug_coord ( type(qg_field), intent(in)  self,
type(unstructured_grid), intent(inout)  ug,
integer, intent(in)  colocated 
)

Definition at line 1422 of file qg_fields.F90.

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

◆ ug_size()

subroutine qg_fields::ug_size ( type(qg_field), intent(in)  self,
type(unstructured_grid), intent(inout)  ug 
)
private

Definition at line 1368 of file qg_fields.F90.

Here is the caller graph for this function:

◆ write_file()

subroutine, public qg_fields::write_file ( type(qg_field), intent(in)  fld,
type(c_ptr), intent(in)  c_conf,
type(datetime), intent(in)  vdate 
)
Parameters
[in]fldFields
[in]c_confConfiguration
[in]vdateDateTime

Definition at line 994 of file qg_fields.F90.

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

◆ zeros()

subroutine, public qg_fields::zeros ( type(qg_field), intent(inout)  self)

Definition at line 151 of file qg_fields.F90.

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

Variable Documentation

◆ netcdfio

logical qg_fields::netcdfio = .true.
private

Definition at line 57 of file qg_fields.F90.

◆ qg_field_registry

type(registry_t), public qg_fields::qg_field_registry

Linked list interface - defines registry_t type.

Global registry

Definition at line 65 of file qg_fields.F90.