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  !
45 
46  !-----------------------------------------------------------------------------
47  !
48  !++ Public parameters & variables
49  !
50  !-----------------------------------------------------------------------------
51  !
52  !++ Private procedure
53  !
54  !-----------------------------------------------------------------------------
55  !
56  !++ Private parameters & variables
57  !
58  real(RP), allocatable :: QTRC_RK_list(:,:,:,:)
59  integer, private :: I_COMM_RK_list(2)
60  integer, private :: nstage = 3
61 
62  !-----------------------------------------------------------------------------
63 contains
64 
65  !-----------------------------------------------------------------------------
68  tinteg_type )
69  use scale_prc, only: &
70  prc_abort
71  use scale_comm_cartesc, only: &
73  implicit none
74 
75  character(len=*) :: tinteg_type
76 
77  !---------------------------------------------------------------------------
78 
79  if ( tinteg_type(1:5) /= 'LINRK' .or. tinteg_type(7:7) /= 's' ) then
80  log_error("ATMOS_DYN_Tinteg_tracer_linrkNs_setup",*) 'TINTEG_TRACER_TYPE is invalid. Check!'
81  call prc_abort
82  end if
83  read(tinteg_type(6:6),*) nstage
84 
85  allocate( qtrc_rk_list(ka,ia,ja,2) )
86  i_comm_rk_list(:) = (/ 1, 2 /)
87  call comm_vars8_init( 'QTRC_RK1', qtrc_rk_list(:,:,:,1), i_comm_rk_list(1) )
88  call comm_vars8_init( 'QTRC_RK2', qtrc_rk_list(:,:,:,2), i_comm_rk_list(2) )
89 
90  return
92 
93  !-----------------------------------------------------------------------------
95  subroutine atmos_dyn_tinteg_tracer_linrk( &
96  QTRC, qflx, & ! (out)
97  qtrc0, rhoq_t, dens0, dens, & ! (in)
98  mflx_hi, num_diff, & ! (in)
99  gsqrt, mapf, & ! (in)
100  cdz, rcdz, rcdx, rcdy, & ! (in)
101  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
102  twod, & ! (in)
103  dtl, & ! (in)
104  flag_fct_tracer, flag_fct_along_stream ) ! (in)
105 
106  use scale_comm_cartesc, 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(out) :: qflx (ka,ia,ja,3)
116  real(rp), intent(in) :: qtrc0 (ka,ia,ja)
117  real(rp), intent(in) :: rhoq_t (ka,ia,ja)
118  real(rp), intent(in) :: dens0 (ka,ia,ja)
119  real(rp), intent(in) :: dens (ka,ia,ja)
120  real(rp), intent(in) :: mflx_hi (ka,ia,ja,3)
121  real(rp), intent(in) :: num_diff(ka,ia,ja,3)
122  real(rp), intent(in) :: gsqrt (ka,ia,ja,7)
123  real(rp), intent(in) :: mapf (ia,ja)
124  real(rp), intent(in) :: cdz(ka)
125  real(rp), intent(in) :: rcdz(ka)
126  real(rp), intent(in) :: rcdx(ia)
127  real(rp), intent(in) :: rcdy(ja)
128  logical, intent(in) :: bnd_w
129  logical, intent(in) :: bnd_e
130  logical, intent(in) :: bnd_s
131  logical, intent(in) :: bnd_n
132  logical, intent(in) :: twod
133  real(rp), intent(in) :: dtl
134  logical, intent(in) :: flag_fct_tracer
135  logical, intent(in) :: flag_fct_along_stream
136 
137  real(rp) :: dens_rk(ka,ia,ja)
138  integer :: k, i, j
139  integer :: nowstage
140  real(rp) :: linrk_coef
141  integer :: i_in, i_out, i_tmp
142 
143  !------------------------------------------------------------------------
144  ! Start RK
145  !------------------------------------------------------------------------
146 
147  i_in = 1; i_out = 2
148 
149  !$omp parallel do collapse(2) private(k)
150  do j = js, je
151  do i = is, ie
152  do k = ks, ke
153  qtrc_rk_list(k,i,j,i_in) = qtrc0(k,i,j)
154  end do
155  end do
156  end do
157 
158  do nowstage=1, nstage
159  linrk_coef = 1.0_rp/real(nstage - nowstage + 1, kind=rp)
160 
161  call atmos_dyn_copy_boundary_tracer( qtrc_rk_list(:,:,:,i_in), & ! [INOUT]
162  qtrc0, & ! [IN]
163  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
164  twod ) ! [IN]
165 
166  call comm_vars8( qtrc_rk_list(:,:,:,i_in), i_comm_rk_list(i_in) )
167  call comm_wait ( qtrc_rk_list(:,:,:,i_in), i_comm_rk_list(i_in), .false. )
168 
169  if (nowstage < nstage) then
170  !$omp parallel do collapse(2) private(k)
171  do j = js-1, je+1
172  do i = max(is-1,1), min(ie+1,ia)
173  do k = ks, ke
174  dens_rk(k,i,j) = dens0(k,i,j) &
175  + ( dens(k,i,j) - dens0(k,i,j) ) * linrk_coef
176  end do
177  end do
178  end do
179 
180  call atmos_dyn_tstep_tracer( &
181  qtrc_rk_list(:,:,:,i_out), qflx, & ! (out)
182  qtrc_rk_list(:,:,:,i_in), qtrc0, rhoq_t, & ! (in)
183  dens0, dens_rk, & ! (in)
184  mflx_hi, num_diff, & ! (in)
185  gsqrt, mapf, & ! (in)
186  cdz, rcdz, rcdx, rcdy, & ! (in)
187  twod, dtl*linrk_coef, & ! (in)
188  .false., flag_fct_along_stream ) ! (in)
189 
190  else ! final stage
191  call atmos_dyn_tstep_tracer( &
192  qtrc, qflx, & ! (out)
193  qtrc_rk_list(:,:,:,i_in), qtrc0, rhoq_t, & ! (in)
194  dens0, dens, & ! (in)
195  mflx_hi, num_diff, & ! (in)
196  gsqrt, mapf, & ! (in)
197  cdz, rcdz, rcdx, rcdy, & ! (in)
198  twod, dtl, & ! (in)
199  flag_fct_tracer, flag_fct_along_stream ) ! (in)
200 
201  end if
202 
203  i_tmp = i_out
204  i_out = i_in; i_in = i_tmp
205  end do
206 
207  return
208  end subroutine atmos_dyn_tinteg_tracer_linrk
209 
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:342
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:38
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:105
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:44
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_comm_cartesc::comm_vars8_init
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
Definition: scale_comm_cartesC.F90:294
scale_debug::check
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:56
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:311
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:69
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_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56