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
History
  • 2016-05-17 (S.Nishizawa) [new]

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 68 of file scale_atmos_dyn_tinteg_tracer_rk3.F90.

References scale_comm::comm_vars8_init(), scale_grid_index::ia, scale_grid_index::ja, scale_grid_index::ka, and scale_process::prc_mpistop().

Referenced by scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_setup().

68  use scale_process, only: &
70  use scale_comm, only: &
72  implicit none
73 
74  character(len=*) :: tinteg_type
75 
76  integer :: iv
77  !---------------------------------------------------------------------------
78 
79  if ( tinteg_type /= 'RK3WS2002' ) then
80  write(*,*) 'xxx TINTEG_LARGE_TYPE is not RK3WS2002. Check!'
81  call prc_mpistop
82  end if
83 
84  allocate( qtrc_rk1(ka,ia,ja) )
85  allocate( qtrc_rk2(ka,ia,ja) )
86 
87  call comm_vars8_init( 'QTRC_RK1', qtrc_rk1, i_comm_rk1 )
88  call comm_vars8_init( 'QTRC_RK2', qtrc_rk2, i_comm_rk2 )
89 
90  return
subroutine, public prc_mpistop
Abort MPI.
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
Definition: scale_comm.F90:313
integer, public ia
of whole cells: x, local, with HALO
integer, public ka
of whole cells: z, local, with HALO
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
integer, public ja
of whole cells: y, 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 106 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_grid_index::ie, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ke, and scale_grid_index::ks.

Referenced by scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_setup().

106  use scale_comm, 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(in) :: QTRC0 (KA,IA,JA)
116  real(RP), intent(in) :: RHOQ_t (KA,IA,JA)
117  real(RP), intent(in) :: DENS0 (KA,IA,JA)
118  real(RP), intent(in) :: DENS (KA,IA,JA)
119  real(RP), intent(in) :: mflx_hi (KA,IA,JA,3)
120  real(RP), intent(in) :: num_diff(KA,IA,JA,3)
121  real(RP), intent(in) :: GSQRT (KA,IA,JA,7)
122  real(RP), intent(in) :: MAPF (IA,JA)
123  real(RP), intent(in) :: CDZ(KA)
124  real(RP), intent(in) :: RCDZ(KA)
125  real(RP), intent(in) :: RCDX(IA)
126  real(RP), intent(in) :: RCDY(JA)
127  logical, intent(in) :: BND_W
128  logical, intent(in) :: BND_E
129  logical, intent(in) :: BND_S
130  logical, intent(in) :: BND_N
131  real(RP), intent(in) :: dtl
132  logical, intent(in) :: FLAG_FCT_TRACER
133  logical, intent(in) :: FLAG_FCT_ALONG_STREAM
134 
135  real(RP) :: DENS_RK(KA,IA,JA)
136  real(RP) :: dtrk
137  integer :: k, i, j
138 
139  !------------------------------------------------------------------------
140  ! Start RK
141  !------------------------------------------------------------------------
142 
143  do j = js-1, je+1
144  do i = is-1, ie+1
145  do k = ks, ke
146  dens_rk(k,i,j) = dens0(k,i,j) &
147  + ( dens(k,i,j) - dens0(k,i,j) ) / 3.0_rp
148  end do
149  end do
150  end do
151 
152  dtrk = dtl / 3.0_rp
153  call atmos_dyn_tstep_tracer( &
154  qtrc_rk1, & ! (out)
155  qtrc, qtrc0, rhoq_t, &! (in)
156  dens0, dens_rk, & ! (in)
157  mflx_hi, num_diff, & ! (in)
158  gsqrt, mapf, & ! (in)
159  cdz, rcdz, rcdx, rcdy, & ! (in)
160  dtrk, & ! (in)
161  .false., flag_fct_along_stream ) ! (in)
162 
163  call atmos_dyn_copy_boundary_tracer( qtrc_rk1, & ! [INOUT]
164  qtrc0, & ! [IN]
165  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
166 
167  call comm_vars8( qtrc_rk1(:,:,:), i_comm_rk1 )
168  call comm_wait ( qtrc_rk1(:,:,:), i_comm_rk1, .false. )
169 
170 
171  do j = js-1, je+1
172  do i = is-1, ie+1
173  do k = ks, ke
174  dens_rk(k,i,j) = dens0(k,i,j) &
175  + ( dens(k,i,j) - dens0(k,i,j) ) * 0.5_rp
176  end do
177  end do
178  end do
179 
180  dtrk = dtl / 2.0_rp
181  call atmos_dyn_tstep_tracer( &
182  qtrc_rk2, & ! (out)
183  qtrc_rk1, qtrc0, rhoq_t, &! (in)
184  dens0, dens_rk, & ! (in)
185  mflx_hi, num_diff, & ! (in)
186  gsqrt, mapf, & ! (in)
187  cdz, rcdz, rcdx, rcdy, & ! (in)
188  dtrk, & ! (in)
189  .false., flag_fct_along_stream ) ! (in)
190 
191  call atmos_dyn_copy_boundary_tracer( qtrc_rk2, & ! [INOUT]
192  qtrc0, & ! [IN]
193  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
194 
195  call comm_vars8( qtrc_rk2(:,:,:), i_comm_rk2 )
196  call comm_wait ( qtrc_rk2(:,:,:), i_comm_rk2, .false. )
197 
198 
199  dtrk = dtl
200  call atmos_dyn_tstep_tracer( &
201  qtrc, & ! (out)
202  qtrc_rk2, qtrc0, rhoq_t, &! (in)
203  dens0, dens, & ! (in)
204  mflx_hi, num_diff, & ! (in)
205  gsqrt, mapf, & ! (in)
206  cdz, rcdz, rcdx, rcdy, & ! (in)
207  dtrk, & ! (in)
208  flag_fct_tracer, flag_fct_along_stream ) ! (in)
209 
210  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
integer, public ke
end point of inner domain: z, local
module COMMUNICATION
Definition: scale_comm.F90:23
integer, public js
start point of inner domain: y, local
module Atmosphere / Dynamics common
integer, public ks
start point of inner domain: z, local
integer, public ie
end point of inner domain: x, local
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: