SCALE-RM
scale_atmos_dyn_tinteg_tracer_euler.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
11 !-------------------------------------------------------------------------------
12 #include "scalelib.h"
14  !-----------------------------------------------------------------------------
15  !
16  !++ used modules
17  !
18  use scale_precision
19  use scale_io
20  use scale_prof
22  use scale_index
23  use scale_tracer
24 
25 #ifdef DEBUG
26  use scale_debug, only: &
27  check
28  use scale_const, only: &
29  undef => const_undef, &
30  iundef => const_undef2
31 #endif
32  !-----------------------------------------------------------------------------
33  implicit none
34  private
35  !-----------------------------------------------------------------------------
36  !
37  !++ Public procedure
38  !
41 
42  !-----------------------------------------------------------------------------
43  !
44  !++ Public parameters & variables
45  !
46  !-----------------------------------------------------------------------------
47  !
48  !++ Private procedure
49  !
50  !-----------------------------------------------------------------------------
51  !
52  !++ Private parameters & variables
53  !
54  !-----------------------------------------------------------------------------
55 contains
56 
57  !-----------------------------------------------------------------------------
60  tinteg_type )
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
77 
78  !-----------------------------------------------------------------------------
80  subroutine atmos_dyn_tinteg_tracer_euler( &
81  QTRC, & ! (out)
82  qtrc0, rhoq_t, &! (in)
83  dens0, dens, & ! (in)
84  mflx_hi, num_diff, & ! (in)
85  gsqrt, mapf, & ! (in)
86  cdz, rcdz, rcdx, rcdy, & ! (in)
87  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
88  dtl, & ! (in)
89  flag_fct_tracer, & ! (in)
90  flag_fct_along_stream ) ! (in)
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
133  end subroutine atmos_dyn_tinteg_tracer_euler
134 
module DEBUG
Definition: scale_debug.F90:11
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 check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:56
real(rp), public const_undef
Definition: scale_const.F90:41
subroutine, public atmos_dyn_tinteg_tracer_euler_setup(tinteg_type)
Setup.
module COMMUNICATION
module TRACER
module Index
Definition: scale_index.F90:11
module atmosphere / grid / cartesC index
module PROCESS
Definition: scale_prc.F90:11
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:38
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
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
procedure(step), pointer, public atmos_dyn_tstep_tracer
module Atmosphere / Dynamical scheme