SCALE-RM
Functions/Subroutines
scale_atmos_dyn_tinteg_tracer_euler Module Reference

module Atmosphere / Dyn Tinteg More...

Functions/Subroutines

subroutine, public atmos_dyn_tinteg_tracer_euler_setup (tinteg_type)
 Setup. More...
 
subroutine, public atmos_dyn_tinteg_tracer_euler (QTRC, QTRC0, RHOQ_t, DENS0, DENS, mflx_hi, num_diff, GSQRT, MAPF, CDZ, RCDZ, RCDX, RCDY, BND_W, BND_E, BND_S, BND_N, dtl, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM)
 EULER. More...
 

Detailed Description

module Atmosphere / Dyn Tinteg

Description
Temporal integration for tracer advection for Atmospheric process Euler scheme
Author
Team SCALE

Function/Subroutine Documentation

◆ atmos_dyn_tinteg_tracer_euler_setup()

subroutine, public scale_atmos_dyn_tinteg_tracer_euler::atmos_dyn_tinteg_tracer_euler_setup ( character(len=*)  tinteg_type)

Setup.

Definition at line 61 of file scale_atmos_dyn_tinteg_tracer_euler.F90.

References scale_prc::prc_abort().

Referenced by scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_setup().

61  use scale_prc, only: &
62  prc_abort
63  implicit none
64 
65  character(len=*) :: tinteg_type
66 
67  integer :: iv
68  !---------------------------------------------------------------------------
69 
70  if ( tinteg_type /= 'EULER' ) then
71  log_error("ATMOS_DYN_Tinteg_tracer_euler_setup",*) 'TINTEG_LARGE_TYPE is not EULER. Check!'
72  call prc_abort
73  end if
74 
75  return
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_dyn_tinteg_tracer_euler()

subroutine, public scale_atmos_dyn_tinteg_tracer_euler::atmos_dyn_tinteg_tracer_euler ( real(rp), dimension (ka,ia,ja), intent(inout)  QTRC,
real(rp), dimension (ka,ia,ja), intent(in)  QTRC0,
real(rp), dimension (ka,ia,ja), intent(in)  RHOQ_t,
real(rp), dimension (ka,ia,ja), intent(in)  DENS0,
real(rp), dimension (ka,ia,ja), intent(in)  DENS,
real(rp), dimension (ka,ia,ja,3), intent(in)  mflx_hi,
real(rp), dimension(ka,ia,ja,3), intent(in)  num_diff,
real(rp), dimension (ka,ia,ja,7), intent(in)  GSQRT,
real(rp), dimension (ia,ja), intent(in)  MAPF,
real(rp), dimension(ka), intent(in)  CDZ,
real(rp), dimension(ka), intent(in)  RCDZ,
real(rp), dimension(ia), intent(in)  RCDX,
real(rp), dimension(ja), intent(in)  RCDY,
logical, intent(in)  BND_W,
logical, intent(in)  BND_E,
logical, intent(in)  BND_S,
logical, intent(in)  BND_N,
real(rp), intent(in)  dtl,
logical, intent(in)  FLAG_FCT_TRACER,
logical, intent(in)  FLAG_FCT_ALONG_STREAM 
)

EULER.

Definition at line 91 of file scale_atmos_dyn_tinteg_tracer_euler.F90.

References scale_atmos_dyn_tstep_tracer::atmos_dyn_tstep_tracer.

Referenced by scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_setup().

91  use scale_comm_cartesc, only: &
92  comm_vars8, &
93  comm_wait
94  use scale_atmos_dyn_tstep_tracer, only: &
96  implicit none
97  real(RP), intent(inout) :: qtrc (ka,ia,ja)
98  real(RP), intent(in) :: qtrc0 (ka,ia,ja)
99  real(RP), intent(in) :: rhoq_t (ka,ia,ja)
100  real(RP), intent(in) :: dens0 (ka,ia,ja)
101  real(RP), intent(in) :: dens (ka,ia,ja)
102  real(RP), intent(in) :: mflx_hi (ka,ia,ja,3)
103  real(RP), intent(in) :: num_diff(ka,ia,ja,3)
104  real(RP), intent(in) :: gsqrt (ka,ia,ja,7)
105  real(RP), intent(in) :: mapf (ia,ja)
106  real(RP), intent(in) :: cdz(ka)
107  real(RP), intent(in) :: rcdz(ka)
108  real(RP), intent(in) :: rcdx(ia)
109  real(RP), intent(in) :: rcdy(ja)
110  logical, intent(in) :: bnd_w
111  logical, intent(in) :: bnd_e
112  logical, intent(in) :: bnd_s
113  logical, intent(in) :: bnd_n
114  real(RP), intent(in) :: dtl
115  logical, intent(in) :: flag_fct_tracer
116  logical, intent(in) :: flag_fct_along_stream
117 
118  !------------------------------------------------------------------------
119  ! Start RK
120  !------------------------------------------------------------------------
121 
122  call atmos_dyn_tstep_tracer( &
123  qtrc, & ! (out)
124  qtrc, qtrc0, rhoq_t, &! (in)
125  dens0, dens, & ! (in)
126  mflx_hi, num_diff, & ! (in)
127  gsqrt, mapf, & ! (in)
128  cdz, rcdz, rcdx, rcdy, & ! (in)
129  dtl, & ! (in)
130  flag_fct_tracer, flag_fct_along_stream ) ! (in)
131 
132  return
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
module COMMUNICATION
integer, public ka
of whole cells: z, local, with HALO
procedure(step), pointer, public atmos_dyn_tstep_tracer
module Atmosphere / Dynamical scheme
Here is the caller graph for this function: