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

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 59 of file scale_atmos_dyn_tstep_tracer_fvm_heve.F90.

References scale_prc::prc_abort().

Referenced by scale_atmos_dyn_tstep_tracer::atmos_dyn_tstep_tracer_setup().

59  use scale_prc, only: &
60  prc_abort
61  implicit none
62  character(len=*), intent(in) :: type
63 
64  if ( type /= 'FVM-HEVE' ) then
65  log_error("ATMOS_DYN_Tstep_tracer_fvm_heve_setup",*) 'Tstep_tracer_type is not "FVM-HEVE"!'
66  call prc_abort
67  end if
68 
69  return
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
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 82 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_atmos_grid_cartesc_index::i_uyz, scale_atmos_grid_cartesc_index::i_xvz, scale_atmos_grid_cartesc_index::i_xyw, scale_atmos_grid_cartesc_index::i_xyz, scale_atmos_grid_cartesc_index::iblock, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::jblock, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, scale_atmos_grid_cartesc_index::xdir, scale_atmos_grid_cartesc_index::ydir, and scale_atmos_grid_cartesc_index::zdir.

Referenced by scale_atmos_dyn_tstep_tracer::atmos_dyn_tstep_tracer_setup().

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