SCALE-RM
scale_atmos_dyn_tinteg_tracer_euler.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14 !-------------------------------------------------------------------------------
15 #include "inc_openmp.h"
17  !-----------------------------------------------------------------------------
18  !
19  !++ used modules
20  !
21  use scale_precision
22  use scale_stdio
23  use scale_prof
25  use scale_index
26  use scale_tracer
27 
28 #ifdef DEBUG
29  use scale_debug, only: &
30  check
31  use scale_const, only: &
32  undef => const_undef, &
33  iundef => const_undef2
34 #endif
35  !-----------------------------------------------------------------------------
36  implicit none
37  private
38  !-----------------------------------------------------------------------------
39  !
40  !++ Public procedure
41  !
44 
45  !-----------------------------------------------------------------------------
46  !
47  !++ Public parameters & variables
48  !
49  !-----------------------------------------------------------------------------
50  !
51  !++ Private procedure
52  !
53  !-----------------------------------------------------------------------------
54  !
55  !++ Private parameters & variables
56  !
57  !-----------------------------------------------------------------------------
58 contains
59 
60  !-----------------------------------------------------------------------------
63  tinteg_type )
64  use scale_process, only: &
66  implicit none
67 
68  character(len=*) :: tinteg_type
69 
70  integer :: iv
71  !---------------------------------------------------------------------------
72 
73  if ( tinteg_type /= 'EULER' ) then
74  write(*,*) 'xxx TINTEG_LARGE_TYPE is not EULER. Check!'
75  call prc_mpistop
76  end if
77 
78  return
80 
81  !-----------------------------------------------------------------------------
83  subroutine atmos_dyn_tinteg_tracer_euler( &
84  QTRC, & ! (out)
85  qtrc0, rhoq_t, &! (in)
86  dens0, dens, & ! (in)
87  mflx_hi, num_diff, & ! (in)
88  gsqrt, mapf, & ! (in)
89  cdz, rcdz, rcdx, rcdy, & ! (in)
90  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
91  dtl, & ! (in)
92  flag_fct_tracer, & ! (in)
93  flag_fct_along_stream ) ! (in)
94  use scale_comm, only: &
95  comm_vars8, &
96  comm_wait
97  use scale_atmos_dyn_tstep_tracer, only: &
99  implicit none
100  real(RP), intent(inout) :: qtrc (ka,ia,ja)
101  real(RP), intent(in) :: qtrc0 (ka,ia,ja)
102  real(RP), intent(in) :: rhoq_t (ka,ia,ja)
103  real(RP), intent(in) :: dens0 (ka,ia,ja)
104  real(RP), intent(in) :: dens (ka,ia,ja)
105  real(RP), intent(in) :: mflx_hi (ka,ia,ja,3)
106  real(RP), intent(in) :: num_diff(ka,ia,ja,3)
107  real(RP), intent(in) :: gsqrt (ka,ia,ja,7)
108  real(RP), intent(in) :: mapf (ia,ja)
109  real(RP), intent(in) :: cdz(ka)
110  real(RP), intent(in) :: rcdz(ka)
111  real(RP), intent(in) :: rcdx(ia)
112  real(RP), intent(in) :: rcdy(ja)
113  logical, intent(in) :: bnd_w
114  logical, intent(in) :: bnd_e
115  logical, intent(in) :: bnd_s
116  logical, intent(in) :: bnd_n
117  real(RP), intent(in) :: dtl
118  logical, intent(in) :: flag_fct_tracer
119  logical, intent(in) :: flag_fct_along_stream
120 
121  !------------------------------------------------------------------------
122  ! Start RK
123  !------------------------------------------------------------------------
124 
125  call atmos_dyn_tstep_tracer( &
126  qtrc, & ! (out)
127  qtrc, qtrc0, rhoq_t, &! (in)
128  dens0, dens, & ! (in)
129  mflx_hi, num_diff, & ! (in)
130  gsqrt, mapf, & ! (in)
131  cdz, rcdz, rcdx, rcdy, & ! (in)
132  dtl, & ! (in)
133  flag_fct_tracer, flag_fct_along_stream ) ! (in)
134 
135  return
136  end subroutine atmos_dyn_tinteg_tracer_euler
137 
module DEBUG
Definition: scale_debug.F90:13
subroutine, public prc_mpistop
Abort MPI.
module STDIO
Definition: scale_stdio.F90:12
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.
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:58
real(rp), public const_undef
Definition: scale_const.F90:43
subroutine, public atmos_dyn_tinteg_tracer_euler_setup(tinteg_type)
Setup.
module grid index
module TRACER
module Index
Definition: scale_index.F90:14
integer, public ia
of whole cells: x, local, with HALO
integer, public ka
of whole cells: z, local, with HALO
module COMMUNICATION
Definition: scale_comm.F90:23
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:40
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
module profiler
Definition: scale_prof.F90:10
module PRECISION
procedure(step), pointer, public atmos_dyn_tstep_tracer
module Atmosphere / Dynamical scheme
integer, public ja
of whole cells: y, local, with HALO