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_TRACER_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  qflx, & ! (out)
83  qtrc0, rhoq_t, &! (in)
84  dens0, dens, & ! (in)
85  mflx_hi, num_diff, & ! (in)
86  gsqrt, mapf, & ! (in)
87  cdz, rcdz, rcdx, rcdy, & ! (in)
88  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
89  twod, & ! (in)
90  dtl, & ! (in)
91  flag_fct_tracer, & ! (in)
92  flag_fct_along_stream ) ! (in)
93  use scale_comm_cartesc, only: &
94  comm_vars8, &
95  comm_wait
96  use scale_atmos_dyn_tstep_tracer, only: &
98  implicit none
99  real(rp), intent(inout) :: qtrc (ka,ia,ja)
100  real(rp), intent(out) :: qflx (ka,ia,ja,3)
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  logical, intent(in) :: twod
118  real(rp), intent(in) :: dtl
119  logical, intent(in) :: flag_fct_tracer
120  logical, intent(in) :: flag_fct_along_stream
121 
122  !------------------------------------------------------------------------
123  ! Start RK
124  !------------------------------------------------------------------------
125 
126  call atmos_dyn_tstep_tracer( &
127  qtrc, & ! (out)
128  qflx, & ! (out)
129  qtrc, qtrc0, rhoq_t, &! (in)
130  dens0, dens, & ! (in)
131  mflx_hi, num_diff, & ! (in)
132  gsqrt, mapf, & ! (in)
133  cdz, rcdz, rcdx, rcdy, & ! (in)
134  twod, & ! (in)
135  dtl, & ! (in)
136  flag_fct_tracer, flag_fct_along_stream ) ! (in)
137 
138  return
139  end subroutine atmos_dyn_tinteg_tracer_euler
140 
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_index
module Index
Definition: scale_index.F90:11
scale_const::const_undef2
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:38
scale_atmos_dyn_tinteg_tracer_euler
module Atmosphere / Dyn Tinteg
Definition: scale_atmos_dyn_tinteg_tracer_euler.F90:13
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_io
module STDIO
Definition: scale_io.F90:10
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_debug::check
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:56
scale_atmos_dyn_tinteg_tracer_euler::atmos_dyn_tinteg_tracer_euler
subroutine, public atmos_dyn_tinteg_tracer_euler(QTRC, qflx, QTRC0, RHOQ_t, DENS0, DENS, mflx_hi, num_diff, GSQRT, MAPF, CDZ, RCDZ, RCDX, RCDY, BND_W, BND_E, BND_S, BND_N, TwoD, dtl, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM)
EULER.
Definition: scale_atmos_dyn_tinteg_tracer_euler.F90:93
scale_atmos_dyn_tstep_tracer
module Atmosphere / Dynamical scheme
Definition: scale_atmos_dyn_tstep_tracer.F90:12
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_atmos_dyn_tinteg_tracer_euler::atmos_dyn_tinteg_tracer_euler_setup
subroutine, public atmos_dyn_tinteg_tracer_euler_setup(tinteg_type)
Setup.
Definition: scale_atmos_dyn_tinteg_tracer_euler.F90:61
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_debug
module DEBUG
Definition: scale_debug.F90:11
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
scale_atmos_dyn_tstep_tracer::atmos_dyn_tstep_tracer
procedure(step), pointer, public atmos_dyn_tstep_tracer
Definition: scale_atmos_dyn_tstep_tracer.F90:72
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41