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, qflx_hi, QTRC, QTRC0, RHOQ_t, DENS0, DENS, mflx_hi, num_diff, GSQRT, MAPF, CDZ, RCDZ, RCDX, RCDY, TwoD, 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.

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

References scale_prc::prc_abort().

Referenced by scale_atmos_dyn_tstep_tracer::atmos_dyn_tstep_tracer_setup().

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,3), intent(out)  qflx_hi,
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,
logical, intent(in)  TwoD,
real(rp), intent(in)  dtl,
logical, intent(in)  FLAG_FCT_TRACER,
logical, intent(in)  FLAG_FCT_ALONG_STREAM 
)

Definition at line 84 of file scale_atmos_dyn_tstep_tracer_fvm_heve.F90.

84  use scale_atmos_dyn_fvm_fct, only: &
86  use scale_atmos_dyn_fvm_flux, only: &
90  use scale_atmos_dyn_fvm_flux_ud1, only: &
94  implicit none
95  real(RP), intent(inout) :: QTRCo (KA,IA,JA) ! could be identical to QTRC0
96  real(RP), intent(out) :: qflx_hi (KA,IA,JA,3) ! rho * vel(x,y,z) * phi @ (u,v,w)-face high order
97  real(RP), intent(in) :: QTRC (KA,IA,JA)
98  real(RP), intent(in) :: QTRC0 (KA,IA,JA)
99  real(RP), intent(in) :: RHOQ_t (KA,IA,JA)
100  real(RP), intent(in) :: DENS0 (KA,IA,JA)
101  real(RP), intent(in) :: DENS (KA,IA,JA)
102  real(RP), intent(in) :: mflx_hi (KA,IA,JA,3)
103  real(RP), intent(in) :: num_diff(KA,IA,JA,3)
104  real(RP), intent(in) :: GSQRT (KA,IA,JA,7)
105  real(RP), intent(in) :: MAPF (IA,JA,2)
106  real(RP), intent(in) :: CDZ(KA)
107  real(RP), intent(in) :: RCDZ(KA)
108  real(RP), intent(in) :: RCDX(IA)
109  real(RP), intent(in) :: RCDY(JA)
110  logical, intent(in) :: TwoD
111  real(RP), intent(in) :: dtl
112  logical, intent(in) :: FLAG_FCT_TRACER
113  logical, intent(in) :: FLAG_FCT_ALONG_STREAM
114 
115 
116  ! For tracer advection
117  real(RP) :: qflx_lo (KA,IA,JA,3) ! rho * vel(x,y,z) * phi, monotone flux
118  real(RP) :: qflx_anti(KA,IA,JA,3) ! anti-diffusive flux
119 
120  integer :: IIS, IIE
121  integer :: JJS, JJE
122  integer :: IIS0, JJS0
123  integer :: i, j, k
124  !---------------------------------------------------------------------------
125 
126 #ifdef DEBUG
127  qflx_hi(:,:,:,:) = undef
128  qflx_lo(:,:,:,:) = undef
129 #endif
130 
131  do jjs = js, je, jblock
132  jje = jjs+jblock-1
133  do iis = is, ie, iblock
134  iie = iis+iblock-1
135 
136  ! at (x, y, w)
137  call atmos_dyn_fvm_fluxz_xyz_tracer( qflx_hi(:,:,:,zdir), & ! (out)
138  mflx_hi(:,:,:,zdir), qtrc, gsqrt(:,:,:,i_xyw), & ! (in)
139  num_diff(:,:,:,zdir), & ! (in)
140  cdz, & ! (in)
141  iis, iie, jjs, jje ) ! (in)
142 
143  ! at (u, y, z)
144  if ( .not. twod ) &
145  call atmos_dyn_fvm_fluxx_xyz_tracer( qflx_hi(:,:,:,xdir), & ! (out)
146  mflx_hi(:,:,:,xdir), qtrc, gsqrt(:,:,:,i_uyz), & ! (in)
147  num_diff(:,:,:,xdir), & ! (in)
148  cdz, & ! (in)
149  iis, iie, jjs, jje ) ! (in)
150 
151  ! at (x, v, z)
152  call atmos_dyn_fvm_fluxy_xyz_tracer( qflx_hi(:,:,:,ydir), & ! (out)
153  mflx_hi(:,:,:,ydir), qtrc, gsqrt(:,:,:,i_xvz), & ! (in)
154  num_diff(:,:,:,ydir), & ! (in)
155  cdz, & ! (in)
156  iis, iie, jjs, jje ) ! (in)
157 
158  if ( flag_fct_tracer ) then
159 
160  call atmos_dyn_fvm_fluxz_xyz_ud1( qflx_lo(:,:,:,zdir), & ! (out)
161  mflx_hi(:,:,:,zdir), qtrc0, gsqrt(:,:,:,i_xyw), & ! (in)
162  num_diff(:,:,:,zdir), & ! (in)
163  cdz, & ! (in)
164  iis-1, iie+1, jjs-1, jje+1 ) ! (in)
165 
166  if ( .not. twod ) &
167  call atmos_dyn_fvm_fluxx_xyz_ud1( qflx_lo(:,:,:,xdir), & ! (out)
168  mflx_hi(:,:,:,xdir), qtrc0, gsqrt(:,:,:,i_uyz), & ! (in)
169  num_diff(:,:,:,xdir), & ! (in)
170  cdz, & ! (in)
171  iis-1, iie+1, jjs-1, jje+1 ) ! (in)
172 
173  call atmos_dyn_fvm_fluxy_xyz_ud1( qflx_lo(:,:,:,ydir), & ! (out)
174  mflx_hi(:,:,:,ydir), qtrc0, gsqrt(:,:,:,i_xvz), & ! (in)
175  num_diff(:,:,:,ydir), & ! (in)
176  cdz, & ! (in)
177  iis-1, iie+1, jjs-1, jje+1 ) ! (in)
178  end if
179 
180  enddo
181  enddo
182 
183  if ( flag_fct_tracer ) then
184 
185  call atmos_dyn_fvm_fct( qflx_anti, & ! (out)
186  qtrc0, dens0, dens, & ! (in)
187  qflx_hi, qflx_lo, & ! (in)
188  mflx_hi, & ! (in)
189  rcdz, rcdx, rcdy, & ! (in)
190  gsqrt(:,:,:,i_xyz), & ! (in)
191  mapf, twod, dtl, & ! (in)
192  flag_fct_along_stream ) ! (in)
193 
194  do jjs = js, je, jblock
195  jje = jjs+jblock-1
196  do iis = is, ie, iblock
197  iie = iis+iblock-1
198 
199  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
200  do j = jjs, jje
201  do i = iis, iie
202  do k = ks-1, ke
203  qflx_hi(k,i,j,zdir) = qflx_hi(k,i,j,zdir) - qflx_anti(k,i,j,zdir)
204  end do
205  end do
206  end do
207 
208  if ( .not. twod ) then
209  if ( iis == is ) then
210  iis0 = iis-1
211  else
212  iis0 = iis
213  end if
214  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
215  do j = jjs, jje
216  do i = iis0, iie
217  do k = ks, ke
218  qflx_hi(k,i,j,xdir) = qflx_hi(k,i,j,xdir) - qflx_anti(k,i,j,xdir)
219  end do
220  end do
221  end do
222  end if
223 
224  if ( jjs == js ) then
225  jjs0 = jjs-1
226  else
227  jjs0 = jjs
228  end if
229  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
230  do j = jjs0, jje
231  do i = iis, iie
232  do k = ks, ke
233  qflx_hi(k,i,j,ydir) = qflx_hi(k,i,j,ydir) - qflx_anti(k,i,j,ydir)
234  end do
235  end do
236  end do
237 
238 
239  if ( twod ) then
240  !$omp parallel do private(i,j,k) OMP_SCHEDULE_
241  do j = jjs, jje
242  do k = ks, ke
243  qtrco(k,is,j) = ( qtrc0(k,is,j) * dens0(k,is,j) &
244  + dtl * ( - ( ( qflx_hi(k,is,j,zdir) - qflx_hi(k-1,is,j ,zdir) ) * rcdz(k) &
245  + ( qflx_hi(k,is,j,ydir) - qflx_hi(k ,is,j-1,ydir) ) * rcdy(j) &
246  ) * mapf(is,j,2) / gsqrt(k,is,j,i_xyz) &
247  + rhoq_t(k,is,j) ) ) / dens(k,is,j)
248  enddo
249  enddo
250  else
251  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
252  do j = jjs, jje
253  do i = iis, iie
254  do k = ks, ke
255  qtrco(k,i,j) = ( qtrc0(k,i,j) * dens0(k,i,j) &
256  + dtl * ( - ( ( qflx_hi(k,i,j,zdir) - qflx_hi(k-1,i ,j ,zdir) ) * rcdz(k) &
257  + ( qflx_hi(k,i,j,xdir) - qflx_hi(k ,i-1,j ,xdir) ) * rcdx(i) &
258  + ( qflx_hi(k,i,j,ydir) - qflx_hi(k ,i ,j-1,ydir) ) * rcdy(j) &
259  ) * mapf(i,j,1) * mapf(i,j,2) / gsqrt(k,i,j,i_xyz) &
260  + rhoq_t(k,i,j) ) ) / dens(k,i,j)
261  enddo
262  enddo
263  enddo
264  end if
265 
266  enddo
267  enddo
268 
269  else ! skip FCT
270 
271  do jjs = js, je, jblock
272  jje = jjs+jblock-1
273  do iis = is, ie, iblock
274  iie = iis+iblock-1
275 
276  if ( twod ) then
277  !$omp parallel do default(none) private(j,k) OMP_SCHEDULE_ &
278  !$omp shared(JJS,JJE,IS,KS,KE,QTRCo,QTRC0,DENS0,dtl,qflx_hi,RCDZ,RCDY,MAPF) &
279  !$omp shared(GSQRT,RHOQ_t,DENS,I_XYZ)
280  do j = jjs, jje
281  do k = ks, ke
282  qtrco(k,is,j) = ( qtrc0(k,is,j) * dens0(k,is,j) &
283  + dtl * ( - ( ( qflx_hi(k,is,j,zdir) - qflx_hi(k-1,is,j ,zdir) ) * rcdz(k) &
284  + ( qflx_hi(k,is,j,ydir) - qflx_hi(k ,is,j-1,ydir) ) * rcdy(j) &
285  ) * mapf(is,j,2) / gsqrt(k,is,j,i_xyz) &
286  + rhoq_t(k,is,j) ) ) / dens(k,is,j)
287  enddo
288  enddo
289  else
290  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
291  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,QTRCo,QTRC0,DENS0,dtl,qflx_hi,RCDZ,RCDX,RCDY,MAPF) &
292  !$omp shared(GSQRT,RHOQ_t,DENS,I_XYZ)
293  do j = jjs, jje
294  do i = iis, iie
295  do k = ks, ke
296  qtrco(k,i,j) = ( qtrc0(k,i,j) * dens0(k,i,j) &
297  + dtl * ( - ( ( qflx_hi(k,i,j,zdir) - qflx_hi(k-1,i ,j ,zdir) ) * rcdz(k) &
298  + ( qflx_hi(k,i,j,xdir) - qflx_hi(k ,i-1,j ,xdir) ) * rcdx(i) &
299  + ( qflx_hi(k,i,j,ydir) - qflx_hi(k ,i ,j-1,ydir) ) * rcdy(j) &
300  ) * mapf(i,j,1) * mapf(i,j,2) / gsqrt(k,i,j,i_xyz) &
301  + rhoq_t(k,i,j) ) ) / dens(k,i,j)
302  enddo
303  enddo
304  enddo
305  end if
306 
307  enddo
308  enddo
309 
310  end if
311 
312  return

References scale_atmos_dyn_fvm_fct::atmos_dyn_fvm_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_tracer::k, 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().

Here is the call graph for this function:
Here is the caller graph for this function:
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_atmos_dyn_fvm_fct::atmos_dyn_fvm_fct
subroutine, public atmos_dyn_fvm_fct(qflx_anti, phi_in, DENS0, DENS, qflx_hi, qflx_lo, mflx_hi, rdz, rdx, rdy, GSQRT, MAPF, TwoD, dt, flag_vect)
Setup.
Definition: scale_atmos_dyn_fvm_fct.F90:72
scale_atmos_dyn_fvm_flux
module scale_atmos_dyn_fvm_flux
Definition: scale_atmos_dyn_fvm_flux.F90:13
scale_atmos_grid_cartesc_index::jblock
integer, public jblock
block size for cache blocking: y
Definition: scale_atmos_grid_cartesC_index.F90:41
scale_atmos_grid_cartesc_index::iblock
integer, public iblock
block size for cache blocking: x
Definition: scale_atmos_grid_cartesC_index.F90:40
scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxx_xyz_ud1
subroutine, public atmos_dyn_fvm_fluxx_xyz_ud1(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XYZ
Definition: scale_atmos_dyn_fvm_flux_ud1.F90:212
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxy_xyz_tracer
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxy_xyz_tracer
Definition: scale_atmos_dyn_fvm_flux.F90:175
scale_atmos_grid_cartesc_index::i_xyz
integer, public i_xyz
Definition: scale_atmos_grid_cartesC_index.F90:90
scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxy_xyz_ud1
subroutine, public atmos_dyn_fvm_fluxy_xyz_ud1(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XYZ
Definition: scale_atmos_dyn_fvm_flux_ud1.F90:260
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
scale_atmos_grid_cartesc_index::i_uyz
integer, public i_uyz
Definition: scale_atmos_grid_cartesC_index.F90:94
scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxz_xyz_ud1
subroutine, public atmos_dyn_fvm_fluxz_xyz_ud1(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at XYZ
Definition: scale_atmos_dyn_fvm_flux_ud1.F90:142
scale_atmos_dyn_fvm_flux_ud1
module scale_atmos_dyn_fvm_flux_ud1
Definition: scale_atmos_dyn_fvm_flux_ud1.F90:16
scale_atmos_grid_cartesc_index::zdir
integer, parameter, public zdir
Definition: scale_atmos_grid_cartesC_index.F90:32
scale_atmos_grid_cartesc_index::i_xvz
integer, public i_xvz
Definition: scale_atmos_grid_cartesC_index.F90:95
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxz_xyz_tracer
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxz_xyz_tracer
Definition: scale_atmos_dyn_fvm_flux.F90:165
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxx_xyz_tracer
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxx_xyz_tracer
Definition: scale_atmos_dyn_fvm_flux.F90:170
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
scale_atmos_grid_cartesc_index::xdir
integer, parameter, public xdir
Definition: scale_atmos_grid_cartesC_index.F90:33
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
scale_atmos_grid_cartesc_index::ydir
integer, parameter, public ydir
Definition: scale_atmos_grid_cartesC_index.F90:34
scale_atmos_dyn_fvm_fct
module Atmosphere / Dynamics common
Definition: scale_atmos_dyn_fvm_fct.F90:12
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56
scale_atmos_grid_cartesc_index::i_xyw
integer, public i_xyw
Definition: scale_atmos_grid_cartesC_index.F90:91