54 real(RP),
private,
allocatable :: DENS_t(:,:,:)
55 real(RP),
private,
allocatable :: MOMZ_t(:,:,:)
56 real(RP),
private,
allocatable :: MOMX_t(:,:,:)
57 real(RP),
private,
allocatable :: MOMY_t(:,:,:)
58 real(RP),
private,
allocatable :: RHOT_t(:,:,:)
59 real(RP),
private,
allocatable,
target :: RHOQ_t(:,:,:,:)
60 real(RP),
private,
allocatable,
target :: ZERO(:,:,:)
61 real(RP),
private,
pointer :: RHOQ_tn(:,:,:)
63 real(RP),
private,
allocatable :: DENS_damp(:,:,:)
66 real(RP),
private,
allocatable,
target :: mflx(:,:,:,:)
69 integer :: I_COMM_DENS = 1
70 integer :: I_COMM_MOMZ = 2
71 integer :: I_COMM_MOMX = 3
72 integer :: I_COMM_MOMY = 4
73 integer :: I_COMM_RHOT = 5
74 integer,
allocatable :: I_COMM_PROG(:)
76 integer :: I_COMM_DENS_t = 1
77 integer :: I_COMM_MOMZ_t = 2
78 integer :: I_COMM_MOMX_t = 3
79 integer :: I_COMM_MOMY_t = 4
80 integer :: I_COMM_RHOT_t = 5
82 integer :: I_COMM_DENS_damp = 6
84 integer,
allocatable :: I_COMM_RHOQ_t(:)
85 integer,
allocatable :: I_COMM_QTRC(:)
87 integer :: I_COMM_mflx_z = 1
88 integer :: I_COMM_mflx_x = 2
89 integer :: I_COMM_mflx_y = 3
92 integer :: HIST_mflx(3)
93 integer :: HIST_tflx(3)
94 integer :: HIST_phys(5)
95 integer :: HIST_damp(5)
96 integer,
allocatable :: HIST_qflx(:,:)
97 integer,
allocatable :: HIST_phys_QTRC(:)
98 integer,
allocatable :: HIST_damp_QTRC(:)
101 real(RP),
allocatable,
target :: zero_x(:,:)
102 real(RP),
allocatable,
target :: zero_y(:,:)
103 integer :: MONIT_damp_mass
104 integer :: MONIT_damp_qtot
105 integer :: MONIT_mflx_west
106 integer :: MONIT_mflx_east
107 integer :: MONIT_mflx_south
108 integer :: MONIT_mflx_north
109 integer :: MONIT_qflx_west
110 integer :: MONIT_qflx_east
111 integer :: MONIT_qflx_south
112 integer :: MONIT_qflx_north
114 character(len=H_SHORT) :: EVAL_TYPE_NUMFILTER =
'TENDENCY'
122 DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG )
145 real(rp),
intent(inout) :: dens(
ka,
ia,
ja)
146 real(rp),
intent(inout) :: momz(
ka,
ia,
ja)
147 real(rp),
intent(inout) :: momx(
ka,
ia,
ja)
148 real(rp),
intent(inout) :: momy(
ka,
ia,
ja)
149 real(rp),
intent(inout) :: rhot(
ka,
ia,
ja)
150 real(rp),
intent(inout) :: qtrc(
ka,
ia,
ja,
qa)
151 real(rp),
intent(inout) :: prog(
ka,
ia,
ja,
va)
165 log_info(
"ATMOS_DYN_Tstep_large_fvm_heve_setup",*)
'Not found namelist. Default used.'
166 elseif( ierr > 0 )
then
167 log_error(
"ATMOS_DYN_Tstep_large_fvm_heve_setup",*)
'Not appropriate names in namelist ATMOS_DYN_TSTEP_LARGE_FVM_HEVE. Check!'
172 select case( eval_type_numfilter )
173 case(
'TENDENCY',
'FILTER' )
175 log_error(
"ATMOS_DYN_Tstep_large_fvm_heve_setup",*)
'The specfied value of EVAL_TYPE_NUMFILTER is not appropriate. Check!'
180 allocate( dens_t(
ka,
ia,
ja) )
181 allocate( momz_t(
ka,
ia,
ja) )
182 allocate( momx_t(
ka,
ia,
ja) )
183 allocate( momy_t(
ka,
ia,
ja) )
184 allocate( rhot_t(
ka,
ia,
ja) )
187 allocate( dens_damp(
ka,
ia,
ja) )
189 allocate( mflx(
ka,
ia,
ja,3) )
191 allocate( i_comm_prog(max(
va,1)) )
192 allocate( i_comm_qtrc(
qa) )
193 allocate( i_comm_rhoq_t(
qa) )
201 i_comm_prog(iv) = 5 + iv
212 i_comm_rhoq_t(iq) = 5 +
va + iq
213 i_comm_qtrc(iq) = 5 +
va + iq
225 allocate( zero(
ka,
ia,
ja) )
228 mflx(:,:,:,:) = undef
230 mflx(:,:,:,
xdir) = 0.0_rp
235 call file_history_reg(
'ZFLX_MOM',
'momentum flux of z-direction',
'kg/m2/s', &
238 call file_history_reg(
'XFLX_MOM',
'momentum flux of x-direction',
'kg/m2/s', &
241 call file_history_reg(
'YFLX_MOM',
'momentum flux of y-direction',
'kg/m2/s', &
245 call file_history_reg(
'ZFLX_RHOT',
'potential temperature flux of z-direction',
'K*kg/m2/s', &
248 call file_history_reg(
'XFLX_RHOT',
'potential temperature flux of x-direction',
'K*kg/m2/s', &
251 call file_history_reg(
'YFLX_RHOT',
'potential temperature flux of y-direction',
'K*kg/m2/s', &
255 call file_history_reg(
'DENS_t_phys',
'tendency of dencity due to physics',
'kg/m3/s', &
257 call file_history_reg(
'MOMZ_t_phys',
'tendency of momentum z due to physics',
'kg/m2/s2', &
260 call file_history_reg(
'MOMX_t_phys',
'tendency of momentum x due to physics',
'kg/m2/s2', &
263 call file_history_reg(
'MOMY_t_phys',
'tendency of momentum y due to physics',
'kg/m2/s2', &
266 call file_history_reg(
'RHOT_t_phys',
'tendency of rho*theta temperature due to physics',
'K*kg/m3/s', &
269 call file_history_reg(
'DENS_t_damp',
'tendency of dencity due to damping',
'kg/m3/s', &
271 call file_history_reg(
'MOMZ_t_damp',
'tendency of momentum z due to damping',
'kg/m2/s2', &
274 call file_history_reg(
'MOMX_t_damp',
'tendency of momentum x due to damping',
'kg/m2/s2', &
277 call file_history_reg(
'MOMY_t_damp',
'tendency of momentum y due to damping',
'kg/m2/s2', &
280 call file_history_reg(
'RHOT_t_damp',
'tendency of rho*theta temperature due to damping',
'K kg/m3/s', &
283 allocate( hist_qflx(3,
qa) )
284 allocate( hist_phys_qtrc(
qa) )
285 allocate( hist_damp_qtrc(
qa) )
304 call file_history_put( hist_mflx(iv), zero(:,:,:) )
305 call file_history_put( hist_tflx(iv), zero(:,:,:) )
308 call file_history_put( hist_phys(iv), zero(:,:,:) )
309 call file_history_put( hist_damp(iv), zero(:,:,:) )
313 call file_history_put( hist_qflx(iv,iq), zero(:,:,:) )
315 call file_history_put( hist_phys_qtrc(iq), zero(:,:,:) )
316 call file_history_put( hist_damp_qtrc(iq), zero(:,:,:) )
321 allocate( zero_x(
ka,
ja) )
322 allocate( zero_y(
ka,
ia) )
326 call monitor_reg(
"MASSTND_DAMP",
"mass tendency by the damping",
"kg", &
330 call monitor_reg(
"MASSFLX_WEST",
"mass flux at the western boundary",
"kg", &
332 dim_type=
"ZY-W", is_tendency=.true. )
333 call monitor_reg(
"MASSFLX_EAST",
"mass flux at the eastern boundary",
"kg", &
335 dim_type=
"ZY-E", is_tendency=.true. )
336 call monitor_reg(
"MASSFLX_SOUTH",
"mass flux at the southern boundary",
"kg", &
338 dim_type=
"ZX-S", is_tendency=.true. )
339 call monitor_reg(
"MASSFLX_NORTH",
"mass flux at the northern boundary",
"kg", &
341 dim_type=
"ZX-N", is_tendency=.true. )
343 call monitor_reg(
"QTOTTND_DAMP",
"water mass tendency by the damping",
"kg", &
347 call monitor_reg(
"QTOTFLX_WEST",
"water mass flux at the western boundary",
"kg", &
349 dim_type=
"ZY-W", is_tendency=.true. )
350 call monitor_reg(
"QTOTFLX_EAST",
"water mass flux at the eastern boundary",
"kg", &
352 dim_type=
"ZY-E", is_tendency=.true. )
353 call monitor_reg(
"QTOTFLX_SOUTH",
"water mass flux at the southern boundary",
"kg", &
355 dim_type=
"ZX-S", is_tendency=.true. )
356 call monitor_reg(
"QTOTFLX_NORTH",
"water mass flux at the northern boundary",
"kg", &
358 dim_type=
"ZX-N", is_tendency=.true. )
361 call monitor_put( monit_damp_mass, zero(:,:,:) )
362 call monitor_put( monit_mflx_west, zero_x(:,:) )
363 call monitor_put( monit_mflx_east, zero_x(:,:) )
364 call monitor_put( monit_mflx_south, zero_y(:,:) )
365 call monitor_put( monit_mflx_north, zero_y(:,:) )
367 call monitor_put( monit_damp_qtot, zero(:,:,:) )
368 call monitor_put( monit_qflx_west, zero_x(:,:) )
369 call monitor_put( monit_qflx_east, zero_x(:,:) )
370 call monitor_put( monit_qflx_south, zero_y(:,:) )
371 call monitor_put( monit_qflx_north, zero_y(:,:) )
379 DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, &
380 DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, &
381 num_diff, num_diff_q, &
383 DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, &
385 CDZ, CDX, CDY, FDZ, FDX, FDY, &
386 RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
388 J13G, J23G, J33G, MAPF, &
389 AQ_R, AQ_CV, AQ_CP, AQ_MASS, &
390 REF_dens, REF_pott, REF_qv, REF_pres, &
391 BND_W, BND_E, BND_S, BND_N, TwoD, &
392 ND_COEF, ND_COEF_Q, ND_LAPLACIAN_NUM, &
393 ND_SFC_FACT, ND_USE_RS, &
394 BND_QA, BND_IQ, BND_SMOOTHER_FACT, &
395 DAMP_DENS, DAMP_VELZ, DAMP_VELX, &
396 DAMP_VELY, DAMP_POTT, DAMP_QTRC, &
397 DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, &
398 DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, &
399 MFLUX_OFFSET_X, MFLUX_OFFSET_Y, &
400 wdamp_coef, divdmp_coef, &
401 FLAG_TRACER_SPLIT_TEND, &
402 FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, &
403 FLAG_FCT_ALONG_STREAM, &
421 file_history_query, &
450 real(rp),
intent(inout) :: dens(
ka,
ia,
ja)
451 real(rp),
intent(inout) :: momz(
ka,
ia,
ja)
452 real(rp),
intent(inout) :: momx(
ka,
ia,
ja)
453 real(rp),
intent(inout) :: momy(
ka,
ia,
ja)
454 real(rp),
intent(inout) :: rhot(
ka,
ia,
ja)
455 real(rp),
intent(inout) :: qtrc(
ka,
ia,
ja,
qa)
456 real(rp),
intent(inout) :: prog(
ka,
ia,
ja,
va)
458 real(rp),
intent(inout) :: dens_av(
ka,
ia,
ja)
459 real(rp),
intent(inout) :: momz_av(
ka,
ia,
ja)
460 real(rp),
intent(inout) :: momx_av(
ka,
ia,
ja)
461 real(rp),
intent(inout) :: momy_av(
ka,
ia,
ja)
462 real(rp),
intent(inout) :: rhot_av(
ka,
ia,
ja)
463 real(rp),
intent(inout) :: qtrc_av(
ka,
ia,
ja,
qa)
465 real(rp),
intent(out) :: num_diff(
ka,
ia,
ja,5,3)
466 real(rp),
intent(out) :: num_diff_q(
ka,
ia,
ja,3)
468 real(rp),
intent(in) :: qtrc0(
ka,
ia,
ja,
qa)
470 real(rp),
intent(in) :: dens_tp(
ka,
ia,
ja)
471 real(rp),
intent(in) :: momz_tp(
ka,
ia,
ja)
472 real(rp),
intent(in) :: momx_tp(
ka,
ia,
ja)
473 real(rp),
intent(in) :: momy_tp(
ka,
ia,
ja)
474 real(rp),
intent(in) :: rhot_tp(
ka,
ia,
ja)
475 real(rp),
intent(in) :: rhoq_tp(
ka,
ia,
ja,
qa)
477 real(rp),
intent(in) :: corioli(
ia,
ja)
479 real(rp),
intent(in) :: cdz (
ka)
480 real(rp),
intent(in) :: cdx (
ia)
481 real(rp),
intent(in) :: cdy (
ja)
482 real(rp),
intent(in) :: fdz (
ka-1)
483 real(rp),
intent(in) :: fdx (
ia-1)
484 real(rp),
intent(in) :: fdy (
ja-1)
485 real(rp),
intent(in) :: rcdz(
ka)
486 real(rp),
intent(in) :: rcdx(
ia)
487 real(rp),
intent(in) :: rcdy(
ja)
488 real(rp),
intent(in) :: rfdz(
ka-1)
489 real(rp),
intent(in) :: rfdx(
ia-1)
490 real(rp),
intent(in) :: rfdy(
ja-1)
492 real(rp),
intent(in) :: phi (
ka,
ia,
ja)
493 real(rp),
intent(in) :: gsqrt(
ka,
ia,
ja,7)
494 real(rp),
intent(in) :: j13g (
ka,
ia,
ja,7)
495 real(rp),
intent(in) :: j23g (
ka,
ia,
ja,7)
496 real(rp),
intent(in) :: j33g
497 real(rp),
intent(in) :: mapf (
ia,
ja,2,4)
499 real(rp),
intent(in) :: aq_r (
qa)
500 real(rp),
intent(in) :: aq_cv (
qa)
501 real(rp),
intent(in) :: aq_cp (
qa)
502 real(rp),
intent(in) :: aq_mass(
qa)
504 real(rp),
intent(in) :: ref_dens(
ka,
ia,
ja)
505 real(rp),
intent(in) :: ref_pott(
ka,
ia,
ja)
506 real(rp),
intent(in) :: ref_qv (
ka,
ia,
ja)
507 real(rp),
intent(in) :: ref_pres(
ka,
ia,
ja)
509 logical,
intent(in) :: bnd_w
510 logical,
intent(in) :: bnd_e
511 logical,
intent(in) :: bnd_s
512 logical,
intent(in) :: bnd_n
513 logical,
intent(in) :: twod
515 real(rp),
intent(in) :: nd_coef
516 real(rp),
intent(in) :: nd_coef_q
517 integer,
intent(in) :: nd_laplacian_num
518 real(rp),
intent(in) :: nd_sfc_fact
519 logical,
intent(in) :: nd_use_rs
521 integer,
intent(in) :: bnd_qa
522 integer,
intent(in) :: bnd_iq(
qa)
523 real(rp),
intent(in) :: bnd_smoother_fact
525 real(rp),
intent(in) :: damp_dens(
ka,
ia,
ja)
526 real(rp),
intent(in) :: damp_velz(
ka,
ia,
ja)
527 real(rp),
intent(in) :: damp_velx(
ka,
ia,
ja)
528 real(rp),
intent(in) :: damp_vely(
ka,
ia,
ja)
529 real(rp),
intent(in) :: damp_pott(
ka,
ia,
ja)
530 real(rp),
intent(in) :: damp_qtrc(
ka,
ia,
ja,bnd_qa)
532 real(rp),
intent(in) :: damp_alpha_dens(
ka,
ia,
ja)
533 real(rp),
intent(in) :: damp_alpha_velz(
ka,
ia,
ja)
534 real(rp),
intent(in) :: damp_alpha_velx(
ka,
ia,
ja)
535 real(rp),
intent(in) :: damp_alpha_vely(
ka,
ia,
ja)
536 real(rp),
intent(in) :: damp_alpha_pott(
ka,
ia,
ja)
537 real(rp),
intent(in) :: damp_alpha_qtrc(
ka,
ia,
ja,bnd_qa)
538 real(rp),
intent(in) :: mflux_offset_x(
ka,
ja,2)
539 real(rp),
intent(in) :: mflux_offset_y(
ka,
ia,2)
541 real(rp),
intent(in) :: wdamp_coef(
ka)
542 real(rp),
intent(in) :: divdmp_coef
544 logical,
intent(in) :: flag_tracer_split_tend
545 logical,
intent(in) :: flag_fct_momentum
546 logical,
intent(in) :: flag_fct_t
547 logical,
intent(in) :: flag_fct_tracer
548 logical,
intent(in) :: flag_fct_along_stream
550 logical,
intent(in) :: use_average
552 integer,
intent(in) :: i_qv
554 real(
dp),
intent(in) :: dtls
555 real(
dp),
intent(in) :: dtss
556 logical ,
intent(in) :: llast
559 real(rp) :: dens00 (
ka,
ia,
ja)
560 real(rp) :: qflx (
ka,
ia,
ja,3)
563 real(rp) :: ddiv (
ka,
ia,
ja)
564 real(rp) :: dpres0 (
ka,
ia,
ja)
565 real(rp) :: rt2p (
ka,
ia,
ja)
566 real(rp) :: ref_rhot(
ka,
ia,
ja)
568 real(rp) :: dens_tq(
ka,
ia,
ja)
571 real(rp) :: damp_t_dens(
ka,
ia,
ja)
572 real(rp) :: damp_t_momz(
ka,
ia,
ja)
573 real(rp) :: damp_t_momx(
ka,
ia,
ja)
574 real(rp) :: damp_t_momy(
ka,
ia,
ja)
575 real(rp) :: damp_t_rhot(
ka,
ia,
ja)
576 real(rp) :: damp_t_qtrc(
ka,
ia,
ja)
578 real(rp) :: tflx(
ka,
ia,
ja,3)
581 real(rp) :: mflx_av(
ka,
ia,
ja,3)
591 real(rp),
target :: qflx_west (
ka,
ja)
592 real(rp),
target :: qflx_east (
ka,
ja)
593 real(rp),
target :: qflx_south(
ka,
ia)
594 real(rp),
target :: qflx_north(
ka,
ia)
595 logical :: monit_lateral_flag(3)
597 integer :: i, j,
k, iq, iqb, step
604 dts = real(dtss, kind=rp)
605 dtl = real(dtls, kind=rp)
606 nstep = ceiling( ( dtl - eps ) / dts )
609 monit_lateral_flag(
zdir) = .false.
610 monit_lateral_flag(
xdir) = monit_mflx_west > 0 .or. monit_mflx_east > 0
611 monit_lateral_flag(
ydir) = monit_mflx_south > 0 .or. monit_mflx_north > 0
614 log_info(
"ATMOS_DYN_Tstep_large_fvm_heve",*)
'Dynamics large time step'
615 log_info_cont(
'(1x,A,F0.2,A,F0.2,A,I0)') &
616 '-> DT_large, DT_small, DT_large/DT_small : ', dtl,
', ', dts,
', ', nstep
618 dens00(:,:,:) = undef
619 num_diff(:,:,:,:,:) = undef
623 dens00(:,:,:) = dens(:,:,:)
625 if ( use_average )
then
631 dens_av(
k,i,j) = 0.0_rp
632 momz_av(
k,i,j) = 0.0_rp
633 momx_av(
k,i,j) = 0.0_rp
634 momy_av(
k,i,j) = 0.0_rp
635 rhot_av(
k,i,j) = 0.0_rp
642 mflx_av(:,:,:,:) = 0.0_rp
649 qflx_west(
k,j) = 0.0_rp
650 qflx_east(
k,j) = 0.0_rp
657 qflx_south(
k,i) = 0.0_rp
658 qflx_north(
k,i) = 0.0_rp
664 dpres0, rt2p, ref_rhot, &
665 rhot, qtrc, ref_pres, aq_r, aq_cv, aq_cp, aq_mass )
676 dens_tq(:,:,:) = 0.0_rp
687 do i = max(
is-1,1), min(
ie+1,
ia)
689 diff(
k,i,j) = qtrc(
k,i,j,iq) - damp_qtrc(
k,i,j,iqb)
694 call file_history_query( hist_damp_qtrc(iq), do_put )
704 damp = - damp_alpha_qtrc(
k,
is,j,iqb) &
706 - ( diff(
k,
is,j-1) + diff(
k,
is,j+1) - diff(
k,
is,j)*2.0_rp ) &
707 * 0.25_rp * bnd_smoother_fact )
708 damp = damp * dens00(
k,
is,j)
709 if ( do_put ) damp_t_qtrc(
k,
is,j) = damp
710 rhoq_t(
k,
is,j,iq) = rhoq_tp(
k,
is,j,iq) + damp
722 diff2(
k,i,j) = diff(
k,i,j)
727 call dft_g2g(
ka,
ks,
ke,
ia,
is,
ie,
ja,
js,
je,
spnudge_qv_lm,
spnudge_qv_mm,diff2)
752 damp = - damp_alpha_qtrc(
k,i,j,iqb) &
754 - ( diff(
k,i-1,j) + diff(
k,i+1,j) + diff(
k,i,j-1) + diff(
k,i,j+1) - diff(
k,i,j)*4.0_rp ) &
755 * 0.125_rp * bnd_smoother_fact ) &
757 damp = damp * dens00(
k,i,j)
758 if ( do_put ) damp_t_qtrc(
k,i,j) = damp
759 rhoq_t(
k,i,j,iq) = rhoq_tp(
k,i,j,iq) + damp
767 if ( do_put )
call file_history_put( hist_damp_qtrc(iq), damp_t_qtrc(:,:,:) )
768 call file_history_put( hist_phys_qtrc(iq), rhoq_tp(:,:,:,iq) )
772 call comm_vars8( rhoq_t(:,:,:,iq), i_comm_rhoq_t(iq) )
773 call comm_wait ( rhoq_t(:,:,:,iq), i_comm_rhoq_t(iq), .false. )
782 rhoq_t(
k,i,j,iq) = rhoq_tp(
k,i,j,iq)
795 if ( bnd_w .and. (.not. twod) )
then
799 mflx(
k,
is-1,j,
xdir) = ( momx(
k,
is-1,j) + mflux_offset_x(
k,j,1) ) &
804 if ( bnd_e .and. (.not. twod) )
then
808 mflx(
k,
ie,j,
xdir) = ( momx(
k,
ie,j) + mflux_offset_x(
k,j,2) ) &
817 mflx(
k,i,
js-1,
ydir) = ( momy(
k,i,
js-1) + mflux_offset_y(
k,i,1) ) &
826 mflx(
k,i,
je,
ydir) = ( momy(
k,i,
je) + mflux_offset_y(
k,i,2) ) &
839 damp_t_dens(
k,i,j) = 0.0_rp
840 damp_t_momz(
k,i,j) = 0.0_rp
841 damp_t_momx(
k,i,j) = 0.0_rp
842 damp_t_momy(
k,i,j) = 0.0_rp
843 damp_t_rhot(
k,i,j) = 0.0_rp
860 do i = max(
is-1,1), min(
ie+1,
ia)
862 diff(
k,i,j) = dens(
k,i,j) - damp_dens(
k,i,j)
867 call file_history_query( hist_damp(1), do_put )
877 dens_damp(
k,
is,j) = - damp_alpha_dens(
k,
is,j) &
879 - ( diff(
k,
is,j-1) + diff(
k,
is,j+1) - diff(
k,
is,j)*2.0_rp ) &
880 * 0.25_rp * bnd_smoother_fact ) &
881 + dens_tq(
k,
is,j) * ( 0.5_rp - sign( 0.5_rp, damp_alpha_dens(
k,
is,j)-eps ) )
882 dens_t(
k,
is,j) = dens_tp(
k,
is,j) &
884 if ( do_put ) damp_t_dens(
k,
is,j) = damp_t_dens(
k,
is,j) + dens_damp(
k,
is,j) / nstep
897 dens_damp(
k,i,j) = - damp_alpha_dens(
k,i,j) &
899 - ( diff(
k,i-1,j) + diff(
k,i+1,j) + diff(
k,i,j-1) + diff(
k,i,j+1) - diff(
k,i,j)*4.0_rp ) &
900 * 0.125_rp * bnd_smoother_fact ) &
901 + dens_tq(
k,i,j) * ( 0.5_rp - sign( 0.5_rp, damp_alpha_dens(
k,i,j)-eps ) )
902 dens_t(
k,i,j) = dens_tp(
k,i,j) &
904 if ( do_put ) damp_t_dens(
k,i,j) = damp_t_dens(
k,i,j) + dens_damp(
k,i,j) / nstep
910 call comm_vars8( dens_damp(:,:,:), i_comm_dens_damp )
911 call comm_vars8( dens_t(:,:,:), i_comm_dens_t )
916 do i = max(
is-1,1), min(
ie+1,
ia)
918 diff(
k,i,j) = momz(
k,i,j) - damp_velz(
k,i,j) * ( dens(
k,i,j)+dens(
k+1,i,j) ) * 0.5_rp
923 call file_history_query( hist_damp(2), do_put )
934 damp = - damp_alpha_velz(
k,
is,j) &
936 - ( diff(
k,
is,j-1) + diff(
k,
is,j+1) - diff(
k,
is,j)*2.0_rp ) &
937 * 0.25_rp * bnd_smoother_fact )
939 momz_t(
k,
is,j) = momz_tp(
k,
is,j) &
941 + ( dens_damp(
k,
is,j) + dens_damp(
k+1,
is,j) ) * momz(
k,
is,j) / ( dens(
k,
is,j) + dens(
k,
is,j) )
942 if ( do_put ) damp_t_momz(
k,
is,j) = damp_t_momz(
k,
is,j) + damp / nstep
956 damp = - damp_alpha_velz(
k,i,j) &
958 - ( diff(
k,i-1,j) + diff(
k,i+1,j) + diff(
k,i,j-1) + diff(
k,i,j+1) - diff(
k,i,j)*4.0_rp ) &
959 * 0.125_rp * bnd_smoother_fact )
961 momz_t(
k,i,j) = momz_tp(
k,i,j) &
963 + ( dens_damp(
k,i,j) + dens_damp(
k+1,i,j) ) * momz(
k,i,j) / ( dens(
k,i,j) + dens(
k,i,j) )
964 if ( do_put ) damp_t_momz(
k,i,j) = damp_t_momz(
k,i,j) + damp / nstep
972 momz_t( 1:
ks-1,i,j) = 0.0_rp
973 momz_t(
ke:
ka ,i,j) = 0.0_rp
976 call comm_vars8( momz_t(:,:,:), i_comm_momz_t )
978 call comm_wait( dens_damp(:,:,:), i_comm_dens_damp )
985 diff(
k,
is,j) = momx(
k,
is,j) - damp_velx(
k,
is,j) * dens(
k,
is,j)
994 diff(
k,i,j) = momx(
k,i,j) - damp_velx(
k,i,j) * ( dens(
k,i,j)+dens(
k,i+1,j) ) * 0.5_rp
1001 call file_history_query( hist_damp(3), do_put )
1012 damp = - damp_alpha_velx(
k,
is,j) &
1014 - ( diff(
k,
is,j-1) + diff(
k,
is,j+1) - diff(
k,
is,j)*2.0_rp ) &
1015 * 0.25_rp * bnd_smoother_fact )
1016 momx_t(
k,
is,j) = momx_tp(
k,
is,j) &
1018 + dens_damp(
k,
is,j) * momx(
k,
is,j) / dens(
k,
is,j)
1019 if ( do_put ) damp_t_momx(
k,
is,j) = damp_t_momx(
k,
is,j) + damp / nstep
1031 diff3(
k,i,j) = momy(
k,i,j) - damp_vely(
k,i,j) * ( dens(
k,i,j)+dens(
k,i,j+1) ) * 0.5_rp
1041 diff2(
k,i,j) = diff(
k,i,j) / ( ( dens(
k,i,j)+dens(
k,i+1,j) ) * 0.5_rp )
1051 diff3(
k,i,j) = diff3(
k,i,j) / ( ( dens(
k,i,j)+dens(
k,i,j+1) ) * 0.5_rp )
1056 call dft_g2g_divfree(
ka,
ks,
ke,
ia,
is,
ie,
ja,
js,
je,
spnudge_uv_lm,
spnudge_uv_mm,diff2,diff3)
1063 diff2(
k,i,j) = diff2(
k,i,j) * ( dens(
k,i,j)+dens(
k,i+1,j) ) * 0.5_rp
1075 diff2(
k,i,j) = diff(
k,i,j) / ( ( dens(
k,i,j)+dens(
k,i+1,j) ) * 0.5_rp )
1080 call dft_g2g(
ka,
ks,
ke,
ia,
is,
ie,
ja,
js,
je,
spnudge_uv_lm,
spnudge_uv_mm,diff2)
1086 diff2(
k,i,j) = diff2(
k,i,j) * ( dens(
k,i,j)+dens(
k,i+1,j) ) * 0.5_rp
1116 damp = - damp_alpha_velx(
k,i,j) &
1118 - ( diff(
k,i-1,j) + diff(
k,i+1,j) + diff(
k,i,j-1) + diff(
k,i,j+1) - diff(
k,i,j)*4.0_rp ) &
1119 * 0.125_rp * bnd_smoother_fact ) &
1121 momx_t(
k,i,j) = momx_tp(
k,i,j) &
1123 + ( dens_damp(
k,i,j) + dens_damp(
k,i+1,j) ) * momx(
k,i,j) / ( dens(
k,i,j) + dens(
k,i+1,j) )
1124 if ( do_put ) damp_t_momx(
k,i,j) = damp_t_momx(
k,i,j) + damp / nstep
1131 call comm_vars8( momx_t(:,:,:), i_comm_momx_t )
1136 do i = max(
is-1,1), min(
ie+1,
ia)
1138 diff(
k,i,j) = momy(
k,i,j) - damp_vely(
k,i,j) * ( dens(
k,i,j)+dens(
k,i,j+1) ) * 0.5_rp
1143 call file_history_query( hist_damp(4), do_put )
1154 damp = - damp_alpha_vely(
k,
is,j) &
1156 - ( diff(
k,
is,j-1) + diff(
k,
is,j+1) - diff(
k,
is,j)*2.0_rp ) &
1157 * 0.25_rp * bnd_smoother_fact )
1158 momy_t(
k,
is,j) = momy_tp(
k,
is,j) &
1160 + ( dens_damp(
k,
is,j) + dens_damp(
k,
is,j+1) ) * momy(
k,
is,j) / ( dens(
k,
is,j) + dens(
k,
is,j+1) )
1161 if ( do_put ) damp_t_momy(
k,
is,j) = damp_t_momy(
k,
is,j) + damp / nstep
1173 diff2(
k,i,j) = diff3(
k,i,j) * ( dens(
k,i,j)+dens(
k,i,j+1) ) * 0.5_rp
1184 diff2(
k,i,j) = diff(
k,i,j) / ( ( dens(
k,i,j)+dens(
k,i,j+1) ) * 0.5_rp )
1189 call dft_g2g(
ka,
ks,
ke,
ia,
is,
ie,
ja,
js,
je,
spnudge_uv_lm,
spnudge_uv_mm,diff2)
1196 diff2(
k,i,j) = diff2(
k,i,j) * ( dens(
k,i,j)+dens(
k,i,j+1) ) * 0.5_rp
1213 damp = - damp_alpha_vely(
k,i,j) &
1215 - ( diff(
k,i-1,j) + diff(
k,i+1,j) + diff(
k,i,j-1) + diff(
k,i,j+1) - diff(
k,i,j)*4.0_rp ) &
1216 * 0.125_rp * bnd_smoother_fact ) &
1218 momy_t(
k,i,j) = momy_tp(
k,i,j) &
1220 + ( dens_damp(
k,i,j) + dens_damp(
k,i,j+1) ) * momy(
k,i,j) / ( dens(
k,i,j) + dens(
k,i,j+1) )
1221 if ( do_put ) damp_t_momy(
k,i,j) = damp_t_momy(
k,i,j) + damp / nstep
1229 call comm_vars8( momy_t(:,:,:), i_comm_momy_t )
1234 do i = max(
is-1,1), min(
ie+2,
ia)
1236 diff(
k,i,j) = rhot(
k,i,j) - damp_pott(
k,i,j) * dens(
k,i,j)
1241 call file_history_query( hist_damp(5), do_put )
1252 damp = - damp_alpha_pott(
k,
is,j) &
1254 - ( diff(
k,
is,j-1) + diff(
k,
is,j+1) - diff(
k,
is,j)*2.0_rp ) &
1255 * 0.25_rp * bnd_smoother_fact )
1256 rhot_t(
k,
is,j) = rhot_tp(
k,
is,j) &
1258 + dens_damp(
k,
is,j) * rhot(
k,
is,j) / dens(
k,
is,j)
1259 if ( do_put ) damp_t_rhot(
k,
is,j) = damp_t_rhot(
k,
is,j) + damp / nstep
1270 diff2(
k,i,j) = diff(
k,i,j) / dens(
k,i,j)
1275 call dft_g2g(
ka,
ks,
ke,
ia,
is,
ie,
ja,
js,
je,
spnudge_pt_lm,
spnudge_pt_mm,diff2)
1282 diff2(
k,i,j) = diff2(
k,i,j) * dens(
k,i,j)
1311 damp = - damp_alpha_pott(
k,i,j) &
1313 - ( diff(
k,i-1,j) + diff(
k,i+1,j) + diff(
k,i,j-1) + diff(
k,i,j+1) - diff(
k,i,j)*4.0_rp ) &
1314 * 0.125_rp * bnd_smoother_fact ) &
1316 rhot_t(
k,i,j) = rhot_tp(
k,i,j) &
1318 + dens_damp(
k,i,j) * rhot(
k,i,j) / dens(
k,i,j)
1319 if ( do_put ) damp_t_rhot(
k,i,j) = damp_t_rhot(
k,i,j) + damp / nstep
1327 call comm_vars8( rhot_t(:,:,:), i_comm_rhot_t )
1329 call comm_wait ( dens_t(:,:,:), i_comm_dens_t, .false. )
1330 call comm_wait ( momz_t(:,:,:), i_comm_momz_t, .false. )
1331 call comm_wait ( momx_t(:,:,:), i_comm_momx_t, .false. )
1332 call comm_wait ( momy_t(:,:,:), i_comm_momy_t, .false. )
1333 call comm_wait ( rhot_t(:,:,:), i_comm_rhot_t, .false. )
1340 if ( nd_coef == 0.0_rp .or. eval_type_numfilter ==
'FILTER' )
then
1342 num_diff(:,:,:,:,:) = 0.0_rp
1345 dens, momz, momx, momy, rhot, &
1346 cdz, cdx, cdy, fdz, fdx, fdy, twod, dts, &
1347 ref_dens, ref_pott, &
1348 nd_coef, nd_laplacian_num, nd_sfc_fact, nd_use_rs )
1355 if ( divdmp_coef > 0.0_rp )
then
1359 gsqrt, j13g, j23g, j33g, mapf, &
1361 rcdz, rcdx, rcdy, rfdz, fdz )
1380 dens_t, momz_t, momx_t, momy_t, rhot_t, &
1381 dpres0, rt2p, corioli, &
1382 num_diff, wdamp_coef, divdmp_coef, ddiv, &
1383 flag_fct_momentum, flag_fct_t, &
1384 flag_fct_along_stream, &
1385 cdz, fdz, fdx, fdy, &
1386 rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, &
1387 phi, gsqrt, j13g, j23g, j33g, mapf, &
1388 ref_dens, ref_rhot, &
1389 bnd_w, bnd_e, bnd_s, bnd_n, twod, &
1397 if ( nd_coef > 0.0_rp .and. eval_type_numfilter ==
'FILTER' )
then
1398 call comm_vars8( dens(:,:,:), i_comm_dens )
1399 call comm_vars8( momz(:,:,:), i_comm_momz )
1400 call comm_vars8( momx(:,:,:), i_comm_momx )
1401 call comm_vars8( momy(:,:,:), i_comm_momy )
1402 call comm_vars8( rhot(:,:,:), i_comm_rhot )
1403 call comm_wait ( dens(:,:,:), i_comm_dens, .false. )
1404 call comm_wait ( momz(:,:,:), i_comm_momz, .false. )
1405 call comm_wait ( momx(:,:,:), i_comm_momx, .false. )
1406 call comm_wait ( momy(:,:,:), i_comm_momy, .false. )
1407 call comm_wait ( rhot(:,:,:), i_comm_rhot, .false. )
1411 dens, momz, momx, momy, rhot, &
1412 cdz, cdx, cdy, fdz, fdx, fdy, &
1413 rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, &
1415 gsqrt, mapf, ref_dens, ref_pott, &
1416 nd_coef, nd_laplacian_num, nd_sfc_fact, nd_use_rs )
1423 dens( 1:
ks-1,i,j) = dens(
ks,i,j)
1424 momz( 1:
ks-1,i,j) = 0.0_rp
1425 momx( 1:
ks-1,i,j) = momx(
ks,i,j)
1426 momy( 1:
ks-1,i,j) = momy(
ks,i,j)
1427 rhot( 1:
ks-1,i,j) = rhot(
ks,i,j)
1429 prog( 1:
ks-1,i,j,iv) = prog(
ks,i,j,iv)
1431 dens(
ke+1:
ka, i,j) = dens(
ke,i,j)
1432 momz(
ke+1:
ka, i,j) = 0.0_rp
1433 momx(
ke+1:
ka, i,j) = momx(
ke,i,j)
1434 momy(
ke+1:
ka, i,j) = momy(
ke,i,j)
1435 rhot(
ke+1:
ka, i,j) = rhot(
ke,i,j)
1437 prog(
ke+1:
ka, i,j,iv) = prog(
ke,i,j,iv)
1442 call comm_vars8( dens(:,:,:), i_comm_dens )
1443 call comm_vars8( momz(:,:,:), i_comm_momz )
1444 call comm_vars8( momx(:,:,:), i_comm_momx )
1445 call comm_vars8( momy(:,:,:), i_comm_momy )
1446 call comm_vars8( rhot(:,:,:), i_comm_rhot )
1448 call comm_vars8( prog(:,:,:,iv), i_comm_prog(iv) )
1450 call comm_wait ( dens(:,:,:), i_comm_dens, .false. )
1451 call comm_wait ( momz(:,:,:), i_comm_momz, .false. )
1452 call comm_wait ( momx(:,:,:), i_comm_momx, .false. )
1453 call comm_wait ( momy(:,:,:), i_comm_momy, .false. )
1454 call comm_wait ( rhot(:,:,:), i_comm_rhot, .false. )
1456 call comm_wait ( prog(:,:,:,iv), i_comm_prog(iv), .false. )
1459 if ( use_average )
then
1464 dens_av(
k,i,j) = dens_av(
k,i,j) + dens(
k,i,j) / nstep
1465 momz_av(
k,i,j) = momz_av(
k,i,j) + momz(
k,i,j) / nstep
1466 momx_av(
k,i,j) = momx_av(
k,i,j) + momx(
k,i,j) / nstep
1467 momy_av(
k,i,j) = momy_av(
k,i,j) + momy(
k,i,j) / nstep
1468 rhot_av(
k,i,j) = rhot_av(
k,i,j) + rhot(
k,i,j) / nstep
1480 mflx_av(
k,i,j,n) = mflx_av(
k,i,j,n) + mflx(
k,i,j,n)
1502 mflx(
k,i,j,n) = mflx_av(
k,i,j,n) / nstep
1510 call comm_vars8( mflx(:,:,:,
zdir), i_comm_mflx_z )
1511 if ( .not. twod )
call comm_vars8( mflx(:,:,:,
xdir), i_comm_mflx_x )
1512 call comm_vars8( mflx(:,:,:,
ydir), i_comm_mflx_y )
1513 call comm_wait ( mflx(:,:,:,
zdir), i_comm_mflx_z, .false. )
1514 if ( .not. twod )
call comm_wait ( mflx(:,:,:,
xdir), i_comm_mflx_x, .false. )
1515 call comm_wait ( mflx(:,:,:,
ydir), i_comm_mflx_y, .false. )
1529 if ( nd_coef_q == 0.0_rp )
then
1531 num_diff_q(:,:,:,:) = 0.0_rp
1534 dens00, qtrc(:,:,:,iq), iq==i_qv, &
1535 cdz, cdx, cdy, twod, dtl, &
1537 nd_coef_q, nd_laplacian_num, nd_sfc_fact, nd_use_rs )
1544 if ( flag_tracer_split_tend )
then
1550 qtrc(
k,i,j,iq) = qtrc(
k,i,j,iq) &
1551 + rhoq_t(
k,i,j,iq) * dtl / dens(
k,i,j)
1557 rhoq_tn => rhoq_t(:,:,:,iq)
1563 qtrc0(:,:,:,iq), rhoq_tn, &
1566 gsqrt, mapf(:,:,:,
i_xy), &
1567 cdz, rcdz, rcdx, rcdy, &
1568 bnd_w, bnd_e, bnd_s, bnd_n, &
1571 llast .AND. flag_fct_tracer, &
1572 flag_fct_along_stream )
1578 call file_history_query( hist_qflx(iv,iq), do_put )
1579 if ( do_put .or. monit_lateral_flag(iv) )
then
1581 call file_history_put( hist_qflx(iv,iq), qflx(:,:,:,iv) )
1586 if ( bnd_w .and. monit_qflx_west > 0 )
then
1590 qflx_west(
k,j) = qflx_west(
k,j) + qflx(
k,
is-1,j,
xdir)
1594 if ( bnd_e .and. monit_qflx_east > 0 )
then
1598 qflx_east(
k,j) = qflx_east(
k,j) + qflx(
k,
ie,j,
xdir)
1602 if ( bnd_s .and. monit_qflx_south > 0 )
then
1606 qflx_south(
k,i) = qflx_south(
k,i) + qflx(
k,i,
js-1,
ydir)
1610 if ( bnd_n .and. monit_qflx_north > 0 )
then
1614 qflx_north(
k,i) = qflx_north(
k,i) + qflx(
k,i,
je,
ydir)
1629 qtrc(
k,i,j,iq) = ( qtrc0(
k,i,j,iq) * dens00(
k,i,j) &
1630 + rhoq_t(
k,i,j,iq) * dtl ) / dens(
k,i,j)
1637 call comm_vars8( qtrc(:,:,:,iq), i_comm_qtrc(iq) )
1639 if ( use_average )
then
1643 qtrc_av(
k,i,j,iq) = qtrc(
k,i,j,iq)
1652 call comm_wait ( qtrc(:,:,:,iq), i_comm_qtrc(iq), .false. )
1658 call file_history_put( hist_phys(1), dens_tp )
1659 call file_history_put( hist_phys(2), momz_tp )
1660 call file_history_put( hist_phys(3), momx_tp )
1661 call file_history_put( hist_phys(4), momy_tp )
1662 call file_history_put( hist_phys(5), rhot_tp )
1664 call file_history_put( hist_damp(1), damp_t_dens )
1665 call file_history_put( hist_damp(2), damp_t_momz )
1666 call file_history_put( hist_damp(3), damp_t_momx )
1667 call file_history_put( hist_damp(4), damp_t_momy )
1668 call file_history_put( hist_damp(5), damp_t_rhot )
1671 call file_history_query( hist_mflx(iv), do_put )
1672 if ( do_put .or. monit_lateral_flag(iv) )
then
1674 call file_history_put( hist_mflx(iv), mflx(:,:,:,iv) )
1679 call file_history_query( hist_tflx(iv), do_put )
1682 call file_history_put( hist_tflx(iv), tflx(:,:,:,iv) )
1687 call monitor_put( monit_damp_mass, damp_t_dens(:,:,:) )
1693 call monitor_put( monit_damp_qtot, dens_tq(:,:,:) )
1710 integer,
intent(in) :: MONIT_ID
1711 logical,
intent(in) :: BND_flag
1712 real(RP),
target,
intent(in) :: flx(:,:)
1713 real(RP),
target,
intent(in) :: flx_zero(:,:)
1715 real(RP),
pointer :: flx_ptr(:,:)
1717 if ( bnd_flag )
then
1722 call monitor_put( monit_id, flx_ptr(:,:) )
1729 integer,
intent(in) :: I_DIR
1730 real(RP),
intent(inout) :: flx(KA,IA,JA,3)
1731 real(RP),
intent(in) :: GSQRT(KA,IA,JA,7)
1732 real(RP),
intent(in) :: MAPF (IA,JA,2,4)
1737 if (i_dir ==
zdir)
then
1742 flx(k,i,j,
zdir) = flx(k,i,j,
zdir) * mapf(i,j,1,
i_xy) * mapf(i,j,2,
i_xy) / gsqrt(k,i,j,
i_xyw)
1748 if (i_dir ==
xdir)
then
1753 flx(k,i,j,
xdir) = flx(k,i,j,
xdir) * mapf(i,j,2,
i_uy) / gsqrt(k,i,j,
i_uyz)
1759 if (i_dir ==
ydir)
then
1764 flx(k,i,j,
ydir) = flx(k,i,j,
ydir) * mapf(i,j,1,
i_xv) / gsqrt(k,i,j,
i_xvz)