SCALE-RM
scale_atmos_dyn_tinteg_tracer_linrk.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
15 !-------------------------------------------------------------------------------
16 #include "scalelib.h"
18  !-----------------------------------------------------------------------------
19  !
20  !++ used modules
21  !
22  use scale_precision
23  use scale_io
24  use scale_prof
26  use scale_index
27  use scale_tracer
28 
29 #ifdef DEBUG
30  use scale_debug, only: &
31  check
32  use scale_const, only: &
33  undef => const_undef, &
34  iundef => const_undef2
35 #endif
36  !-----------------------------------------------------------------------------
37  implicit none
38  private
39  !-----------------------------------------------------------------------------
40  !
41  !++ Public procedure
42  !
46 
47  !-----------------------------------------------------------------------------
48  !
49  !++ Public parameters & variables
50  !
51  !-----------------------------------------------------------------------------
52  !
53  !++ Private procedure
54  !
55  !-----------------------------------------------------------------------------
56  !
57  !++ Private parameters & variables
58  !
59  real(RP), allocatable :: QTRC_RK_list(:,:,:,:)
60  integer, private :: I_COMM_RK_list(2)
61  integer, private :: nstage = 3
62 
63  !-----------------------------------------------------------------------------
64 contains
65 
66  !-----------------------------------------------------------------------------
69  tinteg_type )
70  use scale_prc, only: &
71  prc_abort
72  use scale_comm_cartesc, only: &
74  implicit none
75 
76  character(len=*) :: tinteg_type
77 
78  !---------------------------------------------------------------------------
79 
80  if ( tinteg_type(1:5) /= 'LINRK' .or. tinteg_type(7:7) /= 's' ) then
81  log_error("ATMOS_DYN_Tinteg_tracer_linrkNs_setup",*) 'TINTEG_TRACER_TYPE is invalid. Check!'
82  call prc_abort
83  end if
84  read(tinteg_type(6:6),*) nstage
85 
86  allocate( qtrc_rk_list(ka,ia,ja,2) )
87  i_comm_rk_list(:) = (/ 1, 2 /)
88  call comm_vars8_init( 'QTRC_RK1', qtrc_rk_list(:,:,:,1), i_comm_rk_list(1) )
89  call comm_vars8_init( 'QTRC_RK2', qtrc_rk_list(:,:,:,2), i_comm_rk_list(2) )
90 
91  return
93 
94  !-----------------------------------------------------------------------------
97 
98  deallocate( qtrc_rk_list )
99 
100  return
102 
103  !-----------------------------------------------------------------------------
105  subroutine atmos_dyn_tinteg_tracer_linrk( &
106  QTRC, qflx, & ! (out)
107  qtrc0, rhoq_t, dens0, dens, & ! (in)
108  mflx_hi, num_diff, & ! (in)
109  gsqrt, mapf, & ! (in)
110  cdz, rcdz, rcdx, rcdy, & ! (in)
111  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
112  twod, & ! (in)
113  dtl, & ! (in)
114  flag_fct_tracer, flag_fct_along_stream ) ! (in)
115 
116  use scale_comm_cartesc, only: &
117  comm_vars8, &
118  comm_wait
119  use scale_atmos_dyn_tstep_tracer, only: &
121  use scale_atmos_dyn_common, only: &
123  implicit none
124  real(rp), intent(inout) :: qtrc (ka,ia,ja)
125  real(rp), intent(out) :: qflx (ka,ia,ja,3)
126  real(rp), intent(in) :: qtrc0 (ka,ia,ja)
127  real(rp), intent(in) :: rhoq_t (ka,ia,ja)
128  real(rp), intent(in) :: dens0 (ka,ia,ja)
129  real(rp), intent(in) :: dens (ka,ia,ja)
130  real(rp), intent(in) :: mflx_hi (ka,ia,ja,3)
131  real(rp), intent(in) :: num_diff(ka,ia,ja,3)
132  real(rp), intent(in) :: gsqrt (ka,ia,ja,7)
133  real(rp), intent(in) :: mapf (ia,ja)
134  real(rp), intent(in) :: cdz(ka)
135  real(rp), intent(in) :: rcdz(ka)
136  real(rp), intent(in) :: rcdx(ia)
137  real(rp), intent(in) :: rcdy(ja)
138  logical, intent(in) :: bnd_w
139  logical, intent(in) :: bnd_e
140  logical, intent(in) :: bnd_s
141  logical, intent(in) :: bnd_n
142  logical, intent(in) :: twod
143  real(rp), intent(in) :: dtl
144  logical, intent(in) :: flag_fct_tracer
145  logical, intent(in) :: flag_fct_along_stream
146 
147  real(rp) :: dens_rk(ka,ia,ja)
148  integer :: k, i, j
149  integer :: nowstage
150  real(rp) :: linrk_coef
151  integer :: i_in, i_out, i_tmp
152 
153  !------------------------------------------------------------------------
154  ! Start RK
155  !------------------------------------------------------------------------
156 
157  i_in = 1; i_out = 2
158 
159  !$omp parallel do collapse(2) private(k)
160  do j = js, je
161  do i = is, ie
162  do k = ks, ke
163  qtrc_rk_list(k,i,j,i_in) = qtrc0(k,i,j)
164  end do
165  end do
166  end do
167 
168  do nowstage=1, nstage
169  linrk_coef = 1.0_rp/real(nstage - nowstage + 1, kind=rp)
170 
171  call atmos_dyn_copy_boundary_tracer( qtrc_rk_list(:,:,:,i_in), & ! [INOUT]
172  qtrc0, & ! [IN]
173  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
174  twod ) ! [IN]
175 
176  call comm_vars8( qtrc_rk_list(:,:,:,i_in), i_comm_rk_list(i_in) )
177  call comm_wait ( qtrc_rk_list(:,:,:,i_in), i_comm_rk_list(i_in), .false. )
178 
179  if (nowstage < nstage) then
180  !$omp parallel do collapse(2) private(k)
181  do j = js-1, je+1
182  do i = max(is-1,1), min(ie+1,ia)
183  do k = ks, ke
184  dens_rk(k,i,j) = dens0(k,i,j) &
185  + ( dens(k,i,j) - dens0(k,i,j) ) * linrk_coef
186  end do
187  end do
188  end do
189 
190  call atmos_dyn_tstep_tracer( &
191  qtrc_rk_list(:,:,:,i_out), qflx, & ! (out)
192  qtrc_rk_list(:,:,:,i_in), qtrc0, rhoq_t, & ! (in)
193  dens0, dens_rk, & ! (in)
194  mflx_hi, num_diff, & ! (in)
195  gsqrt, mapf, & ! (in)
196  cdz, rcdz, rcdx, rcdy, & ! (in)
197  twod, dtl*linrk_coef, & ! (in)
198  .false., flag_fct_along_stream ) ! (in)
199 
200  else ! final stage
201  call atmos_dyn_tstep_tracer( &
202  qtrc, qflx, & ! (out)
203  qtrc_rk_list(:,:,:,i_in), qtrc0, rhoq_t, & ! (in)
204  dens0, dens, & ! (in)
205  mflx_hi, num_diff, & ! (in)
206  gsqrt, mapf, & ! (in)
207  cdz, rcdz, rcdx, rcdy, & ! (in)
208  twod, dtl, & ! (in)
209  flag_fct_tracer, flag_fct_along_stream ) ! (in)
210 
211  end if
212 
213  i_tmp = i_out
214  i_out = i_in; i_in = i_tmp
215  end do
216 
217  return
218  end subroutine atmos_dyn_tinteg_tracer_linrk
219 
scale_atmos_dyn_tinteg_tracer_linrk
module Atmosphere / Dyn Tinteg
Definition: scale_atmos_dyn_tinteg_tracer_linrk.F90:17
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_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_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_linrk::atmos_dyn_tinteg_tracer_linrk
subroutine, public atmos_dyn_tinteg_tracer_linrk(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)
linear case RK
Definition: scale_atmos_dyn_tinteg_tracer_linrk.F90:115
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_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_dyn_tinteg_tracer_linrk::atmos_dyn_tinteg_tracer_linrk_setup
subroutine, public atmos_dyn_tinteg_tracer_linrk_setup(tinteg_type)
Setup.
Definition: scale_atmos_dyn_tinteg_tracer_linrk.F90:70
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
scale_atmos_dyn_tinteg_tracer_linrk::atmos_dyn_tinteg_tracer_linrk_finalize
subroutine, public atmos_dyn_tinteg_tracer_linrk_finalize
finalize
Definition: scale_atmos_dyn_tinteg_tracer_linrk.F90:97