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_comm, only: &
72  implicit none
73 
74  character(len=*) :: tinteg_type
75 
76  integer :: iv
77  !---------------------------------------------------------------------------
78 
79  if ( tinteg_type /= 'RK3WS2002' ) then
80  write(*,*) 'xxx TINTEG_LARGE_TYPE is not RK3WS2002. Check!'
81  call prc_mpistop
82  end if
83 
84  allocate( qtrc_rk1(ka,ia,ja) )
85  allocate( qtrc_rk2(ka,ia,ja) )
86 
87  call comm_vars8_init( 'QTRC_RK1', qtrc_rk1, i_comm_rk1 )
88  call comm_vars8_init( 'QTRC_RK2', qtrc_rk2, i_comm_rk2 )
89 
90  return
92 
93  !-----------------------------------------------------------------------------
95  subroutine atmos_dyn_tinteg_tracer_rk3( &
96  QTRC, & ! (out)
97  qtrc0, rhoq_t, &! (in)
98  dens0, dens, & ! (in)
99  mflx_hi, num_diff, & ! (in)
100  gsqrt, mapf, & ! (in)
101  cdz, rcdz, rcdx, rcdy, & ! (in)
102  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
103  dtl, & ! (in)
104  flag_fct_tracer, & ! (in)
105  flag_fct_along_stream ) ! (in)
106  use scale_comm, only: &
107  comm_vars8, &
108  comm_wait
109  use scale_atmos_dyn_tstep_tracer, only: &
111  use scale_atmos_dyn_common, only: &
113  implicit none
114  real(RP), intent(inout) :: qtrc (ka,ia,ja)
115  real(RP), intent(in) :: qtrc0 (ka,ia,ja)
116  real(RP), intent(in) :: rhoq_t (ka,ia,ja)
117  real(RP), intent(in) :: dens0 (ka,ia,ja)
118  real(RP), intent(in) :: dens (ka,ia,ja)
119  real(RP), intent(in) :: mflx_hi (ka,ia,ja,3)
120  real(RP), intent(in) :: num_diff(ka,ia,ja,3)
121  real(RP), intent(in) :: gsqrt (ka,ia,ja,7)
122  real(RP), intent(in) :: mapf (ia,ja)
123  real(RP), intent(in) :: cdz(ka)
124  real(RP), intent(in) :: rcdz(ka)
125  real(RP), intent(in) :: rcdx(ia)
126  real(RP), intent(in) :: rcdy(ja)
127  logical, intent(in) :: bnd_w
128  logical, intent(in) :: bnd_e
129  logical, intent(in) :: bnd_s
130  logical, intent(in) :: bnd_n
131  real(RP), intent(in) :: dtl
132  logical, intent(in) :: flag_fct_tracer
133  logical, intent(in) :: flag_fct_along_stream
134 
135  real(RP) :: dens_rk(ka,ia,ja)
136  real(RP) :: dtrk
137  integer :: k, i, j
138 
139  !------------------------------------------------------------------------
140  ! Start RK
141  !------------------------------------------------------------------------
142 
143  do j = js-1, je+1
144  do i = is-1, ie+1
145  do k = ks, ke
146  dens_rk(k,i,j) = dens0(k,i,j) &
147  + ( dens(k,i,j) - dens0(k,i,j) ) / 3.0_rp
148  end do
149  end do
150  end do
151 
152  dtrk = dtl / 3.0_rp
153  call atmos_dyn_tstep_tracer( &
154  qtrc_rk1, & ! (out)
155  qtrc, qtrc0, rhoq_t, &! (in)
156  dens0, dens_rk, & ! (in)
157  mflx_hi, num_diff, & ! (in)
158  gsqrt, mapf, & ! (in)
159  cdz, rcdz, rcdx, rcdy, & ! (in)
160  dtrk, & ! (in)
161  .false., flag_fct_along_stream ) ! (in)
162 
163  call atmos_dyn_copy_boundary_tracer( qtrc_rk1, & ! [INOUT]
164  qtrc0, & ! [IN]
165  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
166 
167  call comm_vars8( qtrc_rk1(:,:,:), i_comm_rk1 )
168  call comm_wait ( qtrc_rk1(:,:,:), i_comm_rk1, .false. )
169 
170 
171  do j = js-1, je+1
172  do i = is-1, ie+1
173  do k = ks, ke
174  dens_rk(k,i,j) = dens0(k,i,j) &
175  + ( dens(k,i,j) - dens0(k,i,j) ) * 0.5_rp
176  end do
177  end do
178  end do
179 
180  dtrk = dtl / 2.0_rp
181  call atmos_dyn_tstep_tracer( &
182  qtrc_rk2, & ! (out)
183  qtrc_rk1, qtrc0, rhoq_t, &! (in)
184  dens0, dens_rk, & ! (in)
185  mflx_hi, num_diff, & ! (in)
186  gsqrt, mapf, & ! (in)
187  cdz, rcdz, rcdx, rcdy, & ! (in)
188  dtrk, & ! (in)
189  .false., flag_fct_along_stream ) ! (in)
190 
191  call atmos_dyn_copy_boundary_tracer( qtrc_rk2, & ! [INOUT]
192  qtrc0, & ! [IN]
193  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
194 
195  call comm_vars8( qtrc_rk2(:,:,:), i_comm_rk2 )
196  call comm_wait ( qtrc_rk2(:,:,:), i_comm_rk2, .false. )
197 
198 
199  dtrk = dtl
200  call atmos_dyn_tstep_tracer( &
201  qtrc, & ! (out)
202  qtrc_rk2, qtrc0, rhoq_t, &! (in)
203  dens0, dens, & ! (in)
204  mflx_hi, num_diff, & ! (in)
205  gsqrt, mapf, & ! (in)
206  cdz, rcdz, rcdx, rcdy, & ! (in)
207  dtrk, & ! (in)
208  flag_fct_tracer, flag_fct_along_stream ) ! (in)
209 
210  return
211  end subroutine atmos_dyn_tinteg_tracer_rk3
212 
integer, public is
start point of inner domain: x, local
module DEBUG
Definition: scale_debug.F90:13
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
Definition: scale_comm.F90:313
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 whole cells: x, local, with HALO
integer, public ka
of whole cells: z, 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 Atmosphere / Dynamics common
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_copy_boundary_tracer(QTRC, QTRC0, BND_W, BND_E, BND_S, BND_N)
subroutine, public atmos_dyn_tinteg_tracer_rk3_setup(tinteg_type)
Setup.
module Atmosphere / Dynamical scheme
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.
integer, public ja
of whole cells: y, local, with HALO