FV3DYCORE  Version1.0.0
dyn_core_mod Module Reference

The module 'dyn_core' peforms the Lagrangian acoustic dynamics described by [lin2004vertically]. More...

Functions/Subroutines

subroutine, public dyn_core (npx, npy, npz, ng, sphum, nq, bdt, n_split, zvir, cp, akap, cappa, grav, hydrostatic, u, v, w, delz, pt, q, delp, pe, pk, phis, ws, omga, ptop, pfull, ua, va, uc, vc, mfx, mfy, cx, cy, pkz, peln, q_con, ak, bk, ks, gridstruct, flagstruct, neststruct, idiag, bd, domain, init_step, i_pack, end_step, diss_est, time_total)
 
subroutine pk3_halo (is, ie, js, je, isd, ied, jsd, jed, npz, ptop, akap, pk3, delp)
 
subroutine pln_halo (is, ie, js, je, isd, ied, jsd, jed, npz, ptop, pk3, delp)
 
subroutine pe_halo (is, ie, js, je, isd, ied, jsd, jed, npz, ptop, pe, delp)
 
subroutine adv_pe (ua, va, pem, om, gridstruct, bd, npx, npy, npz, ng)
 
subroutine p_grad_c (dt2, npz, delpc, pkc, gz, uc, vc, bd, rdxc, rdyc, hydrostatic)
 
subroutine nh_p_grad (u, v, pp, gz, delp, pk, dt, ng, gridstruct, bd, npx, npy, npz, use_logp)
 
subroutine split_p_grad (u, v, pp, gz, delp, pk, beta, dt, ng, gridstruct, bd, npx, npy, npz, use_logp)
 
subroutine one_grad_p (u, v, pk, gz, divg2, delp, dt, ng, gridstruct, bd, npx, npy, npz, ptop, hydrostatic, a2b_ord, d_ext)
 
subroutine grad1_p_update (divg2, u, v, pk, gz, dt, ng, gridstruct, bd, npx, npy, npz, ptop, beta, a2b_ord)
 
subroutine mix_dp (hydrostatic, w, delp, pt, km, ak, bk, CG, fv_debug, bd)
 
subroutine geopk (ptop, pe, peln, delp, pk, gz, hs, pt, q_con, pkz, km, akap, CG, nested, computehalo, npx, npy, a2b_ord, bd)
 The subroutine 'geopk' calculates geopotential and pressure to the kappa. More...
 
subroutine, public del2_cubed (q, cd, gridstruct, domain, npx, npy, km, nmax, bd)
 The subroutine 'del2-cubed' filters the omega field for the physics. More...
 
subroutine, public init_ijk_mem (i1, i2, j1, j2, km, array, var)
 
subroutine ray_fast (dt, npx, npy, npz, pfull, tau, u, v, w, ks, dp, ptop, hydrostatic, rf_cutoff, bd)
 The subroutine 'Ray_fast' computes a simple "inline" version of the Rayleigh friction (EXPERIMENTAL - NOT FOR GENERAL USE). More...
 

Variables

real ptk
 
real peln1
 
real rgrav
 
real d3_damp
 
real, dimension(:,:,:), allocatable ut
 
real, dimension(:,:,:), allocatable vt
 
real, dimension(:,:,:), allocatable crx
 
real, dimension(:,:,:), allocatable cry
 
real, dimension(:,:,:), allocatable xfx
 
real, dimension(:,:,:), allocatable yfx
 
real, dimension(:,:,:), allocatable divgd
 
real, dimension(:,:,:), allocatable zh
 
real, dimension(:,:,:), allocatable du
 
real, dimension(:,:,:), allocatable dv
 
real, dimension(:,:,:), allocatable pkc
 
real, dimension(:,:,:), allocatable delpc
 
real, dimension(:,:,:), allocatable pk3
 
real, dimension(:,:,:), allocatable ptc
 
real, dimension(:,:,:), allocatable gz
 
real(kind=r_grid), parameter cnst_0p20 =0.20d0
 
real, dimension(:), allocatable rf
 
integer k_rf = 0
 
logical rff_initialized = .false.
 
logical first_call = .true.
 
integer kmax =1
 

Detailed Description

The module 'dyn_core' peforms the Lagrangian acoustic dynamics described by [lin2004vertically].

The forward timestep is handled by routines in 'sw_core.F90'. The backwards-in-time PGF is evaluated in one_grad_p or split_p_grad (hydrostatic) and nh_p_grad (nonhydrostatic) see [lin1997explicit]. The nonhydrostatic components are handled by 'nh_core.F90'.

Function/Subroutine Documentation

◆ adv_pe()

subroutine dyn_core_mod::adv_pe ( real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed,npz), intent(in)  ua,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed,npz), intent(in)  va,
real, dimension(bd%is-1:bd%ie+1,1:npz+1,bd%js-1:bd%je+1), intent(in)  pem,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed,npz), intent(inout)  om,
type(fv_grid_type), intent(inout), target  gridstruct,
type(fv_grid_bounds_type), intent(in)  bd,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
integer, intent(in)  ng 
)
private

Definition at line 1554 of file dyn_core.F90.

◆ del2_cubed()

subroutine, public dyn_core_mod::del2_cubed ( real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed,km), intent(inout)  q,
real(kind=r_grid), intent(in)  cd,
type(fv_grid_type), intent(in), target  gridstruct,
type(domain2d), intent(inout)  domain,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  km,
integer, intent(in)  nmax,
type(fv_grid_bounds_type), intent(in)  bd 
)

The subroutine 'del2-cubed' filters the omega field for the physics.

Parameters
[in]cdcd = K * da_min; 0 < K < 0.25

Definition at line 2421 of file dyn_core.F90.

◆ dyn_core()

subroutine, public dyn_core_mod::dyn_core ( integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
integer, intent(in)  ng,
integer, intent(in)  sphum,
integer, intent(in)  nq,
real, intent(in)  bdt,
integer, intent(in)  n_split,
real, intent(in)  zvir,
real, intent(in)  cp,
real, intent(in)  akap,
real, dimension(bd%isd:,bd%jsd:,1:), intent(inout)  cappa,
real, intent(in)  grav,
logical, intent(in)  hydrostatic,
real, dimension(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz), intent(inout)  u,
real, dimension(bd%isd:bd%ied+1,bd%jsd:bd%jed ,npz), intent(inout)  v,
real, dimension( bd%isd:,bd%jsd:,1:), intent(inout)  w,
real, dimension(bd%isd:,bd%jsd:,1:), intent(inout)  delz,
real, dimension( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz), intent(inout)  pt,
real, dimension( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz, nq), intent(inout)  q,
real, dimension(bd%isd:bd%ied ,bd%jsd:bd%jed ,npz), intent(inout)  delp,
real, dimension(bd%is-1:bd%ie+1, npz+1,bd%js-1:bd%je+1), intent(inout)  pe,
real, dimension(bd%is:bd%ie,bd%js:bd%je, npz+1), intent(inout)  pk,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed), intent(inout)  phis,
real, dimension(bd%is:bd%ie,bd%js:bd%je), intent(out)  ws,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed,npz), intent(inout)  omga,
real, intent(in)  ptop,
real, dimension(npz), intent(in)  pfull,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed,npz), intent(inout)  ua,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed,npz), intent(inout)  va,
real, dimension(bd%isd:bd%ied+1,bd%jsd:bd%jed ,npz), intent(inout)  uc,
real, dimension(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz), intent(inout)  vc,
real, dimension(bd%is:bd%ie+1, bd%js:bd%je, npz), intent(inout)  mfx,
real, dimension(bd%is:bd%ie , bd%js:bd%je+1, npz), intent(inout)  mfy,
real, dimension(bd%is:bd%ie+1, bd%jsd:bd%jed, npz), intent(inout)  cx,
real, dimension(bd%isd:bd%ied ,bd%js:bd%je+1, npz), intent(inout)  cy,
real, dimension(bd%is:bd%ie,bd%js:bd%je,npz), intent(inout)  pkz,
real, dimension(bd%is:bd%ie,npz+1,bd%js:bd%je), intent(inout)  peln,
real, dimension(bd%isd:, bd%jsd:, 1:), intent(inout)  q_con,
real, dimension(npz+1), intent(in)  ak,
real, dimension(npz+1), intent(in)  bk,
integer, intent(in)  ks,
type(fv_grid_type), intent(inout), target  gridstruct,
type(fv_flags_type), intent(in), target  flagstruct,
type(fv_nest_type), intent(inout)  neststruct,
type(fv_diag_type), intent(in)  idiag,
type(fv_grid_bounds_type), intent(in)  bd,
type(domain2d), intent(inout)  domain,
logical, intent(in)  init_step,
type(group_halo_update_type), dimension(*), intent(inout)  i_pack,
logical, intent(in)  end_step,
real, dimension(bd%isd:bd%ied ,bd%jsd:bd%jed ,npz), intent(inout)  diss_est,
real, intent(in), optional  time_total 
)
Parameters
[in,out]uD grid zonal wind (m/s)
[in,out]vD grid meridional wind (m/s)
[in,out]wvertical vel. (m/s)
[in,out]delzdelta-height (m, negative)
[in,out]cappamoist kappa
[in,out]ptpotential temperature (K)
[in,out]delppressure thickness (pascal)
[in]time_totaltotal time (seconds) since start
[in,out]diss_estskeb dissipation estimate
[in,out]phisSurface geopotential (g*Z_surf)
[in,out]peedge pressure (pascal)
[in,out]pelnln(pe)
[in,out]pkpe**kappa
[out]wsw at surface
[in,out]omgaVertical pressure velocity (pa/s)
[in,out]uc(uc, vc) are mostly used as the C grid winds

Definition at line 180 of file dyn_core.F90.

◆ geopk()

subroutine dyn_core_mod::geopk ( real, intent(in)  ptop,
real, dimension(bd%is-1:bd%ie+1,km+1,bd%js-1:bd%je+1), intent(out)  pe,
real, dimension(bd%is:bd%ie,km+1,bd%js:bd%je), intent(out)  peln,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed,km), intent(in)  delp,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed,km+1), intent(out)  pk,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed,km+1), intent(out)  gz,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed), intent(in)  hs,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed,km), intent(in)  pt,
real, dimension(bd%isd:,bd%jsd:,1:), intent(in)  q_con,
real, dimension(bd%is:bd%ie,bd%js:bd%je,km), intent(out)  pkz,
integer, intent(in)  km,
real, intent(in)  akap,
logical, intent(in)  CG,
logical, intent(in)  nested,
logical, intent(in)  computehalo,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  a2b_ord,
type(fv_grid_bounds_type), intent(in)  bd 
)
private

The subroutine 'geopk' calculates geopotential and pressure to the kappa.

Parameters
[out]pelnln(pe)

Definition at line 2247 of file dyn_core.F90.

◆ grad1_p_update()

subroutine dyn_core_mod::grad1_p_update ( real, dimension(bd%is:bd%ie+1,bd%js:bd%je+1), intent(in)  divg2,
real, dimension(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz), intent(inout)  u,
real, dimension(bd%isd:bd%ied+1,bd%jsd:bd%jed ,npz), intent(inout)  v,
real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed ,npz+1), intent(inout)  pk,
real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed ,npz+1), intent(inout)  gz,
real, intent(in)  dt,
integer, intent(in)  ng,
type(fv_grid_type), intent(inout), target  gridstruct,
type(fv_grid_bounds_type), intent(in)  bd,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
real, intent(in)  ptop,
real, intent(in)  beta,
integer, intent(in)  a2b_ord 
)
private

Definition at line 2070 of file dyn_core.F90.

◆ init_ijk_mem()

subroutine, public dyn_core_mod::init_ijk_mem ( integer, intent(in)  i1,
integer, intent(in)  i2,
integer, intent(in)  j1,
integer, intent(in)  j2,
integer, intent(in)  km,
real, dimension(i1:i2,j1:j2,km), intent(inout)  array,
real, intent(in)  var 
)

Definition at line 2532 of file dyn_core.F90.

◆ mix_dp()

subroutine dyn_core_mod::mix_dp ( logical, intent(in)  hydrostatic,
real, dimension(bd%isd:,bd%jsd:,1:), intent(inout)  w,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed,km), intent(inout)  delp,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed,km), intent(inout)  pt,
integer, intent(in)  km,
real, dimension(km+1), intent(in)  ak,
real, dimension(km+1), intent(in)  bk,
logical, intent(in)  CG,
logical, intent(in)  fv_debug,
type(fv_grid_bounds_type), intent(in)  bd 
)
private

Definition at line 2164 of file dyn_core.F90.

◆ nh_p_grad()

subroutine dyn_core_mod::nh_p_grad ( real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed+1,npz), intent(inout)  u,
real, dimension(bd%isd:bd%ied+1,bd%jsd:bd%jed, npz), intent(inout)  v,
real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed, npz+1), intent(inout)  pp,
real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed, npz+1), intent(inout)  gz,
real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed, npz), intent(inout)  delp,
real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed, npz+1), intent(inout)  pk,
real, intent(in)  dt,
integer, intent(in)  ng,
type(fv_grid_type), intent(inout), target  gridstruct,
type(fv_grid_bounds_type), intent(in)  bd,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
logical, intent(in)  use_logp 
)
private
Parameters
[in,out]ppperturbation pressure
[in,out]pkp**kappa
[in,out]gzg * h

Definition at line 1722 of file dyn_core.F90.

◆ one_grad_p()

subroutine dyn_core_mod::one_grad_p ( real, dimension(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz), intent(inout)  u,
real, dimension(bd%isd:bd%ied+1,bd%jsd:bd%jed ,npz), intent(inout)  v,
real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed ,npz+1), intent(inout)  pk,
real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed ,npz+1), intent(inout)  gz,
real, dimension(bd%is:bd%ie+1,bd%js:bd%je+1), intent(in)  divg2,
real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed ,npz), intent(inout)  delp,
real, intent(in)  dt,
integer, intent(in)  ng,
type(fv_grid_type), intent(inout), target  gridstruct,
type(fv_grid_bounds_type), intent(in)  bd,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
real, intent(in)  ptop,
logical, intent(in)  hydrostatic,
integer, intent(in)  a2b_ord,
real, intent(in)  d_ext 
)
private

Definition at line 1935 of file dyn_core.F90.

◆ p_grad_c()

subroutine dyn_core_mod::p_grad_c ( real, intent(in)  dt2,
integer, intent(in)  npz,
real, dimension(bd%isd:, bd%jsd: ,: ), intent(in)  delpc,
real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed ,npz+1), intent(in)  pkc,
real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed ,npz+1), intent(in)  gz,
real, dimension(bd%isd:bd%ied+1,bd%jsd:bd%jed ,npz), intent(inout)  uc,
real, dimension(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz), intent(inout)  vc,
type(fv_grid_bounds_type), intent(in)  bd,
real, dimension(bd%isd:bd%ied+1,bd%jsd:bd%jed+1), intent(in)  rdxc,
real, dimension(bd%isd:bd%ied ,bd%jsd:bd%jed), intent(in)  rdyc,
logical, intent(in)  hydrostatic 
)
private
Parameters
[in]pkcpkc is pe**cappa if hydrostatic pkc is full pressure if non-hydrostatic
[in]gzpkc is pe**cappa if hydrostatic pkc is full pressure if non-hydrostatic

Definition at line 1660 of file dyn_core.F90.

◆ pe_halo()

subroutine dyn_core_mod::pe_halo ( integer, intent(in)  is,
integer, intent(in)  ie,
integer, intent(in)  js,
integer, intent(in)  je,
integer, intent(in)  isd,
integer, intent(in)  ied,
integer, intent(in)  jsd,
integer, intent(in)  jed,
integer, intent(in)  npz,
real, intent(in)  ptop,
real, dimension(is-1:ie+1,npz+1,js-1:je+1), intent(inout)  pe,
real, dimension(isd:ied,jsd:jed,npz), intent(in)  delp 
)
private

Definition at line 1523 of file dyn_core.F90.

◆ pk3_halo()

subroutine dyn_core_mod::pk3_halo ( integer, intent(in)  is,
integer, intent(in)  ie,
integer, intent(in)  js,
integer, intent(in)  je,
integer, intent(in)  isd,
integer, intent(in)  ied,
integer, intent(in)  jsd,
integer, intent(in)  jed,
integer, intent(in)  npz,
real, intent(in)  ptop,
real, intent(in)  akap,
real, dimension(isd:ied,jsd:jed,npz+1), intent(inout)  pk3,
real, dimension(isd:ied,jsd:jed,npz), intent(in)  delp 
)
private

Definition at line 1420 of file dyn_core.F90.

◆ pln_halo()

subroutine dyn_core_mod::pln_halo ( integer, intent(in)  is,
integer, intent(in)  ie,
integer, intent(in)  js,
integer, intent(in)  je,
integer, intent(in)  isd,
integer, intent(in)  ied,
integer, intent(in)  jsd,
integer, intent(in)  jed,
integer, intent(in)  npz,
real, intent(in)  ptop,
real, dimension(isd:ied,jsd:jed,npz+1), intent(inout)  pk3,
real, dimension(isd:ied,jsd:jed,npz), intent(in)  delp 
)
private

Definition at line 1474 of file dyn_core.F90.

◆ ray_fast()

subroutine dyn_core_mod::ray_fast ( real, intent(in)  dt,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
real, dimension(npz), intent(in)  pfull,
real, intent(in)  tau,
real, dimension(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz), intent(inout)  u,
real, dimension(bd%isd:bd%ied+1,bd%jsd:bd%jed,npz), intent(inout)  v,
real, dimension(bd%isd: ,bd%jsd: ,1: ), intent(inout)  w,
integer, intent(in)  ks,
real, dimension(npz), intent(in)  dp,
real, intent(in)  ptop,
logical, intent(in)  hydrostatic,
real, intent(in)  rf_cutoff,
type(fv_grid_bounds_type), intent(in)  bd 
)
private

The subroutine 'Ray_fast' computes a simple "inline" version of the Rayleigh friction (EXPERIMENTAL - NOT FOR GENERAL USE).

Parameters
[in]tautime scale (days)
[in,out]uD grid zonal wind (m/s)
[in,out]vD grid meridional wind (m/s)
[in,out]wcell center vertical wind (m/s)

Definition at line 2551 of file dyn_core.F90.

◆ split_p_grad()

subroutine dyn_core_mod::split_p_grad ( real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed+1,npz), intent(inout)  u,
real, dimension(bd%isd:bd%ied+1,bd%jsd:bd%jed, npz), intent(inout)  v,
real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed, npz+1), intent(inout)  pp,
real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed, npz+1), intent(inout)  gz,
real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed, npz), intent(inout)  delp,
real, dimension(bd%isd:bd%ied, bd%jsd:bd%jed, npz+1), intent(inout)  pk,
real, intent(in)  beta,
real, intent(in)  dt,
integer, intent(in)  ng,
type(fv_grid_type), intent(inout), target  gridstruct,
type(fv_grid_bounds_type), intent(in)  bd,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
logical, intent(in)  use_logp 
)
private
Parameters
[in,out]ppperturbation pressure
[in,out]pkp**kappa
[in,out]gzg * h

Definition at line 1820 of file dyn_core.F90.

Variable Documentation

◆ cnst_0p20

real(kind=r_grid), parameter dyn_core_mod::cnst_0p20 =0.20d0
private

Definition at line 157 of file dyn_core.F90.

◆ crx

real, dimension(:,:,:), allocatable dyn_core_mod::crx
private

Definition at line 153 of file dyn_core.F90.

◆ cry

real, dimension(:,:,:), allocatable dyn_core_mod::cry
private

Definition at line 153 of file dyn_core.F90.

◆ d3_damp

real dyn_core_mod::d3_damp
private

Definition at line 152 of file dyn_core.F90.

◆ delpc

real, dimension(:,:,:), allocatable dyn_core_mod::delpc
private

Definition at line 153 of file dyn_core.F90.

◆ divgd

real, dimension(:,:,:), allocatable dyn_core_mod::divgd
private

Definition at line 153 of file dyn_core.F90.

◆ du

real, dimension(:,:,:), allocatable dyn_core_mod::du
private

Definition at line 153 of file dyn_core.F90.

◆ dv

real, dimension(:,:,:), allocatable dyn_core_mod::dv
private

Definition at line 153 of file dyn_core.F90.

◆ first_call

logical dyn_core_mod::first_call = .true.
private

Definition at line 162 of file dyn_core.F90.

◆ gz

real, dimension(:,:,:), allocatable dyn_core_mod::gz
private

Definition at line 153 of file dyn_core.F90.

◆ k_rf

integer dyn_core_mod::k_rf = 0
private

Definition at line 160 of file dyn_core.F90.

◆ kmax

integer dyn_core_mod::kmax =1
private

Definition at line 163 of file dyn_core.F90.

◆ peln1

real dyn_core_mod::peln1
private

Definition at line 151 of file dyn_core.F90.

◆ pk3

real, dimension(:,:,:), allocatable dyn_core_mod::pk3
private

Definition at line 153 of file dyn_core.F90.

◆ pkc

real, dimension(:,:,:), allocatable dyn_core_mod::pkc
private

Definition at line 153 of file dyn_core.F90.

◆ ptc

real, dimension(:,:,:), allocatable dyn_core_mod::ptc
private

Definition at line 153 of file dyn_core.F90.

◆ ptk

real dyn_core_mod::ptk
private

Definition at line 151 of file dyn_core.F90.

◆ rf

real, dimension(:), allocatable dyn_core_mod::rf
private

Definition at line 159 of file dyn_core.F90.

◆ rff_initialized

logical dyn_core_mod::rff_initialized = .false.
private

Definition at line 161 of file dyn_core.F90.

◆ rgrav

real dyn_core_mod::rgrav
private

Definition at line 151 of file dyn_core.F90.

◆ ut

real, dimension(:,:,:), allocatable dyn_core_mod::ut
private

Definition at line 153 of file dyn_core.F90.

◆ vt

real, dimension(:,:,:), allocatable dyn_core_mod::vt
private

Definition at line 153 of file dyn_core.F90.

◆ xfx

real, dimension(:,:,:), allocatable dyn_core_mod::xfx
private

Definition at line 153 of file dyn_core.F90.

◆ yfx

real, dimension(:,:,:), allocatable dyn_core_mod::yfx
private

Definition at line 153 of file dyn_core.F90.

◆ zh

real, dimension(:,:,:), allocatable dyn_core_mod::zh
private

Definition at line 153 of file dyn_core.F90.