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  real(RP), public, allocatable :: coriolis(:,:) ! coriolis term
46 
47  !-----------------------------------------------------------------------------
48  !
49  !++ Private procedure
50  !
51  !-----------------------------------------------------------------------------
52  !
53  !++ Private parameters & variables
54  !
55  ! for communication
56  integer, private :: i_comm_dens = 1
57  integer, private :: i_comm_momz = 2
58  integer, private :: i_comm_momx = 3
59  integer, private :: i_comm_momy = 4
60  integer, private :: i_comm_rhot = 5
61  integer, private, allocatable :: i_comm_prog(:)
62  integer, private, allocatable :: i_comm_qtrc(:)
63 
64  logical, private :: bnd_w
65  logical, private :: bnd_e
66  logical, private :: bnd_s
67  logical, private :: bnd_n
68 
69  real(RP), private, allocatable :: mflx_hi (:,:,:,:) ! rho * vel(x,y,z) @ (u,v,w)-face high order
70 
71  real(RP), private, allocatable :: num_diff (:,:,:,:,:)
72  real(RP), private, allocatable :: num_diff_q(:,:,:,:)
73  real(RP), private, allocatable :: wdamp_coef(:) ! coefficient for Rayleigh damping of w
74 
75  logical, private :: dyn_none
76 
77  !-----------------------------------------------------------------------------
78 contains
79  !-----------------------------------------------------------------------------
81  subroutine atmos_dyn_setup( &
82  DYN_Tinteg_Short_TYPE, &
83  DYN_Tinteg_Tracer_TYPE, &
84  DYN_Tinteg_Large_TYPE, &
85  DYN_Tstep_Tracer_TYPE, &
86  DYN_Tstep_Large_TYPE, &
87  DYN_FVM_FLUX_TYPE, &
88  DYN_FVM_FLUX_TYPE_TRACER, &
89  DENS, MOMZ, MOMX, MOMY, RHOT, &
90  QTRC, PROG, &
91  CDZ, CDX, CDY, &
92  FDZ, FDX, FDY, &
93  wdamp_tau, &
94  wdamp_height, &
95  FZ, &
96  coriolis_type, &
97  coriolis_f0, &
98  coriolis_beta, &
99  coriolis_y0, &
100  CY, &
101  lat, &
102  none )
103  use scale_prc, only: &
104  prc_abort
105  use scale_prc_cartesc, only: &
106  prc_has_e, &
107  prc_has_w, &
108  prc_has_n, &
109  prc_has_s
110  use scale_const, only: &
111  ohm => const_ohm, &
112  undef => const_undef
113  use scale_comm_cartesc, only: &
115  use scale_atmos_dyn_common, only: &
118  use scale_atmos_dyn_tinteg_short, only: &
120  use scale_atmos_dyn_tinteg_tracer, only: &
122  use scale_atmos_dyn_tinteg_large, only: &
124  use scale_atmos_dyn_tstep_short, only: &
126  use scale_atmos_dyn_tstep_tracer, only: &
128  use scale_atmos_dyn_tstep_large, only: &
130  use scale_atmos_dyn_fvm_flux, only: &
132  implicit none
133 
134  character(len=*), intent(in) :: DYN_Tinteg_Short_TYPE
135  character(len=*), intent(in) :: DYN_Tinteg_Tracer_TYPE
136  character(len=*), intent(in) :: DYN_Tinteg_Large_TYPE
137  character(len=*), intent(in) :: DYN_Tstep_Tracer_TYPE
138  character(len=*), intent(in) :: DYN_Tstep_Large_TYPE
139  character(len=*), intent(in) :: DYN_FVM_FLUX_TYPE
140  character(len=*), intent(in) :: DYN_FVM_FLUX_TYPE_TRACER
141  real(RP), intent(inout) :: DENS(ka,ia,ja) ! MPI_RECV_INIT requires intent(inout)
142  real(RP), intent(inout) :: MOMZ(ka,ia,ja)
143  real(RP), intent(inout) :: MOMX(ka,ia,ja)
144  real(RP), intent(inout) :: MOMY(ka,ia,ja)
145  real(RP), intent(inout) :: RHOT(ka,ia,ja)
146  real(RP), intent(inout) :: QTRC(ka,ia,ja,qa)
147  real(RP), intent(inout) :: PROG(ka,ia,ja,va)
148  real(RP), intent(in) :: CDZ(ka)
149  real(RP), intent(in) :: CDX(ia)
150  real(RP), intent(in) :: CDY(ja)
151  real(RP), intent(in) :: FDZ(ka-1)
152  real(RP), intent(in) :: FDX(ia-1)
153  real(RP), intent(in) :: FDY(ja-1)
154  real(RP), intent(in) :: wdamp_tau
155  real(RP), intent(in) :: wdamp_height
156  real(RP), intent(in) :: FZ(0:ka)
157  character(len=*), intent(in) :: coriolis_type
158  real(RP), intent(in) :: coriolis_f0
159  real(RP), intent(in) :: coriolis_beta
160  real(RP), intent(in) :: coriolis_y0
161  real(RP), intent(in) :: CY(ja)
162  real(RP), intent(in) :: lat(ia,ja)
163  logical, optional, intent(in) :: none
164 
165  integer :: j
166  integer :: iv, iq
167  !---------------------------------------------------------------------------
168 
169  dyn_none = .false.
170 
171  if ( present(none) ) then
172  dyn_none = none
173  endif
174 
175  allocate( i_comm_prog(max(va,1)) )
176  allocate( i_comm_qtrc(qa) )
177 
178  if ( .NOT. dyn_none ) then
179  bnd_w = .NOT. prc_has_w
180  bnd_e = .NOT. prc_has_e
181  bnd_s = .NOT. prc_has_s
182  bnd_n = .NOT. prc_has_n
183 
184  allocate( coriolis(ia,ja) )
185  allocate( mflx_hi(ka,ia,ja,3) )
186  allocate( num_diff(ka,ia,ja,5,3) )
187  allocate( num_diff_q(ka,ia,ja,3) )
188  allocate( wdamp_coef(ka) )
189  mflx_hi(:,:,:,:) = undef
190  num_diff(:,:,:,:,:) = undef
191  num_diff_q(:,:,:,:) = undef
192 
193  call atmos_dyn_fvm_flux_setup ( dyn_fvm_flux_type, & ! [IN]
194  dyn_fvm_flux_type_tracer ) ! [IN]
195 
197 
198  call atmos_dyn_tstep_tracer_setup ( dyn_tstep_tracer_type ) ! [IN]
199 
200  call atmos_dyn_tstep_large_setup ( dyn_tstep_large_type, & ! [IN]
201  dens, momz, momx, momy, rhot, & ! [INOUT]
202  qtrc, prog, & ! [INOUT]
203  mflx_hi ) ! [INOUT]
204 
205  call atmos_dyn_tinteg_short_setup ( dyn_tinteg_short_type ) ! [IN]
206 
207  call atmos_dyn_tinteg_tracer_setup( dyn_tinteg_tracer_type ) ! [IN]
208 
209  call atmos_dyn_tinteg_large_setup ( dyn_tinteg_large_type ) ! [IN]
210 
211  ! numerical diffusion
212  call atmos_dyn_filter_setup( num_diff, num_diff_q, & ! [INOUT]
213  cdz, cdx, cdy, fdz, fdx, fdy ) ! [IN]
214 
215  ! numerical diffusion
216  call atmos_dyn_wdamp_setup( wdamp_coef(:), & ! [INOUT]
217  wdamp_tau, wdamp_height, & ! [IN]
218  fz(:) ) ! [IN]
219 
220  ! coriolis parameter
221  select case ( coriolis_type )
222  case ( 'PLANE' )
223  do j = 1, ja
224  coriolis(:,j) = coriolis_f0 + coriolis_beta * ( cy(j) - coriolis_y0 )
225  end do
226  case ( 'SPHERE' )
227  coriolis(:,:) = 2.0_rp * ohm * sin( lat(:,:) )
228  case default
229  log_error("ATMOS_DYN_setup",*) 'Coriolis type is invalid: ', trim(coriolis_type)
230  log_error_cont(*) 'The type must be PLANE or SPHERE'
231  call prc_abort
232  end select
233 
234  else
235 
236  call comm_vars8_init( 'DENS', dens, i_comm_dens )
237  call comm_vars8_init( 'MOMZ', momz, i_comm_momz )
238  call comm_vars8_init( 'MOMX', momx, i_comm_momx )
239  call comm_vars8_init( 'MOMY', momy, i_comm_momy )
240  call comm_vars8_init( 'RHOT', rhot, i_comm_rhot )
241 
242  do iv = 1, va
243  i_comm_prog(iv) = 5 + iv
244  call comm_vars8_init( 'PROG', prog(:,:,:,iv), i_comm_prog(iv) )
245  enddo
246 
247  do iq = 1, qa
248  i_comm_qtrc(iq) = 5 + va + iq
249  call comm_vars8_init( 'QTRC', qtrc(:,:,:,iq), i_comm_qtrc(iq) )
250  enddo
251 
252  endif
253 
254  return
255  end subroutine atmos_dyn_setup
256 
257  !-----------------------------------------------------------------------------
259  subroutine atmos_dyn( &
260  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, &
261  PROG, &
262  DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, &
263  DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, &
264  CDZ, CDX, CDY, FDZ, FDX, FDY, &
265  RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
266  PHI, GSQRT, &
267  J13G, J23G, J33G, MAPF, &
268  AQ_R, AQ_CV, AQ_CP, AQ_MASS, &
269  REF_dens, REF_pott, REF_qv, REF_pres, &
270  ND_COEF, ND_COEF_Q, ND_ORDER, ND_SFC_FACT, ND_USE_RS, &
271  BND_QA, BND_SMOOTHER_FACT, &
272  DAMP_DENS, DAMP_VELZ, DAMP_VELX, &
273  DAMP_VELY, DAMP_POTT, DAMP_QTRC, &
274  DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, &
275  DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, &
276  divdmp_coef, &
277  FLAG_TRACER_SPLIT_TEND, &
278  FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, &
279  FLAG_FCT_ALONG_STREAM, &
280  USE_AVERAGE, &
281  I_QV, &
282  DTSEC, DTSEC_DYN )
283  use scale_comm_cartesc, only: &
284  comm_vars8, &
285  comm_wait
286  use scale_atmos_dyn_tinteg_large, only: &
288  implicit none
289 
290  real(RP), intent(inout) :: DENS(ka,ia,ja)
291  real(RP), intent(inout) :: MOMZ(ka,ia,ja)
292  real(RP), intent(inout) :: MOMX(ka,ia,ja)
293  real(RP), intent(inout) :: MOMY(ka,ia,ja)
294  real(RP), intent(inout) :: RHOT(ka,ia,ja)
295  real(RP), intent(inout) :: QTRC(ka,ia,ja,qa)
296  real(RP), intent(inout) :: PROG(ka,ia,ja,va)
297 
298  real(RP), intent(inout) :: DENS_av(ka,ia,ja)
299  real(RP), intent(inout) :: MOMZ_av(ka,ia,ja)
300  real(RP), intent(inout) :: MOMX_av(ka,ia,ja)
301  real(RP), intent(inout) :: MOMY_av(ka,ia,ja)
302  real(RP), intent(inout) :: RHOT_av(ka,ia,ja)
303  real(RP), intent(inout) :: QTRC_av(ka,ia,ja,qa)
304 
305  real(RP), intent(in) :: DENS_tp(ka,ia,ja)
306  real(RP), intent(in) :: MOMZ_tp(ka,ia,ja)
307  real(RP), intent(in) :: MOMX_tp(ka,ia,ja)
308  real(RP), intent(in) :: MOMY_tp(ka,ia,ja)
309  real(RP), intent(in) :: RHOT_tp(ka,ia,ja)
310  real(RP), intent(in) :: RHOQ_tp(ka,ia,ja,qa)
311 
312  real(RP), intent(in) :: CDZ (ka)
313  real(RP), intent(in) :: CDX (ia)
314  real(RP), intent(in) :: CDY (ja)
315  real(RP), intent(in) :: FDZ (ka-1)
316  real(RP), intent(in) :: FDX (ia-1)
317  real(RP), intent(in) :: FDY (ja-1)
318  real(RP), intent(in) :: RCDZ(ka)
319  real(RP), intent(in) :: RCDX(ia)
320  real(RP), intent(in) :: RCDY(ja)
321  real(RP), intent(in) :: RFDZ(ka-1)
322  real(RP), intent(in) :: RFDX(ia-1)
323  real(RP), intent(in) :: RFDY(ja-1)
324 
325  real(RP), intent(in) :: PHI (ka,ia,ja)
326  real(RP), intent(in) :: GSQRT(ka,ia,ja,7)
327  real(RP), intent(in) :: J13G (ka,ia,ja,7)
328  real(RP), intent(in) :: J23G (ka,ia,ja,7)
329  real(RP), intent(in) :: J33G
330  real(RP), intent(in) :: MAPF (ia,ja,2,4)
331 
332  real(RP), intent(in) :: AQ_R (qa)
333  real(RP), intent(in) :: AQ_CV (qa)
334  real(RP), intent(in) :: AQ_CP (qa)
335  real(RP), intent(in) :: AQ_MASS(qa)
336 
337  real(RP), intent(in) :: REF_dens(ka,ia,ja)
338  real(RP), intent(in) :: REF_pott(ka,ia,ja)
339  real(RP), intent(in) :: REF_qv (ka,ia,ja)
340  real(RP), intent(in) :: REF_pres(ka,ia,ja)
341  real(RP), intent(in) :: ND_COEF
342  real(RP), intent(in) :: ND_COEF_Q
343  integer, intent(in) :: ND_ORDER
344  real(RP), intent(in) :: ND_SFC_FACT
345  logical, intent(in) :: ND_USE_RS
346 
347  integer, intent(in) :: BND_QA
348  real(RP), intent(in) :: BND_SMOOTHER_FACT
349 
350  real(RP), intent(in) :: DAMP_DENS(ka,ia,ja)
351  real(RP), intent(in) :: DAMP_VELZ(ka,ia,ja)
352  real(RP), intent(in) :: DAMP_VELX(ka,ia,ja)
353  real(RP), intent(in) :: DAMP_VELY(ka,ia,ja)
354  real(RP), intent(in) :: DAMP_POTT(ka,ia,ja)
355  real(RP), intent(in) :: DAMP_QTRC(ka,ia,ja,bnd_qa)
356 
357  real(RP), intent(in) :: DAMP_alpha_DENS(ka,ia,ja)
358  real(RP), intent(in) :: DAMP_alpha_VELZ(ka,ia,ja)
359  real(RP), intent(in) :: DAMP_alpha_VELX(ka,ia,ja)
360  real(RP), intent(in) :: DAMP_alpha_VELY(ka,ia,ja)
361  real(RP), intent(in) :: DAMP_alpha_POTT(ka,ia,ja)
362  real(RP), intent(in) :: DAMP_alpha_QTRC(ka,ia,ja,bnd_qa)
363 
364  real(RP), intent(in) :: divdmp_coef
365 
366  logical, intent(in) :: FLAG_TRACER_SPLIT_TEND
367  logical, intent(in) :: FLAG_FCT_MOMENTUM
368  logical, intent(in) :: FLAG_FCT_T
369  logical, intent(in) :: FLAG_FCT_TRACER
370  logical, intent(in) :: FLAG_FCT_ALONG_STREAM
371 
372  logical, intent(in) :: USE_AVERAGE
373 
374  integer, intent(in) :: I_QV
375 
376  real(DP), intent(in) :: DTSEC
377  real(DP), intent(in) :: DTSEC_DYN
378 
379  ! for time integartion
380  real(RP) :: DENS00 (ka,ia,ja) ! saved density before update
381 
382  ! For tracer advection
383  real(RP) :: tflx_hi(ka,ia,ja,3) ! rho * theta * vel(x,y,z) @ (u,v,w)-face high order
384  real(RP) :: dt
385 
386  integer :: i, j, k, iq
387  !---------------------------------------------------------------------------
388 
389  log_progress(*) 'atmosphere / dynamics'
390 
391  dt = real(dtsec, kind=rp)
392 
393  if ( dyn_none ) then
394 
395  call prof_rapstart("DYN_Tinteg", 2)
396 
397  do j = js, je
398  do i = is, ie
399  do k = ks, ke
400  dens00(k,i,j) = dens(k,i,j) ! save previous value
401 
402  dens(k,i,j) = dens(k,i,j) + dens_tp(k,i,j) * dt
403  enddo
404  enddo
405  enddo
406 
407  do j = js, je
408  do i = is, ie
409  do k = ks, ke
410  momz(k,i,j) = momz(k,i,j) + momz_tp(k,i,j) * dt
411  enddo
412  enddo
413  enddo
414 
415  do j = js, je
416  do i = is, ie
417  do k = ks, ke
418  momx(k,i,j) = momx(k,i,j) + momx_tp(k,i,j) * dt
419  enddo
420  enddo
421  enddo
422 
423  do j = js, je
424  do i = is, ie
425  do k = ks, ke
426  momy(k,i,j) = momy(k,i,j) + momy_tp(k,i,j) * dt
427  enddo
428  enddo
429  enddo
430 
431  do j = js, je
432  do i = is, ie
433  do k = ks, ke
434  rhot(k,i,j) = rhot(k,i,j) + rhot_tp(k,i,j) * dt
435  enddo
436  enddo
437  enddo
438 
439  do iq = 1, qa
440  do j = js, je
441  do i = is, ie
442  do k = ks, ke
443  qtrc(k,i,j,iq) = max( 0.0_rp, &
444  qtrc(k,i,j,iq) * dens00(k,i,j) + rhoq_tp(k,i,j,iq) * dt &
445  ) / dens(k,i,j)
446  enddo
447  enddo
448  enddo
449  enddo
450 
451  call comm_vars8( dens(:,:,:), i_comm_dens )
452  call comm_vars8( momz(:,:,:), i_comm_momz )
453  call comm_vars8( momx(:,:,:), i_comm_momx )
454  call comm_vars8( momy(:,:,:), i_comm_momy )
455  call comm_vars8( rhot(:,:,:), i_comm_rhot )
456  do iq = 1, qa
457  call comm_vars8( qtrc(:,:,:,iq), i_comm_qtrc(iq) )
458  enddo
459 
460  call comm_wait ( dens(:,:,:), i_comm_dens, .false. )
461  call comm_wait ( momz(:,:,:), i_comm_momz, .false. )
462  call comm_wait ( momx(:,:,:), i_comm_momx, .false. )
463  call comm_wait ( momy(:,:,:), i_comm_momy, .false. )
464  call comm_wait ( rhot(:,:,:), i_comm_rhot, .false. )
465  do iq = 1, qa
466  call comm_wait ( qtrc(:,:,:,iq), i_comm_qtrc(iq), .false. )
467  enddo
468 
469  if ( use_average ) then
470  dens_av(:,:,:) = dens(:,:,:)
471  momz_av(:,:,:) = momz(:,:,:)
472  momx_av(:,:,:) = momx(:,:,:)
473  momy_av(:,:,:) = momy(:,:,:)
474  rhot_av(:,:,:) = rhot(:,:,:)
475  qtrc_av(:,:,:,:) = qtrc(:,:,:,:)
476  endif
477 
478  call prof_rapend("DYN_Tinteg", 2)
479 
480  return
481  endif ! if DYN_NONE == .true.
482 
483  call prof_rapstart("DYN_Tinteg", 2)
484 
485  call atmos_dyn_tinteg_large( dens, momz, momx, momy, rhot, qtrc, & ! [INOUT]
486  prog, & ! [INOUT]
487  dens_av, momz_av, momx_av, momy_av, rhot_av, qtrc_av, & ! [INOUT]
488  mflx_hi, tflx_hi, & ! [OUT]
489  num_diff, num_diff_q, & ! [OUT;WORK]
490  dens_tp, momz_tp, momx_tp, momy_tp, rhot_tp, rhoq_tp, & ! [IN]
491  coriolis, & ! [IN]
492  cdz, cdx, cdy, fdz, fdx, fdy, & ! [IN]
493  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
494  phi, gsqrt, & ! [IN]
495  j13g, j23g, j33g, mapf, & ! [IN]
496  aq_r, aq_cv, aq_cp, aq_mass, & ! [IN]
497  ref_dens, ref_pott, ref_qv, ref_pres, & ! [IN]
498  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
499  nd_coef, nd_coef_q, nd_order, nd_sfc_fact, nd_use_rs, & ! [IN]
500  bnd_qa, bnd_smoother_fact, & ! [IN]
501  damp_dens, damp_velz, damp_velx, & ! [IN]
502  damp_vely, damp_pott, damp_qtrc, & ! [IN]
503  damp_alpha_dens, damp_alpha_velz, damp_alpha_velx, & ! [IN]
504  damp_alpha_vely, damp_alpha_pott, damp_alpha_qtrc, & ! [IN]
505  wdamp_coef, & ! [IN]
506  divdmp_coef, & ! [IN]
507  flag_tracer_split_tend, & ! [IN]
508  flag_fct_momentum, flag_fct_t, flag_fct_tracer, & ! [IN]
509  flag_fct_along_stream, & ! [IN]
510  use_average, & ! [IN]
511  i_qv, & ! [IN]
512  dtsec, dtsec_dyn ) ! [IN]
513 
514  call prof_rapend ("DYN_Tinteg", 2)
515 
516 #ifdef CHECK_MASS
517  call check_mass( dens, damp_dens, & ! [IN]
518  mflx_hi, tflx_hi, & ! [IN]
519  gsqrt, mapf, & ! [IN]
520  rcdx, rcdy, & ! [IN]
521  dt, & ! [IN]
522  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
523 #endif
524 
525  return
526  end subroutine atmos_dyn
527 
528 #ifdef CHECK_MASS
529  !-----------------------------------------------------------------------------
530  subroutine check_mass( &
531  DENS, DAMP_DENS, &
532  mflx_hi, tflx_hi, &
533  GSQRT, MAPF, &
534  RCDX, RCDY, &
535  dt, &
536  BND_W, BND_E, BND_S, BND_N )
537  use mpi
538  use scale_atmos_grid_cartesc_real, only: &
540  use scale_comm_cartesc, only: &
541  comm_datatype, &
542  comm_world
543  use scale_file_history, only: &
544  file_history_in
545  implicit none
546 
547  real(RP), intent(in) :: dens (ka,ia,ja)
548  real(RP), intent(in) :: DAMP_DENS(ka,ia,ja)
549  real(RP), intent(in) :: mflx_hi (ka,ia,ja,3)
550  real(RP), intent(in) :: tflx_hi (ka,ia,ja,3)
551  real(RP), intent(in) :: GSQRT (ka,ia,ja,7)
552  real(RP), intent(in) :: MAPF ( ia,ja,2,7)
553  real(RP), intent(in) :: RCDX(ia)
554  real(RP), intent(in) :: RCDY(ja)
555  real(RP), intent(in) :: dt
556  logical, intent(in) :: BND_W
557  logical, intent(in) :: BND_E
558  logical, intent(in) :: BND_S
559  logical, intent(in) :: BND_N
560 
561  ! lateral boundary flux
562  real(RP) :: mflx_lb_horizontal(ka)
563  real(RP) :: allmflx_lb_horizontal(ka)
564  real(RP) :: mflx_lb_total
565  real(RP) :: mass_total
566  real(RP) :: mass_total2
567  real(RP) :: allmflx_lb_total
568  real(RP) :: allmass_total
569  real(RP) :: allmass_total2
570 
571  integer :: k, i, j
572  integer :: ierr
573  !---------------------------------------------------------------------------
574 
575  call file_history_in(mflx_hi(:,:,:,zdir), 'MFLXZ', 'momentum flux of z-direction (w/ CHECK_MASS)', 'kg/m2/s', dim_type='ZHXY' )
576  call file_history_in(mflx_hi(:,:,:,xdir), 'MFLXX', 'momentum flux of x-direction (w/ CHECK_MASS)', 'kg/m2/s', dim_type='ZXHY' )
577  call file_history_in(mflx_hi(:,:,:,ydir), 'MFLXY', 'momentum flux of y-direction (w/ CHECK_MASS)', 'kg/m2/s', dim_type='ZXYH' )
578 
579  call file_history_in(tflx_hi(:,:,:,zdir), 'TFLXZ', 'potential temperature flux of z-direction (w/ CHECK_MASS)', 'K*kg/m2/s', dim_type='ZHXY' )
580  call file_history_in(tflx_hi(:,:,:,xdir), 'TFLXX', 'potential temperature flux of x-direction (w/ CHECK_MASS)', 'K*kg/m2/s', dim_type='ZXHY' )
581  call file_history_in(tflx_hi(:,:,:,ydir), 'TFLXY', 'potential temperature flux of y-direction (w/ CHECK_MASS)', 'K*kg/m2/s', dim_type='ZXYH' )
582 
583  mflx_lb_total = 0.0_rp
584  mflx_lb_horizontal(:) = 0.0_rp
585  allmflx_lb_horizontal(:) = 0.0_rp
586 
587  if ( bnd_w ) then ! for western boundary
588  i = is
589  do j = js, je
590  do k = ks, ke
591  mflx_lb_total = mflx_lb_total + mflx_hi(k,i-1,j,xdir) * rcdx(i) * vol(k,i,j) &
592  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
593  mflx_lb_horizontal(k) = mflx_lb_horizontal(k) + mflx_hi(k,i-1,j,xdir) * rcdx(i) * vol(k,i,j) &
594  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
595 
596  enddo
597  enddo
598  endif
599  if ( bnd_e ) then ! for eastern boundary
600  i = ie
601  do j = js, je
602  do k = ks, ke
603  mflx_lb_total = mflx_lb_total - mflx_hi(k,i,j,xdir) * rcdx(i) * vol(k,i,j) &
604  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
605  mflx_lb_horizontal(k) = mflx_lb_horizontal(k) - mflx_hi(k,i,j,xdir) * rcdx(i) * vol(k,i,j) &
606  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
607  enddo
608  enddo
609  endif
610  if ( bnd_s ) then ! for sourthern boundary
611  j = js
612  do i = is, ie
613  do k = ks, ke
614  mflx_lb_total = mflx_lb_total + mflx_hi(k,i,j-1,ydir) * rcdy(j) * vol(k,i,j) &
615  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
616  mflx_lb_horizontal(k) = mflx_lb_horizontal(k) + mflx_hi(k,i,j-1,ydir) * rcdy(j) * vol(k,i,j) &
617  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
618  enddo
619  enddo
620  endif
621  if ( bnd_n ) then ! for northern boundary
622  j = je
623  do i = is, ie
624  do k = ks, ke
625  mflx_lb_total = mflx_lb_total - mflx_hi(k,i,j,ydir) * rcdy(j) * vol(k,i,j) &
626  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
627  mflx_lb_horizontal(k) = mflx_lb_horizontal(k) - mflx_hi(k,i,j,ydir) * rcdy(j) * vol(k,i,j) &
628  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
629  enddo
630  enddo
631  endif
632 
633  mass_total = 0.0_rp
634  mass_total2 = 0.0_rp
635 
636  ! check total mass in the inner region
637  do j = js, je
638  do i = is, ie
639  do k = ks, ke
640  mass_total = mass_total + dens(k,i,j) * vol(k,i,j)
641  mass_total2 = mass_total2 + damp_dens(k,i,j) * vol(k,i,j)
642  enddo
643  enddo
644  enddo
645 
646  call mpi_allreduce( mflx_lb_total, &
647  allmflx_lb_total, &
648  1, &
649  comm_datatype, &
650  mpi_sum, &
651  comm_world, &
652  ierr )
653 
654  log_info("check_mass",'(A,1x,ES24.17)') 'total mflx_lb:', allmflx_lb_total
655 
656  call mpi_allreduce( mass_total, &
657  allmass_total, &
658  1, &
659  comm_datatype, &
660  mpi_sum, &
661  comm_world, &
662  ierr )
663 
664  log_info("check_mass",'(A,1x,ES24.17)') 'total mass :', allmass_total
665 
666  call mpi_allreduce( mass_total2, &
667  allmass_total2, &
668  1, &
669  comm_datatype, &
670  mpi_sum, &
671  comm_world, &
672  ierr )
673 
674  log_info("check_mass",'(A,1x,ES24.17)') 'total mass2 :', allmass_total2
675 
676  call mpi_allreduce( mflx_lb_horizontal(ks:ke), &
677  allmflx_lb_horizontal(ks:ke), &
678  kmax, &
679  comm_datatype, &
680  mpi_sum, &
681  comm_world, &
682  ierr )
683 
684  call file_history_in(allmflx_lb_horizontal(:), 'ALLMOM_lb_hz', &
685  'horizontally total momentum flux from lateral boundary', 'kg/m2/s' )
686 
687  return
688  end subroutine check_mass
689 #endif
690 
691 end module scale_atmos_dyn
module DEBUG
Definition: scale_debug.F90:11
integer, public comm_world
communication world ID
integer, public va
Definition: scale_index.F90:35
subroutine, public atmos_dyn_tinteg_large_setup(ATMOS_DYN_Tinteg_large_TYPE)
Register.
module Atmosphere / Dynamical scheme
module Atmosphere / Dynamics Temporal integration
integer, public ia
of whole cells: x, local, with HALO
module Atmosphere / Dynamics Temporal integration
subroutine, public atmos_dyn_tstep_tracer_setup(ATMOS_DYN_TSTEP_TRACER_TYPE)
Register.
integer, public comm_datatype
datatype of variable
module Atmosphere / Dynamics FENT + FCT
integer, public qa
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_vol
control volume (zxy) [m3]
integer, public ja
of whole cells: y, local, with HALO
module process / cartesC
module Atmosphere / Dynamical scheme
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:56
module Atmosphere / Dynamics Temporal integration
logical, public prc_has_s
logical, public prc_has_n
logical, public prc_has_e
subroutine, public atmos_dyn_wdamp_setup(wdamp_coef, wdamp_tau, wdamp_height, FZ)
Setup.
real(rp), public const_undef
Definition: scale_const.F90:41
module COMMUNICATION
procedure(large), pointer, public atmos_dyn_tinteg_large
integer, public is
start point of inner domain: x, local
real(rp), public const_ohm
angular velocity of the planet [1/s]
Definition: scale_const.F90:45
integer, public ie
end point of inner domain: x, local
module TRACER
module Index
Definition: scale_index.F90:11
module atmosphere / grid / cartesC index
integer, public ke
end point of inner domain: z, local
module PROCESS
Definition: scale_prc.F90:11
integer, public je
end point of inner domain: y, local
subroutine, public atmos_dyn_filter_setup(num_diff, num_diff_q, CDZ, CDX, CDY, FDZ, FDX, FDY)
Setup.
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:38
module Atmosphere / Dynamics common
integer, public ks
start point of inner domain: z, local
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, 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_ORDER, ND_SFC_FACT, ND_USE_RS, BND_QA, 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, 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.
integer, public kmax
of computational cells: z, local
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
integer, public js
start point of inner domain: y, local
real(rp), dimension(:,:), allocatable, public coriolis
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:157
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_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, coriolis_type, coriolis_f0, coriolis_beta, coriolis_y0, CY, lat, none)
Setup.
module profiler
Definition: scale_prof.F90:11
subroutine, public atmos_dyn_tstep_large_setup(Tstep_large_type, DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, mflx_hi)
Register.
subroutine, public atmos_dyn_fvm_flux_setup(scheme, scheme_tracer)
setup
module Atmosphere GRID CartesC Real(real space)
module scale_atmos_dyn_fvm_flux
module PRECISION
integer, public ka
of whole cells: z, local, with HALO
subroutine, public atmos_dyn_tinteg_short_setup(ATMOS_DYN_Tinteg_short_TYPE)
Register.
module STDIO
Definition: scale_io.F90:10
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:210
integer, parameter, public rp
procedure(short_setup), pointer, public atmos_dyn_tstep_short_setup
module file_history
subroutine, public atmos_dyn_tinteg_tracer_setup(ATMOS_DYN_Tinteg_tracer_TYPE)
Register.
module Atmosphere / Dynamical scheme
logical, public prc_has_w