FV3DYCORE  Version 2.0.0
fv_nesting_mod Module Reference

The module 'fv_nesting' is a collection of routines pertaining to grid nesting [harris2013two]. More...

Functions/Subroutines

subroutine, public setup_nested_grid_bcs (npx, npy, npz, zvir, ncnst, u, v, w, pt, delp, delz, q, uc, vc, nested, inline_q, make_nh, ng, gridstruct, flagstruct, neststruct, nest_timestep, tracer_nest_timestep, domain, parent_grid, bd, nwat, ak, bk)
 The subroutine 'setup_nested_grid_BCs' fetches data from the coarse grid to set up the nested-grid boundary conditions. More...
 
subroutine, public set_physics_bcs (ps, u_dt, v_dt, flagstruct, gridstruct, neststruct, npx, npy, npz, ng, ak, bk, bd)
 
subroutine set_bc_direct (pe_src_BC, pe_dst_BC, buf, var, neststruct, npx, npy, npz, npz_coarse, ng, bd, istag, jstag, iv, kord)
 
subroutine setup_pt_bc (pt_BC, pe_eul_BC, sphum_BC, npx, npy, npz, zvir, bd)
 
subroutine setup_pt_bc_k (ptBC, sphumBC, peBC, zvir, isd_BC, ied_BC, istart, iend, jstart, jend, npz)
 
subroutine setup_eul_delp_bc (delp_lag_BC, delp_eul_BC, pe_lag_BC, pe_eul_BC, ak_dst, bk_dst, npx, npy, npz, npz_coarse, ptop_src, bd)
 
subroutine setup_eul_delp_bc_k (delplagBC, delpeulBC, pelagBC, peeulBC, ptop_src, ak_dst, bk_dst, isd_BC, ied_BC, istart, iend, jstart, jend, npz, npz_coarse)
 
subroutine copy_ps_bc (ps, pe_BC, npx, npy, npz, istag, jstag, bd)
 
subroutine setup_eul_pe_bc (pe_src_BC, pe_eul_BC, ak_dst, bk_dst, npx, npy, npz, npz_src, istag, jstag, bd, make_src_in, ak_src, bk_src)
 
subroutine setup_eul_pe_bc_k (pesrcBC, peeulBC, ak_dst, bk_dst, isd_BC, ied_BC, istart, iend, jstart, jend, npz, npz_src, make_src, ak_src, bk_src)
 
subroutine remap_bc (pe_lag_BC, pe_eul_BC, var_lag_BC, var_eul_BC, npx, npy, npz, npz_coarse, bd, istag, jstag, iv, kord, varname, do_log_pe)
 
subroutine remap_bc_direct (pe_lag_BC, pe_eul_BC, var_lag_BC, var, npx, npy, npz, npz_coarse, bd, istag, jstag, iv, kord, do_log_pe)
 
subroutine remap_bc_k (pe_lagBC, pe_eulBC, var_lagBC, var_eulBC, isd_BC, ied_BC, istart, iend, jstart, jend, npz, npz_coarse, iv, kord, log_pe)
 
subroutine remap_delz_bc (pe_lag_BC, pe_eul_BC, delp_lag_BC, delz_lag_BC, delp_eul_BC, delz_eul_BC, npx, npy, npz, npz_coarse, bd, istag, jstag, iv, kord)
 
subroutine compute_specific_volume_bc_k (delpBC, delzBC, isd_BC, ied_BC, istart, iend, jstart, jend, npz)
 
subroutine compute_delz_bc_k (delpBC, delzBC, isd_BC, ied_BC, istart, iend, jstart, jend, npz)
 
subroutine setup_pt_nh_bc (pt_BC, delp_BC, delz_BC, sphum_BC, q_BC, nq, npx, npy, npz, zvir, bd)
 
subroutine setup_pt_nh_bc_k (ptBC, sphumBC, delpBC, delzBC, liq_watBC, rainwatBC, ice_watBC, snowwatBC, graupelBC, zvir, isd_BC, ied_BC, istart, iend, jstart, jend, npz)
 
subroutine set_nh_bcs_t0 (neststruct)
 
subroutine set_bcs_t0 (ncnst, hydrostatic, neststruct)
 
subroutine d2c_setup (u, v, ua, va, uc, vc, dord4, isd, ied, jsd, jed, is, ie, js, je, npx, npy, grid_type, bounded_domain, se_corner, sw_corner, ne_corner, nw_corner, rsin_u, rsin_v, cosa_s, rsin2)
 
subroutine d2a_setup (u, v, ua, va, dord4, isd, ied, jsd, jed, is, ie, js, je, npx, npy, grid_type, bounded_domain, cosa_s, rsin2)
 
subroutine, public twoway_nesting (Atm, ngrids, grids_on_this_pe, zvir, Time, this_grid)
 The subroutine'twoway_nesting' performs a two-way update of nested-grid data onto the parent grid. More...
 
subroutine twoway_nest_update (npx, npy, npz, zvir, ncnst, sphum, u, v, w, pt, delp, q, pe, pkz, delz, ps, ptop, ak, bk, gridstruct, flagstruct, neststruct, domain, parent_grid, bd, grid_number, conv_theta_in)
 
subroutine level_sum (q, area, domain, bd, npz, L_sum)
 
subroutine remap_up_k (ps_src, ps_dst, ak_src, bk_src, ak_dst, bk_dst, var_src, var_dst, bd, istart, iend, jstart, jend, istag, jstag, npz_src, npz_dst, iv, kord, blend_wt, log_pe)
 
subroutine after_twoway_nest_update (npx, npy, npz, ng, ncnst, u, v, w, delz, pt, delp, q, ps, pe, pk, peln, pkz, phis, ua, va, ptop, gridstruct, flagstruct, domain, bd, Time)
 
subroutine update_remap_tqw (npz, ak_dst, bk_dst, ps_dst, t_dst, q_dst, w_dst, hydrostatic, kmd, ps_src, ak_src, bk_src, t_src, w_src, zvir, ptop, nq, kord_tm, kord_tr, kord_wz, is, ie, js, je, isd, ied, jsd, jed, do_q, istart, iend, jstart, jend, blend_wt)
 The subroutine 'update_remap_tqw' remaps (interpolated) nested-grid data to the coarse-grid's vertical coordinate. More...
 
subroutine update_remap_uv (npz, ak_dst, bk_dst, ps_dst, u_dst, v_dst, kmd, ak_src, bk_src, ps_src, u_src, v_src, kord_mt, is, ie, js, je, isd, ied, jsd, jed, ptop, istart, iend, jstart, jend, blend_wt)
 

Variables

logical rf_initialized = .false.
 
logical bad_range
 
real, dimension(:), allocatable rf
 
real, dimension(:), allocatable rw
 
integer kmax =1
 
real, dimension(:,:), allocatable te_2d_coarse
 
real, dimension(:,:,:), allocatable dp1_coarse
 
type(fv_nest_bc_type_3d) u_buf
 
type(fv_nest_bc_type_3d) v_buf
 
type(fv_nest_bc_type_3d) uc_buf
 
type(fv_nest_bc_type_3d) vc_buf
 
type(fv_nest_bc_type_3d) delp_buf
 
type(fv_nest_bc_type_3d) delz_buf
 
type(fv_nest_bc_type_3d) pt_buf
 
type(fv_nest_bc_type_3d) w_buf
 
type(fv_nest_bc_type_3d) divg_buf
 
type(fv_nest_bc_type_3d) pe_u_buf
 
type(fv_nest_bc_type_3d) pe_v_buf
 
type(fv_nest_bc_type_3d) pe_b_buf
 
type(fv_nest_bc_type_3d), dimension(:), allocatable q_buf
 
real, dimension(:,:,:), allocatable, target dum_west
 
real, dimension(:,:,:), allocatable, target dum_east
 
real, dimension(:,:,:), allocatable, target dum_north
 
real, dimension(:,:,:), allocatable, target dum_south
 

Detailed Description

The module 'fv_nesting' is a collection of routines pertaining to grid nesting [harris2013two].

Function/Subroutine Documentation

◆ after_twoway_nest_update()

subroutine fv_nesting_mod::after_twoway_nest_update ( integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
integer, intent(in)  ng,
integer, intent(in)  ncnst,
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%is: ,bd%js: ,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), intent(inout)  delp,
real, dimension( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz, ncnst), intent(inout)  q,
real, dimension (bd%isd:bd%ied ,bd%jsd:bd%jed), intent(inout)  ps,
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%is:bd%ie,npz+1,bd%js:bd%je), intent(inout)  peln,
real, dimension (bd%is:bd%ie,bd%js:bd%je,npz), intent(inout)  pkz,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed), intent(inout)  phis,
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, intent(in)  ptop,
type(fv_grid_type), intent(in)  gridstruct,
type(fv_flags_type), intent(in)  flagstruct,
type(domain2d), intent(inout)  domain,
type(fv_grid_bounds_type), intent(in)  bd,
type(time_type), intent(in)  Time 
)
private
Parameters
[in,out]uD grid zonal wind (m/s)
[in,out]vD grid meridional wind (m/s)
[in,out]wW (m/s)
[in,out]pttemperature (K)
[in,out]delppressure thickness (pascal)
[in,out]qspecific humidity and constituents
[in,out]delzdelta-height (m); non-hydrostatic only
[in,out]psSurface pressure (pascal)
[in,out]peedge pressure (pascal)
[in,out]pkpe**cappa
[in,out]pelnln(pe)
[in,out]pkzfinite-volume mean pk
[in,out]phisSurface geopotential (g*Z_surf)

Definition at line 2892 of file fv_nesting.F90.

◆ compute_delz_bc_k()

subroutine fv_nesting_mod::compute_delz_bc_k ( real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(in)  delpBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(inout)  delzBC,
integer, intent(in)  isd_BC,
integer, intent(in)  ied_BC,
integer, intent(in)  istart,
integer, intent(in)  iend,
integer, intent(in)  jstart,
integer, intent(in)  jend,
integer, intent(in)  npz 
)
private

Definition at line 1433 of file fv_nesting.F90.

◆ compute_specific_volume_bc_k()

subroutine fv_nesting_mod::compute_specific_volume_bc_k ( real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(in)  delpBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(inout)  delzBC,
integer, intent(in)  isd_BC,
integer, intent(in)  ied_BC,
integer, intent(in)  istart,
integer, intent(in)  iend,
integer, intent(in)  jstart,
integer, intent(in)  jend,
integer, intent(in)  npz 
)
private

Definition at line 1407 of file fv_nesting.F90.

◆ copy_ps_bc()

subroutine fv_nesting_mod::copy_ps_bc ( real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag), intent(in)  ps,
type(fv_nest_bc_type_3d), intent(inout)  pe_BC,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
integer, intent(in)  istag,
integer, intent(in)  jstag,
type(fv_grid_bounds_type), intent(in)  bd 
)
private

Definition at line 1000 of file fv_nesting.F90.

◆ d2a_setup()

subroutine fv_nesting_mod::d2a_setup ( real, dimension(isd:ied,jsd:jed+1), intent(in)  u,
real, dimension(isd:ied+1,jsd:jed), intent(in)  v,
real, dimension(isd:ied ,jsd:jed ), intent(out)  ua,
real, dimension(isd:ied ,jsd:jed ), intent(out)  va,
logical, intent(in)  dord4,
integer, intent(in)  isd,
integer, intent(in)  ied,
integer, intent(in)  jsd,
integer, intent(in)  jed,
integer, intent(in)  is,
integer, intent(in)  ie,
integer, intent(in)  js,
integer, intent(in)  je,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  grid_type,
logical, intent(in)  bounded_domain,
real, dimension(isd:ied,jsd:jed), intent(in)  cosa_s,
real, dimension(isd:ied,jsd:jed), intent(in)  rsin2 
)
private

Definition at line 2153 of file fv_nesting.F90.

◆ d2c_setup()

subroutine fv_nesting_mod::d2c_setup ( real, dimension(isd:ied,jsd:jed+1), intent(in)  u,
real, dimension(isd:ied+1,jsd:jed), intent(in)  v,
real, dimension(isd:ied ,jsd:jed ), intent(out)  ua,
real, dimension(isd:ied ,jsd:jed ), intent(out)  va,
real, dimension(isd:ied+1,jsd:jed ), intent(out)  uc,
real, dimension(isd:ied ,jsd:jed+1), intent(out)  vc,
logical, intent(in)  dord4,
integer, intent(in)  isd,
integer, intent(in)  ied,
integer, intent(in)  jsd,
integer, intent(in)  jed,
integer, intent(in)  is,
integer, intent(in)  ie,
integer, intent(in)  js,
integer, intent(in)  je,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  grid_type,
logical, intent(in)  bounded_domain,
logical, intent(in)  se_corner,
logical, intent(in)  sw_corner,
logical, intent(in)  ne_corner,
logical, intent(in)  nw_corner,
real, dimension(isd:ied+1,jsd:jed), intent(in)  rsin_u,
real, dimension(isd:ied,jsd:jed+1), intent(in)  rsin_v,
real, dimension(isd:ied,jsd:jed), intent(in)  cosa_s,
real, dimension(isd:ied,jsd:jed), intent(in)  rsin2 
)
private

Definition at line 1878 of file fv_nesting.F90.

◆ level_sum()

subroutine fv_nesting_mod::level_sum ( real, dimension( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz), intent(in)  q,
real, dimension( bd%isd:bd%ied ,bd%jsd:bd%jed), intent(in)  area,
type(domain2d), intent(in)  domain,
type(fv_grid_bounds_type), intent(in)  bd,
integer, intent(in)  npz,
real, dimension( npz ), intent(out)  L_sum 
)
private

Definition at line 2726 of file fv_nesting.F90.

◆ remap_bc()

subroutine fv_nesting_mod::remap_bc ( type(fv_nest_bc_type_3d), intent(inout), target  pe_lag_BC,
type(fv_nest_bc_type_3d), intent(inout), target  pe_eul_BC,
type(fv_nest_bc_type_3d), intent(inout), target  var_lag_BC,
type(fv_nest_bc_type_3d), intent(inout), target  var_eul_BC,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
integer, intent(in)  npz_coarse,
type(fv_grid_bounds_type), intent(in)  bd,
integer, intent(in)  istag,
integer, intent(in)  jstag,
integer, intent(in)  iv,
integer, intent(in)  kord,
character(len=*), intent(in)  varname,
logical, intent(in), optional  do_log_pe 
)
private

Definition at line 1168 of file fv_nesting.F90.

◆ remap_bc_direct()

subroutine fv_nesting_mod::remap_bc_direct ( type(fv_nest_bc_type_3d), intent(inout), target  pe_lag_BC,
type(fv_nest_bc_type_3d), intent(inout), target  pe_eul_BC,
type(fv_nest_bc_type_3d), intent(inout), target  var_lag_BC,
real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,npz), intent(inout)  var,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
integer, intent(in)  npz_coarse,
type(fv_grid_bounds_type), intent(in)  bd,
integer, intent(in)  istag,
integer, intent(in)  jstag,
integer, intent(in)  iv,
integer, intent(in)  kord,
logical, intent(in), optional  do_log_pe 
)
private

Definition at line 1224 of file fv_nesting.F90.

◆ remap_bc_k()

subroutine fv_nesting_mod::remap_bc_k ( real, dimension(isd_bc:ied_bc,jstart:jend,npz_coarse+1), intent(inout)  pe_lagBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz+1), intent(inout)  pe_eulBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz_coarse), intent(inout)  var_lagBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(inout)  var_eulBC,
integer, intent(in)  isd_BC,
integer, intent(in)  ied_BC,
integer, intent(in)  istart,
integer, intent(in)  iend,
integer, intent(in)  jstart,
integer, intent(in)  jend,
integer, intent(in)  npz,
integer, intent(in)  npz_coarse,
integer, intent(in)  iv,
integer, intent(in)  kord,
logical, intent(in)  log_pe 
)
private

Definition at line 1282 of file fv_nesting.F90.

◆ remap_delz_bc()

subroutine fv_nesting_mod::remap_delz_bc ( type(fv_nest_bc_type_3d), intent(inout), target  pe_lag_BC,
type(fv_nest_bc_type_3d), intent(inout), target  pe_eul_BC,
type(fv_nest_bc_type_3d), intent(inout), target  delp_lag_BC,
type(fv_nest_bc_type_3d), intent(inout), target  delz_lag_BC,
type(fv_nest_bc_type_3d), intent(inout), target  delp_eul_BC,
type(fv_nest_bc_type_3d), intent(inout), target  delz_eul_BC,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
integer, intent(in)  npz_coarse,
type(fv_grid_bounds_type), intent(in)  bd,
integer, intent(in)  istag,
integer, intent(in)  jstag,
integer, intent(in)  iv,
integer, intent(in)  kord 
)
private

Definition at line 1345 of file fv_nesting.F90.

◆ remap_up_k()

subroutine fv_nesting_mod::remap_up_k ( real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed), intent(inout)  ps_src,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed), intent(inout)  ps_dst,
real, dimension(npz_src+1), intent(in)  ak_src,
real, dimension(npz_src+1), intent(in)  bk_src,
real, dimension(npz_dst+1), intent(in)  ak_dst,
real, dimension(npz_dst+1), intent(in)  bk_dst,
real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,npz_src), intent(inout)  var_src,
real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,npz_dst), intent(inout)  var_dst,
type(fv_grid_bounds_type), intent(in)  bd,
integer, intent(in)  istart,
integer, intent(in)  iend,
integer, intent(in)  jstart,
integer, intent(in)  jend,
integer, intent(in)  istag,
integer, intent(in)  jstag,
integer, intent(in)  npz_src,
integer, intent(in)  npz_dst,
integer, intent(in)  iv,
integer, intent(in)  kord,
real, dimension(npz_dst), intent(in)  blend_wt,
logical, intent(in)  log_pe 
)
private

Definition at line 2758 of file fv_nesting.F90.

◆ set_bc_direct()

subroutine fv_nesting_mod::set_bc_direct ( type(fv_nest_bc_type_3d), intent(inout)  pe_src_BC,
type(fv_nest_bc_type_3d), intent(inout)  pe_dst_BC,
type(fv_nest_bc_type_3d), intent(inout)  buf,
real, dimension(bd%isd:bd%ied+istag,bd%jsd:bd%jed+jstag,npz), intent(inout)  var,
type(fv_nest_type), intent(inout)  neststruct,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
integer, intent(in)  npz_coarse,
integer, intent(in)  ng,
type(fv_grid_bounds_type), intent(in)  bd,
integer, intent(in)  istag,
integer, intent(in)  jstag,
integer, intent(in)  iv,
integer, intent(in)  kord 
)
private

Definition at line 761 of file fv_nesting.F90.

◆ set_bcs_t0()

subroutine fv_nesting_mod::set_bcs_t0 ( integer, intent(in)  ncnst,
logical, intent(in)  hydrostatic,
type(fv_nest_type), intent(inout)  neststruct 
)
private

Definition at line 1801 of file fv_nesting.F90.

◆ set_nh_bcs_t0()

subroutine fv_nesting_mod::set_nh_bcs_t0 ( type(fv_nest_type), intent(inout)  neststruct)
private

Definition at line 1783 of file fv_nesting.F90.

◆ set_physics_bcs()

subroutine, public fv_nesting_mod::set_physics_bcs ( real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed), intent(inout)  ps,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed,npz), intent(inout)  u_dt,
real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed,npz), intent(inout)  v_dt,
type(fv_flags_type), intent(in)  flagstruct,
type(fv_grid_type gridstruct,
type(fv_nest_type), intent(inout), target  neststruct,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
integer, intent(in)  ng,
real, dimension(npz+1), intent(in)  ak,
real, dimension(npz+1), intent(in)  bk,
type(fv_grid_bounds_type), intent(in)  bd 
)

Definition at line 687 of file fv_nesting.F90.

◆ setup_eul_delp_bc()

subroutine fv_nesting_mod::setup_eul_delp_bc ( type(fv_nest_bc_type_3d), intent(inout), target  delp_lag_BC,
type(fv_nest_bc_type_3d), intent(inout), target  delp_eul_BC,
type(fv_nest_bc_type_3d), intent(inout), target  pe_lag_BC,
type(fv_nest_bc_type_3d), intent(inout), target  pe_eul_BC,
real, dimension(npz+1), intent(in)  ak_dst,
real, dimension(npz+1), intent(in)  bk_dst,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
integer, intent(in)  npz_coarse,
real, intent(in)  ptop_src,
type(fv_grid_bounds_type), intent(in)  bd 
)
private

Definition at line 882 of file fv_nesting.F90.

◆ setup_eul_delp_bc_k()

subroutine fv_nesting_mod::setup_eul_delp_bc_k ( real, dimension(isd_bc:ied_bc,jstart:jend,npz_coarse), intent(inout)  delplagBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(inout)  delpeulBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz_coarse+1), intent(inout)  pelagBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz+1), intent(inout)  peeulBC,
real, intent(in)  ptop_src,
real, dimension(npz+1), intent(in)  ak_dst,
real, dimension(npz+1), intent(in)  bk_dst,
integer, intent(in)  isd_BC,
integer, intent(in)  ied_BC,
integer, intent(in)  istart,
integer, intent(in)  iend,
integer, intent(in)  jstart,
integer, intent(in)  jend,
integer, intent(in)  npz,
integer, intent(in)  npz_coarse 
)
private

Definition at line 938 of file fv_nesting.F90.

◆ setup_eul_pe_bc()

subroutine fv_nesting_mod::setup_eul_pe_bc ( type(fv_nest_bc_type_3d), intent(inout), target  pe_src_BC,
type(fv_nest_bc_type_3d), intent(inout), target  pe_eul_BC,
real, dimension(npz+1), intent(in)  ak_dst,
real, dimension(npz+1), intent(in)  bk_dst,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
integer, intent(in)  npz_src,
integer, intent(in)  istag,
integer, intent(in)  jstag,
type(fv_grid_bounds_type), intent(in)  bd,
logical, intent(in), optional  make_src_in,
real, dimension(npz_src), intent(in), optional  ak_src,
real, dimension(npz_src), intent(in), optional  bk_src 
)
private

Definition at line 1071 of file fv_nesting.F90.

◆ setup_eul_pe_bc_k()

subroutine fv_nesting_mod::setup_eul_pe_bc_k ( real, dimension(isd_bc:ied_bc,jstart:jend,npz_src+1), intent(inout)  pesrcBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz+1), intent(inout)  peeulBC,
real, dimension(npz+1), intent(in)  ak_dst,
real, dimension(npz+1), intent(in)  bk_dst,
integer, intent(in)  isd_BC,
integer, intent(in)  ied_BC,
integer, intent(in)  istart,
integer, intent(in)  iend,
integer, intent(in)  jstart,
integer, intent(in)  jend,
integer, intent(in)  npz,
integer, intent(in)  npz_src,
logical, intent(in)  make_src,
real, dimension(npz_src+1), intent(in)  ak_src,
real, dimension(npz_src+1), intent(in)  bk_src 
)
private

Definition at line 1132 of file fv_nesting.F90.

◆ setup_nested_grid_bcs()

subroutine, public fv_nesting_mod::setup_nested_grid_bcs ( integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
real, intent(in)  zvir,
integer, intent(in)  ncnst,
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%ied ,bd%jsd:bd%jed ,npz), intent(inout)  pt,
real, dimension(bd%isd:bd%ied ,bd%jsd:bd%jed ,npz), intent(inout)  delp,
real, dimension(bd%is: ,bd%js: ,1:), intent(inout)  delz,
real, dimension( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz, ncnst), intent(inout)  q,
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,
logical, intent(in)  nested,
logical, intent(in)  inline_q,
logical, intent(in)  make_nh,
integer, intent(in)  ng,
type(fv_grid_type), intent(inout)  gridstruct,
type(fv_flags_type), intent(inout)  flagstruct,
type(fv_nest_type), intent(inout), target  neststruct,
integer, intent(inout)  nest_timestep,
integer, intent(inout)  tracer_nest_timestep,
type(domain2d), intent(inout)  domain,
type(fv_atmos_type), intent(in), pointer  parent_grid,
type(fv_grid_bounds_type), intent(in)  bd,
integer, intent(in)  nwat,
real, dimension(npz), intent(in)  ak,
real, dimension(npz), intent(in)  bk 
)

The subroutine 'setup_nested_grid_BCs' fetches data from the coarse grid to set up the nested-grid boundary conditions.

Parameters
[in,out]uD grid zonal wind (m/s)
[in,out]vD grid meridional wind (m/s)
[in,out]wW (m/s)
[in,out]pttemperature (K)
[in,out]delppressure thickness (pascal)
[in,out]delzheight thickness (m)
[in,out]qspecific humidity and constituents
[in,out]uc(uc,vc) mostly used as the C grid winds

Definition at line 165 of file fv_nesting.F90.

◆ setup_pt_bc()

subroutine fv_nesting_mod::setup_pt_bc ( type(fv_nest_bc_type_3d), intent(inout)  pt_BC,
type(fv_nest_bc_type_3d), intent(in)  pe_eul_BC,
type(fv_nest_bc_type_3d), intent(in)  sphum_BC,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
real, intent(in)  zvir,
type(fv_grid_bounds_type), intent(in)  bd 
)
private

Definition at line 782 of file fv_nesting.F90.

◆ setup_pt_bc_k()

subroutine fv_nesting_mod::setup_pt_bc_k ( real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(inout)  ptBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(in)  sphumBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz+1), intent(in)  peBC,
real, intent(in)  zvir,
integer, intent(in)  isd_BC,
integer, intent(in)  ied_BC,
integer, intent(in)  istart,
integer, intent(in)  iend,
integer, intent(in)  jstart,
integer, intent(in)  jend,
integer, intent(in)  npz 
)
private

Definition at line 851 of file fv_nesting.F90.

◆ setup_pt_nh_bc()

subroutine fv_nesting_mod::setup_pt_nh_bc ( type(fv_nest_bc_type_3d), intent(inout), target  pt_BC,
type(fv_nest_bc_type_3d), intent(in), target  delp_BC,
type(fv_nest_bc_type_3d), intent(in), target  delz_BC,
type(fv_nest_bc_type_3d), intent(in), target  sphum_BC,
type(fv_nest_bc_type_3d), dimension(nq), intent(in), target  q_BC,
integer, intent(in)  nq,
integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
real, intent(in)  zvir,
type(fv_grid_bounds_type), intent(in)  bd 
)
private

Definition at line 1467 of file fv_nesting.F90.

◆ setup_pt_nh_bc_k()

subroutine fv_nesting_mod::setup_pt_nh_bc_k ( real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(out)  ptBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(in)  sphumBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(in)  delpBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(in)  delzBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(in)  liq_watBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(in)  rainwatBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(in)  ice_watBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(in)  snowwatBC,
real, dimension(isd_bc:ied_bc,jstart:jend,npz), intent(in)  graupelBC,
real, intent(in)  zvir,
integer, intent(in)  isd_BC,
integer, intent(in)  ied_BC,
integer, intent(in)  istart,
integer, intent(in)  iend,
integer, intent(in)  jstart,
integer, intent(in)  jend,
integer, intent(in)  npz 
)
private

Definition at line 1714 of file fv_nesting.F90.

◆ twoway_nest_update()

subroutine fv_nesting_mod::twoway_nest_update ( integer, intent(in)  npx,
integer, intent(in)  npy,
integer, intent(in)  npz,
real, intent(in)  zvir,
integer, intent(in)  ncnst,
integer, intent(in)  sphum,
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%ied ,bd%jsd:bd%jed ,npz), intent(inout)  pt,
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, ncnst), intent(inout)  q,
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), intent(inout)  pkz,
real, dimension(bd%isd: ,bd%jsd: ,1: ), intent(inout)  delz,
real, dimension (bd%isd:bd%ied ,bd%jsd:bd%jed), intent(inout)  ps,
real, intent(in)  ptop,
real, dimension(npz+1), intent(in)  ak,
real, dimension(npz+1), intent(in)  bk,
type(fv_grid_type), intent(inout)  gridstruct,
type(fv_flags_type), intent(inout)  flagstruct,
type(fv_nest_type), intent(inout)  neststruct,
type(domain2d), intent(inout)  domain,
type(fv_atmos_type), intent(in), pointer  parent_grid,
type(fv_grid_bounds_type), intent(in)  bd,
integer, intent(in)  grid_number,
logical, intent(in), optional  conv_theta_in 
)
private
Parameters
[in,out]uD grid zonal wind (m/s)
[in,out]vD grid meridional wind (m/s)
[in,out]wW (m/s)
[in,out]pttemperature (K)
[in,out]delppressure thickness (pascal)
[in,out]qspecific humidity and constituents
[in,out]pefinite-volume interface p ! NOTE TRANSPOSITION NEEDED
[in,out]pkzfinite-volume mean pk
[in,out]delzdelta-height (m); non-hydrostatic only
[in,out]psSurface pressure (pascal)

Definition at line 2379 of file fv_nesting.F90.

◆ twoway_nesting()

subroutine, public fv_nesting_mod::twoway_nesting ( type(fv_atmos_type), dimension(ngrids), intent(inout)  Atm,
integer, intent(in)  ngrids,
logical, dimension(ngrids), intent(in)  grids_on_this_pe,
real, intent(in)  zvir,
type(time_type), intent(in)  Time,
integer, intent(in)  this_grid 
)

The subroutine'twoway_nesting' performs a two-way update of nested-grid data onto the parent grid.

Definition at line 2312 of file fv_nesting.F90.

◆ update_remap_tqw()

subroutine fv_nesting_mod::update_remap_tqw ( integer, intent(in)  npz,
real, dimension(npz+1), intent(in)  ak_dst,
real, dimension(npz+1), intent(in)  bk_dst,
real, dimension(isd:ied,jsd:jed), intent(in)  ps_dst,
real, dimension(isd:ied,jsd:jed,npz), intent(inout)  t_dst,
real, dimension(isd:ied,jsd:jed,npz,nq), intent(inout)  q_dst,
real, dimension(isd:ied,jsd:jed,npz), intent(inout)  w_dst,
logical, intent(in)  hydrostatic,
integer, intent(in)  kmd,
real, dimension(isd:ied,jsd:jed), intent(in)  ps_src,
real, dimension(kmd+1), intent(in)  ak_src,
real, dimension(kmd+1), intent(in)  bk_src,
real, dimension(isd:ied,jsd:jed,kmd), intent(in)  t_src,
real, dimension(isd:ied,jsd:jed,kmd), intent(in)  w_src,
real, intent(in)  zvir,
real, intent(in)  ptop,
integer, intent(in)  nq,
integer, intent(in)  kord_tm,
integer, intent(in)  kord_tr,
integer, intent(in)  kord_wz,
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,
logical, intent(in)  do_q,
integer, intent(in)  istart,
integer, intent(in)  iend,
integer, intent(in)  jstart,
integer, intent(in)  jend,
real, dimension(npz), intent(in)  blend_wt 
)
private

The subroutine 'update_remap_tqw' remaps (interpolated) nested-grid data to the coarse-grid's vertical coordinate.

Definition at line 2990 of file fv_nesting.F90.

◆ update_remap_uv()

subroutine fv_nesting_mod::update_remap_uv ( integer, intent(in)  npz,
real, dimension(npz+1), intent(in)  ak_dst,
real, dimension(npz+1), intent(in)  bk_dst,
real, dimension(isd:ied,jsd:jed), intent(in)  ps_dst,
real, dimension(isd:ied,jsd:jed+1,npz), intent(inout)  u_dst,
real, dimension(isd:ied+1,jsd:jed,npz), intent(inout)  v_dst,
integer, intent(in)  kmd,
real, dimension(kmd+1), intent(in)  ak_src,
real, dimension(kmd+1), intent(in)  bk_src,
real, dimension(isd:ied,jsd:jed), intent(in)  ps_src,
real, dimension(isd:ied,jsd:jed+1,kmd), intent(inout)  u_src,
real, dimension(isd:ied+1,jsd:jed,kmd), intent(inout)  v_src,
integer, intent(in)  kord_mt,
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,
real, intent(in)  ptop,
integer, intent(in)  istart,
integer, intent(in)  iend,
integer, intent(in)  jstart,
integer, intent(in)  jend,
real, dimension(npz), intent(in)  blend_wt 
)
private

Definition at line 3095 of file fv_nesting.F90.

Variable Documentation

◆ bad_range

logical fv_nesting_mod::bad_range

Definition at line 132 of file fv_nesting.F90.

◆ delp_buf

type(fv_nest_bc_type_3d) fv_nesting_mod::delp_buf

Definition at line 141 of file fv_nesting.F90.

◆ delz_buf

type(fv_nest_bc_type_3d) fv_nesting_mod::delz_buf

Definition at line 141 of file fv_nesting.F90.

◆ divg_buf

type(fv_nest_bc_type_3d) fv_nesting_mod::divg_buf

Definition at line 141 of file fv_nesting.F90.

◆ dp1_coarse

real, dimension(:,:,:), allocatable fv_nesting_mod::dp1_coarse

Definition at line 137 of file fv_nesting.F90.

◆ dum_east

real, dimension(:,:,:), allocatable, target fv_nesting_mod::dum_east

Definition at line 144 of file fv_nesting.F90.

◆ dum_north

real, dimension(:,:,:), allocatable, target fv_nesting_mod::dum_north

Definition at line 144 of file fv_nesting.F90.

◆ dum_south

real, dimension(:,:,:), allocatable, target fv_nesting_mod::dum_south

Definition at line 144 of file fv_nesting.F90.

◆ dum_west

real, dimension(:,:,:), allocatable, target fv_nesting_mod::dum_west

Definition at line 144 of file fv_nesting.F90.

◆ kmax

integer fv_nesting_mod::kmax =1

Definition at line 134 of file fv_nesting.F90.

◆ pe_b_buf

type(fv_nest_bc_type_3d) fv_nesting_mod::pe_b_buf

Definition at line 141 of file fv_nesting.F90.

◆ pe_u_buf

type(fv_nest_bc_type_3d) fv_nesting_mod::pe_u_buf

Definition at line 141 of file fv_nesting.F90.

◆ pe_v_buf

type(fv_nest_bc_type_3d) fv_nesting_mod::pe_v_buf

Definition at line 141 of file fv_nesting.F90.

◆ pt_buf

type(fv_nest_bc_type_3d) fv_nesting_mod::pt_buf

Definition at line 141 of file fv_nesting.F90.

◆ q_buf

type(fv_nest_bc_type_3d), dimension(:), allocatable fv_nesting_mod::q_buf

Definition at line 142 of file fv_nesting.F90.

◆ rf

real, dimension(:), allocatable fv_nesting_mod::rf

Definition at line 133 of file fv_nesting.F90.

◆ rf_initialized

logical fv_nesting_mod::rf_initialized = .false.

Definition at line 131 of file fv_nesting.F90.

◆ rw

real, dimension(:), allocatable fv_nesting_mod::rw

Definition at line 133 of file fv_nesting.F90.

◆ te_2d_coarse

real, dimension(:,:), allocatable fv_nesting_mod::te_2d_coarse

Definition at line 136 of file fv_nesting.F90.

◆ u_buf

type(fv_nest_bc_type_3d) fv_nesting_mod::u_buf

Definition at line 141 of file fv_nesting.F90.

◆ uc_buf

type(fv_nest_bc_type_3d) fv_nesting_mod::uc_buf

Definition at line 141 of file fv_nesting.F90.

◆ v_buf

type(fv_nest_bc_type_3d) fv_nesting_mod::v_buf

Definition at line 141 of file fv_nesting.F90.

◆ vc_buf

type(fv_nest_bc_type_3d) fv_nesting_mod::vc_buf

Definition at line 141 of file fv_nesting.F90.

◆ w_buf

type(fv_nest_bc_type_3d) fv_nesting_mod::w_buf

Definition at line 141 of file fv_nesting.F90.