SCALE-RM
scale_atmos_dyn_tinteg_tracer_rk3.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  real(RP), allocatable :: qtrc_rk1(:,:,:)
58  real(RP), allocatable :: qtrc_rk2(:,:,:)
59  integer :: i_comm_rk1 = 1
60  integer :: i_comm_rk2 = 1
61  !-----------------------------------------------------------------------------
62 contains
63 
64  !-----------------------------------------------------------------------------
67  tinteg_type )
68  use scale_process, only: &
70  use scale_const, only: &
71  undef => const_undef
72  use scale_comm, only: &
74  implicit none
75 
76  character(len=*) :: tinteg_type
77 
78  integer :: iv
79  !---------------------------------------------------------------------------
80 
81  if ( tinteg_type /= 'RK3WS2002' ) then
82  write(*,*) 'xxx TINTEG_LARGE_TYPE is not RK3WS2002. Check!'
83  call prc_mpistop
84  end if
85 
86  allocate( qtrc_rk1(ka,ia,ja) )
87  allocate( qtrc_rk2(ka,ia,ja) )
88 
89  call comm_vars8_init( qtrc_rk1, i_comm_rk1 )
90  call comm_vars8_init( qtrc_rk2, i_comm_rk2 )
91 
92  return
94 
95  !-----------------------------------------------------------------------------
97  subroutine atmos_dyn_tinteg_tracer_rk3( &
98  QTRC, & ! (out)
99  qtrc0, rhoq_t, &! (in)
100  dens0, dens, & ! (in)
101  mflx_hi, num_diff, & ! (in)
102  gsqrt, mapf, & ! (in)
103  cdz, rcdz, rcdx, rcdy, & ! (in)
104  dtl, & ! (in)
105  flag_fct_tracer, & ! (in)
106  flag_fct_along_stream ) ! (in)
107  use scale_comm, only: &
108  comm_vars8, &
109  comm_wait
110  use scale_atmos_dyn_tstep_tracer, only: &
112  implicit none
113  real(RP), intent(inout) :: QTRC (ka,ia,ja)
114  real(RP), intent(in) :: QTRC0 (ka,ia,ja)
115  real(RP), intent(in) :: RHOQ_t (ka,ia,ja)
116  real(RP), intent(in) :: DENS0 (ka,ia,ja)
117  real(RP), intent(in) :: DENS (ka,ia,ja)
118  real(RP), intent(in) :: mflx_hi (ka,ia,ja,3)
119  real(RP), intent(in) :: num_diff(ka,ia,ja,3)
120  real(RP), intent(in) :: GSQRT (ka,ia,ja,7)
121  real(RP), intent(in) :: MAPF (ia,ja)
122  real(RP), intent(in) :: CDZ(ka)
123  real(RP), intent(in) :: RCDZ(ka)
124  real(RP), intent(in) :: RCDX(ia)
125  real(RP), intent(in) :: RCDY(ja)
126  real(RP), intent(in) :: dtl
127  logical, intent(in) :: FLAG_FCT_TRACER
128  logical, intent(in) :: FLAG_FCT_ALONG_STREAM
129 
130  real(RP) :: DENS_RK(ka,ia,ja)
131  real(RP) :: dtrk
132  integer :: k, i, j
133 
134  !------------------------------------------------------------------------
135  ! Start RK
136  !------------------------------------------------------------------------
137 
138  do j = js-1, je+1
139  do i = is-1, ie+1
140  do k = ks, ke
141  dens_rk(k,i,j) = dens0(k,i,j) &
142  + ( dens(k,i,j) - dens0(k,i,j) ) / 3.0_rp
143  end do
144  end do
145  end do
146 
147  dtrk = dtl / 3.0_rp
148  call atmos_dyn_tstep_tracer( &
149  qtrc_rk1, & ! (out)
150  qtrc, qtrc0, rhoq_t, &! (in)
151  dens0, dens_rk, & ! (in)
152  mflx_hi, num_diff, & ! (in)
153  gsqrt, mapf, & ! (in)
154  cdz, rcdz, rcdx, rcdy, & ! (in)
155  dtrk, & ! (in)
156  .false., flag_fct_along_stream ) ! (in)
157 
158  call comm_vars8( qtrc_rk1(:,:,:), i_comm_rk1 )
159  call comm_wait ( qtrc_rk1(:,:,:), i_comm_rk1, .false. )
160 
161 
162  do j = js-1, je+1
163  do i = is-1, ie+1
164  do k = ks, ke
165  dens_rk(k,i,j) = dens0(k,i,j) &
166  + ( dens(k,i,j) - dens0(k,i,j) ) * 0.5_rp
167  end do
168  end do
169  end do
170 
171  dtrk = dtl / 2.0_rp
172  call atmos_dyn_tstep_tracer( &
173  qtrc_rk2, & ! (out)
174  qtrc_rk1, qtrc0, rhoq_t, &! (in)
175  dens0, dens_rk, & ! (in)
176  mflx_hi, num_diff, & ! (in)
177  gsqrt, mapf, & ! (in)
178  cdz, rcdz, rcdx, rcdy, & ! (in)
179  dtrk, & ! (in)
180  .false., flag_fct_along_stream ) ! (in)
181 
182  call comm_vars8( qtrc_rk2(:,:,:), i_comm_rk2 )
183  call comm_wait ( qtrc_rk2(:,:,:), i_comm_rk2, .false. )
184 
185 
186  dtrk = dtl
187  call atmos_dyn_tstep_tracer( &
188  qtrc, & ! (out)
189  qtrc_rk2, qtrc0, rhoq_t, &! (in)
190  dens0, dens, & ! (in)
191  mflx_hi, num_diff, & ! (in)
192  gsqrt, mapf, & ! (in)
193  cdz, rcdz, rcdx, rcdy, & ! (in)
194  dtrk, & ! (in)
195  flag_fct_tracer, flag_fct_along_stream ) ! (in)
196 
197  return
198  end subroutine atmos_dyn_tinteg_tracer_rk3
199 
integer, public is
start point of inner domain: x, local
module DEBUG
Definition: scale_debug.F90:13
subroutine, public atmos_dyn_tinteg_tracer_rk3(QTRC, QTRC0, RHOQ_t, DENS0, DENS, mflx_hi, num_diff, GSQRT, MAPF, CDZ, RCDZ, RCDX, RCDY, dtl, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM)
RK3.
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:58
real(rp), public const_undef
Definition: scale_const.F90:43
module grid index
module TRACER
module Index
Definition: scale_index.F90:14
integer, public ia
of x whole cells (local, with HALO)
subroutine, public comm_vars8_init(var, vid)
Register variables.
Definition: scale_comm.F90:328
integer, public ka
of z whole cells (local, with HALO)
module COMMUNICATION
Definition: scale_comm.F90:23
integer, public js
start point of inner domain: y, local
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:40
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
module PRECISION
procedure(step), pointer, public atmos_dyn_tstep_tracer
subroutine, public atmos_dyn_tinteg_tracer_rk3_setup(tinteg_type)
Setup.
module Atmosphere / Dynamical scheme
integer, public ja
of y whole cells (local, with HALO)