SCALE-RM
Functions/Subroutines
scale_atmos_dyn_tinteg_tracer_rk3 Module Reference

module Atmosphere / Dyn Tinteg More...

Functions/Subroutines

subroutine, public atmos_dyn_tinteg_tracer_rk3_setup (tinteg_type)
 Setup. More...
 
subroutine, public atmos_dyn_tinteg_tracer_rk3 (QTRC, QTRC0, RHOQ_t, DENS0, DENS, mflx_hi, num_diff, GSQRT, MAPF, CDZ, RCDZ, RCDX, RCDY, BND_W, BND_E, BND_S, BND_N, dtl, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM)
 RK3. More...
 

Detailed Description

module Atmosphere / Dyn Tinteg

Description
Temporal integration for tracer advection for Atmospheric process three step Runge-Kutta scheme
Author
Team SCALE

Function/Subroutine Documentation

◆ atmos_dyn_tinteg_tracer_rk3_setup()

subroutine, public scale_atmos_dyn_tinteg_tracer_rk3::atmos_dyn_tinteg_tracer_rk3_setup ( character(len=*)  tinteg_type)

Setup.

Definition at line 65 of file scale_atmos_dyn_tinteg_tracer_rk3.F90.

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

65  use scale_prc, only: &
66  prc_abort
67  use scale_comm_cartesc, only: &
69  implicit none
70 
71  character(len=*) :: tinteg_type
72 
73  integer :: iv
74  !---------------------------------------------------------------------------
75 
76  if ( tinteg_type /= 'RK3WS2002' ) then
77  log_error("ATMOS_DYN_Tinteg_tracer_rk3_setup",*) 'TINTEG_LARGE_TYPE is not RK3WS2002. Check!'
78  call prc_abort
79  end if
80 
81  allocate( qtrc_rk1(ka,ia,ja) )
82  allocate( qtrc_rk2(ka,ia,ja) )
83 
84  call comm_vars8_init( 'QTRC_RK1', qtrc_rk1, i_comm_rk1 )
85  call comm_vars8_init( 'QTRC_RK2', qtrc_rk2, i_comm_rk2 )
86 
87  return
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
module COMMUNICATION
module PROCESS
Definition: scale_prc.F90:11
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
integer, public ka
of whole cells: z, local, with HALO
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_dyn_tinteg_tracer_rk3()

subroutine, public scale_atmos_dyn_tinteg_tracer_rk3::atmos_dyn_tinteg_tracer_rk3 ( real(rp), dimension (ka,ia,ja), intent(inout)  QTRC,
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,
real(rp), intent(in)  dtl,
logical, intent(in)  FLAG_FCT_TRACER,
logical, intent(in)  FLAG_FCT_ALONG_STREAM 
)

RK3.

Definition at line 103 of file scale_atmos_dyn_tinteg_tracer_rk3.F90.

References scale_atmos_dyn_common::atmos_dyn_copy_boundary_tracer(), scale_atmos_dyn_tstep_tracer::atmos_dyn_tstep_tracer, 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_atmos_grid_cartesc_index::ke, and scale_atmos_grid_cartesc_index::ks.

Referenced by scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_setup().

103  use scale_comm_cartesc, only: &
104  comm_vars8, &
105  comm_wait
106  use scale_atmos_dyn_tstep_tracer, only: &
108  use scale_atmos_dyn_common, only: &
110  implicit none
111  real(RP), intent(inout) :: qtrc (ka,ia,ja)
112  real(RP), intent(in) :: qtrc0 (ka,ia,ja)
113  real(RP), intent(in) :: rhoq_t (ka,ia,ja)
114  real(RP), intent(in) :: dens0 (ka,ia,ja)
115  real(RP), intent(in) :: dens (ka,ia,ja)
116  real(RP), intent(in) :: mflx_hi (ka,ia,ja,3)
117  real(RP), intent(in) :: num_diff(ka,ia,ja,3)
118  real(RP), intent(in) :: gsqrt (ka,ia,ja,7)
119  real(RP), intent(in) :: mapf (ia,ja)
120  real(RP), intent(in) :: cdz(ka)
121  real(RP), intent(in) :: rcdz(ka)
122  real(RP), intent(in) :: rcdx(ia)
123  real(RP), intent(in) :: rcdy(ja)
124  logical, intent(in) :: bnd_w
125  logical, intent(in) :: bnd_e
126  logical, intent(in) :: bnd_s
127  logical, intent(in) :: bnd_n
128  real(RP), intent(in) :: dtl
129  logical, intent(in) :: flag_fct_tracer
130  logical, intent(in) :: flag_fct_along_stream
131 
132  real(RP) :: dens_rk(ka,ia,ja)
133  real(RP) :: dtrk
134  integer :: k, i, j
135 
136  !------------------------------------------------------------------------
137  ! Start RK
138  !------------------------------------------------------------------------
139 
140  do j = js-1, je+1
141  do i = is-1, ie+1
142  do k = ks, ke
143  dens_rk(k,i,j) = dens0(k,i,j) &
144  + ( dens(k,i,j) - dens0(k,i,j) ) / 3.0_rp
145  end do
146  end do
147  end do
148 
149  dtrk = dtl / 3.0_rp
150  call atmos_dyn_tstep_tracer( &
151  qtrc_rk1, & ! (out)
152  qtrc, qtrc0, rhoq_t, &! (in)
153  dens0, dens_rk, & ! (in)
154  mflx_hi, num_diff, & ! (in)
155  gsqrt, mapf, & ! (in)
156  cdz, rcdz, rcdx, rcdy, & ! (in)
157  dtrk, & ! (in)
158  .false., flag_fct_along_stream ) ! (in)
159 
160  call atmos_dyn_copy_boundary_tracer( qtrc_rk1, & ! [INOUT]
161  qtrc0, & ! [IN]
162  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
163 
164  call comm_vars8( qtrc_rk1(:,:,:), i_comm_rk1 )
165  call comm_wait ( qtrc_rk1(:,:,:), i_comm_rk1, .false. )
166 
167 
168  do j = js-1, je+1
169  do i = is-1, ie+1
170  do k = ks, ke
171  dens_rk(k,i,j) = dens0(k,i,j) &
172  + ( dens(k,i,j) - dens0(k,i,j) ) * 0.5_rp
173  end do
174  end do
175  end do
176 
177  dtrk = dtl / 2.0_rp
178  call atmos_dyn_tstep_tracer( &
179  qtrc_rk2, & ! (out)
180  qtrc_rk1, qtrc0, rhoq_t, &! (in)
181  dens0, dens_rk, & ! (in)
182  mflx_hi, num_diff, & ! (in)
183  gsqrt, mapf, & ! (in)
184  cdz, rcdz, rcdx, rcdy, & ! (in)
185  dtrk, & ! (in)
186  .false., flag_fct_along_stream ) ! (in)
187 
188  call atmos_dyn_copy_boundary_tracer( qtrc_rk2, & ! [INOUT]
189  qtrc0, & ! [IN]
190  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
191 
192  call comm_vars8( qtrc_rk2(:,:,:), i_comm_rk2 )
193  call comm_wait ( qtrc_rk2(:,:,:), i_comm_rk2, .false. )
194 
195 
196  dtrk = dtl
197  call atmos_dyn_tstep_tracer( &
198  qtrc, & ! (out)
199  qtrc_rk2, qtrc0, rhoq_t, &! (in)
200  dens0, dens, & ! (in)
201  mflx_hi, num_diff, & ! (in)
202  gsqrt, mapf, & ! (in)
203  cdz, rcdz, rcdx, rcdy, & ! (in)
204  dtrk, & ! (in)
205  flag_fct_tracer, flag_fct_along_stream ) ! (in)
206 
207  return
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
module COMMUNICATION
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
integer, public ke
end point of inner domain: z, local
integer, public je
end point of inner domain: y, local
module Atmosphere / Dynamics common
integer, public ks
start point of inner domain: z, local
integer, public js
start point of inner domain: y, local
integer, public ka
of whole cells: z, local, with HALO
procedure(step), pointer, public atmos_dyn_tstep_tracer
subroutine, public atmos_dyn_copy_boundary_tracer(QTRC, QTRC0, BND_W, BND_E, BND_S, BND_N)
module Atmosphere / Dynamical scheme
Here is the call graph for this function:
Here is the caller graph for this function: