SCALE-RM
scale_atmos_dyn_tstep_tracer_fvm_heve.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
21  use scale_index
22  use scale_tracer
23 
24 #ifdef DEBUG
25  use scale_debug, only: &
26  check
27  use scale_const, only: &
28  undef => const_undef, &
29  iundef => const_undef2
30 #endif
31  !-----------------------------------------------------------------------------
32  implicit none
33  private
34  !-----------------------------------------------------------------------------
35  !
36  !++ Public procedure
37  !
40 
41  !-----------------------------------------------------------------------------
42  !
43  !++ Public parameters & variables
44  !
45  !-----------------------------------------------------------------------------
46  !
47  !++ Private procedure
48  !
49  !-----------------------------------------------------------------------------
50  !
51  !++ Private parameters & variables
52  !
53  !-----------------------------------------------------------------------------
54 contains
55 
56  !-----------------------------------------------------------------------------
58  subroutine atmos_dyn_tstep_tracer_fvm_heve_setup( type )
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
71 
73  QTRCo, & ! (out)
74  qflx_hi, & ! (out)
75  qtrc, qtrc0, rhoq_t, &! (in)
76  dens0, dens, & ! (in)
77  mflx_hi, num_diff, & ! (in)
78  gsqrt, mapf, & ! (in)
79  cdz, rcdz, rcdx, rcdy, & ! (in)
80  twod, & ! (in)
81  dtl, & ! (in)
82  flag_fct_tracer, & ! (in)
83  flag_fct_along_stream ) ! (in)
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  !$acc data copy(QTRCo) &
127  !$acc copyout(qflx_hi) &
128  !$acc copyin(QTRC, QTRC0, RHOQ_t, DENS0, DENS, mflx_hi, num_diff, &
129  !$acc GSQRT, MAPF, CDZ, RCDZ, RCDX, RCDY) &
130  !$acc create(qflx_lo, qflx_anti)
131 
132 #ifdef DEBUG
133  !$acc kernels
134  qflx_hi(:,:,:,:) = undef
135  qflx_lo(:,:,:,:) = undef
136  !$acc end kernels
137 #endif
138 
139  do jjs = js, je, jblock
140  jje = jjs+jblock-1
141  do iis = is, ie, iblock
142  iie = iis+iblock-1
143 
144  ! at (x, y, w)
145  call atmos_dyn_fvm_fluxz_xyz_tracer( qflx_hi(:,:,:,zdir), & ! (out)
146  mflx_hi(:,:,:,zdir), qtrc, gsqrt(:,:,:,i_xyw), & ! (in)
147  num_diff(:,:,:,zdir), & ! (in)
148  cdz, & ! (in)
149  iis, iie, jjs, jje ) ! (in)
150 
151  ! at (u, y, z)
152  if ( .not. twod ) &
153  call atmos_dyn_fvm_fluxx_xyz_tracer( qflx_hi(:,:,:,xdir), & ! (out)
154  mflx_hi(:,:,:,xdir), qtrc, gsqrt(:,:,:,i_uyz), & ! (in)
155  num_diff(:,:,:,xdir), & ! (in)
156  cdz, & ! (in)
157  iis, iie, jjs, jje ) ! (in)
158 
159  ! at (x, v, z)
160  call atmos_dyn_fvm_fluxy_xyz_tracer( qflx_hi(:,:,:,ydir), & ! (out)
161  mflx_hi(:,:,:,ydir), qtrc, gsqrt(:,:,:,i_xvz), & ! (in)
162  num_diff(:,:,:,ydir), & ! (in)
163  cdz, & ! (in)
164  iis, iie, jjs, jje ) ! (in)
165 
166  if ( flag_fct_tracer ) then
167 
168  call atmos_dyn_fvm_fluxz_xyz_ud1( qflx_lo(:,:,:,zdir), & ! (out)
169  mflx_hi(:,:,:,zdir), qtrc0, gsqrt(:,:,:,i_xyw), & ! (in)
170  num_diff(:,:,:,zdir), & ! (in)
171  cdz, & ! (in)
172  iis-1, iie+1, jjs-1, jje+1 ) ! (in)
173 
174  if ( .not. twod ) &
175  call atmos_dyn_fvm_fluxx_xyz_ud1( qflx_lo(:,:,:,xdir), & ! (out)
176  mflx_hi(:,:,:,xdir), qtrc0, gsqrt(:,:,:,i_uyz), & ! (in)
177  num_diff(:,:,:,xdir), & ! (in)
178  cdz, & ! (in)
179  iis-1, iie+1, jjs-1, jje+1 ) ! (in)
180 
181  call atmos_dyn_fvm_fluxy_xyz_ud1( qflx_lo(:,:,:,ydir), & ! (out)
182  mflx_hi(:,:,:,ydir), qtrc0, gsqrt(:,:,:,i_xvz), & ! (in)
183  num_diff(:,:,:,ydir), & ! (in)
184  cdz, & ! (in)
185  iis-1, iie+1, jjs-1, jje+1 ) ! (in)
186  end if
187 
188  enddo
189  enddo
190 
191  if ( flag_fct_tracer ) then
192 
193  call atmos_dyn_fvm_fct( qflx_anti, & ! (out)
194  qtrc0, dens0, dens, & ! (in)
195  qflx_hi, qflx_lo, & ! (in)
196  mflx_hi, & ! (in)
197  rcdz, rcdx, rcdy, & ! (in)
198  gsqrt(:,:,:,i_xyz), & ! (in)
199  mapf, twod, dtl, & ! (in)
200  flag_fct_along_stream ) ! (in)
201 
202  do jjs = js, je, jblock
203  jje = jjs+jblock-1
204  do iis = is, ie, iblock
205  iie = iis+iblock-1
206 
207  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
208  do j = jjs, jje
209  do i = iis, iie
210  do k = ks-1, ke
211  qflx_hi(k,i,j,zdir) = qflx_hi(k,i,j,zdir) - qflx_anti(k,i,j,zdir)
212  end do
213  end do
214  end do
215 
216  if ( .not. twod ) then
217  if ( iis == is ) then
218  iis0 = iis-1
219  else
220  iis0 = iis
221  end if
222  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
223  do j = jjs, jje
224  do i = iis0, iie
225  do k = ks, ke
226  qflx_hi(k,i,j,xdir) = qflx_hi(k,i,j,xdir) - qflx_anti(k,i,j,xdir)
227  end do
228  end do
229  end do
230  end if
231 
232  if ( jjs == js ) then
233  jjs0 = jjs-1
234  else
235  jjs0 = jjs
236  end if
237  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
238  do j = jjs0, jje
239  do i = iis, iie
240  do k = ks, ke
241  qflx_hi(k,i,j,ydir) = qflx_hi(k,i,j,ydir) - qflx_anti(k,i,j,ydir)
242  end do
243  end do
244  end do
245 
246 
247  if ( twod ) then
248  !$omp parallel do private(i,j,k) OMP_SCHEDULE_
249  do j = jjs, jje
250  do k = ks, ke
251  qtrco(k,is,j) = ( qtrc0(k,is,j) * dens0(k,is,j) &
252  + dtl * ( - ( ( qflx_hi(k,is,j,zdir) - qflx_hi(k-1,is,j ,zdir) ) * rcdz(k) &
253  + ( qflx_hi(k,is,j,ydir) - qflx_hi(k ,is,j-1,ydir) ) * rcdy(j) &
254  ) * mapf(is,j,2) / gsqrt(k,is,j,i_xyz) &
255  + rhoq_t(k,is,j) ) ) / dens(k,is,j)
256  enddo
257  enddo
258  else
259  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
260  do j = jjs, jje
261  do i = iis, iie
262  do k = ks, ke
263  qtrco(k,i,j) = ( qtrc0(k,i,j) * dens0(k,i,j) &
264  + dtl * ( - ( ( qflx_hi(k,i,j,zdir) - qflx_hi(k-1,i ,j ,zdir) ) * rcdz(k) &
265  + ( qflx_hi(k,i,j,xdir) - qflx_hi(k ,i-1,j ,xdir) ) * rcdx(i) &
266  + ( qflx_hi(k,i,j,ydir) - qflx_hi(k ,i ,j-1,ydir) ) * rcdy(j) &
267  ) * mapf(i,j,1) * mapf(i,j,2) / gsqrt(k,i,j,i_xyz) &
268  + rhoq_t(k,i,j) ) ) / dens(k,i,j)
269  enddo
270  enddo
271  enddo
272  end if
273 
274  enddo
275  enddo
276 
277  else ! skip FCT
278 
279  do jjs = js, je, jblock
280  jje = jjs+jblock-1
281  do iis = is, ie, iblock
282  iie = iis+iblock-1
283 
284  if ( twod ) then
285  !$omp parallel do default(none) private(j,k) OMP_SCHEDULE_ &
286  !$omp shared(JJS,JJE,IS,KS,KE,QTRCo,QTRC0,DENS0,dtl,qflx_hi,RCDZ,RCDY,MAPF) &
287  !$omp shared(GSQRT,RHOQ_t,DENS,I_XYZ)
288  !$acc kernels
289  do j = jjs, jje
290  do k = ks, ke
291  qtrco(k,is,j) = ( qtrc0(k,is,j) * dens0(k,is,j) &
292  + dtl * ( - ( ( qflx_hi(k,is,j,zdir) - qflx_hi(k-1,is,j ,zdir) ) * rcdz(k) &
293  + ( qflx_hi(k,is,j,ydir) - qflx_hi(k ,is,j-1,ydir) ) * rcdy(j) &
294  ) * mapf(is,j,2) / gsqrt(k,is,j,i_xyz) &
295  + rhoq_t(k,is,j) ) ) / dens(k,is,j)
296  enddo
297  enddo
298  !$acc end kernels
299  else
300  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
301  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,QTRCo,QTRC0,DENS0,dtl,qflx_hi,RCDZ,RCDX,RCDY,MAPF) &
302  !$omp shared(GSQRT,RHOQ_t,DENS,I_XYZ)
303  !$acc kernels
304  do j = jjs, jje
305  do i = iis, iie
306  do k = ks, ke
307  qtrco(k,i,j) = ( qtrc0(k,i,j) * dens0(k,i,j) &
308  + dtl * ( - ( ( qflx_hi(k,i,j,zdir) - qflx_hi(k-1,i ,j ,zdir) ) * rcdz(k) &
309  + ( qflx_hi(k,i,j,xdir) - qflx_hi(k ,i-1,j ,xdir) ) * rcdx(i) &
310  + ( qflx_hi(k,i,j,ydir) - qflx_hi(k ,i ,j-1,ydir) ) * rcdy(j) &
311  ) * mapf(i,j,1) * mapf(i,j,2) / gsqrt(k,i,j,i_xyz) &
312  + rhoq_t(k,i,j) ) ) / dens(k,i,j)
313  enddo
314  enddo
315  enddo
316  !$acc end kernels
317  end if
318 
319  enddo
320  enddo
321 
322  end if
323  !$acc end data
324 
325  return
326  end subroutine atmos_dyn_tstep_tracer_fvm_heve
327 
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:350
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_index
module Index
Definition: scale_index.F90:11
scale_const::const_undef2
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:40
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
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:221
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:177
scale_atmos_grid_cartesc_index::i_xyz
integer, public i_xyz
Definition: scale_atmos_grid_cartesC_index.F90:91
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:271
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
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:95
scale_io
module STDIO
Definition: scale_io.F90:10
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:143
scale_tracer::k
real(rp), public k
Definition: scale_tracer.F90:45
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_debug::check
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:59
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:96
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:167
scale_atmos_dyn_tstep_tracer_fvm_heve::atmos_dyn_tstep_tracer_fvm_heve
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)
Definition: scale_atmos_dyn_tstep_tracer_fvm_heve.F90:84
scale_prof
module profiler
Definition: scale_prof.F90:11
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:172
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::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
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_dyn_tstep_tracer_fvm_heve::atmos_dyn_tstep_tracer_fvm_heve_setup
subroutine, public atmos_dyn_tstep_tracer_fvm_heve_setup(type)
Setup.
Definition: scale_atmos_dyn_tstep_tracer_fvm_heve.F90:59
scale_debug
module DEBUG
Definition: scale_debug.F90:11
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_dyn_tstep_tracer_fvm_heve
module Atmosphere / Dynamics
Definition: scale_atmos_dyn_tstep_tracer_fvm_heve.F90:12
scale_atmos_grid_cartesc_index::ydir
integer, parameter, public ydir
Definition: scale_atmos_grid_cartesC_index.F90:34
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
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:92