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, 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_const::const_undef, 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_const, only: &
71  undef => const_undef
72  use scale_comm, only: &
74  implicit none
75 
76  character(len=*) :: tinteg_type
77 
78  integer :: iv
79  !---------------------------------------------------------------------------
80 
81  if ( tinteg_type /= 'RK3WS2002' ) then
82  write(*,*) 'xxx TINTEG_LARGE_TYPE is not RK3WS2002. Check!'
83  call prc_mpistop
84  end if
85 
86  allocate( qtrc_rk1(ka,ia,ja) )
87  allocate( qtrc_rk2(ka,ia,ja) )
88 
89  call comm_vars8_init( qtrc_rk1, i_comm_rk1 )
90  call comm_vars8_init( qtrc_rk2, i_comm_rk2 )
91 
92  return
subroutine, public prc_mpistop
Abort MPI.
real(rp), public const_undef
Definition: scale_const.F90:43
integer, public ia
of x whole cells (local, with HALO)
subroutine, public comm_vars8_init(var, vid)
Register variables.
Definition: scale_comm.F90:328
integer, public ka
of z whole cells (local, with HALO)
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, public ja
of y whole cells (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,
real(rp), intent(in)  dtl,
logical, intent(in)  FLAG_FCT_TRACER,
logical, intent(in)  FLAG_FCT_ALONG_STREAM 
)

RK3.

Definition at line 107 of file scale_atmos_dyn_tinteg_tracer_rk3.F90.

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

107  use scale_comm, only: &
108  comm_vars8, &
109  comm_wait
110  use scale_atmos_dyn_tstep_tracer, only: &
112  implicit none
113  real(RP), intent(inout) :: qtrc (ka,ia,ja)
114  real(RP), intent(in) :: qtrc0 (ka,ia,ja)
115  real(RP), intent(in) :: rhoq_t (ka,ia,ja)
116  real(RP), intent(in) :: dens0 (ka,ia,ja)
117  real(RP), intent(in) :: dens (ka,ia,ja)
118  real(RP), intent(in) :: mflx_hi (ka,ia,ja,3)
119  real(RP), intent(in) :: num_diff(ka,ia,ja,3)
120  real(RP), intent(in) :: gsqrt (ka,ia,ja,7)
121  real(RP), intent(in) :: mapf (ia,ja)
122  real(RP), intent(in) :: cdz(ka)
123  real(RP), intent(in) :: rcdz(ka)
124  real(RP), intent(in) :: rcdx(ia)
125  real(RP), intent(in) :: rcdy(ja)
126  real(RP), intent(in) :: dtl
127  logical, intent(in) :: flag_fct_tracer
128  logical, intent(in) :: flag_fct_along_stream
129 
130  real(RP) :: dens_rk(ka,ia,ja)
131  real(RP) :: dtrk
132  integer :: k, i, j
133 
134  !------------------------------------------------------------------------
135  ! Start RK
136  !------------------------------------------------------------------------
137 
138  do j = js-1, je+1
139  do i = is-1, ie+1
140  do k = ks, ke
141  dens_rk(k,i,j) = dens0(k,i,j) &
142  + ( dens(k,i,j) - dens0(k,i,j) ) / 3.0_rp
143  end do
144  end do
145  end do
146 
147  dtrk = dtl / 3.0_rp
148  call atmos_dyn_tstep_tracer( &
149  qtrc_rk1, & ! (out)
150  qtrc, qtrc0, rhoq_t, &! (in)
151  dens0, dens_rk, & ! (in)
152  mflx_hi, num_diff, & ! (in)
153  gsqrt, mapf, & ! (in)
154  cdz, rcdz, rcdx, rcdy, & ! (in)
155  dtrk, & ! (in)
156  .false., flag_fct_along_stream ) ! (in)
157 
158  call comm_vars8( qtrc_rk1(:,:,:), i_comm_rk1 )
159  call comm_wait ( qtrc_rk1(:,:,:), i_comm_rk1, .false. )
160 
161 
162  do j = js-1, je+1
163  do i = is-1, ie+1
164  do k = ks, ke
165  dens_rk(k,i,j) = dens0(k,i,j) &
166  + ( dens(k,i,j) - dens0(k,i,j) ) * 0.5_rp
167  end do
168  end do
169  end do
170 
171  dtrk = dtl / 2.0_rp
172  call atmos_dyn_tstep_tracer( &
173  qtrc_rk2, & ! (out)
174  qtrc_rk1, qtrc0, rhoq_t, &! (in)
175  dens0, dens_rk, & ! (in)
176  mflx_hi, num_diff, & ! (in)
177  gsqrt, mapf, & ! (in)
178  cdz, rcdz, rcdx, rcdy, & ! (in)
179  dtrk, & ! (in)
180  .false., flag_fct_along_stream ) ! (in)
181 
182  call comm_vars8( qtrc_rk2(:,:,:), i_comm_rk2 )
183  call comm_wait ( qtrc_rk2(:,:,:), i_comm_rk2, .false. )
184 
185 
186  dtrk = dtl
187  call atmos_dyn_tstep_tracer( &
188  qtrc, & ! (out)
189  qtrc_rk2, qtrc0, rhoq_t, &! (in)
190  dens0, dens, & ! (in)
191  mflx_hi, num_diff, & ! (in)
192  gsqrt, mapf, & ! (in)
193  cdz, rcdz, rcdx, rcdy, & ! (in)
194  dtrk, & ! (in)
195  flag_fct_tracer, flag_fct_along_stream ) ! (in)
196 
197  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
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
module COMMUNICATION
Definition: scale_comm.F90:23
integer, public js
start point of inner domain: y, local
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
module Atmosphere / Dynamical scheme
integer, public ja
of y whole cells (local, with HALO)
Here is the caller graph for this function: