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, 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
History
  • 2016-05-17 (S.Nishizawa) [new]

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 64 of file scale_atmos_dyn_tinteg_tracer_euler.F90.

References scale_comm::comm_vars8_init(), scale_const::const_undef, and scale_process::prc_mpistop().

Referenced by scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_setup().

64  use scale_process, only: &
66  use scale_const, only: &
67  undef => const_undef
68  use scale_comm, only: &
70  implicit none
71 
72  character(len=*) :: tinteg_type
73 
74  integer :: iv
75  !---------------------------------------------------------------------------
76 
77  if ( tinteg_type /= 'EULER' ) then
78  write(*,*) 'xxx TINTEG_LARGE_TYPE is not EULER. Check!'
79  call prc_mpistop
80  end if
81 
82  return
subroutine, public prc_mpistop
Abort MPI.
real(rp), public const_undef
Definition: scale_const.F90:43
subroutine, public comm_vars8_init(var, vid)
Register variables.
Definition: scale_comm.F90:328
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
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,
real(rp), intent(in)  dtl,
logical, intent(in)  FLAG_FCT_TRACER,
logical, intent(in)  FLAG_FCT_ALONG_STREAM 
)

EULER.

Definition at line 97 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().

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