SCALE-RM
scale_atmos_dyn_tinteg_tracer_rk3.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  real(RP), allocatable :: qtrc_rk1(:,:,:)
55  real(RP), allocatable :: qtrc_rk2(:,:,:)
56  integer :: i_comm_rk1 = 1
57  integer :: i_comm_rk2 = 1
58  !-----------------------------------------------------------------------------
59 contains
60 
61  !-----------------------------------------------------------------------------
64  tinteg_type )
65  use scale_prc, only: &
66  prc_abort
67  use scale_comm_cartesc, only: &
69  implicit none
70 
71  character(len=*) :: tinteg_type
72 
73  integer :: iv
74  !---------------------------------------------------------------------------
75 
76  if ( tinteg_type /= 'RK3WS2002' ) then
77  log_error("ATMOS_DYN_Tinteg_tracer_rk3_setup",*) 'TINTEG_LARGE_TYPE is not RK3WS2002. Check!'
78  call prc_abort
79  end if
80 
81  allocate( qtrc_rk1(ka,ia,ja) )
82  allocate( qtrc_rk2(ka,ia,ja) )
83 
84  call comm_vars8_init( 'QTRC_RK1', qtrc_rk1, i_comm_rk1 )
85  call comm_vars8_init( 'QTRC_RK2', qtrc_rk2, i_comm_rk2 )
86 
87  return
89 
90  !-----------------------------------------------------------------------------
92  subroutine atmos_dyn_tinteg_tracer_rk3( &
93  QTRC, & ! (out)
94  qtrc0, rhoq_t, &! (in)
95  dens0, dens, & ! (in)
96  mflx_hi, num_diff, & ! (in)
97  gsqrt, mapf, & ! (in)
98  cdz, rcdz, rcdx, rcdy, & ! (in)
99  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
100  dtl, & ! (in)
101  flag_fct_tracer, & ! (in)
102  flag_fct_along_stream ) ! (in)
103  use scale_comm_cartesc, only: &
104  comm_vars8, &
105  comm_wait
106  use scale_atmos_dyn_tstep_tracer, only: &
108  use scale_atmos_dyn_common, only: &
110  implicit none
111  real(RP), intent(inout) :: QTRC (ka,ia,ja)
112  real(RP), intent(in) :: QTRC0 (ka,ia,ja)
113  real(RP), intent(in) :: RHOQ_t (ka,ia,ja)
114  real(RP), intent(in) :: DENS0 (ka,ia,ja)
115  real(RP), intent(in) :: DENS (ka,ia,ja)
116  real(RP), intent(in) :: mflx_hi (ka,ia,ja,3)
117  real(RP), intent(in) :: num_diff(ka,ia,ja,3)
118  real(RP), intent(in) :: GSQRT (ka,ia,ja,7)
119  real(RP), intent(in) :: MAPF (ia,ja)
120  real(RP), intent(in) :: CDZ(ka)
121  real(RP), intent(in) :: RCDZ(ka)
122  real(RP), intent(in) :: RCDX(ia)
123  real(RP), intent(in) :: RCDY(ja)
124  logical, intent(in) :: BND_W
125  logical, intent(in) :: BND_E
126  logical, intent(in) :: BND_S
127  logical, intent(in) :: BND_N
128  real(RP), intent(in) :: dtl
129  logical, intent(in) :: FLAG_FCT_TRACER
130  logical, intent(in) :: FLAG_FCT_ALONG_STREAM
131 
132  real(RP) :: DENS_RK(ka,ia,ja)
133  real(RP) :: dtrk
134  integer :: k, i, j
135 
136  !------------------------------------------------------------------------
137  ! Start RK
138  !------------------------------------------------------------------------
139 
140  do j = js-1, je+1
141  do i = is-1, ie+1
142  do k = ks, ke
143  dens_rk(k,i,j) = dens0(k,i,j) &
144  + ( dens(k,i,j) - dens0(k,i,j) ) / 3.0_rp
145  end do
146  end do
147  end do
148 
149  dtrk = dtl / 3.0_rp
150  call atmos_dyn_tstep_tracer( &
151  qtrc_rk1, & ! (out)
152  qtrc, qtrc0, rhoq_t, &! (in)
153  dens0, dens_rk, & ! (in)
154  mflx_hi, num_diff, & ! (in)
155  gsqrt, mapf, & ! (in)
156  cdz, rcdz, rcdx, rcdy, & ! (in)
157  dtrk, & ! (in)
158  .false., flag_fct_along_stream ) ! (in)
159 
160  call atmos_dyn_copy_boundary_tracer( qtrc_rk1, & ! [INOUT]
161  qtrc0, & ! [IN]
162  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
163 
164  call comm_vars8( qtrc_rk1(:,:,:), i_comm_rk1 )
165  call comm_wait ( qtrc_rk1(:,:,:), i_comm_rk1, .false. )
166 
167 
168  do j = js-1, je+1
169  do i = is-1, ie+1
170  do k = ks, ke
171  dens_rk(k,i,j) = dens0(k,i,j) &
172  + ( dens(k,i,j) - dens0(k,i,j) ) * 0.5_rp
173  end do
174  end do
175  end do
176 
177  dtrk = dtl / 2.0_rp
178  call atmos_dyn_tstep_tracer( &
179  qtrc_rk2, & ! (out)
180  qtrc_rk1, qtrc0, rhoq_t, &! (in)
181  dens0, dens_rk, & ! (in)
182  mflx_hi, num_diff, & ! (in)
183  gsqrt, mapf, & ! (in)
184  cdz, rcdz, rcdx, rcdy, & ! (in)
185  dtrk, & ! (in)
186  .false., flag_fct_along_stream ) ! (in)
187 
188  call atmos_dyn_copy_boundary_tracer( qtrc_rk2, & ! [INOUT]
189  qtrc0, & ! [IN]
190  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
191 
192  call comm_vars8( qtrc_rk2(:,:,:), i_comm_rk2 )
193  call comm_wait ( qtrc_rk2(:,:,:), i_comm_rk2, .false. )
194 
195 
196  dtrk = dtl
197  call atmos_dyn_tstep_tracer( &
198  qtrc, & ! (out)
199  qtrc_rk2, qtrc0, rhoq_t, &! (in)
200  dens0, dens, & ! (in)
201  mflx_hi, num_diff, & ! (in)
202  gsqrt, mapf, & ! (in)
203  cdz, rcdz, rcdx, rcdy, & ! (in)
204  dtrk, & ! (in)
205  flag_fct_tracer, flag_fct_along_stream ) ! (in)
206 
207  return
208  end subroutine atmos_dyn_tinteg_tracer_rk3
209 
module DEBUG
Definition: scale_debug.F90:11
integer, public ia
of whole cells: x, local, with HALO
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
module COMMUNICATION
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
module TRACER
module Index
Definition: scale_index.F90:11
module atmosphere / grid / cartesC index
integer, public ke
end point of inner domain: z, local
module PROCESS
Definition: scale_prc.F90:11
integer, public je
end point of inner domain: y, local
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:38
module Atmosphere / Dynamics common
integer, public ks
start point of inner domain: z, local
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
integer, public js
start point of inner domain: y, local
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
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.