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
40 
41  !-----------------------------------------------------------------------------
42  !
43  !++ Public parameters & variables
44  !
45  !-----------------------------------------------------------------------------
46  !
47  !++ Private procedure
48  !
49  !-----------------------------------------------------------------------------
50  !
51  !++ Private parameters & variables
52  !
53  ! for communication
54  integer, private :: I_COMM_DENS = 1
55  integer, private :: I_COMM_MOMZ = 2
56  integer, private :: I_COMM_MOMX = 3
57  integer, private :: I_COMM_MOMY = 4
58  integer, private :: I_COMM_RHOT = 5
59  integer, private, allocatable :: I_COMM_PROG(:)
60  integer, private, allocatable :: I_COMM_QTRC(:)
61 
62  logical, private :: BND_W
63  logical, private :: BND_E
64  logical, private :: BND_S
65  logical, private :: BND_N
66 
67  real(RP), private, allocatable :: num_diff (:,:,:,:,:)
68  real(RP), private, allocatable :: num_diff_q(:,:,:,:)
69  real(RP), private, allocatable :: wdamp_coef(:) ! coefficient for Rayleigh damping of w
70 
71  logical, private :: DYN_NONE
72 
73  !-----------------------------------------------------------------------------
74 contains
75  !-----------------------------------------------------------------------------
77  subroutine atmos_dyn_setup( &
78  DYN_Tinteg_Short_TYPE, &
79  DYN_Tinteg_Tracer_TYPE, &
80  DYN_Tinteg_Large_TYPE, &
81  DYN_Tstep_Tracer_TYPE, &
82  DYN_Tstep_Large_TYPE, &
83  DYN_Tstep_Short_TYPE, &
84  DYN_FVM_FLUX_TYPE, &
85  DYN_FVM_FLUX_TYPE_TRACER, &
86  DENS, MOMZ, MOMX, MOMY, RHOT, &
87  QTRC, PROG, &
88  CDZ, CDX, CDY, &
89  FDZ, FDX, FDY, &
90  wdamp_tau, &
91  wdamp_height, &
92  FZ, &
93  none )
94  use scale_prc, only: &
95  prc_abort
96  use scale_prc_cartesc, only: &
97  prc_twod, &
98  prc_has_e, &
99  prc_has_w, &
100  prc_has_n, &
101  prc_has_s
102  use scale_const, only: &
103  undef => const_undef
104  use scale_comm_cartesc, only: &
106  use scale_atmos_dyn_common, only: &
108  use scale_atmos_dyn_fvm_numfilter, only: &
110  use scale_atmos_dyn_tinteg_short, only: &
112  use scale_atmos_dyn_tinteg_tracer, only: &
114  use scale_atmos_dyn_tinteg_large, only: &
116  use scale_atmos_dyn_tstep_short, only: &
118  use scale_atmos_dyn_tstep_tracer, only: &
120  use scale_atmos_dyn_tstep_large, only: &
122  use scale_atmos_dyn_fvm_flux, only: &
124  use scale_spnudge, only: &
126  implicit none
127 
128  character(len=*), intent(in) :: dyn_tinteg_short_type
129  character(len=*), intent(in) :: dyn_tinteg_tracer_type
130  character(len=*), intent(in) :: dyn_tinteg_large_type
131  character(len=*), intent(in) :: dyn_tstep_tracer_type
132  character(len=*), intent(in) :: dyn_tstep_large_type
133  character(len=*), intent(in) :: dyn_tstep_short_type
134  character(len=*), intent(in) :: dyn_fvm_flux_type
135  character(len=*), intent(in) :: dyn_fvm_flux_type_tracer
136  real(rp), intent(inout) :: dens(ka,ia,ja) ! MPI_RECV_INIT requires intent(inout)
137  real(rp), intent(inout) :: momz(ka,ia,ja)
138  real(rp), intent(inout) :: momx(ka,ia,ja)
139  real(rp), intent(inout) :: momy(ka,ia,ja)
140  real(rp), intent(inout) :: rhot(ka,ia,ja)
141  real(rp), intent(inout) :: qtrc(ka,ia,ja,qa)
142  real(rp), intent(inout) :: prog(ka,ia,ja,va)
143  real(rp), intent(in) :: cdz(ka)
144  real(rp), intent(in) :: cdx(ia)
145  real(rp), intent(in) :: cdy(ja)
146  real(rp), intent(in) :: fdz(ka-1)
147  real(rp), intent(in) :: fdx(ia-1)
148  real(rp), intent(in) :: fdy(ja-1)
149  real(rp), intent(in) :: wdamp_tau
150  real(rp), intent(in) :: wdamp_height
151  real(rp), intent(in) :: fz(0:ka)
152  logical, optional, intent(in) :: none
153 
154  integer :: iv, iq
155  !---------------------------------------------------------------------------
156 
157  dyn_none = .false.
158 
159  if ( present(none) ) then
160  dyn_none = none
161  endif
162 
163  allocate( i_comm_prog(max(va,1)) )
164  allocate( i_comm_qtrc(qa) )
165 
166  if ( .NOT. dyn_none ) then
167  bnd_w = ( .NOT. prc_has_w ) .and. ( .NOT. prc_twod )
168  bnd_e = ( .NOT. prc_has_e ) .and. ( .NOT. prc_twod )
169  bnd_s = .NOT. prc_has_s
170  bnd_n = .NOT. prc_has_n
171 
172  allocate( num_diff(ka,ia,ja,5,3) )
173  allocate( num_diff_q(ka,ia,ja,3) )
174  allocate( wdamp_coef(ka) )
175  num_diff(:,:,:,:,:) = undef
176  num_diff_q(:,:,:,:) = undef
177 
178  call atmos_dyn_fvm_flux_setup ( dyn_fvm_flux_type, & ! [IN]
179  dyn_fvm_flux_type_tracer ) ! [IN]
180 
182 
183  call atmos_dyn_tstep_tracer_setup ( dyn_tstep_tracer_type ) ! [IN]
184 
185  call atmos_dyn_tstep_large_setup ( dyn_tstep_large_type, & ! [IN]
186  dens, momz, momx, momy, rhot, & ! [INOUT]
187  qtrc, prog ) ! [INOUT]
188 
189  call atmos_dyn_tinteg_short_setup ( dyn_tinteg_short_type, dyn_tstep_short_type ) ! [IN]
190 
191  call atmos_dyn_tinteg_tracer_setup( dyn_tinteg_tracer_type ) ! [IN]
192 
193  call atmos_dyn_tinteg_large_setup ( dyn_tinteg_large_type ) ! [IN]
194 
195  ! numerical diffusion
196  call atmos_dyn_fvm_numfilter_setup( num_diff, num_diff_q, & ! [INOUT]
197  cdz, cdx, cdy, fdz, fdx, fdy ) ! [IN]
198 
199  ! sponge layer
200  call atmos_dyn_wdamp_setup( wdamp_coef(:), & ! [INOUT]
201  wdamp_tau, wdamp_height, & ! [IN]
202  fz(:) ) ! [IN]
203 
204  ! setup for spectral nudging
205  call spnudge_setup( ka, ks, ke, ia, is, ie, ja, js, je )
206 
207 
208  else
209 
210  call comm_vars8_init( 'DENS', dens, i_comm_dens )
211  call comm_vars8_init( 'MOMZ', momz, i_comm_momz )
212  call comm_vars8_init( 'MOMX', momx, i_comm_momx )
213  call comm_vars8_init( 'MOMY', momy, i_comm_momy )
214  call comm_vars8_init( 'RHOT', rhot, i_comm_rhot )
215 
216  do iv = 1, va
217  i_comm_prog(iv) = 5 + iv
218  call comm_vars8_init( 'PROG', prog(:,:,:,iv), i_comm_prog(iv) )
219  enddo
220 
221  do iq = 1, qa
222  i_comm_qtrc(iq) = 5 + va + iq
223  call comm_vars8_init( 'QTRC', qtrc(:,:,:,iq), i_comm_qtrc(iq) )
224  enddo
225 
226  endif
227 
228  return
229  end subroutine atmos_dyn_setup
230 
231  !-----------------------------------------------------------------------------
233  subroutine atmos_dyn( &
234  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, &
235  PROG, &
236  DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, &
237  DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, &
238  CORIOLIS, &
239  CDZ, CDX, CDY, FDZ, FDX, FDY, &
240  RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
241  PHI, GSQRT, &
242  J13G, J23G, J33G, MAPF, &
243  AQ_R, AQ_CV, AQ_CP, AQ_MASS, &
244  REF_dens, REF_pott, REF_qv, REF_pres, &
245  ND_COEF, ND_COEF_Q, ND_LAPLACIAN_NUM, &
246  ND_SFC_FACT, ND_USE_RS, &
247  BND_QA, BND_IQ, BND_SMOOTHER_FACT, &
248  DAMP_DENS, DAMP_VELZ, DAMP_VELX, &
249  DAMP_VELY, DAMP_POTT, DAMP_QTRC, &
250  DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, &
251  DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, &
252  MFLUX_OFFSET_X, MFLUX_OFFSET_Y, &
253  divdmp_coef, &
254  FLAG_TRACER_SPLIT_TEND, &
255  FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, &
256  FLAG_FCT_ALONG_STREAM, &
257  USE_AVERAGE, &
258  I_QV, &
259  DTSEC, DTSEC_DYN )
260  use scale_prc_cartesc, only: &
261  prc_twod
262  use scale_comm_cartesc, only: &
263  comm_vars8, &
264  comm_wait
265  use scale_atmos_dyn_tinteg_large, only: &
267  implicit none
268 
269  real(rp), intent(inout) :: dens(ka,ia,ja)
270  real(rp), intent(inout) :: momz(ka,ia,ja)
271  real(rp), intent(inout) :: momx(ka,ia,ja)
272  real(rp), intent(inout) :: momy(ka,ia,ja)
273  real(rp), intent(inout) :: rhot(ka,ia,ja)
274  real(rp), intent(inout) :: qtrc(ka,ia,ja,qa)
275  real(rp), intent(inout) :: prog(ka,ia,ja,va)
276 
277  real(rp), intent(inout) :: dens_av(ka,ia,ja)
278  real(rp), intent(inout) :: momz_av(ka,ia,ja)
279  real(rp), intent(inout) :: momx_av(ka,ia,ja)
280  real(rp), intent(inout) :: momy_av(ka,ia,ja)
281  real(rp), intent(inout) :: rhot_av(ka,ia,ja)
282  real(rp), intent(inout) :: qtrc_av(ka,ia,ja,qa)
283 
284  real(rp), intent(in) :: dens_tp(ka,ia,ja)
285  real(rp), intent(in) :: momz_tp(ka,ia,ja)
286  real(rp), intent(in) :: momx_tp(ka,ia,ja)
287  real(rp), intent(in) :: momy_tp(ka,ia,ja)
288  real(rp), intent(in) :: rhot_tp(ka,ia,ja)
289  real(rp), intent(in) :: rhoq_tp(ka,ia,ja,qa)
290 
291  real(rp), intent(in) :: coriolis(ia,ja)
292 
293  real(rp), intent(in) :: cdz (ka)
294  real(rp), intent(in) :: cdx (ia)
295  real(rp), intent(in) :: cdy (ja)
296  real(rp), intent(in) :: fdz (ka-1)
297  real(rp), intent(in) :: fdx (ia-1)
298  real(rp), intent(in) :: fdy (ja-1)
299  real(rp), intent(in) :: rcdz(ka)
300  real(rp), intent(in) :: rcdx(ia)
301  real(rp), intent(in) :: rcdy(ja)
302  real(rp), intent(in) :: rfdz(ka-1)
303  real(rp), intent(in) :: rfdx(ia-1)
304  real(rp), intent(in) :: rfdy(ja-1)
305 
306  real(rp), intent(in) :: phi (ka,ia,ja)
307  real(rp), intent(in) :: gsqrt(ka,ia,ja,7)
308  real(rp), intent(in) :: j13g (ka,ia,ja,7)
309  real(rp), intent(in) :: j23g (ka,ia,ja,7)
310  real(rp), intent(in) :: j33g
311  real(rp), intent(in) :: mapf (ia,ja,2,4)
312 
313  real(rp), intent(in) :: aq_r (qa)
314  real(rp), intent(in) :: aq_cv (qa)
315  real(rp), intent(in) :: aq_cp (qa)
316  real(rp), intent(in) :: aq_mass(qa)
317 
318  real(rp), intent(in) :: ref_dens(ka,ia,ja)
319  real(rp), intent(in) :: ref_pott(ka,ia,ja)
320  real(rp), intent(in) :: ref_qv (ka,ia,ja)
321  real(rp), intent(in) :: ref_pres(ka,ia,ja)
322  real(rp), intent(in) :: nd_coef
323  real(rp), intent(in) :: nd_coef_q
324  integer, intent(in) :: nd_laplacian_num
325  real(rp), intent(in) :: nd_sfc_fact
326  logical, intent(in) :: nd_use_rs
327 
328  integer, intent(in) :: bnd_qa
329  integer, intent(in) :: bnd_iq(qa)
330  real(rp), intent(in) :: bnd_smoother_fact
331 
332  real(rp), intent(in) :: damp_dens(ka,ia,ja)
333  real(rp), intent(in) :: damp_velz(ka,ia,ja)
334  real(rp), intent(in) :: damp_velx(ka,ia,ja)
335  real(rp), intent(in) :: damp_vely(ka,ia,ja)
336  real(rp), intent(in) :: damp_pott(ka,ia,ja)
337  real(rp), intent(in) :: damp_qtrc(ka,ia,ja,bnd_qa)
338 
339  real(rp), intent(in) :: damp_alpha_dens(ka,ia,ja)
340  real(rp), intent(in) :: damp_alpha_velz(ka,ia,ja)
341  real(rp), intent(in) :: damp_alpha_velx(ka,ia,ja)
342  real(rp), intent(in) :: damp_alpha_vely(ka,ia,ja)
343  real(rp), intent(in) :: damp_alpha_pott(ka,ia,ja)
344  real(rp), intent(in) :: damp_alpha_qtrc(ka,ia,ja,bnd_qa)
345  real(rp), intent(in) :: mflux_offset_x(ka,ja,2)
346  real(rp), intent(in) :: mflux_offset_y(ka,ja,2)
347 
348  real(rp), intent(in) :: divdmp_coef
349 
350  logical, intent(in) :: flag_tracer_split_tend
351  logical, intent(in) :: flag_fct_momentum
352  logical, intent(in) :: flag_fct_t
353  logical, intent(in) :: flag_fct_tracer
354  logical, intent(in) :: flag_fct_along_stream
355 
356  logical, intent(in) :: use_average
357 
358  integer, intent(in) :: i_qv
359 
360  real(dp), intent(in) :: dtsec
361  real(dp), intent(in) :: dtsec_dyn
362 
363  ! for time integartion
364  real(rp) :: dens00 (ka,ia,ja) ! saved density before update
365 
366  ! For tracer advection
367  real(rp) :: dt
368 
369  integer :: i, j, k, iq
370  !---------------------------------------------------------------------------
371 
372  log_progress(*) 'atmosphere / dynamics'
373 
374  dt = real(dtsec, kind=rp)
375 
376  if ( dyn_none ) then
377 
378  call prof_rapstart("DYN_Tinteg", 2)
379 
380  do j = js, je
381  do i = is, ie
382  do k = ks, ke
383  dens00(k,i,j) = dens(k,i,j) ! save previous value
384 
385  dens(k,i,j) = dens(k,i,j) + dens_tp(k,i,j) * dt
386  enddo
387  enddo
388  enddo
389 
390  do j = js, je
391  do i = is, ie
392  do k = ks, ke
393  momz(k,i,j) = momz(k,i,j) + momz_tp(k,i,j) * dt
394  enddo
395  enddo
396  enddo
397 
398  do j = js, je
399  do i = is, ie
400  do k = ks, ke
401  momx(k,i,j) = momx(k,i,j) + momx_tp(k,i,j) * dt
402  enddo
403  enddo
404  enddo
405 
406  do j = js, je
407  do i = is, ie
408  do k = ks, ke
409  momy(k,i,j) = momy(k,i,j) + momy_tp(k,i,j) * dt
410  enddo
411  enddo
412  enddo
413 
414  do j = js, je
415  do i = is, ie
416  do k = ks, ke
417  rhot(k,i,j) = rhot(k,i,j) + rhot_tp(k,i,j) * dt
418  enddo
419  enddo
420  enddo
421 
422  do iq = 1, qa
423  do j = js, je
424  do i = is, ie
425  do k = ks, ke
426  qtrc(k,i,j,iq) = max( 0.0_rp, &
427  qtrc(k,i,j,iq) * dens00(k,i,j) + rhoq_tp(k,i,j,iq) * dt &
428  ) / dens(k,i,j)
429  enddo
430  enddo
431  enddo
432  enddo
433 
434  call comm_vars8( dens(:,:,:), i_comm_dens )
435  call comm_vars8( momz(:,:,:), i_comm_momz )
436  call comm_vars8( momx(:,:,:), i_comm_momx )
437  call comm_vars8( momy(:,:,:), i_comm_momy )
438  call comm_vars8( rhot(:,:,:), i_comm_rhot )
439  do iq = 1, qa
440  call comm_vars8( qtrc(:,:,:,iq), i_comm_qtrc(iq) )
441  enddo
442 
443  call comm_wait ( dens(:,:,:), i_comm_dens, .false. )
444  call comm_wait ( momz(:,:,:), i_comm_momz, .false. )
445  call comm_wait ( momx(:,:,:), i_comm_momx, .false. )
446  call comm_wait ( momy(:,:,:), i_comm_momy, .false. )
447  call comm_wait ( rhot(:,:,:), i_comm_rhot, .false. )
448  do iq = 1, qa
449  call comm_wait ( qtrc(:,:,:,iq), i_comm_qtrc(iq), .false. )
450  enddo
451 
452  if ( use_average ) then
453  dens_av(:,:,:) = dens(:,:,:)
454  momz_av(:,:,:) = momz(:,:,:)
455  momx_av(:,:,:) = momx(:,:,:)
456  momy_av(:,:,:) = momy(:,:,:)
457  rhot_av(:,:,:) = rhot(:,:,:)
458  qtrc_av(:,:,:,:) = qtrc(:,:,:,:)
459  endif
460 
461  call prof_rapend("DYN_Tinteg", 2)
462 
463  return
464  endif ! if DYN_NONE == .true.
465 
466  call prof_rapstart("DYN_Tinteg", 2)
467 
468  call atmos_dyn_tinteg_large( dens, momz, momx, momy, rhot, qtrc, & ! [INOUT]
469  prog, & ! [INOUT]
470  dens_av, momz_av, momx_av, momy_av, rhot_av, qtrc_av, & ! [INOUT]
471  num_diff, num_diff_q, & ! [OUT;WORK]
472  dens_tp, momz_tp, momx_tp, momy_tp, rhot_tp, rhoq_tp, & ! [IN]
473  coriolis, & ! [IN]
474  cdz, cdx, cdy, fdz, fdx, fdy, & ! [IN]
475  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
476  phi, gsqrt, & ! [IN]
477  j13g, j23g, j33g, mapf, & ! [IN]
478  aq_r, aq_cv, aq_cp, aq_mass, & ! [IN]
479  ref_dens, ref_pott, ref_qv, ref_pres, & ! [IN]
480  bnd_w, bnd_e, bnd_s, bnd_n, prc_twod, & ! [IN]
481  nd_coef, nd_coef_q, nd_laplacian_num, & ! [IN]
482  nd_sfc_fact, nd_use_rs, & ! [IN]
483  bnd_qa, bnd_iq, bnd_smoother_fact, & ! [IN]
484  damp_dens, damp_velz, damp_velx, & ! [IN]
485  damp_vely, damp_pott, damp_qtrc, & ! [IN]
486  damp_alpha_dens, damp_alpha_velz, damp_alpha_velx, & ! [IN]
487  damp_alpha_vely, damp_alpha_pott, damp_alpha_qtrc, & ! [IN]
488  mflux_offset_x, mflux_offset_y, & ! [IN]
489  wdamp_coef, divdmp_coef, & ! [IN]
490  flag_tracer_split_tend, & ! [IN]
491  flag_fct_momentum, flag_fct_t, flag_fct_tracer, & ! [IN]
492  flag_fct_along_stream, & ! [IN]
493  use_average, & ! [IN]
494  i_qv, & ! [IN]
495  dtsec, dtsec_dyn ) ! [IN]
496 
497  call prof_rapend ("DYN_Tinteg", 2)
498 
499  return
500  end subroutine atmos_dyn
501 
502 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:50
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_tracer::qa
integer, public qa
Definition: scale_tracer.F90:34
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:38
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:50
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:159
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:130
scale_prc_cartesc::prc_has_n
logical, public prc_has_n
Definition: scale_prc_cartesC.F90:48
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:201
scale_prc_cartesc::prc_has_e
logical, public prc_has_e
Definition: scale_prc_cartesC.F90:49
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:94
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:44
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_comm_cartesc::comm_vars8_init
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
Definition: scale_comm_cartesC.F90:294
scale_debug::check
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:56
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:103
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:217
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_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:281
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41
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:95
scale_prc_cartesc::prc_has_w
logical, public prc_has_w
Definition: scale_prc_cartesC.F90:47
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:260
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:55