SCALE-RM
scale_atmos_dyn_tinteg_large_rk3.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14 !-------------------------------------------------------------------------------
15 #include "inc_openmp.h"
17  !-----------------------------------------------------------------------------
18  !
19  !++ used modules
20  !
21  use scale_precision
22  use scale_stdio
23  use scale_prof
25  use scale_index
26  use scale_tracer
27 
28 #ifdef DEBUG
29  use scale_debug, only: &
30  check
31  use scale_const, only: &
32  undef => const_undef, &
33  iundef => const_undef2
34 #endif
35  !-----------------------------------------------------------------------------
36  implicit none
37  private
38  !-----------------------------------------------------------------------------
39  !
40  !++ Public procedure
41  !
44 
45  !-----------------------------------------------------------------------------
46  !
47  !++ Public parameters & variables
48  !
49  !-----------------------------------------------------------------------------
50  !
51  !++ Private procedure
52  !
53  !-----------------------------------------------------------------------------
54  !
55  !++ Private parameters & variables
56  !
57  !-----------------------------------------------------------------------------
58 contains
59 
60  !-----------------------------------------------------------------------------
63  tinteg_type )
64  use scale_process, only: &
66  implicit none
67 
68  character(len=*) :: tinteg_type
69 
70  integer :: iv
71  !---------------------------------------------------------------------------
72 
73  if ( tinteg_type /= 'RK3' ) then
74  write(*,*) 'xxx TINTEG_LARGE_TYPE is not RK3. Check!'
75  call prc_mpistop
76  end if
77 
78  return
80 
81  !-----------------------------------------------------------------------------
83  subroutine atmos_dyn_tinteg_large_rk3( &
84  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, &
85  PROG, &
86  DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, &
87  mflx_hi, tflx_hi, &
88  num_diff, num_diff_q, &
89  DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, &
90  CORIOLI, &
91  CDZ, CDX, CDY, FDZ, FDX, FDY, &
92  RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
93  PHI, GSQRT, &
94  J13G, J23G, J33G, MAPF, &
95  AQ_R, AQ_CV, AQ_CP, AQ_MASS, &
96  REF_dens, REF_pott, REF_qv, REF_pres, &
97  BND_W, BND_E, BND_S, BND_N, &
98  ND_COEF, ND_COEF_Q, ND_ORDER, ND_SFC_FACT, ND_USE_RS, &
99  DAMP_DENS, DAMP_VELZ, DAMP_VELX, &
100  DAMP_VELY, DAMP_POTT, DAMP_QTRC, &
101  DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, &
102  DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, &
103  wdamp_coef, &
104  divdmp_coef, &
105  FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, &
106  FLAG_FCT_ALONG_STREAM, &
107  USE_AVERAGE, &
108  DTL, DTS )
109  use scale_const, only: &
110  rdry => const_rdry, &
111  rvap => const_rvap, &
112  cvdry => const_cvdry
113  use scale_comm, only: &
114  comm_vars8, &
115  comm_wait
116  use scale_gridtrans, only: &
117  i_xyz, &
118  i_xyw, &
119  i_uyz, &
120  i_xvz, &
121  i_xy, &
122  i_uy, &
123  i_xv
124  use scale_atmos_dyn_common, only: &
128  use scale_atmos_dyn_fvm_flux_ud1, only: &
132  use scale_atmos_dyn_fvm_flux, only: &
136  use scale_atmos_dyn_tstep_large, only: &
138  use scale_atmos_boundary, only: &
139  bnd_qa, &
140  bnd_smoother_fact => atmos_boundary_smoother_fact
141 #ifdef HIST_TEND
142  use scale_history, only: &
143  hist_in
144 #endif
145  implicit none
146 
147  real(RP), intent(inout) :: dens(ka,ia,ja)
148  real(RP), intent(inout) :: momz(ka,ia,ja)
149  real(RP), intent(inout) :: momx(ka,ia,ja)
150  real(RP), intent(inout) :: momy(ka,ia,ja)
151  real(RP), intent(inout) :: rhot(ka,ia,ja)
152  real(RP), intent(inout) :: qtrc(ka,ia,ja,qa)
153 
154  real(RP), intent(inout) :: prog(ka,ia,ja,va)
155 
156  real(RP), intent(inout) :: dens_av(ka,ia,ja)
157  real(RP), intent(inout) :: momz_av(ka,ia,ja)
158  real(RP), intent(inout) :: momx_av(ka,ia,ja)
159  real(RP), intent(inout) :: momy_av(ka,ia,ja)
160  real(RP), intent(inout) :: rhot_av(ka,ia,ja)
161  real(RP), intent(inout) :: qtrc_av(ka,ia,ja,qa)
162 
163  real(RP), intent(out) :: mflx_hi(ka,ia,ja,3)
164  real(RP), intent(out) :: tflx_hi(ka,ia,ja,3)
165 
166  real(RP), intent(out) :: num_diff(ka,ia,ja,5,3)
167  real(RP), intent(out) :: num_diff_q(ka,ia,ja,3)
168 
169  real(RP), intent(in) :: dens_tp(ka,ia,ja)
170  real(RP), intent(in) :: momz_tp(ka,ia,ja)
171  real(RP), intent(in) :: momx_tp(ka,ia,ja)
172  real(RP), intent(in) :: momy_tp(ka,ia,ja)
173  real(RP), intent(in) :: rhot_tp(ka,ia,ja)
174  real(RP), intent(in) :: rhoq_tp(ka,ia,ja,qa)
175 
176  real(RP), intent(in) :: corioli(ia,ja)
177 
178  real(RP), intent(in) :: cdz (ka)
179  real(RP), intent(in) :: cdx (ia)
180  real(RP), intent(in) :: cdy (ja)
181  real(RP), intent(in) :: fdz (ka-1)
182  real(RP), intent(in) :: fdx (ia-1)
183  real(RP), intent(in) :: fdy (ja-1)
184  real(RP), intent(in) :: rcdz(ka)
185  real(RP), intent(in) :: rcdx(ia)
186  real(RP), intent(in) :: rcdy(ja)
187  real(RP), intent(in) :: rfdz(ka-1)
188  real(RP), intent(in) :: rfdx(ia-1)
189  real(RP), intent(in) :: rfdy(ja-1)
190 
191  real(RP), intent(in) :: phi (ka,ia,ja)
192  real(RP), intent(in) :: gsqrt(ka,ia,ja,7)
193  real(RP), intent(in) :: j13g (ka,ia,ja,7)
194  real(RP), intent(in) :: j23g (ka,ia,ja,7)
195  real(RP), intent(in) :: j33g
196  real(RP), intent(in) :: mapf (ia,ja,2,4)
197 
198  real(RP), intent(in) :: aq_r (qa)
199  real(RP), intent(in) :: aq_cv (qa)
200  real(RP), intent(in) :: aq_cp (qa)
201  real(RP), intent(in) :: aq_mass(qa)
202 
203  real(RP), intent(in) :: ref_dens(ka,ia,ja)
204  real(RP), intent(in) :: ref_pott(ka,ia,ja)
205  real(RP), intent(in) :: ref_qv (ka,ia,ja)
206  real(RP), intent(in) :: ref_pres(ka,ia,ja)
207 
208  logical, intent(in) :: bnd_w
209  logical, intent(in) :: bnd_e
210  logical, intent(in) :: bnd_s
211  logical, intent(in) :: bnd_n
212 
213  real(RP), intent(in) :: nd_coef
214  real(RP), intent(in) :: nd_coef_q
215  integer, intent(in) :: nd_order
216  real(RP), intent(in) :: nd_sfc_fact
217  logical, intent(in) :: nd_use_rs
218 
219  real(RP), intent(in) :: damp_dens(ka,ia,ja)
220  real(RP), intent(in) :: damp_velz(ka,ia,ja)
221  real(RP), intent(in) :: damp_velx(ka,ia,ja)
222  real(RP), intent(in) :: damp_vely(ka,ia,ja)
223  real(RP), intent(in) :: damp_pott(ka,ia,ja)
224  real(RP), intent(in) :: damp_qtrc(ka,ia,ja,bnd_qa)
225 
226  real(RP), intent(in) :: damp_alpha_dens(ka,ia,ja)
227  real(RP), intent(in) :: damp_alpha_velz(ka,ia,ja)
228  real(RP), intent(in) :: damp_alpha_velx(ka,ia,ja)
229  real(RP), intent(in) :: damp_alpha_vely(ka,ia,ja)
230  real(RP), intent(in) :: damp_alpha_pott(ka,ia,ja)
231  real(RP), intent(in) :: damp_alpha_qtrc(ka,ia,ja,bnd_qa)
232 
233  real(RP), intent(in) :: wdamp_coef(ka)
234  real(RP), intent(in) :: divdmp_coef
235 
236  logical, intent(in) :: flag_fct_momentum
237  logical, intent(in) :: flag_fct_t
238  logical, intent(in) :: flag_fct_tracer
239  logical, intent(in) :: flag_fct_along_stream
240 
241  logical, intent(in) :: use_average
242 
243  real(DP), intent(in) :: dtl
244  real(DP), intent(in) :: dts
245 
246  ! for time integartion
247  real(RP) :: dens0 (ka,ia,ja)
248  real(RP) :: momz0 (ka,ia,ja)
249  real(RP) :: momx0 (ka,ia,ja)
250  real(RP) :: momy0 (ka,ia,ja)
251  real(RP) :: rhot0 (ka,ia,ja)
252  real(RP) :: qtrc0 (ka,ia,ja,qa)
253  real(RP) :: prog0 (ka,ia,ja,va)
254 
255  real(DP) :: dtrk
256 
257  logical :: last
258 
259  integer :: n, iq
260 
261  !---------------------------------------------------------------------------
262 
263  call prof_rapstart("DYN_Large_RK3_Prep", 3)
264 
265  !##### SAVE #####
266 !OCL XFILL
267  dens0(:,:,:) = dens(:,:,:)
268 !OCL XFILL
269  momz0(:,:,:) = momz(:,:,:)
270 !OCL XFILL
271  momx0(:,:,:) = momx(:,:,:)
272 !OCL XFILL
273  momy0(:,:,:) = momy(:,:,:)
274 !OCL XFILL
275  rhot0(:,:,:) = rhot(:,:,:)
276 !OCL XFILL
277  qtrc0(:,:,:,:) = qtrc(:,:,:,:)
278 
279  if ( va > 0 ) then
280 !OCL XFILL
281  prog0(:,:,:,:) = prog(:,:,:,:)
282  end if
283 
284  call prof_rapend ("DYN_Large_RK3_Prep", 3)
285 
286  !------------------------------------------------------------------------
287  ! Start RK
288  !------------------------------------------------------------------------
289 
290  call prof_rapstart("DYN_Large_RK3", 3)
291 
292  do n = 1, 3
293 
294  dtrk = dtl / real( 4-n, kind=dp )
295  last = n == 3
296 
297  if ( n > 1 ) then
298  dens = dens0
299  momz = momz0
300  momx = momx0
301  momy = momy0
302  rhot = rhot0
303  if ( va > 0 ) prog = prog0
304  end if
305 
306  call atmos_dyn_tstep_large( &
307  dens, momz, momx, momy, rhot, qtrc, prog, & ! (inout)
308  dens_av, momz_av, momx_av, momy_av, rhot_av, qtrc_av, & ! (inout)
309  mflx_hi, tflx_hi, & ! (out)
310  num_diff, num_diff_q, & ! (out)
311  qtrc0, & ! (in)
312  dens_tp, momz_tp, momx_tp, momy_tp, rhot_tp, rhoq_tp, & ! (in)
313  corioli, & ! (in)
314  cdz, cdx, cdy, fdz, fdx, fdy, & ! (in)
315  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! (in)
316  phi, gsqrt, & ! (in)
317  j13g, j23g, j33g, mapf, & ! (in)
318  aq_r, aq_cv, aq_cp, aq_mass, & ! (in)
319  ref_dens, ref_pott, ref_qv, ref_pres, & ! (in)
320  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
321  nd_coef, nd_coef_q, nd_order, nd_sfc_fact, nd_use_rs, & ! (in)
322  damp_dens, damp_velz, damp_velx, & ! (in)
323  damp_vely, damp_pott, damp_qtrc, & ! (in)
324  damp_alpha_dens, damp_alpha_velz, damp_alpha_velx, & ! (in)
325  damp_alpha_vely, damp_alpha_pott, damp_alpha_qtrc, & ! (in)
326  wdamp_coef, & ! (in)
327  divdmp_coef, & ! (in)
328  flag_fct_momentum, flag_fct_t, flag_fct_tracer, & ! (in)
329  flag_fct_along_stream, & ! (in)
330  use_average .AND. last, & ! (in)
331  dtrk, dts, last ) ! (in)
332 
333  end do
334 
335  call prof_rapend ("DYN_Large_RK3", 3)
336 
337  return
338  end subroutine atmos_dyn_tinteg_large_rk3
339 
module DEBUG
Definition: scale_debug.F90:13
integer, public i_xvz
real(rp), public const_cvdry
specific heat (dry air,constant volume) [J/kg/K]
Definition: scale_const.F90:59
subroutine, public atmos_dyn_fvm_fluxy_xyz_ud1(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XYZ
subroutine, public prc_mpistop
Abort MPI.
subroutine, public atmos_dyn_numfilter_coef_q(num_diff_q, DENS, QTRC, CDZ, CDX, CDY, dt, REF_qv, iq, ND_COEF, ND_ORDER, ND_SFC_FACT, ND_USE_RS)
Calc coefficient of numerical filter.
integer, public va
Definition: scale_index.F90:38
subroutine, public atmos_dyn_fvm_fluxz_xyz_ud1(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at XYZ
subroutine, public atmos_dyn_numfilter_coef(num_diff, DENS, MOMZ, MOMX, MOMY, RHOT, CDZ, CDX, CDY, FDZ, FDX, FDY, DT, REF_dens, REF_pott, ND_COEF, ND_ORDER, ND_SFC_FACT, ND_USE_RS)
Calc coefficient of numerical filter.
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxx_xyz
module STDIO
Definition: scale_stdio.F90:12
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
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
real(rp), public const_undef
Definition: scale_const.F90:43
module grid index
module TRACER
module Index
Definition: scale_index.F90:14
integer, public ia
of whole cells: x, local, with HALO
module GRIDTRANS
integer, public ka
of whole cells: z, local, with HALO
integer, public i_uy
integer, public i_xyw
subroutine, public atmos_dyn_tinteg_large_rk3_setup(tinteg_type)
Setup.
module COMMUNICATION
Definition: scale_comm.F90:23
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:40
module Atmosphere / Dynamics common
module PROCESS
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:65
module CONSTANT
Definition: scale_const.F90:14
subroutine, public atmos_dyn_tinteg_large_rk3(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, mflx_hi, tflx_hi, num_diff, num_diff_q, DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, CORIOLI, 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, BND_W, BND_E, BND_S, BND_N, 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, wdamp_coef, divdmp_coef, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM, USE_AVERAGE, DTL, DTS)
RK3.
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:156
integer, public i_xv
module profiler
Definition: scale_prof.F90:10
module scale_atmos_dyn_fvm_flux
integer, public i_uyz
module PRECISION
subroutine, public atmos_dyn_fct(qflx_anti, phi_in, DENS0, DENS, qflx_hi, qflx_lo, mflx_hi, rdz, rdx, rdy, GSQRT, MAPF, dt, flag_vect)
Flux Correction Transport Limiter.
module ATMOSPHERE / Boundary treatment
module HISTORY
integer, public i_xyz
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:204
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxz_xyz
integer, parameter, public dp
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxy_xyz
procedure(large), pointer, public atmos_dyn_tstep_large
subroutine, public atmos_dyn_fvm_fluxx_xyz_ud1(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XYZ
module scale_atmos_dyn_fvm_flux_ud1
real(rp), public atmos_boundary_smoother_fact
integer, public ja
of whole cells: y, local, with HALO