SCALE-RM
Functions/Subroutines
scale_atmos_dyn_tinteg_tracer_linrk Module Reference

module Atmosphere / Dyn Tinteg More...

Functions/Subroutines

subroutine, public atmos_dyn_tinteg_tracer_linrk_setup (tinteg_type)
 Setup. More...
 
subroutine, public atmos_dyn_tinteg_tracer_linrk_finalize
 finalize More...
 
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 More...
 

Detailed Description

module Atmosphere / Dyn Tinteg

Description
Temporal integration for tracer advection for Atmospheric process 'Linear case' Runge-Kutta scheme with N stage

For a linear, homogeneous and time-indepent ODE system defined by equation (9) in Baldauf (2008, JCP), the 'Linear case' Runge-Kutta scheme with N stage has Nth order accuracy.

The 'Linear case' Runge-Kutta scheme with 3 stage corresponds to the 3-stage Runge-Kutta scheme in Wicker and Skamarock (2002).

Author
Team SCALE

Function/Subroutine Documentation

◆ atmos_dyn_tinteg_tracer_linrk_setup()

subroutine, public scale_atmos_dyn_tinteg_tracer_linrk::atmos_dyn_tinteg_tracer_linrk_setup ( character(len=*)  tinteg_type)

Setup.

Definition at line 70 of file scale_atmos_dyn_tinteg_tracer_linrk.F90.

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

References scale_comm_cartesc::comm_vars8_init(), scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::ka, and scale_prc::prc_abort().

Referenced by scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_setup().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_dyn_tinteg_tracer_linrk_finalize()

subroutine, public scale_atmos_dyn_tinteg_tracer_linrk::atmos_dyn_tinteg_tracer_linrk_finalize

finalize

Definition at line 97 of file scale_atmos_dyn_tinteg_tracer_linrk.F90.

97 
98  deallocate( qtrc_rk_list )
99 
100  return

Referenced by scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_setup().

Here is the caller graph for this function:

◆ atmos_dyn_tinteg_tracer_linrk()

subroutine, public scale_atmos_dyn_tinteg_tracer_linrk::atmos_dyn_tinteg_tracer_linrk ( real(rp), dimension (ka,ia,ja), intent(inout)  QTRC,
real(rp), dimension (ka,ia,ja,3), intent(out)  qflx,
real(rp), dimension (ka,ia,ja), intent(in)  QTRC0,
real(rp), dimension (ka,ia,ja), intent(in)  RHOQ_t,
real(rp), dimension (ka,ia,ja), intent(in)  DENS0,
real(rp), dimension (ka,ia,ja), intent(in)  DENS,
real(rp), dimension (ka,ia,ja,3), intent(in)  mflx_hi,
real(rp), dimension(ka,ia,ja,3), intent(in)  num_diff,
real(rp), dimension (ka,ia,ja,7), intent(in)  GSQRT,
real(rp), dimension (ia,ja), intent(in)  MAPF,
real(rp), dimension(ka), intent(in)  CDZ,
real(rp), dimension(ka), intent(in)  RCDZ,
real(rp), dimension(ia), intent(in)  RCDX,
real(rp), dimension(ja), intent(in)  RCDY,
logical, intent(in)  BND_W,
logical, intent(in)  BND_E,
logical, intent(in)  BND_S,
logical, intent(in)  BND_N,
logical, intent(in)  TwoD,
real(rp), intent(in)  dtl,
logical, intent(in)  FLAG_FCT_TRACER,
logical, intent(in)  FLAG_FCT_ALONG_STREAM 
)

linear case RK

Definition at line 115 of file scale_atmos_dyn_tinteg_tracer_linrk.F90.

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

References scale_atmos_dyn_common::atmos_dyn_copy_boundary_tracer(), scale_atmos_dyn_tstep_tracer::atmos_dyn_tstep_tracer, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_tracer::k, scale_atmos_grid_cartesc_index::ke, and scale_atmos_grid_cartesc_index::ks.

Referenced by scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_setup().

Here is the call graph for this function:
Here is the caller graph for this function:
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_common
module Atmosphere / Dynamics common
Definition: scale_atmos_dyn_common.F90:12
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
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_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
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_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_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
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_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56