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