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  !
42 
43  !-----------------------------------------------------------------------------
44  !
45  !++ Public parameters & variables
46  !
47  !-----------------------------------------------------------------------------
48  !
49  !++ Private procedure
50  !
51  !-----------------------------------------------------------------------------
52  !
53  !++ Private parameters & variables
54  !
55  real(RP), allocatable :: QTRC_RK1(:,:,:)
56  real(RP), allocatable :: QTRC_RK2(:,:,:)
57  integer :: I_COMM_RK1
58  integer :: I_COMM_RK2
59  !-----------------------------------------------------------------------------
60 contains
61 
62  !-----------------------------------------------------------------------------
65  tinteg_type )
66  use scale_const, only: &
67  undef => const_undef
68  use scale_prc, only: &
69  prc_abort
70  use scale_comm_cartesc, only: &
72  implicit none
73 
74  character(len=*) :: tinteg_type
75 
76  integer :: iv
77  !---------------------------------------------------------------------------
78 
79  if ( tinteg_type /= 'RK3WS2002' ) then
80  log_error("ATMOS_DYN_Tinteg_tracer_rk3_setup",*) 'TINTEG_TRACER_TYPE is not RK3WS2002. Check!'
81  call prc_abort
82  end if
83 
84  allocate( qtrc_rk1(ka,ia,ja) )
85  allocate( qtrc_rk2(ka,ia,ja) )
86  qtrc_rk1(:,:,:) = undef
87  qtrc_rk2(:,:,:) = undef
88  !$acc enter data create(QTRC_RK1, QTRC_RK2)
89 
90  i_comm_rk1 = 1
91  call comm_vars8_init( 'QTRC_RK1', qtrc_rk1, i_comm_rk1 )
92 
93  i_comm_rk2 = 1
94  call comm_vars8_init( 'QTRC_RK2', qtrc_rk2, i_comm_rk2 )
95 
96  return
98 
99  !-----------------------------------------------------------------------------
102 
103  !$acc exit data delete(QTRC_RK1, QTRC_RK2)
104  deallocate( qtrc_rk1 )
105  deallocate( qtrc_rk2 )
106 
107  return
109 
110  !-----------------------------------------------------------------------------
112  subroutine atmos_dyn_tinteg_tracer_rk3( &
113  QTRC, & ! (out)
114  qflx, & ! (out)
115  qtrc0, rhoq_t, &! (in)
116  dens0, dens, & ! (in)
117  mflx_hi, num_diff, & ! (in)
118  gsqrt, mapf, & ! (in)
119  cdz, rcdz, rcdx, rcdy, & ! (in)
120  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
121  twod, & ! (in)
122  dtl, & ! (in)
123  flag_fct_tracer, & ! (in)
124  flag_fct_along_stream ) ! (in)
125  use scale_comm_cartesc, only: &
126  comm_vars8, &
127  comm_wait
128  use scale_atmos_dyn_tstep_tracer, only: &
130  use scale_atmos_dyn_common, only: &
132  implicit none
133  real(rp), intent(inout) :: qtrc (ka,ia,ja)
134  real(rp), intent(out) :: qflx (ka,ia,ja,3)
135  real(rp), intent(in) :: qtrc0 (ka,ia,ja)
136  real(rp), intent(in) :: rhoq_t (ka,ia,ja)
137  real(rp), intent(in) :: dens0 (ka,ia,ja)
138  real(rp), intent(in) :: dens (ka,ia,ja)
139  real(rp), intent(in) :: mflx_hi (ka,ia,ja,3)
140  real(rp), intent(in) :: num_diff(ka,ia,ja,3)
141  real(rp), intent(in) :: gsqrt (ka,ia,ja,7)
142  real(rp), intent(in) :: mapf (ia,ja)
143  real(rp), intent(in) :: cdz(ka)
144  real(rp), intent(in) :: rcdz(ka)
145  real(rp), intent(in) :: rcdx(ia)
146  real(rp), intent(in) :: rcdy(ja)
147  logical, intent(in) :: bnd_w
148  logical, intent(in) :: bnd_e
149  logical, intent(in) :: bnd_s
150  logical, intent(in) :: bnd_n
151  logical, intent(in) :: twod
152  real(rp), intent(in) :: dtl
153  logical, intent(in) :: flag_fct_tracer
154  logical, intent(in) :: flag_fct_along_stream
155 
156  real(rp) :: dens_rk(ka,ia,ja)
157  real(rp) :: dtrk
158  integer :: k, i, j
159 
160  !$acc data copy(QTRC), copyout(qflx) &
161  !$acc copyin(QTRC0, RHOQ_t, DENS0, DENS, mflx_hi, num_diff, &
162  !$acc GSQRT, MAPF, CDZ, RCDZ, RCDX, RCDY) &
163  !$acc create(DENS_RK)
164  !------------------------------------------------------------------------
165  ! Start RK
166  !------------------------------------------------------------------------
167 
168  !$acc kernels
169  do j = js-1, je+1
170  do i = max(is-1,1), min(ie+1,ia)
171  do k = ks, ke
172  dens_rk(k,i,j) = dens0(k,i,j) &
173  + ( dens(k,i,j) - dens0(k,i,j) ) / 3.0_rp
174  end do
175  end do
176  end do
177  !$acc end kernels
178 
179  dtrk = dtl / 3.0_rp
180  call atmos_dyn_tstep_tracer( &
181  qtrc_rk1, & ! (out)
182  qflx, & ! (out)
183  qtrc, 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  twod, dtrk, & ! (in)
189  .false., flag_fct_along_stream ) ! (in)
190 
191  call atmos_dyn_copy_boundary_tracer( qtrc_rk1, & ! [INOUT]
192  qtrc0, & ! [IN]
193  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
194  twod ) ! [IN]
195 
196  call comm_vars8( qtrc_rk1(:,:,:), i_comm_rk1 )
197  call comm_wait ( qtrc_rk1(:,:,:), i_comm_rk1, .false. )
198 
199 
200  !$acc kernels
201  do j = js-1, je+1
202  do i = max(is-1,1), min(ie+1,ia)
203  do k = ks, ke
204  dens_rk(k,i,j) = dens0(k,i,j) &
205  + ( dens(k,i,j) - dens0(k,i,j) ) * 0.5_rp
206  end do
207  end do
208  end do
209  !$acc end kernels
210 
211  dtrk = dtl / 2.0_rp
212  call atmos_dyn_tstep_tracer( &
213  qtrc_rk2, & ! (out)
214  qflx, & ! (out)
215  qtrc_rk1, qtrc0, rhoq_t, &! (in)
216  dens0, dens_rk, & ! (in)
217  mflx_hi, num_diff, & ! (in)
218  gsqrt, mapf, & ! (in)
219  cdz, rcdz, rcdx, rcdy, & ! (in)
220  twod, dtrk, & ! (in)
221  .false., flag_fct_along_stream ) ! (in)
222 
223  call atmos_dyn_copy_boundary_tracer( qtrc_rk2, & ! [INOUT]
224  qtrc0, & ! [IN]
225  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
226  twod ) ! [IN]
227 
228  call comm_vars8( qtrc_rk2(:,:,:), i_comm_rk2 )
229  call comm_wait ( qtrc_rk2(:,:,:), i_comm_rk2, .false. )
230 
231 
232  dtrk = dtl
233  call atmos_dyn_tstep_tracer( &
234  qtrc, & ! (out)
235  qflx, & ! (out)
236  qtrc_rk2, qtrc0, rhoq_t, &! (in)
237  dens0, dens, & ! (in)
238  mflx_hi, num_diff, & ! (in)
239  gsqrt, mapf, & ! (in)
240  cdz, rcdz, rcdx, rcdy, & ! (in)
241  twod, & ! (in)
242  dtrk, & ! (in)
243  flag_fct_tracer, flag_fct_along_stream ) ! (in)
244 
245  !$acc end data
246 
247  return
248  end subroutine atmos_dyn_tinteg_tracer_rk3
249 
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_atmos_dyn_tinteg_tracer_rk3::atmos_dyn_tinteg_tracer_rk3_setup
subroutine, public atmos_dyn_tinteg_tracer_rk3_setup(tinteg_type)
Setup.
Definition: scale_atmos_dyn_tinteg_tracer_rk3.F90:66
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:40
scale_atmos_dyn_tinteg_tracer_rk3
module Atmosphere / Dyn Tinteg
Definition: scale_atmos_dyn_tinteg_tracer_rk3.F90:13
scale_atmos_dyn_common
module Atmosphere / Dynamics common
Definition: scale_atmos_dyn_common.F90:12
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_atmos_dyn_tinteg_tracer_rk3::atmos_dyn_tinteg_tracer_rk3_finalize
subroutine, public atmos_dyn_tinteg_tracer_rk3_finalize
finalize
Definition: scale_atmos_dyn_tinteg_tracer_rk3.F90:102
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
scale_io
module STDIO
Definition: scale_io.F90:10
scale_tracer::k
real(rp), public k
Definition: scale_tracer.F90:45
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:59
scale_atmos_dyn_tinteg_tracer_rk3::atmos_dyn_tinteg_tracer_rk3
subroutine, public atmos_dyn_tinteg_tracer_rk3(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)
RK3.
Definition: scale_atmos_dyn_tinteg_tracer_rk3.F90:125
scale_atmos_dyn_common::atmos_dyn_copy_boundary_tracer
subroutine, public atmos_dyn_copy_boundary_tracer(QTRC, QTRC0, BND_W, BND_E, BND_S, BND_N, TwoD)
Definition: scale_atmos_dyn_common.F90:359
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::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
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_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
scale_comm_cartesc::comm_vars8_init
subroutine, public comm_vars8_init(varname, var, vid, gid)
Register variables.
Definition: scale_comm_cartesC.F90:811
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56