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 (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 69 of file scale_atmos_dyn_tinteg_tracer_linrk.F90.

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

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()

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 105 of file scale_atmos_dyn_tinteg_tracer_linrk.F90.

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

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:342
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_comm_cartesc::comm_vars8_init
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
Definition: scale_comm_cartesC.F90:294
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_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_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56