SCALE-RM
Functions/Subroutines
scale_atmos_dyn_tstep_tracer_fvm_heve Module Reference

module Atmosphere / Dynamics More...

Functions/Subroutines

subroutine, public atmos_dyn_tstep_tracer_fvm_heve_setup (type)
 Setup. More...
 
subroutine, public atmos_dyn_tstep_tracer_fvm_heve (QTRCo, QTRC, QTRC0, RHOQ_t, DENS0, DENS, mflx_hi, num_diff, GSQRT, MAPF, CDZ, RCDZ, RCDX, RCDY, dtl, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM)
 

Detailed Description

module Atmosphere / Dynamics

Description
HEVE FVM scheme for tracer advection in Atmospheric dynamical process
Author
Team SCALE
History
  • 2016-05-17 (S.Nishizawa) [new]

Function/Subroutine Documentation

◆ atmos_dyn_tstep_tracer_fvm_heve_setup()

subroutine, public scale_atmos_dyn_tstep_tracer_fvm_heve::atmos_dyn_tstep_tracer_fvm_heve_setup ( character(len=*), intent(in)  type)

Setup.

Definition at line 62 of file scale_atmos_dyn_tstep_tracer_fvm_heve.F90.

References scale_process::prc_mpistop().

Referenced by scale_atmos_dyn_tstep_tracer::atmos_dyn_tstep_tracer_setup().

62  use scale_process, only: &
64  implicit none
65  character(len=*), intent(in) :: type
66 
67  if ( type /= 'FVM-HEVE' ) then
68  write(*,*) 'xxx Tstep_tracer_type is not "FVM-HEVE"!'
69  call prc_mpistop
70  end if
71 
72  return
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_dyn_tstep_tracer_fvm_heve()

subroutine, public scale_atmos_dyn_tstep_tracer_fvm_heve::atmos_dyn_tstep_tracer_fvm_heve ( real(rp), dimension (ka,ia,ja), intent(inout)  QTRCo,
real(rp), dimension (ka,ia,ja), intent(in)  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,2), 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 
)

Definition at line 85 of file scale_atmos_dyn_tstep_tracer_fvm_heve.F90.

References scale_atmos_dyn_common::atmos_dyn_fct(), scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxx_xyz_tracer, scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxx_xyz_ud1(), scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxy_xyz_tracer, scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxy_xyz_ud1(), scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxz_xyz_tracer, scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxz_xyz_ud1(), scale_gridtrans::i_uy, scale_gridtrans::i_uyz, scale_gridtrans::i_xv, scale_gridtrans::i_xvz, scale_gridtrans::i_xy, scale_gridtrans::i_xyw, scale_gridtrans::i_xyz, scale_grid_index::iblock, scale_grid_index::ie, scale_grid_index::is, scale_grid_index::jblock, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ke, scale_grid_index::ks, scale_grid_index::xdir, scale_grid_index::ydir, and scale_grid_index::zdir.

Referenced by scale_atmos_dyn_tstep_tracer::atmos_dyn_tstep_tracer_setup().

85  use scale_gridtrans, only: &
86  i_xyz, &
87  i_xyw, &
88  i_uyz, &
89  i_xvz, &
90  i_xy, &
91  i_uy, &
92  i_xv
93  use scale_atmos_dyn_common, only: &
95  use scale_atmos_dyn_fvm_flux, only: &
99  use scale_atmos_dyn_fvm_flux_ud1, only: &
103  implicit none
104  real(RP), intent(inout) :: QTRCo (KA,IA,JA) ! could be identical to QTRC0
105  real(RP), intent(in) :: QTRC (KA,IA,JA)
106  real(RP), intent(in) :: QTRC0 (KA,IA,JA)
107  real(RP), intent(in) :: RHOQ_t (KA,IA,JA)
108  real(RP), intent(in) :: DENS0 (KA,IA,JA)
109  real(RP), intent(in) :: DENS (KA,IA,JA)
110  real(RP), intent(in) :: mflx_hi (KA,IA,JA,3)
111  real(RP), intent(in) :: num_diff(KA,IA,JA,3)
112  real(RP), intent(in) :: GSQRT (KA,IA,JA,7)
113  real(RP), intent(in) :: MAPF (IA,JA,2)
114  real(RP), intent(in) :: CDZ(KA)
115  real(RP), intent(in) :: RCDZ(KA)
116  real(RP), intent(in) :: RCDX(IA)
117  real(RP), intent(in) :: RCDY(JA)
118  real(RP), intent(in) :: dtl
119  logical, intent(in) :: FLAG_FCT_TRACER
120  logical, intent(in) :: FLAG_FCT_ALONG_STREAM
121 
122 
123  ! For tracer advection
124  real(RP) :: qflx_hi (KA,IA,JA,3) ! rho * vel(x,y,z) * phi @ (u,v,w)-face high order
125  real(RP) :: qflx_lo (KA,IA,JA,3) ! rho * vel(x,y,z) * phi, monotone flux
126  real(RP) :: qflx_anti(KA,IA,JA,3) ! anti-diffusive flux
127 
128  integer :: IIS, IIE
129  integer :: JJS, JJE
130  integer :: i, j, k
131  !---------------------------------------------------------------------------
132 
133 #ifdef DEBUG
134  qflx_hi(:,:,:,:) = undef
135  qflx_lo(:,:,:,:) = undef
136 #endif
137 
138  do jjs = js, je, jblock
139  jje = jjs+jblock-1
140  do iis = is, ie, iblock
141  iie = iis+iblock-1
142 
143  ! at (x, y, w)
144  call atmos_dyn_fvm_fluxz_xyz_tracer( qflx_hi(:,:,:,zdir), & ! (out)
145  mflx_hi(:,:,:,zdir), qtrc, gsqrt(:,:,:,i_xyw), & ! (in)
146  num_diff(:,:,:,zdir), & ! (in)
147  cdz, & ! (in)
148  iis, iie, jjs, jje ) ! (in)
149 
150  ! at (u, y, z)
151  call atmos_dyn_fvm_fluxx_xyz_tracer( qflx_hi(:,:,:,xdir), & ! (out)
152  mflx_hi(:,:,:,xdir), qtrc, gsqrt(:,:,:,i_uyz), & ! (in)
153  num_diff(:,:,:,xdir), & ! (in)
154  cdz, & ! (in)
155  iis, iie, jjs, jje ) ! (in)
156 
157  ! at (x, v, z)
158  call atmos_dyn_fvm_fluxy_xyz_tracer( qflx_hi(:,:,:,ydir), & ! (out)
159  mflx_hi(:,:,:,ydir), qtrc, gsqrt(:,:,:,i_xvz), & ! (in)
160  num_diff(:,:,:,ydir), & ! (in)
161  cdz, & ! (in)
162  iis, iie, jjs, jje ) ! (in)
163 
164  if ( flag_fct_tracer ) then
165 
166  call atmos_dyn_fvm_fluxz_xyz_ud1( qflx_lo(:,:,:,zdir), & ! (out)
167  mflx_hi(:,:,:,zdir), qtrc0, gsqrt(:,:,:,i_xyw), & ! (in)
168  num_diff(:,:,:,zdir), & ! (in)
169  cdz, & ! (in)
170  iis-1, iie+1, jjs-1, jje+1 ) ! (in)
171 
172  call atmos_dyn_fvm_fluxx_xyz_ud1( qflx_lo(:,:,:,xdir), & ! (out)
173  mflx_hi(:,:,:,xdir), qtrc0, gsqrt(:,:,:,i_uyz), & ! (in)
174  num_diff(:,:,:,xdir), & ! (in)
175  cdz, & ! (in)
176  iis-1, iie+1, jjs-1, jje+1 ) ! (in)
177 
178  call atmos_dyn_fvm_fluxy_xyz_ud1( qflx_lo(:,:,:,ydir), & ! (out)
179  mflx_hi(:,:,:,ydir), qtrc0, gsqrt(:,:,:,i_xvz), & ! (in)
180  num_diff(:,:,:,ydir), & ! (in)
181  cdz, & ! (in)
182  iis-1, iie+1, jjs-1, jje+1 ) ! (in)
183  end if
184 
185  enddo
186  enddo
187 
188  if ( flag_fct_tracer ) then
189 
190  call atmos_dyn_fct( qflx_anti, & ! (out)
191  qtrc0, dens0, dens, & ! (in)
192  qflx_hi, qflx_lo, & ! (in)
193  mflx_hi, & ! (in)
194  rcdz, rcdx, rcdy, & ! (in)
195  gsqrt(:,:,:,i_xyz), & ! (in)
196  mapf, dtl, & ! (in)
197  flag_fct_along_stream ) ! (in)
198 
199  do jjs = js, je, jblock
200  jje = jjs+jblock-1
201  do iis = is, ie, iblock
202  iie = iis+iblock-1
203 
204  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
205  do j = jjs, jje
206  do i = iis, iie
207  do k = ks, ke
208  qtrco(k,i,j) = ( qtrc0(k,i,j) * dens0(k,i,j) &
209  + dtl * ( - ( ( qflx_hi(k ,i ,j ,zdir) - qflx_anti(k ,i ,j ,zdir) &
210  - qflx_hi(k-1,i ,j ,zdir) + qflx_anti(k-1,i ,j ,zdir) ) * rcdz(k) &
211  + ( qflx_hi(k ,i ,j ,xdir) - qflx_anti(k ,i ,j ,xdir) &
212  - qflx_hi(k ,i-1,j ,xdir) + qflx_anti(k ,i-1,j ,xdir) ) * rcdx(i) &
213  + ( qflx_hi(k ,i ,j ,ydir) - qflx_anti(k ,i ,j ,ydir) &
214  - qflx_hi(k ,i ,j-1,ydir) + qflx_anti(k ,i ,j-1,ydir) ) * rcdy(j) &
215  ) * mapf(i,j,1) * mapf(i,j,2) / gsqrt(k,i,j,i_xyz) &
216  + rhoq_t(k,i,j) ) ) / dens(k,i,j)
217  enddo
218  enddo
219  enddo
220 
221  enddo
222  enddo
223 
224  else ! skip FCT
225 
226  do jjs = js, je, jblock
227  jje = jjs+jblock-1
228  do iis = is, ie, iblock
229  iie = iis+iblock-1
230 
231  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
232  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,QTRCo,QTRC0,DENS0,dtl,qflx_hi,RCDZ,RCDX,RCDY,MAPF) &
233  !$omp shared(GSQRT,RHOQ_t,DENS,I_XYZ)
234  do j = jjs, jje
235  do i = iis, iie
236  do k = ks, ke
237  qtrco(k,i,j) = ( qtrc0(k,i,j) * dens0(k,i,j) &
238  + dtl * ( - ( ( qflx_hi(k,i,j,zdir) - qflx_hi(k-1,i ,j ,zdir) ) * rcdz(k) &
239  + ( qflx_hi(k,i,j,xdir) - qflx_hi(k ,i-1,j ,xdir) ) * rcdx(i) &
240  + ( qflx_hi(k,i,j,ydir) - qflx_hi(k ,i ,j-1,ydir) ) * rcdy(j) &
241  ) * mapf(i,j,1) * mapf(i,j,2) / gsqrt(k,i,j,i_xyz) &
242  + rhoq_t(k,i,j) ) ) / dens(k,i,j)
243  enddo
244  enddo
245  enddo
246 
247  enddo
248  enddo
249 
250  end if
251 
252  return
integer, public is
start point of inner domain: x, local
integer, public i_xvz
integer, public je
end point of inner domain: y, local
subroutine, public atmos_dyn_fvm_fluxy_xyz_ud1(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XYZ
subroutine, public atmos_dyn_fvm_fluxz_xyz_ud1(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at XYZ
integer, public iblock
block size for cache blocking: x
integer, parameter, public zdir
integer, parameter, public ydir
integer, public ke
end point of inner domain: z, local
integer, parameter, public xdir
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxy_xyz_tracer
integer, public i_xy
module GRIDTRANS
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxx_xyz_tracer
integer, public i_uy
integer, public jblock
block size for cache blocking: y
integer, public i_xyw
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxz_xyz_tracer
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 i_xv
integer, public ie
end point of inner domain: x, local
module scale_atmos_dyn_fvm_flux
integer, public i_uyz
subroutine, public atmos_dyn_fct(qflx_anti, phi_in, DENS0, DENS, qflx_hi, qflx_lo, mflx_hi, rdz, rdx, rdy, GSQRT, MAPF, dt, flag_vect)
Flux Correction Transport Limiter.
integer, public i_xyz
subroutine, public atmos_dyn_fvm_fluxx_xyz_ud1(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XYZ
module scale_atmos_dyn_fvm_flux_ud1
Here is the call graph for this function:
Here is the caller graph for this function: