SCALE-RM
scale_atmos_dyn.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
11 !-------------------------------------------------------------------------------
12 #include "scalelib.h"
14  !-----------------------------------------------------------------------------
15  !
16  !++ used modules
17  !
18  use scale_precision
19  use scale_io
20  use scale_prof
22  use scale_index
23  use scale_tracer
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  !
38  public :: atmos_dyn_setup
39  public :: atmos_dyn_finalize
40  public :: atmos_dyn
41 
42  !-----------------------------------------------------------------------------
43  !
44  !++ Public parameters & variables
45  !
46  !-----------------------------------------------------------------------------
47  !
48  !++ Private procedure
49  !
50  !-----------------------------------------------------------------------------
51  !
52  !++ Private parameters & variables
53  !
54  ! for communication
55  integer, private :: I_COMM_DENS = 1
56  integer, private :: I_COMM_MOMZ = 2
57  integer, private :: I_COMM_MOMX = 3
58  integer, private :: I_COMM_MOMY = 4
59  integer, private :: I_COMM_RHOT = 5
60  integer, private, allocatable :: I_COMM_PROG(:)
61  integer, private, allocatable :: I_COMM_QTRC(:)
62 
63  logical, private :: BND_W
64  logical, private :: BND_E
65  logical, private :: BND_S
66  logical, private :: BND_N
67 
68  real(RP), private, allocatable :: num_diff (:,:,:,:,:)
69  real(RP), private, allocatable :: num_diff_q(:,:,:,:)
70  real(RP), private, allocatable :: wdamp_coef(:) ! coefficient for Rayleigh damping of w
71 
72  logical, private :: DYN_NONE
73 
74  !-----------------------------------------------------------------------------
75 contains
76  !-----------------------------------------------------------------------------
78  subroutine atmos_dyn_setup( &
79  DYN_Tinteg_Short_TYPE, &
80  DYN_Tinteg_Tracer_TYPE, &
81  DYN_Tinteg_Large_TYPE, &
82  DYN_Tstep_Tracer_TYPE, &
83  DYN_Tstep_Large_TYPE, &
84  DYN_Tstep_Short_TYPE, &
85  DYN_FVM_FLUX_TYPE, &
86  DYN_FVM_FLUX_TYPE_TRACER, &
87  DENS, MOMZ, MOMX, MOMY, RHOT, &
88  QTRC, PROG, &
89  CDZ, CDX, CDY, &
90  FDZ, FDX, FDY, &
91  wdamp_tau, &
92  wdamp_height, &
93  FZ, &
94  none )
95  use scale_prc, only: &
96  prc_abort
97  use scale_prc_cartesc, only: &
98  prc_twod, &
99  prc_has_e, &
100  prc_has_w, &
101  prc_has_n, &
102  prc_has_s
103  use scale_const, only: &
104  undef => const_undef
105  use scale_comm_cartesc, only: &
107  use scale_atmos_dyn_common, only: &
109  use scale_atmos_dyn_fvm_numfilter, only: &
111  use scale_atmos_dyn_tinteg_short, only: &
113  use scale_atmos_dyn_tinteg_tracer, only: &
115  use scale_atmos_dyn_tinteg_large, only: &
117  use scale_atmos_dyn_tstep_short, only: &
119  use scale_atmos_dyn_tstep_tracer, only: &
121  use scale_atmos_dyn_tstep_large, only: &
123  use scale_atmos_dyn_fvm_flux, only: &
125  use scale_spnudge, only: &
127  implicit none
128 
129  character(len=*), intent(in) :: dyn_tinteg_short_type
130  character(len=*), intent(in) :: dyn_tinteg_tracer_type
131  character(len=*), intent(in) :: dyn_tinteg_large_type
132  character(len=*), intent(in) :: dyn_tstep_tracer_type
133  character(len=*), intent(in) :: dyn_tstep_large_type
134  character(len=*), intent(in) :: dyn_tstep_short_type
135  character(len=*), intent(in) :: dyn_fvm_flux_type
136  character(len=*), intent(in) :: dyn_fvm_flux_type_tracer
137  real(rp), intent(inout) :: dens(ka,ia,ja) ! MPI_RECV_INIT requires intent(inout)
138  real(rp), intent(inout) :: momz(ka,ia,ja)
139  real(rp), intent(inout) :: momx(ka,ia,ja)
140  real(rp), intent(inout) :: momy(ka,ia,ja)
141  real(rp), intent(inout) :: rhot(ka,ia,ja)
142  real(rp), intent(inout) :: qtrc(ka,ia,ja,qa)
143  real(rp), intent(inout) :: prog(ka,ia,ja,va)
144  real(rp), intent(in) :: cdz(ka)
145  real(rp), intent(in) :: cdx(ia)
146  real(rp), intent(in) :: cdy(ja)
147  real(rp), intent(in) :: fdz(ka-1)
148  real(rp), intent(in) :: fdx(ia-1)
149  real(rp), intent(in) :: fdy(ja-1)
150  real(rp), intent(in) :: wdamp_tau
151  real(rp), intent(in) :: wdamp_height
152  real(rp), intent(in) :: fz(0:ka)
153  logical, optional, intent(in) :: none
154 
155  integer :: iv, iq
156  !---------------------------------------------------------------------------
157 
158  dyn_none = .false.
159 
160  if ( present(none) ) then
161  dyn_none = none
162  endif
163 
164  allocate( i_comm_prog(max(va,1)) )
165  allocate( i_comm_qtrc(qa) )
166 
167  if ( .NOT. dyn_none ) then
168  bnd_w = ( .NOT. prc_has_w ) .and. ( .NOT. prc_twod )
169  bnd_e = ( .NOT. prc_has_e ) .and. ( .NOT. prc_twod )
170  bnd_s = .NOT. prc_has_s
171  bnd_n = .NOT. prc_has_n
172 
173  allocate( num_diff(ka,ia,ja,5,3) )
174  allocate( num_diff_q(ka,ia,ja,3) )
175  allocate( wdamp_coef(ka) )
176  num_diff(:,:,:,:,:) = undef
177  num_diff_q(:,:,:,:) = undef
178  !$acc enter data create(num_diff, num_diff_q, wdamp_coef)
179 
180  call atmos_dyn_fvm_flux_setup ( dyn_fvm_flux_type, & ! [IN]
181  dyn_fvm_flux_type_tracer ) ! [IN]
182 
184 
185  call atmos_dyn_tstep_tracer_setup ( dyn_tstep_tracer_type ) ! [IN]
186 
187  call atmos_dyn_tstep_large_setup ( dyn_tstep_large_type, & ! [IN]
188  dens, momz, momx, momy, rhot, & ! [INOUT]
189  qtrc, prog ) ! [INOUT]
190 
191  call atmos_dyn_tinteg_short_setup ( dyn_tinteg_short_type, dyn_tstep_short_type ) ! [IN]
192 
193  call atmos_dyn_tinteg_tracer_setup( dyn_tinteg_tracer_type ) ! [IN]
194 
195  call atmos_dyn_tinteg_large_setup ( dyn_tinteg_large_type ) ! [IN]
196 
197  ! numerical diffusion
198  call atmos_dyn_fvm_numfilter_setup( num_diff, num_diff_q, & ! [INOUT]
199  cdz, cdx, cdy, fdz, fdx, fdy ) ! [IN]
200 
201  ! sponge layer
202  call atmos_dyn_wdamp_setup( wdamp_coef(:), & ! [INOUT]
203  wdamp_tau, wdamp_height, & ! [IN]
204  fz(:) ) ! [IN]
205 
206  ! setup for spectral nudging
207  call spnudge_setup( ka, ks, ke, ia, is, ie, ja, js, je )
208 
209 
210  else
211 
212  call comm_vars8_init( 'DENS', dens, i_comm_dens )
213  call comm_vars8_init( 'MOMZ', momz, i_comm_momz )
214  call comm_vars8_init( 'MOMX', momx, i_comm_momx )
215  call comm_vars8_init( 'MOMY', momy, i_comm_momy )
216  call comm_vars8_init( 'RHOT', rhot, i_comm_rhot )
217 
218  do iv = 1, va
219  i_comm_prog(iv) = 5 + iv
220  call comm_vars8_init( 'PROG', prog(:,:,:,iv), i_comm_prog(iv) )
221  enddo
222 
223  do iq = 1, qa
224  i_comm_qtrc(iq) = 5 + va + iq
225  call comm_vars8_init( 'QTRC', qtrc(:,:,:,iq), i_comm_qtrc(iq) )
226  enddo
227 
228  endif
229 
230  return
231  end subroutine atmos_dyn_setup
232 
233  !-----------------------------------------------------------------------------
235  subroutine atmos_dyn_finalize
238  use scale_atmos_dyn_tinteg_short, only: &
240  use scale_atmos_dyn_tinteg_tracer, only: &
242  use scale_atmos_dyn_fvm_numfilter, only: &
244  use scale_spnudge, only: &
246 
247  deallocate( i_comm_prog )
248  deallocate( i_comm_qtrc )
249 
250  if ( .NOT. dyn_none ) then
251 
253 
255 
257 
259 
260  call spnudge_finalize
261 
262  !$acc exit data delete(num_diff, num_diff_q, wdamp_coef)
263  deallocate( num_diff )
264  deallocate( num_diff_q )
265  deallocate( wdamp_coef )
266  end if
267 
268  return
269  end subroutine atmos_dyn_finalize
270  !-----------------------------------------------------------------------------
272  subroutine atmos_dyn( &
273  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, &
274  PROG, &
275  DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, &
276  DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, &
277  CORIOLIS, &
278  CDZ, CDX, CDY, FDZ, FDX, FDY, &
279  RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
280  PHI, GSQRT, &
281  J13G, J23G, J33G, MAPF, &
282  AQ_R, AQ_CV, AQ_CP, AQ_MASS, &
283  REF_dens, REF_pott, REF_qv, REF_pres, &
284  ND_COEF, ND_COEF_Q, ND_LAPLACIAN_NUM, &
285  ND_SFC_FACT, ND_USE_RS, &
286  BND_QA, BND_IQ, BND_SMOOTHER_FACT, &
287  DAMP_DENS, DAMP_VELZ, DAMP_VELX, &
288  DAMP_VELY, DAMP_POTT, DAMP_QTRC, &
289  DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, &
290  DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, &
291  MFLUX_OFFSET_X, MFLUX_OFFSET_Y, &
292  divdmp_coef, &
293  FLAG_TRACER_SPLIT_TEND, &
294  FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, &
295  FLAG_FCT_ALONG_STREAM, &
296  USE_AVERAGE, &
297  I_QV, &
298  DTSEC, DTSEC_DYN )
299  use scale_prc_cartesc, only: &
300  prc_twod
301  use scale_comm_cartesc, only: &
302  comm_vars8, &
303  comm_wait
304  use scale_atmos_dyn_tinteg_large, only: &
306  implicit none
307 
308  real(rp), intent(inout) :: dens(ka,ia,ja)
309  real(rp), intent(inout) :: momz(ka,ia,ja)
310  real(rp), intent(inout) :: momx(ka,ia,ja)
311  real(rp), intent(inout) :: momy(ka,ia,ja)
312  real(rp), intent(inout) :: rhot(ka,ia,ja)
313  real(rp), intent(inout) :: qtrc(ka,ia,ja,qa)
314  real(rp), intent(inout) :: prog(ka,ia,ja,va)
315 
316  real(rp), intent(inout) :: dens_av(ka,ia,ja)
317  real(rp), intent(inout) :: momz_av(ka,ia,ja)
318  real(rp), intent(inout) :: momx_av(ka,ia,ja)
319  real(rp), intent(inout) :: momy_av(ka,ia,ja)
320  real(rp), intent(inout) :: rhot_av(ka,ia,ja)
321  real(rp), intent(inout) :: qtrc_av(ka,ia,ja,qa)
322 
323  real(rp), intent(in) :: dens_tp(ka,ia,ja)
324  real(rp), intent(in) :: momz_tp(ka,ia,ja)
325  real(rp), intent(in) :: momx_tp(ka,ia,ja)
326  real(rp), intent(in) :: momy_tp(ka,ia,ja)
327  real(rp), intent(in) :: rhot_tp(ka,ia,ja)
328  real(rp), intent(in) :: rhoq_tp(ka,ia,ja,qa)
329 
330  real(rp), intent(in) :: coriolis(ia,ja)
331 
332  real(rp), intent(in) :: cdz (ka)
333  real(rp), intent(in) :: cdx (ia)
334  real(rp), intent(in) :: cdy (ja)
335  real(rp), intent(in) :: fdz (ka-1)
336  real(rp), intent(in) :: fdx (ia-1)
337  real(rp), intent(in) :: fdy (ja-1)
338  real(rp), intent(in) :: rcdz(ka)
339  real(rp), intent(in) :: rcdx(ia)
340  real(rp), intent(in) :: rcdy(ja)
341  real(rp), intent(in) :: rfdz(ka-1)
342  real(rp), intent(in) :: rfdx(ia-1)
343  real(rp), intent(in) :: rfdy(ja-1)
344 
345  real(rp), intent(in) :: phi (ka,ia,ja)
346  real(rp), intent(in) :: gsqrt(ka,ia,ja,7)
347  real(rp), intent(in) :: j13g (ka,ia,ja,7)
348  real(rp), intent(in) :: j23g (ka,ia,ja,7)
349  real(rp), intent(in) :: j33g
350  real(rp), intent(in) :: mapf (ia,ja,2,4)
351 
352  real(rp), intent(in) :: aq_r (qa)
353  real(rp), intent(in) :: aq_cv (qa)
354  real(rp), intent(in) :: aq_cp (qa)
355  real(rp), intent(in) :: aq_mass(qa)
356 
357  real(rp), intent(in) :: ref_dens(ka,ia,ja)
358  real(rp), intent(in) :: ref_pott(ka,ia,ja)
359  real(rp), intent(in) :: ref_qv (ka,ia,ja)
360  real(rp), intent(in) :: ref_pres(ka,ia,ja)
361  real(rp), intent(in) :: nd_coef
362  real(rp), intent(in) :: nd_coef_q
363  integer, intent(in) :: nd_laplacian_num
364  real(rp), intent(in) :: nd_sfc_fact
365  logical, intent(in) :: nd_use_rs
366 
367  integer, intent(in) :: bnd_qa
368  integer, intent(in) :: bnd_iq(qa)
369  real(rp), intent(in) :: bnd_smoother_fact
370 
371  real(rp), intent(in) :: damp_dens(ka,ia,ja)
372  real(rp), intent(in) :: damp_velz(ka,ia,ja)
373  real(rp), intent(in) :: damp_velx(ka,ia,ja)
374  real(rp), intent(in) :: damp_vely(ka,ia,ja)
375  real(rp), intent(in) :: damp_pott(ka,ia,ja)
376  real(rp), intent(in) :: damp_qtrc(ka,ia,ja,bnd_qa)
377 
378  real(rp), intent(in) :: damp_alpha_dens(ka,ia,ja)
379  real(rp), intent(in) :: damp_alpha_velz(ka,ia,ja)
380  real(rp), intent(in) :: damp_alpha_velx(ka,ia,ja)
381  real(rp), intent(in) :: damp_alpha_vely(ka,ia,ja)
382  real(rp), intent(in) :: damp_alpha_pott(ka,ia,ja)
383  real(rp), intent(in) :: damp_alpha_qtrc(ka,ia,ja,bnd_qa)
384  real(rp), intent(in) :: mflux_offset_x(ka,ja,2)
385  real(rp), intent(in) :: mflux_offset_y(ka,ia,2)
386 
387  real(rp), intent(in) :: divdmp_coef
388 
389  logical, intent(in) :: flag_tracer_split_tend
390  logical, intent(in) :: flag_fct_momentum
391  logical, intent(in) :: flag_fct_t
392  logical, intent(in) :: flag_fct_tracer
393  logical, intent(in) :: flag_fct_along_stream
394 
395  logical, intent(in) :: use_average
396 
397  integer, intent(in) :: i_qv
398 
399  real(dp), intent(in) :: dtsec
400  real(dp), intent(in) :: dtsec_dyn
401 
402  ! for time integartion
403  real(rp) :: dens00 (ka,ia,ja) ! saved density before update
404 
405  ! For tracer advection
406  real(rp) :: dt
407 
408  integer :: i, j, k, iq
409  !---------------------------------------------------------------------------
410 
411  log_progress(*) 'atmosphere / dynamics'
412 
413  dt = real(dtsec, kind=rp)
414 
415  if ( dyn_none ) then
416 
417  call prof_rapstart("DYN_Tinteg", 2)
418 
419  !$acc data create(DENS00)
420 
421  !$acc kernels
422  do j = js, je
423  do i = is, ie
424  do k = ks, ke
425  dens00(k,i,j) = dens(k,i,j) ! save previous value
426 
427  dens(k,i,j) = dens(k,i,j) + dens_tp(k,i,j) * dt
428  enddo
429  enddo
430  enddo
431  !$acc end kernels
432 
433  !$acc kernels
434  do j = js, je
435  do i = is, ie
436  do k = ks, ke
437  momz(k,i,j) = momz(k,i,j) + momz_tp(k,i,j) * dt
438  enddo
439  enddo
440  enddo
441  !$acc end kernels
442 
443  !$acc kernels
444  do j = js, je
445  do i = is, ie
446  do k = ks, ke
447  momx(k,i,j) = momx(k,i,j) + momx_tp(k,i,j) * dt
448  enddo
449  enddo
450  enddo
451  !$acc end kernels
452 
453  !$acc kernels
454  do j = js, je
455  do i = is, ie
456  do k = ks, ke
457  momy(k,i,j) = momy(k,i,j) + momy_tp(k,i,j) * dt
458  enddo
459  enddo
460  enddo
461  !$acc end kernels
462 
463  !$acc kernels
464  do j = js, je
465  do i = is, ie
466  do k = ks, ke
467  rhot(k,i,j) = rhot(k,i,j) + rhot_tp(k,i,j) * dt
468  enddo
469  enddo
470  enddo
471  !$acc end kernels
472 
473  !$acc kernels
474  do iq = 1, qa
475  do j = js, je
476  do i = is, ie
477  do k = ks, ke
478  qtrc(k,i,j,iq) = max( 0.0_rp, &
479  qtrc(k,i,j,iq) * dens00(k,i,j) + rhoq_tp(k,i,j,iq) * dt &
480  ) / dens(k,i,j)
481  enddo
482  enddo
483  enddo
484  enddo
485  !$acc end kernels
486 
487  !$acc end data
488 
489  call comm_vars8( dens(:,:,:), i_comm_dens )
490  call comm_vars8( momz(:,:,:), i_comm_momz )
491  call comm_vars8( momx(:,:,:), i_comm_momx )
492  call comm_vars8( momy(:,:,:), i_comm_momy )
493  call comm_vars8( rhot(:,:,:), i_comm_rhot )
494  do iq = 1, qa
495  call comm_vars8( qtrc(:,:,:,iq), i_comm_qtrc(iq) )
496  enddo
497 
498  call comm_wait ( dens(:,:,:), i_comm_dens, .false. )
499  call comm_wait ( momz(:,:,:), i_comm_momz, .false. )
500  call comm_wait ( momx(:,:,:), i_comm_momx, .false. )
501  call comm_wait ( momy(:,:,:), i_comm_momy, .false. )
502  call comm_wait ( rhot(:,:,:), i_comm_rhot, .false. )
503  do iq = 1, qa
504  call comm_wait ( qtrc(:,:,:,iq), i_comm_qtrc(iq), .false. )
505  enddo
506 
507  if ( use_average ) then
508  !$acc kernels
509  dens_av(:,:,:) = dens(:,:,:)
510  !$acc end kernels
511  !$acc kernels
512  momz_av(:,:,:) = momz(:,:,:)
513  !$acc end kernels
514  !$acc kernels
515  momx_av(:,:,:) = momx(:,:,:)
516  !$acc end kernels
517  !$acc kernels
518  momy_av(:,:,:) = momy(:,:,:)
519  !$acc end kernels
520  !$acc kernels
521  rhot_av(:,:,:) = rhot(:,:,:)
522  !$acc end kernels
523  !$acc kernels
524  qtrc_av(:,:,:,:) = qtrc(:,:,:,:)
525  !$acc end kernels
526  endif
527 
528  call prof_rapend("DYN_Tinteg", 2)
529 
530  return
531  endif ! if DYN_NONE == .true.
532 
533  call prof_rapstart("DYN_Tinteg", 2)
534 
535  !$acc data copy(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG) &
536  !$acc copyin(DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, &
537  !$acc CORIOLIS, &
538  !$acc CDZ, CDX, CDY, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
539  !$acc PHI, GSQRT, J13G, J23G, MAPF, &
540  !$acc AQ_R, AQ_CV, AQ_CP, AQ_MASS, &
541  !$acc REF_dens, REF_pott, REF_qv, REF_pres, &
542  !$acc BND_IQ, DAMP_DENS, DAMP_VELZ, DAMP_VELX, DAMP_VELY, DAMP_POTT, DAMP_QTRC, &
543  !$acc DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, &
544  !$acc MFLUX_OFFSET_X, MFLUX_OFFSET_Y)
545 
546  !$acc data copy(DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av) if(USE_AVERAGE)
547 
548 
549  call atmos_dyn_tinteg_large( dens, momz, momx, momy, rhot, qtrc, & ! [INOUT]
550  prog, & ! [INOUT]
551  dens_av, momz_av, momx_av, momy_av, rhot_av, qtrc_av, & ! [INOUT]
552  num_diff, num_diff_q, & ! [OUT;WORK]
553  dens_tp, momz_tp, momx_tp, momy_tp, rhot_tp, rhoq_tp, & ! [IN]
554  coriolis, & ! [IN]
555  cdz, cdx, cdy, fdz, fdx, fdy, & ! [IN]
556  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
557  phi, gsqrt, & ! [IN]
558  j13g, j23g, j33g, mapf, & ! [IN]
559  aq_r, aq_cv, aq_cp, aq_mass, & ! [IN]
560  ref_dens, ref_pott, ref_qv, ref_pres, & ! [IN]
561  bnd_w, bnd_e, bnd_s, bnd_n, prc_twod, & ! [IN]
562  nd_coef, nd_coef_q, nd_laplacian_num, & ! [IN]
563  nd_sfc_fact, nd_use_rs, & ! [IN]
564  bnd_qa, bnd_iq, bnd_smoother_fact, & ! [IN]
565  damp_dens, damp_velz, damp_velx, & ! [IN]
566  damp_vely, damp_pott, damp_qtrc, & ! [IN]
567  damp_alpha_dens, damp_alpha_velz, damp_alpha_velx, & ! [IN]
568  damp_alpha_vely, damp_alpha_pott, damp_alpha_qtrc, & ! [IN]
569  mflux_offset_x, mflux_offset_y, & ! [IN]
570  wdamp_coef, divdmp_coef, & ! [IN]
571  flag_tracer_split_tend, & ! [IN]
572  flag_fct_momentum, flag_fct_t, flag_fct_tracer, & ! [IN]
573  flag_fct_along_stream, & ! [IN]
574  use_average, & ! [IN]
575  i_qv, & ! [IN]
576  dtsec, dtsec_dyn ) ! [IN]
577 
578  call prof_rapend ("DYN_Tinteg", 2)
579 
580  !$acc end data
581  !$acc end data
582 
583  return
584  end subroutine atmos_dyn
585 
586 end module scale_atmos_dyn
scale_spnudge::spnudge_setup
subroutine, public spnudge_setup(KA, KS, KE, IA, IS, IE, JA, JS, JE)
Definition: scale_spnudge.F90:51
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_tracer::qa
integer, public qa
Definition: scale_tracer.F90:35
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_atmos_dyn_tinteg_large::atmos_dyn_tinteg_large_setup
subroutine, public atmos_dyn_tinteg_large_setup(ATMOS_DYN_Tinteg_large_TYPE)
Register.
Definition: scale_atmos_dyn_tinteg_large.F90:190
scale_prc_cartesc::prc_has_s
logical, public prc_has_s
Definition: scale_prc_cartesC.F90:51
scale_atmos_dyn_common
module Atmosphere / Dynamics common
Definition: scale_atmos_dyn_common.F90:12
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_dyn
module Atmosphere / Dynamics FENT + FCT
Definition: scale_atmos_dyn.F90:13
scale_prof::prof_rapstart
subroutine, public prof_rapstart(rapname_base, level, disable_barrier)
Start raptime.
Definition: scale_prof.F90:174
scale_atmos_dyn_tinteg_large::atmos_dyn_tinteg_large
procedure(large), pointer, public atmos_dyn_tinteg_large
Definition: scale_atmos_dyn_tinteg_large.F90:169
scale_atmos_dyn_tinteg_short::atmos_dyn_tinteg_short_setup
subroutine, public atmos_dyn_tinteg_short_setup(ATMOS_DYN_Tinteg_short_TYPE, ATMOS_DYN_Tstep_short_TYPE)
Register.
Definition: scale_atmos_dyn_tinteg_short.F90:135
scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_finalize
procedure(finalize), pointer, public atmos_dyn_tinteg_tracer_finalize
Definition: scale_atmos_dyn_tinteg_tracer.F90:78
scale_prc_cartesc::prc_has_n
logical, public prc_has_n
Definition: scale_prc_cartesC.F90:49
scale_atmos_dyn_tinteg_short
module Atmosphere / Dynamics Temporal integration
Definition: scale_atmos_dyn_tinteg_short.F90:11
scale_atmos_dyn_tstep_large::atmos_dyn_tstep_large_setup
subroutine, public atmos_dyn_tstep_large_setup(Tstep_large_type, DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG)
Register.
Definition: scale_atmos_dyn_tstep_large.F90:205
scale_prc_cartesc::prc_has_e
logical, public prc_has_e
Definition: scale_prc_cartesC.F90:50
scale_index::va
integer, public va
Definition: scale_index.F90:35
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_atmos_dyn::atmos_dyn_setup
subroutine, public atmos_dyn_setup(DYN_Tinteg_Short_TYPE, DYN_Tinteg_Tracer_TYPE, DYN_Tinteg_Large_TYPE, DYN_Tstep_Tracer_TYPE, DYN_Tstep_Large_TYPE, DYN_Tstep_Short_TYPE, DYN_FVM_FLUX_TYPE, DYN_FVM_FLUX_TYPE_TRACER, DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, CDZ, CDX, CDY, FDZ, FDX, FDY, wdamp_tau, wdamp_height, FZ, none)
Setup.
Definition: scale_atmos_dyn.F90:95
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_io
module STDIO
Definition: scale_io.F90:10
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_atmos_dyn_tstep_large::atmos_dyn_tstep_large_finalize
procedure(finalize), pointer, public atmos_dyn_tstep_large_finalize
Definition: scale_atmos_dyn_tstep_large.F90:183
scale_spnudge::spnudge_finalize
subroutine, public spnudge_finalize
Definition: scale_spnudge.F90:268
scale_debug::check
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:59
scale_atmos_dyn_common::atmos_dyn_wdamp_setup
subroutine, public atmos_dyn_wdamp_setup(wdamp_coef, wdamp_tau, wdamp_height, FZ)
Setup.
Definition: scale_atmos_dyn_common.F90:69
scale_prc_cartesc
module process / cartesC
Definition: scale_prc_cartesC.F90:11
scale_atmos_dyn_tstep_tracer
module Atmosphere / Dynamical scheme
Definition: scale_atmos_dyn_tstep_tracer.F90:12
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_dyn_fvm_numfilter::atmos_dyn_fvm_numfilter_setup
subroutine, public atmos_dyn_fvm_numfilter_setup(num_diff, num_diff_q, CDZ, CDX, CDY, FDZ, FDX, FDY)
Setup.
Definition: scale_atmos_dyn_fvm_numfilter.F90:104
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_dyn_fvm_numfilter
module Atmosphere / Dynamics common
Definition: scale_atmos_dyn_fvm_numfilter.F90:12
scale_precision::dp
integer, parameter, public dp
Definition: scale_precision.F90:32
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_atmos_dyn_tstep_tracer::atmos_dyn_tstep_tracer_setup
subroutine, public atmos_dyn_tstep_tracer_setup(ATMOS_DYN_TSTEP_TRACER_TYPE)
Register.
Definition: scale_atmos_dyn_tstep_tracer.F90:93
scale_atmos_dyn_tstep_short
module Atmosphere / Dynamical scheme
Definition: scale_atmos_dyn_tstep_short.F90:12
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_spnudge
Definition: scale_spnudge.F90:2
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_tinteg_tracer
module Atmosphere / Dynamics Temporal integration
Definition: scale_atmos_dyn_tinteg_tracer.F90:12
scale_atmos_dyn_tstep_large
module Atmosphere / Dynamical scheme
Definition: scale_atmos_dyn_tstep_large.F90:12
scale_debug
module DEBUG
Definition: scale_debug.F90:11
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.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_tinteg_large
module Atmosphere / Dynamics Temporal integration
Definition: scale_atmos_dyn_tinteg_large.F90:12
scale_prof::prof_rapend
subroutine, public prof_rapend(rapname_base, level, disable_barrier)
Save raptime.
Definition: scale_prof.F90:246
scale_atmos_dyn_tstep_short::atmos_dyn_tstep_short_setup
procedure(short_setup), pointer, public atmos_dyn_tstep_short_setup
Definition: scale_atmos_dyn_tstep_short.F90:133
scale_atmos_dyn_tinteg_short::atmos_dyn_tinteg_short_finalize
procedure(finalize), pointer, public atmos_dyn_tinteg_short_finalize
Definition: scale_atmos_dyn_tinteg_short.F90:114
scale_comm_cartesc::comm_vars8_init
subroutine, public comm_vars8_init(varname, var, vid, gid)
Register variables.
Definition: scale_comm_cartesC.F90:811
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_flux_setup
subroutine, public atmos_dyn_fvm_flux_setup(scheme, scheme_tracer)
setup
Definition: scale_atmos_dyn_fvm_flux.F90:306
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_setup
subroutine, public atmos_dyn_tinteg_tracer_setup(ATMOS_DYN_Tinteg_tracer_TYPE)
Register.
Definition: scale_atmos_dyn_tinteg_tracer.F90:99
scale_prc_cartesc::prc_has_w
logical, public prc_has_w
Definition: scale_prc_cartesC.F90:48
scale_atmos_dyn::atmos_dyn
subroutine, public atmos_dyn(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, CORIOLIS, CDZ, CDX, CDY, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, PHI, GSQRT, J13G, J23G, J33G, MAPF, AQ_R, AQ_CV, AQ_CP, AQ_MASS, REF_dens, REF_pott, REF_qv, REF_pres, ND_COEF, ND_COEF_Q, ND_LAPLACIAN_NUM, ND_SFC_FACT, ND_USE_RS, BND_QA, BND_IQ, BND_SMOOTHER_FACT, DAMP_DENS, DAMP_VELZ, DAMP_VELX, DAMP_VELY, DAMP_POTT, DAMP_QTRC, DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, MFLUX_OFFSET_X, MFLUX_OFFSET_Y, divdmp_coef, FLAG_TRACER_SPLIT_TEND, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM, USE_AVERAGE, I_QV, DTSEC, DTSEC_DYN)
Dynamical Process.
Definition: scale_atmos_dyn.F90:299
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_prc_cartesc::prc_twod
logical, public prc_twod
2D experiment
Definition: scale_prc_cartesC.F90:56
scale_atmos_dyn_fvm_numfilter::atmos_dyn_fvm_numfilter_finalize
subroutine, public atmos_dyn_fvm_numfilter_finalize
finalize
Definition: scale_atmos_dyn_fvm_numfilter.F90:396
scale_atmos_dyn::atmos_dyn_finalize
subroutine, public atmos_dyn_finalize
finalize
Definition: scale_atmos_dyn.F90:236