SCALE-RM
Functions/Subroutines
scale_atmos_dyn_common Module Reference

module Atmosphere / Dynamics common More...

Functions/Subroutines

subroutine, public atmos_dyn_wdamp_setup (wdamp_coef, wdamp_tau, wdamp_height, FZ)
 Setup. More...
 
subroutine, public atmos_dyn_fill_halo (var, fill_constval, lateral_halo, top_bottom_halo)
 
subroutine, public atmos_dyn_copy_boundary (DENS, MOMZ, MOMX, MOMY, RHOT, PROG, DENS0, MOMZ0, MOMX0, MOMY0, RHOT0, PROG0, BND_W, BND_E, BND_S, BND_N, TwoD)
 
subroutine, public atmos_dyn_copy_boundary_tracer (QTRC, QTRC0, BND_W, BND_E, BND_S, BND_N, TwoD)
 
subroutine, public atmos_dyn_divergence (DDIV, MOMZ, MOMX, MOMY, GSQRT, J13G, J23G, J33G, MAPF, TwoD, RCDZ, RCDX, RCDY, RFDZ, FDZ)
 
subroutine, public atmos_dyn_prep_pres_linearization (DPRES, RT2P, REF_rhot, RHOT, QTRC, REF_pres, AQ_R, AQ_CV, AQ_CP, AQ_MASS)
 

Detailed Description

module Atmosphere / Dynamics common

Description
common subroutines for Atmospheric dynamical process
Author
Team SCALE

Function/Subroutine Documentation

◆ atmos_dyn_wdamp_setup()

subroutine, public scale_atmos_dyn_common::atmos_dyn_wdamp_setup ( real(rp), dimension(ka), intent(out)  wdamp_coef,
real(rp), intent(in)  wdamp_tau,
real(rp), intent(in)  wdamp_height,
real(rp), dimension(0:ka), intent(in)  FZ 
)

Setup.

Definition at line 69 of file scale_atmos_dyn_common.F90.

69  use scale_const, only: &
70  pi => const_pi, &
71  eps => const_eps
72  implicit none
73 
74  real(RP), intent(out) :: wdamp_coef(KA)
75  real(RP), intent(in) :: wdamp_tau
76  real(RP), intent(in) :: wdamp_height
77  real(RP), intent(in) :: FZ(0:KA)
78 
79  real(RP) :: alpha, sw
80 
81  integer :: k
82  !---------------------------------------------------------------------------
83 
84  !$acc data copyout(wdamp_coef) copyin(FZ)
85 
86  if ( wdamp_height < 0.0_rp ) then
87  !$acc kernels
88  wdamp_coef(:) = 0.0_rp
89  !$acc end kernels
90  elseif( fz(ke)-wdamp_height < eps ) then
91  !$acc kernels
92  wdamp_coef(:) = 0.0_rp
93  !$acc end kernels
94  else
95  alpha = 1.0_rp / wdamp_tau
96 
97  !$acc kernels
98  do k = ks, ke
99  sw = 0.5_rp + sign( 0.5_rp, fz(k)-wdamp_height )
100 
101  wdamp_coef(k) = alpha * sw &
102  * 0.5_rp * ( 1.0_rp - cos( pi * (fz(k)-wdamp_height) / (fz(ke)-wdamp_height)) )
103  enddo
104  !$acc end kernels
105  !$acc kernels
106  wdamp_coef( 1:ks-1) = wdamp_coef(ks)
107  !$acc end kernels
108  !$acc kernels
109  wdamp_coef(ke+1:ka ) = wdamp_coef(ke)
110  !$acc end kernels
111 
112  !$acc update host(wdamp_coef)
113  log_newline
114  log_info("ATMOS_DYN_wdamp_setup",*) 'Setup Rayleigh damping coefficient'
115  log_info_cont('(1x,A)') '|=== Rayleigh Damping Coef ===|'
116  log_info_cont('(1x,A)') '| k zh[m] coef[/s] |'
117  do k = ka, ke+1, -1
118  log_info_cont('(1x,A,I5,F10.2,ES12.4,A)') '| ',k, fz(k), wdamp_coef(k),' |'
119  enddo
120  k = ke
121  log_info_cont('(1x,A,I5,F10.2,ES12.4,A)') '| ',k, fz(k), wdamp_coef(k),' | KE = TOA'
122  do k = ke-1, ks, -1
123  log_info_cont('(1x,A,I5,F10.2,ES12.4,A)') '| ',k, fz(k), wdamp_coef(k),' |'
124  enddo
125  k = ks-1
126  log_info_cont('(1x,A,I5,F10.2,ES12.4,A)') '| ',k, fz(k), wdamp_coef(k),' | KS-1 = surface'
127  do k = ks-2, 1, -1
128  log_info_cont('(1x,A,I5,F10.2,ES12.4,A)') '| ',k, fz(k), wdamp_coef(k),' |'
129  enddo
130  k = 0
131  log_info_cont('(1x,A,I5,F10.2,12x,A)') '| ',k, fz(k), ' |'
132  log_info_cont('(1x,A)') '|=============================|'
133  endif
134 
135  !$acc end data
136 
137  return

References scale_const::const_eps, scale_const::const_pi, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, and scale_atmos_grid_cartesc_index::ks.

Referenced by scale_atmos_dyn::atmos_dyn_setup().

Here is the caller graph for this function:

◆ atmos_dyn_fill_halo()

subroutine, public scale_atmos_dyn_common::atmos_dyn_fill_halo ( real(rp), dimension(ka,ia,ja), intent(inout)  var,
real(rp), intent(in)  fill_constval,
logical, intent(in)  lateral_halo,
logical, intent(in)  top_bottom_halo 
)

Definition at line 144 of file scale_atmos_dyn_common.F90.

144  implicit none
145 
146  real(RP), intent(inout) :: var(KA,IA,JA)
147  real(RP), intent(in) :: fill_constval
148  logical, intent(in) :: lateral_halo
149  logical, intent(in) :: top_bottom_halo
150 
151  integer :: i, j, k
152  !----------------------------
153 
154  !$acc data copy(var)
155 
156  if (lateral_halo) then
157  !$acc kernels
158 !OCL XFILL
159  do j = 1, ja
160  do i = 1, isb-1
161  do k = 1, ka
162  var(k,i,j) = fill_constval
163  enddo
164  enddo
165  do i = ieb+1, ia
166  do k = 1, ka
167  var(k,i,j) = fill_constval
168  enddo
169  enddo
170  enddo
171  !$acc end kernels
172  !$acc kernels
173 !OCL XFILL
174  do j = 1, jsb-1
175  do i = 1, ia
176  do k = 1, ka
177  var(k,i,j) = fill_constval
178  enddo
179  enddo
180  enddo
181  !$acc end kernels
182  !$acc kernels
183 !OCL XFILL
184  do j = jeb+1, ja
185  do i = 1, ia
186  do k = 1, ka
187  var(k,i,j) = fill_constval
188  enddo
189  enddo
190  enddo
191  !$acc end kernels
192  end if
193 
194  if (top_bottom_halo) then
195  !$acc kernels
196 !OCL XFILL
197  do j = js, je
198  do i = is, ie
199  var( 1:ks-1,i,j) = fill_constval
200  var(ke+1:ka ,i,j) = fill_constval
201  enddo
202  enddo
203  !$acc end kernels
204  end if
205 
206  !$acc end data
207 
208  return

References scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::ieb, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::isb, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::jeb, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::jsb, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, and scale_atmos_grid_cartesc_index::ks.

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve().

Here is the caller graph for this function:

◆ atmos_dyn_copy_boundary()

subroutine, public scale_atmos_dyn_common::atmos_dyn_copy_boundary ( real(rp), dimension (ka,ia,ja), intent(inout)  DENS,
real(rp), dimension (ka,ia,ja), intent(inout)  MOMZ,
real(rp), dimension (ka,ia,ja), intent(inout)  MOMX,
real(rp), dimension (ka,ia,ja), intent(inout)  MOMY,
real(rp), dimension (ka,ia,ja), intent(inout)  RHOT,
real(rp), dimension (ka,ia,ja,va), intent(inout)  PROG,
real(rp), dimension(ka,ia,ja), intent(in)  DENS0,
real(rp), dimension(ka,ia,ja), intent(in)  MOMZ0,
real(rp), dimension(ka,ia,ja), intent(in)  MOMX0,
real(rp), dimension(ka,ia,ja), intent(in)  MOMY0,
real(rp), dimension(ka,ia,ja), intent(in)  RHOT0,
real(rp), dimension(ka,ia,ja,va), intent(in)  PROG0,
logical, intent(in)  BND_W,
logical, intent(in)  BND_E,
logical, intent(in)  BND_S,
logical, intent(in)  BND_N,
logical, intent(in)  TwoD 
)

Definition at line 215 of file scale_atmos_dyn_common.F90.

215  implicit none
216  real(RP), intent(inout) :: DENS (KA,IA,JA)
217  real(RP), intent(inout) :: MOMZ (KA,IA,JA)
218  real(RP), intent(inout) :: MOMX (KA,IA,JA)
219  real(RP), intent(inout) :: MOMY (KA,IA,JA)
220  real(RP), intent(inout) :: RHOT (KA,IA,JA)
221  real(RP), intent(inout) :: PROG (KA,IA,JA,VA)
222  real(RP), intent(in) :: DENS0(KA,IA,JA)
223  real(RP), intent(in) :: MOMZ0(KA,IA,JA)
224  real(RP), intent(in) :: MOMX0(KA,IA,JA)
225  real(RP), intent(in) :: MOMY0(KA,IA,JA)
226  real(RP), intent(in) :: RHOT0(KA,IA,JA)
227  real(RP), intent(in) :: PROG0(KA,IA,JA,VA)
228  logical, intent(in) :: BND_W
229  logical, intent(in) :: BND_E
230  logical, intent(in) :: BND_S
231  logical, intent(in) :: BND_N
232  logical, intent(in) :: TwoD
233 
234  integer :: k, i, j, iv
235 
236  !$acc data copy(DENS, MOMZ, MOMX, MOMY, RHOT, PROG) &
237  !$acc copyin(DENS0, MOMZ0, MOMX0, MOMY0, RHOT0, PROG0)
238 
239  if ( bnd_w .and. (.not. twod) ) then
240  !$omp parallel do default(none) private(j,k) OMP_SCHEDULE_ collapse(2) &
241  !$omp private(i,iv) &
242  !$omp shared(JA,IS,KS,KE,DENS,DENS0,MOMZ,MOMZ0,MOMX,MOMX0,MOMY,MOMY0,RHOT,RHOT0,VA,PROG,PROG0)
243  !$acc kernels
244 !OCL XFILL
245  do j = 1, ja
246  do i = 1, is-1
247  do k = ks, ke
248  dens(k,i,j) = dens0(k,i,j)
249  momz(k,i,j) = momz0(k,i,j)
250  momx(k,i,j) = momx0(k,i,j)
251  momy(k,i,j) = momy0(k,i,j)
252  rhot(k,i,j) = rhot0(k,i,j)
253  !$acc loop seq
254  do iv = 1, va
255  prog(k,i,j,iv) = prog0(k,i,j,iv)
256  end do
257  enddo
258  enddo
259  enddo
260  !$acc end kernels
261  end if
262  if ( bnd_e .and. (.not. twod) ) then
263  !$omp parallel do default(none) private(j,k) OMP_SCHEDULE_ collapse(2) &
264  !$omp private(i,iv) &
265  !$omp shared(JA,IE,IA,KS,KE,DENS,DENS0,MOMZ,MOMZ0,MOMX,MOMX0,MOMY,MOMY0,RHOT,RHOT0,VA,PROG,PROG0)
266  !$acc kernels
267 !OCL XFILL
268  do j = 1, ja
269  do i = ie+1, ia
270  do k = ks, ke
271  dens(k,i,j) = dens0(k,i,j)
272  momz(k,i,j) = momz0(k,i,j)
273  momx(k,i,j) = momx0(k,i,j)
274  momy(k,i,j) = momy0(k,i,j)
275  rhot(k,i,j) = rhot0(k,i,j)
276  !$acc loop seq
277  do iv = 1, va
278  prog(k,i,j,iv) = prog0(k,i,j,iv)
279  end do
280  enddo
281  enddo
282  enddo
283  !$acc end kernels
284  !$omp parallel do private(j,k) OMP_SCHEDULE_ collapse(2)
285  !$acc kernels
286 !OCL XFILL
287  do j = 1, ja
288  do k = ks, ke
289  momx(k,ie,j) = momx0(k,ie,j)
290  enddo
291  enddo
292  !$acc end kernels
293  end if
294  if ( bnd_s ) then
295  !$omp parallel do default(none) private(j,k) OMP_SCHEDULE_ collapse(2) &
296  !$omp private(i,iv) &
297  !$omp shared(JS,IA,KS,KE,DENS,DENS0,MOMZ,MOMZ0,MOMX,MOMX0,MOMY,MOMY0,RHOT,RHOT0,VA,PROG,PROG0)
298  !$acc kernels
299 !OCL XFILL
300  do j = 1, js-1
301  do i = 1, ia
302  do k = ks, ke
303  dens(k,i,j) = dens0(k,i,j)
304  momz(k,i,j) = momz0(k,i,j)
305  momx(k,i,j) = momx0(k,i,j)
306  momy(k,i,j) = momy0(k,i,j)
307  rhot(k,i,j) = rhot0(k,i,j)
308  !$acc loop seq
309  do iv = 1, va
310  prog(k,i,j,iv) = prog0(k,i,j,iv)
311  end do
312  enddo
313  enddo
314  enddo
315  !$acc end kernels
316  end if
317  if ( bnd_n ) then
318  !$omp parallel do default(none) private(j,k) OMP_SCHEDULE_ collapse(2) &
319  !$omp private(i,iv) &
320  !$omp shared(JA,JE,IA,KS,KE,DENS,DENS0,MOMZ,MOMZ0,MOMX,MOMX0,MOMY,MOMY0,RHOT,RHOT0,VA,PROG,PROG0)
321  !$acc kernels
322 !OCL XFILL
323  do j = je+1, ja
324  do i = 1, ia
325  do k = ks, ke
326  dens(k,i,j) = dens0(k,i,j)
327  momz(k,i,j) = momz0(k,i,j)
328  momx(k,i,j) = momx0(k,i,j)
329  momy(k,i,j) = momy0(k,i,j)
330  rhot(k,i,j) = rhot0(k,i,j)
331  !$acc loop seq
332  do iv = 1, va
333  prog(k,i,j,iv) = prog0(k,i,j,iv)
334  end do
335  enddo
336  enddo
337  enddo
338  !$acc end kernels
339  !$omp parallel do private(i,k) OMP_SCHEDULE_ collapse(2)
340  !$acc kernels
341 !OCL XFILL
342  do i = 1, ia
343  do k = ks, ke
344  momy(k,i,je) = momy0(k,i,je)
345  enddo
346  enddo
347  !$acc end kernels
348  end if
349 
350  !$acc end data
351 
352  return

References scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::ja, 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, and scale_index::va.

Referenced by scale_atmos_dyn_tinteg_short_rk11s8o::atmos_dyn_tinteg_short_rk11s8o(), scale_atmos_dyn_tinteg_short_rk3::atmos_dyn_tinteg_short_rk3(), scale_atmos_dyn_tinteg_short_rk4::atmos_dyn_tinteg_short_rk4(), and scale_atmos_dyn_tinteg_short_rk7s6o::atmos_dyn_tinteg_short_rk7s6o().

Here is the caller graph for this function:

◆ atmos_dyn_copy_boundary_tracer()

subroutine, public scale_atmos_dyn_common::atmos_dyn_copy_boundary_tracer ( real(rp), dimension (ka,ia,ja), intent(inout)  QTRC,
real(rp), dimension(ka,ia,ja), intent(in)  QTRC0,
logical, intent(in)  BND_W,
logical, intent(in)  BND_E,
logical, intent(in)  BND_S,
logical, intent(in)  BND_N,
logical, intent(in)  TwoD 
)

Definition at line 359 of file scale_atmos_dyn_common.F90.

359  implicit none
360  real(RP), intent(inout) :: QTRC (KA,IA,JA)
361  real(RP), intent(in) :: QTRC0(KA,IA,JA)
362  logical, intent(in) :: BND_W
363  logical, intent(in) :: BND_E
364  logical, intent(in) :: BND_S
365  logical, intent(in) :: BND_N
366  logical, intent(in) :: TwoD
367 
368  integer :: k, i, j
369 
370  !$acc data copy(QTRC) copyin(QTRC0)
371 
372  if ( bnd_w .and. (.not. twod) ) then
373  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
374  !$omp shared(JA,IS,KS,KE,QTRC,QTRC0)
375  !$acc kernels
376 !OCL XFILL
377  do j = 1, ja
378  do i = 1, is-1
379  do k = ks, ke
380  qtrc(k,i,j) = qtrc0(k,i,j)
381  enddo
382  enddo
383  enddo
384  !$acc end kernels
385  end if
386  if ( bnd_e .and. (.not. twod) ) then
387  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
388  !$omp shared(JA,IE,IA,KS,KE,QTRC,QTRC0)
389  !$acc kernels
390 !OCL XFILL
391  do j = 1, ja
392  do i = ie+1, ia
393  do k = ks, ke
394  qtrc(k,i,j) = qtrc0(k,i,j)
395  enddo
396  enddo
397  enddo
398  !$acc end kernels
399  end if
400  if ( bnd_s ) then
401  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
402  !$omp shared(JS,IA,KS,KE,QTRC,QTRC0)
403  !$acc kernels
404 !OCL XFILL
405  do j = 1, js-1
406  do i = 1, ia
407  do k = ks, ke
408  qtrc(k,i,j) = qtrc0(k,i,j)
409  enddo
410  enddo
411  enddo
412  !$acc end kernels
413  end if
414  if ( bnd_n ) then
415  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
416  !$omp shared(JA,JE,IA,KS,KE,QTRC,QTRC0)
417  !$acc kernels
418 !OCL XFILL
419  do j = je+1, ja
420  do i = 1, ia
421  do k = ks, ke
422  qtrc(k,i,j) = qtrc0(k,i,j)
423  enddo
424  enddo
425  enddo
426  !$acc end kernels
427  end if
428 
429  !$acc end data
430 
431  return

References scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_tracer::k, scale_atmos_grid_cartesc_index::ke, and scale_atmos_grid_cartesc_index::ks.

Referenced by scale_atmos_dyn_tinteg_tracer_linrk::atmos_dyn_tinteg_tracer_linrk(), and scale_atmos_dyn_tinteg_tracer_rk3::atmos_dyn_tinteg_tracer_rk3().

Here is the caller graph for this function:

◆ atmos_dyn_divergence()

subroutine, public scale_atmos_dyn_common::atmos_dyn_divergence ( real(rp), dimension(ka,ia,ja), intent(out)  DDIV,
real(rp), dimension(ka,ia,ja), intent(in)  MOMZ,
real(rp), dimension(ka,ia,ja), intent(in)  MOMX,
real(rp), dimension(ka,ia,ja), intent(in)  MOMY,
real(rp), dimension(ka,ia,ja,7), intent(in)  GSQRT,
real(rp), dimension(ka,ia,ja,7), intent(in)  J13G,
real(rp), dimension(ka,ia,ja,7), intent(in)  J23G,
real(rp), intent(in)  J33G,
real(rp), dimension(ia,ja,2,7), intent(in)  MAPF,
logical, intent(in)  TwoD,
real(rp), dimension(ka), intent(in)  RCDZ,
real(rp), dimension(ia), intent(in)  RCDX,
real(rp), dimension(ja), intent(in)  RCDY,
real(rp), dimension(ka-1), intent(in)  RFDZ,
real(rp), dimension(ka-1), intent(in)  FDZ 
)

Definition at line 442 of file scale_atmos_dyn_common.F90.

442  implicit none
443  real(RP), intent(out) :: DDIV(KA,IA,JA)
444  real(RP), intent(in) :: MOMZ(KA,IA,JA)
445  real(RP), intent(in) :: MOMX(KA,IA,JA)
446  real(RP), intent(in) :: MOMY(KA,IA,JA)
447  real(RP), intent(in) :: GSQRT(KA,IA,JA,7)
448  real(RP), intent(in) :: J13G(KA,IA,JA,7)
449  real(RP), intent(in) :: J23G(KA,IA,JA,7)
450  real(RP), intent(in) :: J33G
451  real(RP), intent(in) :: MAPF(IA,JA,2,7)
452  logical, intent(in) :: TwoD
453  real(RP), intent(in) :: RCDZ(KA)
454  real(RP), intent(in) :: RCDX(IA)
455  real(RP), intent(in) :: RCDY(JA)
456  real(RP), intent(in) :: RFDZ(KA-1)
457  real(RP), intent(in) :: FDZ(KA-1)
458 
459  integer :: k, i, j
460 
461  call prof_rapstart("DYN_divercence", 2)
462 
463  !$acc data copyout(DDIV) copyin(MOMZ, MOMX, MOMY, GSQRT, J13G, J23G, MAPF, RCDZ, RCDX, RCDY, RFDZ, FDZ)
464 
465  ! 3D divergence
466 
467  if ( twod ) then
468  !$omp parallel do private(j,k) OMP_SCHEDULE_
469  !$acc kernels
470  do j = js, je+1
471  do k = ks-1, ke+1
472  ddiv(k,is,j) = j33g * ( momz(k,is,j) - momz(k-1,is,j) ) * rcdz(k) &
473  + ( ( momy(k+1,is,j) + momy(k+1,is,j-1) ) * j23g(k+1,is,j,i_xyw) &
474  - ( momy(k-1,is,j) + momy(k-1,is,j-1) ) * j23g(k-1,is,j,i_xyw) ) / ( fdz(k)+fdz(k-1) ) &
475  + mapf(is,j,2,i_xy) &
476  * ( momy(k,is,j ) * gsqrt(k,is,j ,i_xvz) / mapf(is,j ,1,i_xv) &
477  - momy(k,is,j-1) * gsqrt(k,is,j-1,i_xvz) / mapf(is,j-1,1,i_xv) ) * rcdy(j)
478  enddo
479  enddo
480  !$acc end kernels
481 #ifdef DEBUG
482  k = iundef; i = iundef; j = iundef
483 #endif
484  !$omp parallel do private(j) OMP_SCHEDULE_
485  !$acc kernels
486  do j = js, je+1
487  ddiv(ks,is,j) = j33g * ( momz(ks,is,j) ) * rcdz(ks) &
488  + ( ( momy(ks+1,is,j) + momy(ks+1,is,j-1) ) * j23g(ks+1,is,j,i_xyw) &
489  - ( momy(ks ,is,j) + momy(ks ,is,j-1) ) * j23g(ks ,is,j,i_xyw) ) * rfdz(ks) &
490  + mapf(is,j,2,i_xy) &
491  * ( momy(ks,is,j ) * gsqrt(ks,is,j ,i_xvz) / mapf(is,j ,1,i_xv) &
492  - momy(ks,is,j-1) * gsqrt(ks,is,j-1,i_xvz) / mapf(is,j-1,1,i_xv) ) * rcdy(j)
493  ddiv(ke,is,j) = j33g * ( - momz(ke-1,is,j ) ) * rcdz(ke) &
494  + ( ( momy(ke ,is,j) + momy(ke ,is,j-1) ) * j23g(ke ,is,j,i_xyw) &
495  - ( momy(ke-1,is,j) + momy(ke-1,is,j-1) ) * j23g(ke-1,is,j,i_xyw) ) * rfdz(ke) &
496  + mapf(is,j,2,i_xy) &
497  * ( momy(ke,is,j ) * gsqrt(ke,is,j ,i_xvz) / mapf(is,j ,1,i_xv) &
498  - momy(ke,is,j-1) * gsqrt(ke,is,j-1,i_xvz) / mapf(is,j-1,1,i_xv) ) * rcdy(j)
499  enddo
500  !$acc end kernels
501 #ifdef DEBUG
502  k = iundef; i = iundef; j = iundef
503 #endif
504  else
505  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
506  !$acc kernels
507  do j = js, je+1
508  do i = is, ie+1
509  do k = ks-1, ke+1
510  ddiv(k,i,j) = j33g * ( momz(k,i,j) - momz(k-1,i ,j ) ) * rcdz(k) &
511  + ( ( momx(k+1,i,j) + momx(k+1,i-1,j ) ) * j13g(k+1,i,j,i_xyw) &
512  - ( momx(k-1,i,j) + momx(k-1,i-1,j ) ) * j13g(k-1,i,j,i_xyw) &
513  + ( momy(k+1,i,j) + momy(k+1,i ,j-1) ) * j23g(k+1,i,j,i_xyw) &
514  - ( momy(k-1,i,j) + momy(k-1,i ,j-1) ) * j23g(k-1,i,j,i_xyw) ) / ( fdz(k)+fdz(k-1) ) &
515  + mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) &
516  * ( ( momx(k,i ,j ) * gsqrt(k,i ,j ,i_uyz) / mapf(i ,j ,2,i_uy) &
517  - momx(k,i-1,j ) * gsqrt(k,i-1,j ,i_uyz) / mapf(i-1,j ,2,i_uy) ) * rcdx(i) &
518  + ( momy(k,i ,j ) * gsqrt(k,i ,j ,i_xvz) / mapf(i ,j ,1,i_xv) &
519  - momy(k,i, j-1) * gsqrt(k,i ,j-1,i_xvz) / mapf(i ,j-1,1,i_xv) ) * rcdy(j) )
520  enddo
521  enddo
522  enddo
523  !$acc end kernels
524 #ifdef DEBUG
525  k = iundef; i = iundef; j = iundef
526 #endif
527  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
528  !$acc kernels
529  do j = js, je+1
530  do i = is, ie+1
531  ddiv(ks,i,j) = j33g * ( momz(ks,i,j) ) * rcdz(ks) &
532  + ( ( momx(ks+1,i,j) + momx(ks+1,i-1,j ) ) * j13g(ks+1,i,j,i_xyw) &
533  - ( momx(ks-1,i,j) + momx(ks ,i-1,j ) ) * j13g(ks ,i,j,i_xyw) &
534  + ( momy(ks+1,i,j) + momy(ks+1,i ,j-1) ) * j23g(ks+1,i,j,i_xyw) &
535  - ( momy(ks ,i,j) + momy(ks ,i ,j-1) ) * j23g(ks ,i,j,i_xyw) ) * rfdz(ks) &
536  + mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) &
537  * ( ( momx(ks,i ,j ) * gsqrt(ks,i ,j ,i_uyz) / mapf(i ,j ,2,i_uy) &
538  - momx(ks,i-1,j ) * gsqrt(ks,i-1,j ,i_uyz) / mapf(i-1,j ,2,i_uy) ) * rcdx(i) &
539  + ( momy(ks,i ,j ) * gsqrt(ks,i ,j ,i_xvz) / mapf(i ,j ,1,i_xv) &
540  - momy(ks,i, j-1) * gsqrt(ks,i ,j-1,i_xvz) / mapf(i ,j-1,1,i_xv) ) * rcdy(j) )
541  ddiv(ke,i,j) = j33g * ( - momz(ke-1,i ,j ) ) * rcdz(ke) &
542  + ( ( momx(ke ,i,j) + momx(ke ,i-1,j ) ) * j13g(ke ,i,j,i_xyw) &
543  - ( momx(ke-1,i,j) + momx(ke-1,i-1,j ) ) * j13g(ke-1,i,j,i_xyw) &
544  + ( momy(ke ,i,j) + momy(ke ,i ,j-1) ) * j23g(ke ,i,j,i_xyw) &
545  - ( momy(ke-1,i,j) + momy(ke-1,i ,j-1) ) * j23g(ke-1,i,j,i_xyw) ) * rfdz(ke) &
546  + mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) &
547  * ( ( momx(ke,i ,j ) * gsqrt(ke,i ,j ,i_uyz) / mapf(i ,j ,2,i_uy) &
548  - momx(ke,i-1,j ) * gsqrt(ke,i-1,j ,i_uyz) / mapf(i-1,j ,2,i_uy) ) * rcdx(i) &
549  + ( momy(ke,i ,j ) * gsqrt(ke,i ,j ,i_xvz) / mapf(i ,j ,1,i_xv) &
550  - momy(ke,i, j-1) * gsqrt(ke,i ,j-1,i_xvz) / mapf(i ,j-1,1,i_xv) ) * rcdy(j) )
551  enddo
552  enddo
553  !$acc end kernels
554 #ifdef DEBUG
555  k = iundef; i = iundef; j = iundef
556 #endif
557  end if
558 
559  !$acc end data
560 
561  call prof_rapend ("DYN_divercence", 2)
562 
563  return

References scale_atmos_grid_cartesc_index::i_uy, scale_atmos_grid_cartesc_index::i_uyz, scale_atmos_grid_cartesc_index::i_xv, scale_atmos_grid_cartesc_index::i_xvz, scale_atmos_grid_cartesc_index::i_xy, scale_atmos_grid_cartesc_index::i_xyw, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, 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_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_dyn_prep_pres_linearization()

subroutine, public scale_atmos_dyn_common::atmos_dyn_prep_pres_linearization ( real(rp), dimension(ka,ia,ja), intent(out)  DPRES,
real(rp), dimension(ka,ia,ja), intent(out)  RT2P,
real(rp), dimension(ka,ia,ja), intent(out)  REF_rhot,
real(rp), dimension(ka,ia,ja), intent(in)  RHOT,
real(rp), dimension(ka,ia,ja,qa), intent(in)  QTRC,
real(rp), dimension(ka,ia,ja), intent(in)  REF_pres,
real(rp), dimension(qa), intent(in)  AQ_R,
real(rp), dimension(qa), intent(in)  AQ_CV,
real(rp), dimension(qa), intent(in)  AQ_CP,
real(rp), dimension(qa), intent(in)  AQ_MASS 
)

Definition at line 583 of file scale_atmos_dyn_common.F90.

583  use scale_const, only: &
584  p0 => const_pre00, &
585  rdry => const_rdry, &
586  cvdry => const_cvdry, &
587  cpdry => const_cpdry
588  implicit none
589 
590  real(RP), intent(out) :: DPRES(KA,IA,JA)
591  real(RP), intent(out) :: RT2P(KA,IA,JA)
592  real(RP), intent(out) :: REF_rhot(KA,IA,JA)
593  real(RP), intent(in) :: RHOT(KA,IA,JA)
594  real(RP), intent(in) :: QTRC(KA,IA,JA,QA)
595  real(RP), intent(in) :: REF_pres(KA,IA,JA)
596  real(RP), intent(in) :: AQ_R(QA)
597  real(RP), intent(in) :: AQ_CV(QA)
598  real(RP), intent(in) :: AQ_CP(QA)
599  real(RP), intent(in) :: AQ_MASS(QA)
600 
601  integer :: i, j, k
602  integer :: iq
603  real(RP) :: QDRY (KA) ! dry air
604  real(RP) :: Rtot (KA) ! total R
605  real(RP) :: CVtot(KA) ! total CV
606  real(RP) :: CPtot(KA) ! total CP
607  real(RP) :: PRES ! pressure
608 
609 #ifdef DRY
610  real(RP) :: CPovCV
611 #endif
612 
613  !--------------------------------------
614 
615 #ifdef DRY
616  cpovcv = cpdry / cvdry
617 #endif
618 
619 !OCL XFILL
620  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
621  !$omp shared(JA,IA,KS,KE) &
622  !$omp shared(P0,Rdry,RHOT,AQ_R,AQ_CV,AQ_CP,QTRC,AQ_MASS,REF_rhot,REF_pres,CPdry,CVdry,QA,RT2P,DPRES) &
623 #ifdef DRY
624  !$omp shared(CPovCV) &
625 #endif
626  !$omp private(i,j,k,iq) &
627  !$omp private(PRES,Rtot,CVtot,CPtot,QDRY)
628  !$acc kernels copyout(DPRES, RT2P, REF_rhot) copyin(RHOT, QTRC, REF_pres, AQ_R, AQ_CV, AQ_CP, AQ_MASS)
629  do j = 1, ja
630  !$acc loop private(Rtot, CVtot, CPtot, QDRY)
631  do i = 1, ia
632  do k = ks, ke
633  rtot(k) = 0.0_rp
634  cvtot(k) = 0.0_rp
635  cptot(k) = 0.0_rp
636  qdry(k) = 1.0_rp
637  end do
638  !$acc loop seq
639  do iq = 1, qa
640  do k = ks, ke
641  rtot(k) = rtot(k) + aq_r(iq) * qtrc(k,i,j,iq)
642  cvtot(k) = cvtot(k) + aq_cv(iq) * qtrc(k,i,j,iq)
643  cptot(k) = cptot(k) + aq_cp(iq) * qtrc(k,i,j,iq)
644  qdry(k) = qdry(k) - qtrc(k,i,j,iq) * aq_mass(iq)
645  enddo
646  end do
647  do k = ks, ke
648  rtot(k) = rtot(k) + rdry * qdry(k)
649  cvtot(k) = cvtot(k) + cvdry * qdry(k)
650  cptot(k) = cptot(k) + cpdry * qdry(k)
651  end do
652  do k = ks, ke
653  pres = p0 * ( rtot(k) * rhot(k,i,j) / p0 )**( cptot(k) / cvtot(k) )
654  rt2p(k,i,j) = cptot(k) / cvtot(k) * pres / rhot(k,i,j)
655  dpres(k,i,j) = pres - ref_pres(k,i,j)
656  ref_rhot(k,i,j) = rhot(k,i,j)
657  end do
658  dpres(ks-1,i,j) = dpres(ks+1,i,j) + ( ref_pres(ks+1,i,j) - ref_pres(ks-1,i,j) )
659  dpres(ke+1,i,j) = dpres(ke-1,i,j) + ( ref_pres(ke-1,i,j) - ref_pres(ke+1,i,j) )
660  end do
661  end do
662  !$acc end kernels
663 
664  return

References scale_const::const_cpdry, scale_const::const_cvdry, scale_const::const_pre00, scale_const::const_rdry, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ja, scale_tracer::k, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, and scale_tracer::qa.

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve().

Here is the caller graph for this function:
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_const::const_cvdry
real(rp), public const_cvdry
specific heat (dry air,constant volume) [J/kg/K]
Definition: scale_const.F90:61
scale_const::const_cpdry
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:60
scale_const::const_pi
real(rp), parameter, public const_pi
pi
Definition: scale_const.F90:32
scale_const::const_rdry
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:59
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:97
scale_atmos_grid_cartesc_index::i_xyw
integer, public i_xyw
Definition: scale_atmos_grid_cartesC_index.F90:92