SCALE-RM
scale_atmos_dyn_tinteg_short_rk11s8o.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
13 !-------------------------------------------------------------------------------
14 #include "scalelib.h"
16  !-----------------------------------------------------------------------------
17  !
18  !++ used modules
19  !
20  use scale_precision
21  use scale_io
22  use scale_prof
24  use scale_index
25  use scale_tracer
26 #if defined DEBUG || defined QUICKDEBUG
27  use scale_debug, only: &
28  check
29  use scale_const, only: &
30  undef => const_undef, &
31  iundef => const_undef2
32 #endif
33 
35  rkinfo, &
36  rkcommon_setup => atmos_dyn_tinteg_rkcommon_setup, &
37  rkcommon_rkwork_alloc => atmos_dyn_tinteg_rkcommon_rkwork_alloc, &
38  rkcommon_rkwork_dealloc => atmos_dyn_tinteg_rkcommon_rkwork_dealloc, &
39  rkcommon_comm => atmos_dyn_tinteg_rkcommon_comm, &
40  rkcommon_comm_wait => atmos_dyn_tinteg_rkcommon_comm_wait, &
41  rkcommon_nextstage => atmos_dyn_tinteg_rkcommon_nextstage, &
42  rkcommon_updatevar => atmos_dyn_tinteg_rkcommon_updatevar, &
43  rkcommon_updateflux => atmos_dyn_tinteg_rkcommon_updateflux
44 
45 
46  !-----------------------------------------------------------------------------
47  implicit none
48  private
49  !-----------------------------------------------------------------------------
50  !
51  !++ Public procedure
52  !
55 
56  !-----------------------------------------------------------------------------
57  !
58  !++ Public parameters & variables
59  !
60  !-----------------------------------------------------------------------------
61 
62  !- coeffecients for 11 stage RK with 8th order accuracy by Cooper and Verner (1972) --------------
63  real(RP), parameter :: sq = 21.0_rp**0.5
64  real(RP), parameter :: RC41 = 1.0_rp/7.0_rp
65  real(RP), parameter :: RC42 = (-7.0_rp+3.0_rp*sq)/98.0_rp
66  real(RP), parameter :: RC43 = (21.0_rp-5.0_rp*sq)/49.0_rp
67  real(RP), parameter :: RC51 = (11.0_rp-sq)/84.0_rp
68  real(RP), parameter :: RC53 = (18.0_rp-4.0_rp*sq)/63.0_rp
69  real(RP), parameter :: RC54 = (21.0_rp+sq)/252.0_rp
70  real(RP), parameter :: RC61 = (5.0_rp-sq)/48.0_rp
71  real(RP), parameter :: RC63 = (9.0_rp-sq)/36.0_rp
72  real(RP), parameter :: RC64 = (-231.0_rp-14.0_rp*sq)/360.0_rp
73  real(RP), parameter :: RC65 = (63.0_rp+7.0_rp*sq)/80.0_rp
74  real(RP), parameter :: RC71 = (10.0_rp+sq)/42.0_rp
75  real(RP), parameter :: RC73 = (-432.0_rp-92.0_rp*sq)/315.0_rp
76  real(RP), parameter :: RC74 = (633.0_rp+145.0_rp*sq)/90.0_rp
77  real(RP), parameter :: RC75 = (-504.0_rp-115.0_rp*sq)/70.0_rp
78  real(RP), parameter :: RC76 = (63.0_rp+13.0_rp*sq)/35.0_rp
79  real(RP), parameter :: RC81 = 1.0_rp/14.0_rp
80  real(RP), parameter :: RC85 = (14.0_rp+3.0_rp*sq)/126.0_rp
81  real(RP), parameter :: RC86 = (13.0_rp+3.0_rp*sq)/63.0_rp
82  real(RP), parameter :: RC87 = 1.0_rp/9.0_rp
83  real(RP), parameter :: RC91 = 1.0_rp/32.0_rp
84  real(RP), parameter :: RC95 = (91.0_rp+21.0_rp*sq)/576.0_rp
85  real(RP), parameter :: RC96 = 11.0_rp/72.0_rp
86  real(RP), parameter :: RC97 = (-385.0_rp+75.0_rp*sq)/1152.0_rp
87  real(RP), parameter :: RC98 = (63.0_rp-13.0_rp*sq)/128.0_rp
88  real(RP), parameter :: RC101 = 1.0_rp/14.0_rp
89  real(RP), parameter :: RC105 = 1.0_rp/9.0_rp
90  real(RP), parameter :: RC106 = (-733.0_rp+147.0_rp*sq)/2205.0_rp
91  real(RP), parameter :: RC107 = (515.0_rp-111.0_rp*sq)/504.0_rp
92  real(RP), parameter :: RC108 = (-51.0_rp+11.0_rp*sq)/56.0_rp
93  real(RP), parameter :: RC109 = (132.0_rp-28.0_rp*sq)/245.0_rp
94  real(RP), parameter :: RC115 = (-42.0_rp-7.0_rp*sq)/18.0_rp
95  real(RP), parameter :: RC116 = (-18.0_rp-28.0_rp*sq)/45.0_rp
96  real(RP), parameter :: RC117 = (-273.0_rp+53.0_rp*sq)/72.0_rp
97  real(RP), parameter :: RC118 = (301.0_rp-53.0_rp*sq)/72.0_rp
98  real(RP), parameter :: RC119 = (28.0_rp+28.0_rp*sq)/45.0_rp
99  real(RP), parameter :: RC1110 = (49.0_rp+7.0_rp*sq)/18.0_rp
100  real(RP), parameter :: ZERO = 0.0_rp
101 
102  real(rp), parameter, public :: rkcoef_a_11s8o_cooperverner1972(11,11) = reshape( &
103  (/ zero, 0.5_rp, 0.25_rp, rc41, rc51, rc61, rc71, rc81, rc91, rc101, zero, &
104  zero, zero, 0.25_rp, rc42, zero, zero, zero, zero, zero, zero, zero, &
105  zero, zero, zero, rc43, rc53, rc63, rc73, zero, zero, zero, zero, &
106  zero, zero, zero, zero, rc54, rc64, rc74, zero, zero, zero, zero, &
107  zero, zero, zero, zero, zero, rc65, rc75, rc85, rc95, rc105, rc115, &
108  zero, zero, zero, zero, zero, zero, rc76, rc86, rc96, rc106, rc116, &
109  zero, zero, zero, zero, zero, zero, zero, rc87, rc97, rc107, rc117, &
110  zero, zero, zero, zero, zero, zero, zero, zero, rc98, rc108, rc118, &
111  zero, zero, zero, zero, zero, zero, zero, zero, zero, rc109, rc119, &
112  zero, zero, zero, zero, zero, zero, zero, zero, zero, zero, rc1110, &
113  zero, zero, zero, zero, zero, zero, zero, zero, zero, zero, zero /), &
115 
116  real(rp), parameter, public :: rkcoef_b_11s8o_cooperverner1972(11) = &
117  (/ 0.05_rp, zero, zero, zero, zero, zero, zero, 49.0_rp/180.0_rp, 16.0_rp/45.0_rp, &
118  49.0_rp/180.0_rp, 0.05_rp /)
119 
120  !-----------------------------------------------------------------------------
121  !
122  !++ Private procedure
123  !
124  !-----------------------------------------------------------------------------
125  !
126  !++ Private parameters & variables
127  !
128  integer, private, parameter :: rk_nstage = 11
129  integer, private, parameter :: rk_nregister = 11
130 
131  type(rkinfo), private :: rk_dynvar
132  integer, private, parameter :: i_rk_dens = 1
133  integer, private, parameter :: i_rk_momz = 2
134  integer, private, parameter :: i_rk_momx = 3
135  integer, private, parameter :: i_rk_momy = 4
136  integer, private, parameter :: i_rk_rhot = 5
137 
138  type(rkinfo), private :: rk_prgvar
139 
140  type(rkinfo), private :: rk_mflx_hi
141  type(rkinfo), private :: rk_tflx_hi
142 
143  !-----------------------------------------------------------------------------
144 contains
145  !-----------------------------------------------------------------------------
148  tinteg_type )
149  use scale_prc, only: &
150  prc_abort
151  use scale_const, only: &
152  undef => const_undef
153  use scale_comm_cartesc, only: &
155  implicit none
156 
157  character(len=*) :: tinteg_type
158 
159  integer :: iv, d
160  character(H_SHORT) :: dynvar_name_list(5)
161  character(H_SHORT) :: prgvar_name_list(va)
162  character(H_SHORT) :: flux_name_list(3)
163 
164  real(rp) :: rkcoef_a(rk_nstage,rk_nstage)
165  real(rp) :: rkcoef_b(rk_nstage)
166 
167  !---------------------------------------------------------------------------
168 
169  select case( trim(tinteg_type) )
170  case ('RK11s8o', 'RK11s8oCooperVerner1972')
171  rkcoef_a(:,:) = rkcoef_a_11s8o_cooperverner1972
172  rkcoef_b(:) = rkcoef_b_11s8o_cooperverner1972
173  case default
174  log_error("ATMOS_DYN_Tinteg_short_rk11s8o_setup",*) 'The specified TINTEG_TYPE is invalid. Check!', tinteg_type
175  call prc_abort
176  end select
177 
178  dynvar_name_list(1) = 'DENS'
179  dynvar_name_list(2) = 'MOMZ'
180  dynvar_name_list(3) = 'MOMX'
181  dynvar_name_list(4) = 'MOMY'
182  dynvar_name_list(5) = 'RHOT'
183  call rkcommon_setup( rk_dynvar, rk_nstage, rk_nregister, rkcoef_a, rkcoef_b, dynvar_name_list )
184 
185  do iv = 1, va
186  flux_name_list(iv) = 'PROG'
187  end do
188  call rkcommon_setup( rk_prgvar, rk_nstage, rk_nregister, rkcoef_a, rkcoef_b,prgvar_name_list, comm_id_offset=5 )
189 
190  do d = 1, 3
191  flux_name_list(iv) = 'mflx_hi'
192  end do
193  call rkcommon_setup( rk_mflx_hi, rk_nstage, rk_nregister, rkcoef_a, rkcoef_b, flux_name_list, is_type_flux=.true. )
194 
195 
196  do d = 1, 3
197  flux_name_list(iv) = 'tflx_hi'
198  end do
199  call rkcommon_setup( rk_tflx_hi, rk_nstage, rk_nregister, rkcoef_a, rkcoef_b, flux_name_list, is_type_flux=.true. )
200 
201  !----------------------------------------------------
202 
203  return
205 
206  !-----------------------------------------------------------------------------
208  subroutine atmos_dyn_tinteg_short_rk11s8o( &
209  DENS, MOMZ, MOMX, MOMY, RHOT, PROG, &
210  mflx_hi, tflx_hi, &
211  DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, &
212  DPRES0, CVtot, CORIOLI, &
213  num_diff, wdamp_coef, divdmp_coef, DDIV, &
214  FLAG_FCT_MOMENTUM, FLAG_FCT_T, &
215  FLAG_FCT_ALONG_STREAM, &
216  CDZ, FDZ, FDX, FDY, &
217  RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
218  PHI, GSQRT, J13G, J23G, J33G, MAPF, &
219  REF_pres, REF_dens, &
220  BND_W, BND_E, BND_S, BND_N, TwoD, &
221  dt )
222  use scale_comm_cartesc, only: &
223  comm_vars8, &
224  comm_wait
225  use scale_atmos_dyn_tstep_short, only: &
226  atmos_dyn_tstep => atmos_dyn_tstep_short
227  use scale_atmos_dyn_common, only: &
229  implicit none
230 
231  real(rp), intent(inout) :: dens(ka,ia,ja)
232  real(rp), intent(inout) :: momz(ka,ia,ja)
233  real(rp), intent(inout) :: momx(ka,ia,ja)
234  real(rp), intent(inout) :: momy(ka,ia,ja)
235  real(rp), intent(inout) :: rhot(ka,ia,ja)
236  real(rp), intent(inout) :: prog(ka,ia,ja,va)
237 
238  real(rp), intent(inout) :: mflx_hi(ka,ia,ja,3)
239  real(rp), intent(out) :: tflx_hi(ka,ia,ja,3)
240 
241  real(rp), intent(in) :: dens_t(ka,ia,ja)
242  real(rp), intent(in) :: momz_t(ka,ia,ja)
243  real(rp), intent(in) :: momx_t(ka,ia,ja)
244  real(rp), intent(in) :: momy_t(ka,ia,ja)
245  real(rp), intent(in) :: rhot_t(ka,ia,ja)
246 
247  real(rp), intent(in) :: dpres0(ka,ia,ja)
248  real(rp), intent(in) :: cvtot(ka,ia,ja)
249  real(rp), intent(in) :: corioli(ia,ja)
250  real(rp), intent(in) :: num_diff(ka,ia,ja,5,3)
251  real(rp), intent(in) :: wdamp_coef(ka)
252  real(rp), intent(in) :: divdmp_coef
253  real(rp), intent(in) :: ddiv(ka,ia,ja)
254 
255  logical, intent(in) :: flag_fct_momentum
256  logical, intent(in) :: flag_fct_t
257  logical, intent(in) :: flag_fct_along_stream
258 
259  real(rp), intent(in) :: cdz (ka)
260  real(rp), intent(in) :: fdz (ka-1)
261  real(rp), intent(in) :: fdx (ia-1)
262  real(rp), intent(in) :: fdy (ja-1)
263  real(rp), intent(in) :: rcdz(ka)
264  real(rp), intent(in) :: rcdx(ia)
265  real(rp), intent(in) :: rcdy(ja)
266  real(rp), intent(in) :: rfdz(ka-1)
267  real(rp), intent(in) :: rfdx(ia-1)
268  real(rp), intent(in) :: rfdy(ja-1)
269 
270  real(rp), intent(in) :: phi (ka,ia,ja)
271  real(rp), intent(in) :: gsqrt(ka,ia,ja,7)
272  real(rp), intent(in) :: j13g (ka,ia,ja,7)
273  real(rp), intent(in) :: j23g (ka,ia,ja,7)
274  real(rp), intent(in) :: j33g
275  real(rp), intent(in) :: mapf (ia,ja,2,4)
276 
277  real(rp), intent(in) :: ref_pres(ka,ia,ja)
278  real(rp), intent(in) :: ref_dens(ka,ia,ja)
279 
280  logical, intent(in) :: bnd_w
281  logical, intent(in) :: bnd_e
282  logical, intent(in) :: bnd_s
283  logical, intent(in) :: bnd_n
284  logical, intent(in) :: twod
285 
286  real(rp), intent(in) :: dt
287 
288  integer :: i, j, k, iv, n, s
289  integer :: stage
290 
291  integer, parameter :: ko_dynvar(5) = (/ 0, 1, 0, 0, 0 /)
292  integer, parameter :: io_dynvar(5) = (/ 0, 0, 1, 0, 0 /)
293  integer, parameter :: jo_dynvar(5) = (/ 0, 0, 0, 1, 0 /)
294  integer :: ko_prgvar(va)
295  integer :: io_prgvar(va)
296  integer :: jo_prgvar(va)
297 
298  !---------------------------------------------------------------------------
299 
300  call prof_rapstart("DYN_RK11s8o_Prep",3)
301 
302  call rkcommon_rkwork_alloc( rk_dynvar )
303  call rkcommon_rkwork_alloc( rk_prgvar )
304 
305 #ifdef QUICKDEBUG
306  rk_mflx_hi%buf( 1:ks-1,:,:,:) = undef
307  rk_mflx_hi%buf(ke+1:ka ,:,:,:) = undef
308 #endif
309 
310  !$omp parallel
311 
312  !$omp workshare
313 !OCL XFILL
314  rk_dynvar%work0(:,:,:,i_rk_dens) = dens(:,:,:)
315 !OCL XFILL
316  rk_dynvar%work0(:,:,:,i_rk_momz) = momz(:,:,:)
317 !OCL XFILL
318  rk_dynvar%work0(:,:,:,i_rk_momx) = momx(:,:,:)
319 !OCL XFILL
320  rk_dynvar%work0(:,:,:,i_rk_momy) = momy(:,:,:)
321 !OCL XFILL
322  rk_dynvar%work0(:,:,:,i_rk_rhot) = rhot(:,:,:)
323 !OCL XFILL
324  rk_dynvar%buf(:,:,:,:) = rk_dynvar%work0(:,:,:,:)
325  !$omp end workshare
326 
327  if (va > 0) then
328  !$omp workshare
329 !OCL XFILL
330  rk_prgvar%work0(:,:,:,:) = prog
331 !OCL XFILL
332  rk_prgvar%buf(:,:,:,:) = rk_prgvar%work0(:,:,:,:)
333  !$omp end workshare
334  end if
335  !$omp end parallel
336 
337  if ( bnd_w ) then
338  do j = js, je
339  do k = ks, ke
340  rk_mflx_hi%buf(k,is-1,j,2) = mflx_hi(k,is-1,j,2)
341  end do
342  end do
343  end if
344  if ( bnd_e ) then
345  do j = js, je
346  do k = ks, ke
347  rk_mflx_hi%buf(k,ie,j,2) = mflx_hi(k,ie,j,2)
348  end do
349  end do
350  end if
351  if ( bnd_s ) then
352  do i = is, ie
353  do k = ks, ke
354  rk_mflx_hi%buf(k,i,js-1,3) = mflx_hi(k,i,js-1,3)
355  end do
356  end do
357  end if
358  if ( bnd_n ) then
359  do i = is, ie
360  do k = ks, ke
361  rk_mflx_hi%buf(k,i,je,3) = mflx_hi(k,i,je,3)
362  end do
363  end do
364  end if
365 
366  io_prgvar(:) = 0; jo_prgvar(:) = 0; ko_prgvar(:) = 0
367 
368  call prof_rapend ("DYN_RK11s8o_Prep",3)
369 
370  !------------------------------------------------------------------------
371  ! Start RK
372  !------------------------------------------------------------------------
373 
374  do stage = 1, rk_nstage
375 
376  if ( stage > 1) then
377  call prof_rapstart("DYN_RK11s8o_BND",3)
378  call atmos_dyn_copy_boundary( rk_dynvar%buf(:,:,:,i_rk_dens), & ! [INOUT]
379  rk_dynvar%buf(:,:,:,i_rk_momz), & ! [INOUT]
380  rk_dynvar%buf(:,:,:,i_rk_momx), & ! [INOUT]
381  rk_dynvar%buf(:,:,:,i_rk_momy), & ! [INOUT]
382  rk_dynvar%buf(:,:,:,i_rk_rhot), & ! [INOUT]
383  rk_prgvar%buf(:,:,:,: ), & ! [INOUT]
384  rk_dynvar%work0(:,:,:,i_rk_dens), & ! [IN]
385  rk_dynvar%work0(:,:,:,i_rk_momz), & ! [IN]
386  rk_dynvar%work0(:,:,:,i_rk_momx), & ! [IN]
387  rk_dynvar%work0(:,:,:,i_rk_momy), & ! [IN]
388  rk_dynvar%work0(:,:,:,i_rk_rhot), & ! [IN]
389  rk_prgvar%work0(:,:,:,:), & ! [IN]
390  bnd_w, bnd_e, bnd_s, bnd_n, twod ) ! [IN]
391  call prof_rapend ("DYN_RK11s8o_BND",3)
392 
393  call rkcommon_comm( rk_dynvar )
394  call rkcommon_comm( rk_prgvar )
395  call rkcommon_comm_wait( rk_dynvar )
396  call rkcommon_comm_wait( rk_prgvar )
397  end if
398  !--
399 
400  call prof_rapstart("DYN_RK11s8o",3)
401 
402  call atmos_dyn_tstep( &
403  rk_dynvar%work(:,:,:,i_rk_dens,stage), & ! [OUT]
404  rk_dynvar%work(:,:,:,i_rk_momz,stage), & ! [OUT]
405  rk_dynvar%work(:,:,:,i_rk_momx,stage), & ! [OUT]
406  rk_dynvar%work(:,:,:,i_rk_momy,stage), & ! [OUT]
407  rk_dynvar%work(:,:,:,i_rk_rhot,stage), & ! [OUT]
408  rk_prgvar%work(:,:,:,: ,stage), & ! [OUT]
409  rk_mflx_hi%buf(:,:,:,:), rk_tflx_hi%buf(:,:,:,:), & ! [INOUT,OUT]
410  rk_dynvar%work0(:,:,:,i_rk_dens), & ! [IN]
411  rk_dynvar%work0(:,:,:,i_rk_momz), & ! [IN]
412  rk_dynvar%work0(:,:,:,i_rk_momx), & ! [IN]
413  rk_dynvar%work0(:,:,:,i_rk_momy), & ! [IN]
414  rk_dynvar%work0(:,:,:,i_rk_rhot), & ! [IN]
415  rk_dynvar%buf(:,:,:,i_rk_dens), & ! [IN]
416  rk_dynvar%buf(:,:,:,i_rk_momz), & ! [IN]
417  rk_dynvar%buf(:,:,:,i_rk_momx), & ! [IN]
418  rk_dynvar%buf(:,:,:,i_rk_momy), & ! [IN]
419  rk_dynvar%buf(:,:,:,i_rk_rhot), & ! [IN]
420  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
421  rk_prgvar%work0, rk_prgvar%buf(:,:,:,:), & ! [IN]
422  dpres0, cvtot, corioli, & ! [IN]
423  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
424  flag_fct_momentum, flag_fct_t, & ! [IN]
425  flag_fct_along_stream, & ! [IN]
426  cdz, fdz, fdx, fdy, & ! [IN]
427  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
428  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
429  ref_pres, ref_dens, & ! [IN]
430  bnd_w, bnd_e, bnd_s, bnd_n, twod, & ! [IN]
431  1.0_rp, .false. ) ! [IN]
432 
433  if ( stage < rk_nstage) then
434  call rkcommon_nextstage( rk_dynvar, stage, io_dynvar, jo_dynvar, ko_dynvar, dt )
435  call rkcommon_nextstage( rk_prgvar, stage, io_prgvar, jo_prgvar, ko_prgvar, dt )
436  else
437  call rkcommon_updatevar( rk_dynvar, io_dynvar, jo_dynvar, ko_dynvar, i_rk_dens, i_rk_dens, dt, dens )
438  call rkcommon_updatevar( rk_dynvar, io_dynvar, jo_dynvar, ko_dynvar, i_rk_momz, i_rk_momz, dt, momz )
439  call rkcommon_updatevar( rk_dynvar, io_dynvar, jo_dynvar, ko_dynvar, i_rk_momx, i_rk_momx, dt, momx )
440  call rkcommon_updatevar( rk_dynvar, io_dynvar, jo_dynvar, ko_dynvar, i_rk_momy, i_rk_momy, dt, momy )
441  call rkcommon_updatevar( rk_dynvar, io_dynvar, jo_dynvar, ko_dynvar, i_rk_rhot, i_rk_rhot, dt, rhot )
442  call rkcommon_updatevar( rk_prgvar, io_prgvar, jo_prgvar, ko_prgvar, 1, va, dt, prog )
443  end if
444  call rkcommon_updateflux( rk_mflx_hi, stage, 0, 0, 0, 3, mflx_hi )
445  call rkcommon_updateflux( rk_tflx_hi, stage, 0, 0, 0, 3, tflx_hi )
446 
447  call prof_rapend ("DYN_RK11s8o",3)
448  end do
449 
450  call rkcommon_rkwork_dealloc( rk_dynvar )
451  call rkcommon_rkwork_dealloc( rk_prgvar )
452 
453  return
454  end subroutine atmos_dyn_tinteg_short_rk11s8o
455 
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
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_common
module Atmosphere / Dynamics common
Definition: scale_atmos_dyn_common.F90:12
scale_atmos_dyn_tinteg_short_rk11s8o::atmos_dyn_tinteg_short_rk11s8o_setup
subroutine, public atmos_dyn_tinteg_short_rk11s8o_setup(tinteg_type)
Setup.
Definition: scale_atmos_dyn_tinteg_short_rk11s8o.F90:149
scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_setup
subroutine, public atmos_dyn_tinteg_rkcommon_setup(this, rk_stage_num, rk_register_num, rkcoef_a, rkcoef_b, varname_list, is_type_flux, alloc_rkwork_flag, comm_id_offset)
Definition: scale_atmos_dyn_tinteg_rkcommon.F90:81
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_atmos_dyn_tinteg_rkcommon::rkinfo
Definition: scale_atmos_dyn_tinteg_rkcommon.F90:33
scale_prof::prof_rapstart
subroutine, public prof_rapstart(rapname_base, level, disable_barrier)
Start raptime.
Definition: scale_prof.F90:159
scale_atmos_dyn_tinteg_short_rk11s8o::rkcoef_a_11s8o_cooperverner1972
real(rp), dimension(11, 11), parameter, public rkcoef_a_11s8o_cooperverner1972
Definition: scale_atmos_dyn_tinteg_short_rk11s8o.F90:102
scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_rkwork_dealloc
subroutine, public atmos_dyn_tinteg_rkcommon_rkwork_dealloc(this)
Definition: scale_atmos_dyn_tinteg_rkcommon.F90:171
scale_atmos_dyn_tstep_short::atmos_dyn_tstep_short
procedure(short), pointer, public atmos_dyn_tstep_short
Definition: scale_atmos_dyn_tstep_short.F90:135
scale_index::va
integer, public va
Definition: scale_index.F90:35
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
scale_io
module STDIO
Definition: scale_io.F90:10
scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_updatevar
subroutine, public atmos_dyn_tinteg_rkcommon_updatevar(this, io, jo, ko, vs, ve, dt, var)
Definition: scale_atmos_dyn_tinteg_rkcommon.F90:270
scale_tracer::k
real(rp), public k
Definition: scale_tracer.F90:44
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_comm_cartesc::comm_vars8_init
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
Definition: scale_comm_cartesC.F90:294
scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_rkwork_alloc
subroutine, public atmos_dyn_tinteg_rkcommon_rkwork_alloc(this)
Definition: scale_atmos_dyn_tinteg_rkcommon.F90:149
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_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_atmos_dyn_tstep_short
module Atmosphere / Dynamical scheme
Definition: scale_atmos_dyn_tstep_short.F90:12
scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_comm_wait
subroutine, public atmos_dyn_tinteg_rkcommon_comm_wait(this)
Definition: scale_atmos_dyn_tinteg_rkcommon.F90:201
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_atmos_dyn_tinteg_short_rk11s8o
module Atmosphere / Dyn Tinteg
Definition: scale_atmos_dyn_tinteg_short_rk11s8o.F90:15
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_comm
subroutine, public atmos_dyn_tinteg_rkcommon_comm(this)
Definition: scale_atmos_dyn_tinteg_rkcommon.F90:185
scale_debug
module DEBUG
Definition: scale_debug.F90:11
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
scale_atmos_dyn_tinteg_rkcommon
module scale_atmos_dyn_tinteg_rkcommon
Definition: scale_atmos_dyn_tinteg_rkcommon.F90:12
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
scale_atmos_dyn_tinteg_short_rk11s8o::atmos_dyn_tinteg_short_rk11s8o
subroutine, public atmos_dyn_tinteg_short_rk11s8o(DENS, MOMZ, MOMX, MOMY, RHOT, PROG, mflx_hi, tflx_hi, DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, DPRES0, CVtot, CORIOLI, num_diff, wdamp_coef, divdmp_coef, DDIV, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_ALONG_STREAM, CDZ, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, PHI, GSQRT, J13G, J23G, J33G, MAPF, REF_pres, REF_dens, BND_W, BND_E, BND_S, BND_N, TwoD, dt)
RK7s6o.
Definition: scale_atmos_dyn_tinteg_short_rk11s8o.F90:222
scale_prof::prof_rapend
subroutine, public prof_rapend(rapname_base, level, disable_barrier)
Save raptime.
Definition: scale_prof.F90:217
scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_updateflux
subroutine, public atmos_dyn_tinteg_rkcommon_updateflux(this, nowstage, io, jo, ko, va_, flux)
Definition: scale_atmos_dyn_tinteg_rkcommon.F90:320
scale_atmos_dyn_tinteg_rkcommon::atmos_dyn_tinteg_rkcommon_nextstage
subroutine, public atmos_dyn_tinteg_rkcommon_nextstage(this, nowstage, io, jo, ko, dt)
Definition: scale_atmos_dyn_tinteg_rkcommon.F90:217
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41
scale_atmos_dyn_common::atmos_dyn_copy_boundary
subroutine, public atmos_dyn_copy_boundary(DENS, MOMZ, MOMX, MOMY, RHOT, PROG, DENS0, MOMZ0, MOMX0, MOMY0, RHOT0, PROG0, BND_W, BND_E, BND_S, BND_N, TwoD)
Definition: scale_atmos_dyn_common.F90:188
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56
scale_atmos_dyn_tinteg_short_rk11s8o::rkcoef_b_11s8o_cooperverner1972
real(rp), dimension(11), parameter, public rkcoef_b_11s8o_cooperverner1972
Definition: scale_atmos_dyn_tinteg_short_rk11s8o.F90:116