SCALE-RM
scale_atmos_dyn_tinteg_tracer.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
21  use scale_index
22  use scale_tracer
23  !-----------------------------------------------------------------------------
24  implicit none
25  private
26  !-----------------------------------------------------------------------------
27  !
28  !++ Public procedure
29  !
31 
33  abstract interface
34  subroutine tinteg( &
35  QTRC, & ! (out)
36  qtrc0, rhoq_t, &! (in)
37  dens0, dens, & ! (in)
38  mflx_hi, num_diff, & ! (in)
39  gsqrt, mapf, & ! (in)
40  cdz, rcdz, rcdx, rcdy, & ! (in)
41  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
42  dtl, & ! (in)
43  flag_fct_tracer, & ! (in)
44  flag_fct_along_stream ) ! (in)
45  use scale_precision
47  use scale_index
48  real(RP), intent(inout) :: qtrc (ka,ia,ja)
49  real(RP), intent(in) :: qtrc0 (ka,ia,ja)
50  real(RP), intent(in) :: rhoq_t (ka,ia,ja)
51  real(RP), intent(in) :: dens0 (ka,ia,ja)
52  real(RP), intent(in) :: dens (ka,ia,ja)
53  real(RP), intent(in) :: mflx_hi (ka,ia,ja,3)
54  real(RP), intent(in) :: num_diff(ka,ia,ja,3)
55  real(RP), intent(in) :: gsqrt (ka,ia,ja,7)
56  real(RP), intent(in) :: mapf (ia,ja)
57  real(RP), intent(in) :: cdz(ka)
58  real(RP), intent(in) :: rcdz(ka)
59  real(RP), intent(in) :: rcdx(ia)
60  real(RP), intent(in) :: rcdy(ja)
61  logical, intent(in) :: bnd_w
62  logical, intent(in) :: bnd_e
63  logical, intent(in) :: bnd_s
64  logical, intent(in) :: bnd_n
65  real(RP), intent(in) :: dtl
66  logical, intent(in) :: flag_fct_tracer
67  logical, intent(in) :: flag_fct_along_stream
68  end subroutine tinteg
69  end interface
70  procedure(tinteg), pointer :: atmos_dyn_tinteg_tracer => null()
71  public :: atmos_dyn_tinteg_tracer
72 
73  !-----------------------------------------------------------------------------
74  !
75  !++ Public parameters & variables
76  !
77  !-----------------------------------------------------------------------------
78  !
79  !++ Private procedure
80  !
81  !-----------------------------------------------------------------------------
82  !
83  !++ Private parameters & variables
84  !
85  !-----------------------------------------------------------------------------
86 contains
87  !-----------------------------------------------------------------------------
89  subroutine atmos_dyn_tinteg_tracer_setup( &
90  ATMOS_DYN_Tinteg_tracer_TYPE )
91 
92  use scale_precision
94  use scale_index
95  use scale_prc, only: &
96  prc_abort
103  implicit none
104  character(len=*), intent(in) :: ATMOS_DYN_Tinteg_tracer_TYPE
105  !---------------------------------------------------------------------------
106 
107  select case( atmos_dyn_tinteg_tracer_type )
108  case( 'EULER' )
110  atmos_dyn_tinteg_tracer_type )
112  case( 'RK3WS2002' )
114  atmos_dyn_tinteg_tracer_type )
116  case( 'OFF', 'NONE' )
117  ! do nothing
118  case default
119  log_error("ATMOS_DYN_Tinteg_tracer_setup",*) 'ATMOS_DYN_TINTEG_TRACER_TYPE is invalid: ', atmos_dyn_tinteg_tracer_type
120  call prc_abort
121  end select
122 
123  return
124  end subroutine atmos_dyn_tinteg_tracer_setup
125 
module Atmosphere / Dynamics Temporal integration
integer, public ia
of whole cells: x, local, with HALO
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.
integer, public ja
of whole cells: y, local, with HALO
subroutine, public atmos_dyn_tinteg_tracer_euler_setup(tinteg_type)
Setup.
module TRACER
module Index
Definition: scale_index.F90:11
module atmosphere / grid / cartesC index
module PROCESS
Definition: scale_prc.F90:11
procedure(tinteg), pointer, public atmos_dyn_tinteg_tracer
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module profiler
Definition: scale_prof.F90:11
module PRECISION
integer, public ka
of whole cells: z, local, with HALO
module STDIO
Definition: scale_io.F90:10
subroutine, public atmos_dyn_tinteg_tracer_rk3_setup(tinteg_type)
Setup.
subroutine, public atmos_dyn_tinteg_tracer_setup(ATMOS_DYN_Tinteg_tracer_TYPE)
Register.
subroutine, public atmos_dyn_tinteg_tracer_rk3(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)
RK3.