SCALE-RM
scale_atmos_dyn_tinteg_large_rk3.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
9 
10 !! This module provides two types of 3rd order and 3 stage Runge=Kutta method in Wichere and Skamarock (2002).
11 !! Note that, it ensures 3rd order accuracy only for the case of linear eqautions, and is generally 2nd order accuracy.
12 !!
13 
14 !-------------------------------------------------------------------------------
15 #include "scalelib.h"
17  !-----------------------------------------------------------------------------
18  !
19  !++ used modules
20  !
21  use scale_precision
22  use scale_io
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_prc, only: &
65  prc_abort
66  implicit none
67 
68  character(len=*) :: tinteg_type
69 
70  integer :: iv
71  !---------------------------------------------------------------------------
72 
73  if ( tinteg_type /= 'RK3' ) then
74  log_error("ATMOS_DYN_Tinteg_large_rk3_setup",*) 'TINTEG_LARGE_TYPE is not RK3. Check!'
75  call prc_abort
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  num_diff, num_diff_q, &
88  DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, &
89  CORIOLI, &
90  CDZ, CDX, CDY, FDZ, FDX, FDY, &
91  RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
92  PHI, GSQRT, &
93  J13G, J23G, J33G, MAPF, &
94  AQ_R, AQ_CV, AQ_CP, AQ_MASS, &
95  REF_dens, REF_pott, REF_qv, REF_pres, &
96  BND_W, BND_E, BND_S, BND_N, TwoD, &
97  ND_COEF, ND_COEF_Q, ND_LAPLACIAN_NUM, &
98  ND_SFC_FACT, ND_USE_RS, &
99  BND_QA, BND_IQ, BND_SMOOTHER_FACT, &
100  DAMP_DENS, DAMP_VELZ, DAMP_VELX, &
101  DAMP_VELY, DAMP_POTT, DAMP_QTRC, &
102  DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, &
103  DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, &
104  MFLUX_OFFSET_X, MFLUX_OFFSET_Y, &
105  wdamp_coef, divdmp_coef, &
106  FLAG_TRACER_SPLIT_TEND, &
107  FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, &
108  FLAG_FCT_ALONG_STREAM, &
109  USE_AVERAGE, &
110  I_QV, &
111  DTL, DTS )
112 
113  use scale_atmos_dyn_tstep_large, only: &
115 
116  implicit none
117 
118  real(rp), intent(inout) :: dens(ka,ia,ja)
119  real(rp), intent(inout) :: momz(ka,ia,ja)
120  real(rp), intent(inout) :: momx(ka,ia,ja)
121  real(rp), intent(inout) :: momy(ka,ia,ja)
122  real(rp), intent(inout) :: rhot(ka,ia,ja)
123  real(rp), intent(inout) :: qtrc(ka,ia,ja,qa)
124 
125  real(rp), intent(inout) :: prog(ka,ia,ja,va)
126 
127  real(rp), intent(inout) :: dens_av(ka,ia,ja)
128  real(rp), intent(inout) :: momz_av(ka,ia,ja)
129  real(rp), intent(inout) :: momx_av(ka,ia,ja)
130  real(rp), intent(inout) :: momy_av(ka,ia,ja)
131  real(rp), intent(inout) :: rhot_av(ka,ia,ja)
132  real(rp), intent(inout) :: qtrc_av(ka,ia,ja,qa)
133 
134  real(rp), intent(out) :: num_diff(ka,ia,ja,5,3)
135  real(rp), intent(out) :: num_diff_q(ka,ia,ja,3)
136 
137  real(rp), intent(in) :: dens_tp(ka,ia,ja)
138  real(rp), intent(in) :: momz_tp(ka,ia,ja)
139  real(rp), intent(in) :: momx_tp(ka,ia,ja)
140  real(rp), intent(in) :: momy_tp(ka,ia,ja)
141  real(rp), intent(in) :: rhot_tp(ka,ia,ja)
142  real(rp), intent(in) :: rhoq_tp(ka,ia,ja,qa)
143 
144  real(rp), intent(in) :: corioli(ia,ja)
145 
146  real(rp), intent(in) :: cdz (ka)
147  real(rp), intent(in) :: cdx (ia)
148  real(rp), intent(in) :: cdy (ja)
149  real(rp), intent(in) :: fdz (ka-1)
150  real(rp), intent(in) :: fdx (ia-1)
151  real(rp), intent(in) :: fdy (ja-1)
152  real(rp), intent(in) :: rcdz(ka)
153  real(rp), intent(in) :: rcdx(ia)
154  real(rp), intent(in) :: rcdy(ja)
155  real(rp), intent(in) :: rfdz(ka-1)
156  real(rp), intent(in) :: rfdx(ia-1)
157  real(rp), intent(in) :: rfdy(ja-1)
158 
159  real(rp), intent(in) :: phi (ka,ia,ja)
160  real(rp), intent(in) :: gsqrt(ka,ia,ja,7)
161  real(rp), intent(in) :: j13g (ka,ia,ja,7)
162  real(rp), intent(in) :: j23g (ka,ia,ja,7)
163  real(rp), intent(in) :: j33g
164  real(rp), intent(in) :: mapf (ia,ja,2,4)
165 
166  real(rp), intent(in) :: aq_r (qa)
167  real(rp), intent(in) :: aq_cv (qa)
168  real(rp), intent(in) :: aq_cp (qa)
169  real(rp), intent(in) :: aq_mass(qa)
170 
171  real(rp), intent(in) :: ref_dens(ka,ia,ja)
172  real(rp), intent(in) :: ref_pott(ka,ia,ja)
173  real(rp), intent(in) :: ref_qv (ka,ia,ja)
174  real(rp), intent(in) :: ref_pres(ka,ia,ja)
175 
176  logical, intent(in) :: bnd_w
177  logical, intent(in) :: bnd_e
178  logical, intent(in) :: bnd_s
179  logical, intent(in) :: bnd_n
180  logical, intent(in) :: twod
181 
182  real(rp), intent(in) :: nd_coef
183  real(rp), intent(in) :: nd_coef_q
184  integer, intent(in) :: nd_laplacian_num
185  real(rp), intent(in) :: nd_sfc_fact
186  logical, intent(in) :: nd_use_rs
187 
188  integer, intent(in) :: bnd_qa
189  integer, intent(in) :: bnd_iq(qa)
190  real(rp), intent(in) :: bnd_smoother_fact
191 
192  real(rp), intent(in) :: damp_dens(ka,ia,ja)
193  real(rp), intent(in) :: damp_velz(ka,ia,ja)
194  real(rp), intent(in) :: damp_velx(ka,ia,ja)
195  real(rp), intent(in) :: damp_vely(ka,ia,ja)
196  real(rp), intent(in) :: damp_pott(ka,ia,ja)
197  real(rp), intent(in) :: damp_qtrc(ka,ia,ja,bnd_qa)
198 
199  real(rp), intent(in) :: damp_alpha_dens(ka,ia,ja)
200  real(rp), intent(in) :: damp_alpha_velz(ka,ia,ja)
201  real(rp), intent(in) :: damp_alpha_velx(ka,ia,ja)
202  real(rp), intent(in) :: damp_alpha_vely(ka,ia,ja)
203  real(rp), intent(in) :: damp_alpha_pott(ka,ia,ja)
204  real(rp), intent(in) :: damp_alpha_qtrc(ka,ia,ja,bnd_qa)
205  real(rp), intent(in) :: mflux_offset_x(ka,ja,2)
206  real(rp), intent(in) :: mflux_offset_y(ka,ia,2)
207 
208  real(rp), intent(in) :: wdamp_coef(ka)
209  real(rp), intent(in) :: divdmp_coef
210 
211  logical, intent(in) :: flag_tracer_split_tend
212  logical, intent(in) :: flag_fct_momentum
213  logical, intent(in) :: flag_fct_t
214  logical, intent(in) :: flag_fct_tracer
215  logical, intent(in) :: flag_fct_along_stream
216 
217  logical, intent(in) :: use_average
218 
219  integer, intent(in) :: i_qv
220 
221  real(dp), intent(in) :: dtl
222  real(dp), intent(in) :: dts
223 
224  ! for time integartion
225  real(rp) :: dens0 (ka,ia,ja)
226  real(rp) :: momz0 (ka,ia,ja)
227  real(rp) :: momx0 (ka,ia,ja)
228  real(rp) :: momy0 (ka,ia,ja)
229  real(rp) :: rhot0 (ka,ia,ja)
230  real(rp) :: qtrc0 (ka,ia,ja,qa)
231  real(rp) :: prog0 (ka,ia,ja,va)
232 
233  real(dp) :: dtrk
234 
235  logical :: last
236 
237  integer :: n, iq
238 
239  !---------------------------------------------------------------------------
240 
241  call prof_rapstart("DYN_Large_RK3_Prep", 3)
242 
243  !##### SAVE #####
244 !OCL XFILL
245  dens0(:,:,:) = dens(:,:,:)
246 !OCL XFILL
247  momz0(:,:,:) = momz(:,:,:)
248 !OCL XFILL
249  momx0(:,:,:) = momx(:,:,:)
250 !OCL XFILL
251  momy0(:,:,:) = momy(:,:,:)
252 !OCL XFILL
253  rhot0(:,:,:) = rhot(:,:,:)
254 !OCL XFILL
255  qtrc0(:,:,:,:) = qtrc(:,:,:,:)
256 
257  if ( va > 0 ) then
258 !OCL XFILL
259  prog0(:,:,:,:) = prog(:,:,:,:)
260  end if
261 
262  call prof_rapend ("DYN_Large_RK3_Prep", 3)
263 
264  !------------------------------------------------------------------------
265  ! Start RK
266  !------------------------------------------------------------------------
267 
268  call prof_rapstart("DYN_Large_RK3", 3)
269 
270  do n = 1, 3
271 
272  dtrk = dtl / real( 4-n, kind=dp )
273  last = n == 3
274 
275  if ( n > 1 ) then
276  dens = dens0
277  momz = momz0
278  momx = momx0
279  momy = momy0
280  rhot = rhot0
281  if ( va > 0 ) prog = prog0
282  end if
283 
284  call atmos_dyn_tstep_large( &
285  dens, momz, momx, momy, rhot, qtrc, prog, & ! (inout)
286  dens_av, momz_av, momx_av, momy_av, rhot_av, qtrc_av, & ! (inout)
287  num_diff, num_diff_q, & ! (out)
288  qtrc0, & ! (in)
289  dens_tp, momz_tp, momx_tp, momy_tp, rhot_tp, rhoq_tp, & ! (in)
290  corioli, & ! (in)
291  cdz, cdx, cdy, fdz, fdx, fdy, & ! (in)
292  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! (in)
293  phi, gsqrt, & ! (in)
294  j13g, j23g, j33g, mapf, & ! (in)
295  aq_r, aq_cv, aq_cp, aq_mass, & ! (in)
296  ref_dens, ref_pott, ref_qv, ref_pres, & ! (in)
297  bnd_w, bnd_e, bnd_s, bnd_n, twod, & ! (in)
298  nd_coef, nd_coef_q, nd_laplacian_num, nd_sfc_fact, nd_use_rs, & ! (in)
299  bnd_qa, bnd_iq, bnd_smoother_fact, & ! (in)
300  damp_dens, damp_velz, damp_velx, & ! (in)
301  damp_vely, damp_pott, damp_qtrc, & ! (in)
302  damp_alpha_dens, damp_alpha_velz, damp_alpha_velx, & ! (in)
303  damp_alpha_vely, damp_alpha_pott, damp_alpha_qtrc, & ! (in)
304  mflux_offset_x, mflux_offset_y, & ! (in)
305  wdamp_coef, divdmp_coef, & ! (in)
306  flag_tracer_split_tend, & ! (in)
307  flag_fct_momentum, flag_fct_t, flag_fct_tracer, & ! (in)
308  flag_fct_along_stream, & ! (in)
309  use_average .AND. last, & ! (in)
310  i_qv, & ! (in)
311  dtrk, dts, last ) ! (in)
312 
313  end do
314 
315  call prof_rapend ("DYN_Large_RK3", 3)
316 
317  return
318  end subroutine atmos_dyn_tinteg_large_rk3
319 
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_tracer::qa
integer, public qa
Definition: scale_tracer.F90:34
scale_index
module Index
Definition: scale_index.F90:11
scale_const::const_undef2
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:38
scale_atmos_dyn_tstep_large::atmos_dyn_tstep_large
procedure(large), pointer, public atmos_dyn_tstep_large
Definition: scale_atmos_dyn_tstep_large.F90:179
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
scale_prof::prof_rapstart
subroutine, public prof_rapstart(rapname_base, level, disable_barrier)
Start raptime.
Definition: scale_prof.F90:159
scale_atmos_dyn_tinteg_large_rk3
module Atmosphere / Dyn Tinteg
Definition: scale_atmos_dyn_tinteg_large_rk3.F90:16
scale_atmos_dyn_tinteg_large_rk3::atmos_dyn_tinteg_large_rk3
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, 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, TwoD, ND_COEF, ND_COEF_Q, ND_LAPLACIAN_NUM, ND_SFC_FACT, ND_USE_RS, BND_QA, BND_IQ, BND_SMOOTHER_FACT, DAMP_DENS, DAMP_VELZ, DAMP_VELX, DAMP_VELY, DAMP_POTT, DAMP_QTRC, DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, MFLUX_OFFSET_X, MFLUX_OFFSET_Y, wdamp_coef, divdmp_coef, FLAG_TRACER_SPLIT_TEND, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM, USE_AVERAGE, I_QV, DTL, DTS)
RK3.
Definition: scale_atmos_dyn_tinteg_large_rk3.F90:112
scale_index::va
integer, public va
Definition: scale_index.F90:35
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_io
module STDIO
Definition: scale_io.F90:10
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_debug::check
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:56
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_dyn_tinteg_large_rk3::atmos_dyn_tinteg_large_rk3_setup
subroutine, public atmos_dyn_tinteg_large_rk3_setup(tinteg_type)
Setup.
Definition: scale_atmos_dyn_tinteg_large_rk3.F90:64
scale_precision::dp
integer, parameter, public dp
Definition: scale_precision.F90:32
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_atmos_dyn_tstep_large
module Atmosphere / Dynamical scheme
Definition: scale_atmos_dyn_tstep_large.F90:12
scale_debug
module DEBUG
Definition: scale_debug.F90:11
scale_prof::prof_rapend
subroutine, public prof_rapend(rapname_base, level, disable_barrier)
Save raptime.
Definition: scale_prof.F90:217
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41