FE-Project
Loading...
Searching...
No Matches
scale_timeint_rk.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
2! Warning: This file was generated from common/scale_timeint_rk.F90.erb.
3! Do not edit this file.
4!-------------------------------------------------------------------------------
5!-------------------------------------------------------------------------------
6!> Module common / Runge-Kutta scheme
7!!
8!! @par Description
9!! Driver module to provide various Runge-Kutta schemes.
10!!
11!! @author Yuta Kawai, Team SCALE
12!<
13#include "scaleFElib.h"
15 !-----------------------------------------------------------------------------
16 !
17 !++ Used modules
18 !
19 use scale_precision
20 use scale_io
21 use scale_prc
22 use scale_prof
23 use scale_prc
24
25 !-----------------------------------------------------------------------------
26 implicit none
27 private
28
29 !-----------------------------------------------------------------------------
30 !
31 !++ Public type, procedures
32 !
33
34 !> Derived type to provide RK scheme
35 type, public :: timeint_rk
36 real(RP), private :: dt !< Timestep
37
38 integer, private :: rk_scheme_id_ex !< ID of explicit RK scheme
39 integer, private :: rk_scheme_id_im !< ID of implicit RK scheme
40
41 integer, private, allocatable :: size_each_var(:) !< Dimension information for each variable
42 integer, public :: nstage !< Number of stage in RK scheme
43 integer, public :: tend_buf_size !< Buffer size to store tendencies
44
45 ! For Butcher representation (explicit part)
46 real(rp), public, allocatable :: coef_a_ex(:,:) !< Coefficients A in the Butcher table (explicit part)
47 real(rp), public, allocatable :: coef_b_ex(:) !< Coefficients b in the Butcher table (explicit part)
48 real(rp), public, allocatable :: coef_c_ex(:) !< Coefficients c in the Butcher table (explicit part)
49
50 ! For Shu-Osher representation (explicit part)
51 real(rp), public, allocatable :: coef_sig_ex(:,:)
52 real(rp), public, allocatable :: coef_gam_ex(:,:)
53
54 ! For Butcher representation (implicit part)
55 real(rp), public, allocatable :: coef_a_im(:,:) !< Coefficients A in the Butcher table (implicit part)
56 real(rp), public, allocatable :: coef_b_im(:) !< Coefficients b in the Butcher table (implicit part)
57 real(rp), public, allocatable :: coef_c_im(:) !< Coefficients c in the Butcher table (implicit part)
58
59 integer, allocatable :: tend_buf_indmap(:)
60 integer :: var_num
61 real(rp), public, allocatable :: var_buf1d_ex(:,:,:)
62 real(rp), public, allocatable :: tend_buf1d_ex(:,:,:)
63 real(rp), public, allocatable :: tend_buf1d_im(:,:,:)
64 real(rp), public, allocatable :: var0_1d(:,:)
65 real(rp), private, allocatable :: vartmp_1d(:,:)
66 real(rp), public, allocatable :: var_buf2d_ex(:,:,:,:)
67 real(rp), public, allocatable :: tend_buf2d_ex(:,:,:,:)
68 real(rp), public, allocatable :: tend_buf2d_im(:,:,:,:)
69 real(rp), public, allocatable :: var0_2d(:,:,:)
70 real(rp), private, allocatable :: vartmp_2d(:,:,:)
71 real(rp), public, allocatable :: var_buf3d_ex(:,:,:,:,:)
72 real(rp), public, allocatable :: tend_buf3d_ex(:,:,:,:,:)
73 real(rp), public, allocatable :: tend_buf3d_im(:,:,:,:,:)
74 real(rp), public, allocatable :: var0_3d(:,:,:,:)
75 real(rp), private, allocatable :: vartmp_3d(:,:,:,:)
76
77 logical, private :: low_storage_flag
78 logical, public :: imex_flag
79 integer, public :: ndim
80 contains
81 procedure, public :: init => timeint_rk_init
82 procedure, public :: final => timeint_rk_final
83 procedure, public :: get_implicit_diagfac => timeint_rk_get_implicit_diagfac
84 procedure, public :: get_deltime => timeint_rk_get_deltime
85 procedure, public :: advance1d => timeint_rk_advance1d
86 procedure, public :: advance1d_varlist => timeint_rk_advance1d_varlist
87 procedure, public :: advance_trcvar_1d => timeint_rk_advance_trcvar1d
88 procedure, public :: storevar0_1d => timeint_rk_store_var0_1d
89 procedure, public :: storeimplicit1d => timeint_rk_storeimpl1d
90 procedure, public :: advance2d => timeint_rk_advance2d
91 procedure, public :: advance2d_varlist => timeint_rk_advance2d_varlist
92 procedure, public :: advance_trcvar_2d => timeint_rk_advance_trcvar2d
93 procedure, public :: storevar0_2d => timeint_rk_store_var0_2d
94 procedure, public :: storeimplicit2d => timeint_rk_storeimpl2d
95 procedure, public :: advance3d => timeint_rk_advance3d
96 procedure, public :: advance3d_varlist => timeint_rk_advance3d_varlist
97 procedure, public :: advance_trcvar_3d => timeint_rk_advance_trcvar3d
98 procedure, public :: storevar0_3d => timeint_rk_store_var0_3d
99 procedure, public :: storeimplicit3d => timeint_rk_storeimpl3d
100 generic, public :: advance => advance1d, advance2d, advance3d
101 generic, public :: advance_varlist => advance1d_varlist, advance2d_varlist, advance3d_varlist
102 generic, public :: advance_trcvar => advance_trcvar_1d, advance_trcvar_2d, advance_trcvar_3d
103 generic, public :: storevar0 => storevar0_1d, storevar0_2d, storevar0_3d
104 generic, public :: storeimplicit => storeimplicit1d, storeimplicit2d, storeimplicit3d
105 end type timeint_rk
106
107 type, public :: timeint_rk_var
108 real(rp), pointer :: var1d(:)
109 real(rp), pointer :: var2d(:,:)
110 real(rp), pointer :: var3d(:,:,:)
111 end type timeint_rk_var
112
113 !-----------------------------------------------------------------------------
114 !
115 !++ Public parameters & variables
116 !
117 !-----------------------------------------------------------------------------
118 !
119 !++ Private procedures
120 !
121 !-------------------
122
123contains
124
125!----------------
126
127!> Initialize a object to provide RK scheme
128!!
129!! @param rk_scheme_name Name of RK scheme
130!! @param dt Timestep
131!! @param ndim Number of spatial Dimension
132!! @param var_num Number of variables
133!! @param size_each_var Array to store size of each dimension
134 subroutine timeint_rk_init( this, &
135 rk_scheme_name, dt, var_num, ndim, size_each_var )
136
140 implicit none
141 class(timeint_rk), intent(inout) :: this
142 character(*), intent(in) :: rk_scheme_name
143 real(RP), intent(in) :: dt
144 integer, intent(in) :: var_num
145 integer, intent(in) :: ndim
146 integer, intent(in) :: size_each_var(ndim)
147 !----------------------------------------
148
149 this%dt = dt
150 this%ndim = ndim
151 this%var_num = var_num
152 allocate( this%size_each_var(ndim) )
153 this%size_each_var(:) = size_each_var(:)
154
155 call timeint_rk_butcher_tab_get_info( rk_scheme_name, & ! (in)
156 this%nstage, this%tend_buf_size, & ! (out)
157 this%low_storage_flag, this%imex_flag ) ! (out)
158
159 allocate ( this%coef_a_ex(this%nstage,this%nstage), this%coef_b_ex(this%nstage), this%coef_c_ex(this%nstage) )
160 allocate ( this%coef_sig_ex(this%nstage+1,this%nstage), this%coef_gam_ex(this%nstage+1,this%nstage) )
161! if (this%imex_flag) then
162 allocate ( this%coef_a_im(this%nstage,this%nstage), this%coef_b_im(this%nstage), this%coef_c_im(this%nstage) )
163! end if
164
165 select case(this%ndim)
166 case(1)
167 allocate( this%tend_buf1D_ex(size_each_var(1), var_num, this%tend_buf_size) )
168 if ( this%imex_flag ) then
169 allocate( this%tend_buf1D_im(size_each_var(1), var_num, this%tend_buf_size) )
170 end if
171 allocate( this%var0_1D(size_each_var(1), var_num) )
172 allocate( this%varTmp_1D(size_each_var(1), var_num) )
173 case(2)
174 allocate( this%tend_buf2D_ex(size_each_var(1),size_each_var(2), var_num, this%tend_buf_size) )
175 if ( this%imex_flag ) then
176 allocate( this%tend_buf2D_im(size_each_var(1),size_each_var(2), var_num, this%tend_buf_size) )
177 end if
178 allocate( this%var0_2D(size_each_var(1),size_each_var(2), var_num) )
179 allocate( this%varTmp_2D(size_each_var(1),size_each_var(2), var_num) )
180 case(3)
181 allocate( this%tend_buf3D_ex(size_each_var(1),size_each_var(2),size_each_var(3), var_num, this%tend_buf_size) )
182 if ( this%imex_flag ) then
183 allocate( this%tend_buf3D_im(size_each_var(1),size_each_var(2),size_each_var(3), var_num, this%tend_buf_size) )
184 end if
185 allocate( this%var0_3D(size_each_var(1),size_each_var(2),size_each_var(3), var_num) )
186 allocate( this%varTmp_3D(size_each_var(1),size_each_var(2),size_each_var(3), var_num) )
187 end select
188 allocate( this%tend_buf_indmap(this%nstage) )
189
191 rk_scheme_name, this%nstage, this%imex_flag, & ! (in)
192 this%coef_a_ex, this%coef_b_ex, this%coef_c_ex, & ! (out)
193 this%coef_sig_ex, this%coef_gam_ex, & ! (out)
194 this%coef_a_im, this%coef_b_im, this%coef_c_im, & ! (out)
195 this%tend_buf_indmap ) ! (out)
196
197 return
198 end subroutine timeint_rk_init
199
200!> Finalize a object to provide RK scheme
201!!
202 subroutine timeint_rk_final( this )
203 implicit none
204 class(timeint_rk), intent(inout) :: this
205 !----------------------------------------
206
207 deallocate( this%coef_a_ex, this%coef_b_ex, this%coef_c_ex )
208 deallocate( this%coef_sig_ex, this%coef_gam_ex )
209! if (this%imex_flag) then
210 deallocate( this%coef_a_im, this%coef_b_im, this%coef_c_im )
211! end if
212
213 deallocate( this%size_each_var, this%tend_buf_indmap )
214
215 select case(this%ndim)
216 case(1)
217 deallocate( this%var0_1d )
218 deallocate( this%tend_buf1D_ex )
219 if (this%imex_flag) deallocate( this%tend_buf1D_im )
220 deallocate( this%varTmp_1d )
221 case(2)
222 deallocate( this%var0_2d )
223 deallocate( this%tend_buf2D_ex )
224 if (this%imex_flag) deallocate( this%tend_buf2D_im )
225 deallocate( this%varTmp_2d )
226 case(3)
227 deallocate( this%var0_3d )
228 deallocate( this%tend_buf3D_ex )
229 if (this%imex_flag) deallocate( this%tend_buf3D_im )
230 deallocate( this%varTmp_3d )
231 end select
232
233 return
234 end subroutine timeint_rk_final
235
236!> Get timestep
237!!
238 elemental function timeint_rk_get_deltime( this ) result(dt)
239 implicit none
240 class(timeint_rk), intent(in) :: this
241 real(RP) :: dt
242 !-----------------------------------------------------
243
244 dt = this%dt
245
246 return
247 end function timeint_rk_get_deltime
248
249!> Get a diagonal component of coefficient matrix at a stage for implicit RK scheme
250!!
251!! @param nowstage Current stage
252 elemental function timeint_rk_get_implicit_diagfac( this, nowstage ) result(fac)
253 implicit none
254 class(timeint_rk), intent(in) :: this
255 integer, intent(in) :: nowstage
256 real(RP) :: fac
257 !-----------------------------------------------------
258
259 fac = this%coef_a_im(nowstage,nowstage) * this%dt
260
261 return
262 end function timeint_rk_get_implicit_diagfac
263 !----------------
264
265
266!> Advance a variable at the current stage with explicit RK part
267!!
268!! @param nowstage Current stage
269!! @param q Array to store variable data
270!! @param varID Index of the targeting variable
271!! @param is Index at which begins the loop for the corresponding direction
272!! @param ie Index at which finishes the loop for the corresponding direction
273 subroutine timeint_rk_advance1d( this, nowstage, q, varID, is, ie )
274 implicit none
275 class(timeint_rk), intent(inout) :: this
276 integer, intent(in) :: nowstage
277 real(RP), intent(inout) :: q(:)
278 integer, intent(in) :: varID
279 integer, intent(in) :: is, ie
280 !----------------------------------------
281 if (this%low_storage_flag) then
282 call rk_advance_low_storage1d( this, nowstage, q, varid, is, ie , this%size_each_var(1), this%var_num, &
283 this%var0_1D, this%varTmp_1d, this%tend_buf1D_ex )
284 else
285 if ( this%imex_flag ) then
286 call rk_advance_general1d( this, nowstage, q, varid, is, ie , this%size_each_var(1), this%var_num, &
287 this%var0_1D, this%varTmp_1d, this%tend_buf1D_ex, this%tend_buf1D_im )
288 else
289 call rk_advance_general1d( this, nowstage, q, varid, is, ie , this%size_each_var(1), this%var_num, &
290 this%var0_1D, this%varTmp_1d, this%tend_buf1D_ex )
291 end if
292 end if
293
294 return
295 end subroutine timeint_rk_advance1d
296
297!> Advance variables at the current stage with explicit RK part
298!!
299!! @param nowstage Current stage
300!! @param var_list Array of timeint_rk_var object which saves pointer to variable data
301!! @param varIDs Indices of the targeting variables
302!! @param is Index at which begins the loop for the corresponding direction
303!! @param ie Index at which finishes the loop for the corresponding direction
304!OCL SERIAL
305 subroutine timeint_rk_advance1d_varlist( this, nowstage, var_list, varIDs, is, ie )
306 implicit none
307 class(timeint_rk), intent(inout) :: this
308 integer, intent(in) :: nowstage
309 integer, intent(in) :: varIDs(:)
310 type(timeint_rk_var), intent(inout) :: var_list(size(varIDs))
311 integer, intent(in) :: is, ie
312
313 integer :: i
314 integer :: n
315 !----------------------------------------
316
317 n = size(varids)
318
319 if (this%low_storage_flag) then
320 i = 1
321 do while(i <= n)
322 if ( i+1 <= n ) then
323 call rk_advance_low_storage1d_var2( this, nowstage, var_list(i)%var1D, var_list(i+1)%var1D, varids(i), varids(i+1), is, ie , this%size_each_var(1), this%var_num, &
324 this%var0_1D, this%varTmp_1d, this%tend_buf1D_ex )
325 i = i + 2
326 else
327 call rk_advance_low_storage1d( this, nowstage, var_list(n)%var1D, varids(n), is, ie , this%size_each_var(1), this%var_num, &
328 this%var0_1D, this%varTmp_1d, this%tend_buf1D_ex )
329 i = i + 1
330 end if
331 end do
332 else
333 do i=1, n
334 call this%Advance1D( nowstage, var_list(i)%var1D, varids(i), is, ie )
335 end do
336 end if
337
338 return
339 end subroutine timeint_rk_advance1d_varlist
340
341!> Advance tracer variable data at the current stage with explicit RK part
342!!
343!! @param nowstage Current stage
344!! @param q Array to store tracer variable data
345!! @param varID Index of the targeting variable
346!! @param is Index at which begins the loop for the corresponding direction
347!! @param ie Index at which finishes the loop for the corresponding direction
348!! @param DDENS Array to store density perturbation data at the time level n+1
349!! @param DDENS0 Array to store density perturbation data at the time level n
350!! @param DENS_hyd Array to hydrostatic part of density
351 subroutine timeint_rk_advance_trcvar1d( this, nowstage, q, varID, is, ie , &
352 DDENS, DDENS0, DENS_hyd )
353 implicit none
354 class(timeint_rk), intent(inout) :: this
355 integer, intent(in) :: nowstage
356 real(RP), intent(inout) :: q(:)
357 integer, intent(in) :: varID
358 integer, intent(in) :: is, ie
359 real(RP), intent(in) :: DDENS (:)
360 real(RP), intent(in) :: DDENS0(:)
361 real(RP), intent(in) :: DENS_hyd(:)
362 !----------------------------------------
363
364 if ( this%imex_flag ) then
365 if ( prc_ismaster ) write(*,*) "timeint_rk_advence_trcvar: IMEX is not supported. Check!"
366 call prc_abort
367 end if
368
369 if (this%low_storage_flag) then
370 call rk_advance_trcvar_low_storage1d( this, nowstage, q, varid, is, ie , this%size_each_var(1), this%var_num, &
371 ddens, ddens0, dens_hyd, this%var0_1D, this%varTmp_1d, this%tend_buf1D_ex )
372 else
373 call rk_advance_trcvar_general1d( this, nowstage, q, varid, is, ie , this%size_each_var(1), this%var_num, &
374 ddens, ddens0, dens_hyd, this%var0_1D, this%varTmp_1d, this%tend_buf1D_ex )
375 end if
376
377 return
378 end subroutine timeint_rk_advance_trcvar1d
379
380!> Store variable data at the time level n for the case of IMEX RK scheme
381!!
382!! @param q Array to store variable data at the time level n
383!! @param varID Index of the targeting variable
384!! @param is Index at which begins the loop for the corresponding direction
385!! @param ie Index at which finishes the loop for the corresponding direction
386 subroutine timeint_rk_store_var0_1d( this, q, varID, is, ie )
387 implicit none
388 class(timeint_rk), intent(inout) :: this
389 real(RP), intent(inout) :: q(:)
390 integer, intent(in) :: varID
391 integer, intent(in) :: is, ie
392
393 integer :: i
394 !----------------------------------------
395
396 !$omp parallel private(i)
397 !$omp do
398 do i=is, ie
399 this%var0_1D(i,varid) = q(i)
400 end do
401 !$omp end parallel
402
403 return
404 end subroutine timeint_rk_store_var0_1d
405
406!> Store variable data after the implicit part of current stage in IMEX RK scheme
407!!
408!! @param nowstage Current stage
409!! @param q Array to store variable data after the implicit part of current stage
410!! @param varID Index of the targeting variable
411!! @param is Index at which begins the loop for the corresponding direction
412!! @param ie Index at which finishes the loop for the corresponding direction
413 subroutine timeint_rk_storeimpl1d( this, nowstage, q, varID, is, ie )
414 implicit none
415 class(timeint_rk), intent(inout) :: this
416 integer, intent(in) :: nowstage
417 real(RP), intent(inout) :: q(:)
418 integer, intent(in) :: varID
419 integer, intent(in) :: is, ie
420
421 !----------------------------------------
422
423 call rk_storeimpl_general1d( this, nowstage, q, varid, is, ie , this%size_each_var(1), this%var_num, &
424 this%var0_1D, this%varTmp_1d, this%tend_buf1D_im )
425
426 return
427 end subroutine timeint_rk_storeimpl1d
428
429!> Advance a variable at the current stage with explicit RK part
430!!
431!! @param nowstage Current stage
432!! @param q Array to store variable data
433!! @param varID Index of the targeting variable
434!! @param is Index at which begins the loop for the corresponding direction
435!! @param ie Index at which finishes the loop for the corresponding direction
436!! @param js Index at which begins the loop for the corresponding direction
437!! @param je Index at which finishes the loop for the corresponding direction
438 subroutine timeint_rk_advance2d( this, nowstage, q, varID, is, ie ,js, je )
439 implicit none
440 class(timeint_rk), intent(inout) :: this
441 integer, intent(in) :: nowstage
442 real(RP), intent(inout) :: q(:,:)
443 integer, intent(in) :: varID
444 integer, intent(in) :: is, ie ,js, je
445 !----------------------------------------
446 if (this%low_storage_flag) then
447 call rk_advance_low_storage2d( this, nowstage, q, varid, is, ie ,js, je , this%size_each_var(1),this%size_each_var(2), this%var_num, &
448 this%var0_2D, this%varTmp_2d, this%tend_buf2D_ex )
449 else
450 if ( this%imex_flag ) then
451 call rk_advance_general2d( this, nowstage, q, varid, is, ie ,js, je , this%size_each_var(1),this%size_each_var(2), this%var_num, &
452 this%var0_2D, this%varTmp_2d, this%tend_buf2D_ex, this%tend_buf2D_im )
453 else
454 call rk_advance_general2d( this, nowstage, q, varid, is, ie ,js, je , this%size_each_var(1),this%size_each_var(2), this%var_num, &
455 this%var0_2D, this%varTmp_2d, this%tend_buf2D_ex )
456 end if
457 end if
458
459 return
460 end subroutine timeint_rk_advance2d
461
462!> Advance variables at the current stage with explicit RK part
463!!
464!! @param nowstage Current stage
465!! @param var_list Array of timeint_rk_var object which saves pointer to variable data
466!! @param varIDs Indices of the targeting variables
467!! @param is Index at which begins the loop for the corresponding direction
468!! @param ie Index at which finishes the loop for the corresponding direction
469!! @param js Index at which begins the loop for the corresponding direction
470!! @param je Index at which finishes the loop for the corresponding direction
471!OCL SERIAL
472 subroutine timeint_rk_advance2d_varlist( this, nowstage, var_list, varIDs, is, ie ,js, je )
473 implicit none
474 class(timeint_rk), intent(inout) :: this
475 integer, intent(in) :: nowstage
476 integer, intent(in) :: varIDs(:)
477 type(timeint_rk_var), intent(inout) :: var_list(size(varIDs))
478 integer, intent(in) :: is, ie ,js, je
479
480 integer :: i
481 integer :: n
482 !----------------------------------------
483
484 n = size(varids)
485
486 if (this%low_storage_flag) then
487 i = 1
488 do while(i <= n)
489 if ( i+1 <= n ) then
490 call rk_advance_low_storage2d_var2( this, nowstage, var_list(i)%var2D, var_list(i+1)%var2D, varids(i), varids(i+1), is, ie ,js, je , this%size_each_var(1),this%size_each_var(2), this%var_num, &
491 this%var0_2D, this%varTmp_2d, this%tend_buf2D_ex )
492 i = i + 2
493 else
494 call rk_advance_low_storage2d( this, nowstage, var_list(n)%var2D, varids(n), is, ie ,js, je , this%size_each_var(1),this%size_each_var(2), this%var_num, &
495 this%var0_2D, this%varTmp_2d, this%tend_buf2D_ex )
496 i = i + 1
497 end if
498 end do
499 else
500 do i=1, n
501 call this%Advance2D( nowstage, var_list(i)%var2D, varids(i), is, ie ,js, je )
502 end do
503 end if
504
505 return
506 end subroutine timeint_rk_advance2d_varlist
507
508!> Advance tracer variable data at the current stage with explicit RK part
509!!
510!! @param nowstage Current stage
511!! @param q Array to store tracer variable data
512!! @param varID Index of the targeting variable
513!! @param is Index at which begins the loop for the corresponding direction
514!! @param ie Index at which finishes the loop for the corresponding direction
515!! @param js Index at which begins the loop for the corresponding direction
516!! @param je Index at which finishes the loop for the corresponding direction
517!! @param DDENS Array to store density perturbation data at the time level n+1
518!! @param DDENS0 Array to store density perturbation data at the time level n
519!! @param DENS_hyd Array to hydrostatic part of density
520 subroutine timeint_rk_advance_trcvar2d( this, nowstage, q, varID, is, ie ,js, je , &
521 DDENS, DDENS0, DENS_hyd )
522 implicit none
523 class(timeint_rk), intent(inout) :: this
524 integer, intent(in) :: nowstage
525 real(RP), intent(inout) :: q(:,:)
526 integer, intent(in) :: varID
527 integer, intent(in) :: is, ie ,js, je
528 real(RP), intent(in) :: DDENS (:,:)
529 real(RP), intent(in) :: DDENS0(:,:)
530 real(RP), intent(in) :: DENS_hyd(:,:)
531 !----------------------------------------
532
533 if ( this%imex_flag ) then
534 if ( prc_ismaster ) write(*,*) "timeint_rk_advence_trcvar: IMEX is not supported. Check!"
535 call prc_abort
536 end if
537
538 if (this%low_storage_flag) then
539 call rk_advance_trcvar_low_storage2d( this, nowstage, q, varid, is, ie ,js, je , this%size_each_var(1),this%size_each_var(2), this%var_num, &
540 ddens, ddens0, dens_hyd, this%var0_2D, this%varTmp_2d, this%tend_buf2D_ex )
541 else
542 call rk_advance_trcvar_general2d( this, nowstage, q, varid, is, ie ,js, je , this%size_each_var(1),this%size_each_var(2), this%var_num, &
543 ddens, ddens0, dens_hyd, this%var0_2D, this%varTmp_2d, this%tend_buf2D_ex )
544 end if
545
546 return
547 end subroutine timeint_rk_advance_trcvar2d
548
549!> Store variable data at the time level n for the case of IMEX RK scheme
550!!
551!! @param q Array to store variable data at the time level n
552!! @param varID Index of the targeting variable
553!! @param is Index at which begins the loop for the corresponding direction
554!! @param ie Index at which finishes the loop for the corresponding direction
555!! @param js Index at which begins the loop for the corresponding direction
556!! @param je Index at which finishes the loop for the corresponding direction
557 subroutine timeint_rk_store_var0_2d( this, q, varID, is, ie ,js, je )
558 implicit none
559 class(timeint_rk), intent(inout) :: this
560 real(RP), intent(inout) :: q(:,:)
561 integer, intent(in) :: varID
562 integer, intent(in) :: is, ie ,js, je
563
564 integer :: i,j
565 !----------------------------------------
566
567 !$omp parallel private(i,j)
568 !$omp do
569 do j=js, je
570 do i=is, ie
571 this%var0_2D(i,j,varid) = q(i,j)
572 end do
573 end do
574 !$omp end parallel
575
576 return
577 end subroutine timeint_rk_store_var0_2d
578
579!> Store variable data after the implicit part of current stage in IMEX RK scheme
580!!
581!! @param nowstage Current stage
582!! @param q Array to store variable data after the implicit part of current stage
583!! @param varID Index of the targeting variable
584!! @param is Index at which begins the loop for the corresponding direction
585!! @param ie Index at which finishes the loop for the corresponding direction
586!! @param js Index at which begins the loop for the corresponding direction
587!! @param je Index at which finishes the loop for the corresponding direction
588 subroutine timeint_rk_storeimpl2d( this, nowstage, q, varID, is, ie ,js, je )
589 implicit none
590 class(timeint_rk), intent(inout) :: this
591 integer, intent(in) :: nowstage
592 real(RP), intent(inout) :: q(:,:)
593 integer, intent(in) :: varID
594 integer, intent(in) :: is, ie ,js, je
595
596 !----------------------------------------
597
598 call rk_storeimpl_general2d( this, nowstage, q, varid, is, ie ,js, je , this%size_each_var(1),this%size_each_var(2), this%var_num, &
599 this%var0_2D, this%varTmp_2d, this%tend_buf2D_im )
600
601 return
602 end subroutine timeint_rk_storeimpl2d
603
604!> Advance a variable at the current stage with explicit RK part
605!!
606!! @param nowstage Current stage
607!! @param q Array to store variable data
608!! @param varID Index of the targeting variable
609!! @param is Index at which begins the loop for the corresponding direction
610!! @param ie Index at which finishes the loop for the corresponding direction
611!! @param js Index at which begins the loop for the corresponding direction
612!! @param je Index at which finishes the loop for the corresponding direction
613!! @param ks Index at which begins the loop for the corresponding direction
614!! @param ke Index at which finishes the loop for the corresponding direction
615 subroutine timeint_rk_advance3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke )
616 implicit none
617 class(timeint_rk), intent(inout) :: this
618 integer, intent(in) :: nowstage
619 real(RP), intent(inout) :: q(:,:,:)
620 integer, intent(in) :: varID
621 integer, intent(in) :: is, ie ,js, je ,ks, ke
622 !----------------------------------------
623 if (this%low_storage_flag) then
624 call rk_advance_low_storage3d( this, nowstage, q, varid, is, ie ,js, je ,ks, ke , this%size_each_var(1),this%size_each_var(2),this%size_each_var(3), this%var_num, &
625 this%var0_3D, this%varTmp_3d, this%tend_buf3D_ex )
626 else
627 if ( this%imex_flag ) then
628 call rk_advance_general3d( this, nowstage, q, varid, is, ie ,js, je ,ks, ke , this%size_each_var(1),this%size_each_var(2),this%size_each_var(3), this%var_num, &
629 this%var0_3D, this%varTmp_3d, this%tend_buf3D_ex, this%tend_buf3D_im )
630 else
631 call rk_advance_general3d( this, nowstage, q, varid, is, ie ,js, je ,ks, ke , this%size_each_var(1),this%size_each_var(2),this%size_each_var(3), this%var_num, &
632 this%var0_3D, this%varTmp_3d, this%tend_buf3D_ex )
633 end if
634 end if
635
636 return
637 end subroutine timeint_rk_advance3d
638
639!> Advance variables at the current stage with explicit RK part
640!!
641!! @param nowstage Current stage
642!! @param var_list Array of timeint_rk_var object which saves pointer to variable data
643!! @param varIDs Indices of the targeting variables
644!! @param is Index at which begins the loop for the corresponding direction
645!! @param ie Index at which finishes the loop for the corresponding direction
646!! @param js Index at which begins the loop for the corresponding direction
647!! @param je Index at which finishes the loop for the corresponding direction
648!! @param ks Index at which begins the loop for the corresponding direction
649!! @param ke Index at which finishes the loop for the corresponding direction
650!OCL SERIAL
651 subroutine timeint_rk_advance3d_varlist( this, nowstage, var_list, varIDs, is, ie ,js, je ,ks, ke )
652 implicit none
653 class(timeint_rk), intent(inout) :: this
654 integer, intent(in) :: nowstage
655 integer, intent(in) :: varIDs(:)
656 type(timeint_rk_var), intent(inout) :: var_list(size(varIDs))
657 integer, intent(in) :: is, ie ,js, je ,ks, ke
658
659 integer :: i
660 integer :: n
661 !----------------------------------------
662
663 n = size(varids)
664
665 if (this%low_storage_flag) then
666 i = 1
667 do while(i <= n)
668 if ( i+1 <= n ) then
669 call rk_advance_low_storage3d_var2( this, nowstage, var_list(i)%var3D, var_list(i+1)%var3D, varids(i), varids(i+1), is, ie ,js, je ,ks, ke , this%size_each_var(1),this%size_each_var(2),this%size_each_var(3), this%var_num, &
670 this%var0_3D, this%varTmp_3d, this%tend_buf3D_ex )
671 i = i + 2
672 else
673 call rk_advance_low_storage3d( this, nowstage, var_list(n)%var3D, varids(n), is, ie ,js, je ,ks, ke , this%size_each_var(1),this%size_each_var(2),this%size_each_var(3), this%var_num, &
674 this%var0_3D, this%varTmp_3d, this%tend_buf3D_ex )
675 i = i + 1
676 end if
677 end do
678 else
679 do i=1, n
680 call this%Advance3D( nowstage, var_list(i)%var3D, varids(i), is, ie ,js, je ,ks, ke )
681 end do
682 end if
683
684 return
685 end subroutine timeint_rk_advance3d_varlist
686
687!> Advance tracer variable data at the current stage with explicit RK part
688!!
689!! @param nowstage Current stage
690!! @param q Array to store tracer variable data
691!! @param varID Index of the targeting variable
692!! @param is Index at which begins the loop for the corresponding direction
693!! @param ie Index at which finishes the loop for the corresponding direction
694!! @param js Index at which begins the loop for the corresponding direction
695!! @param je Index at which finishes the loop for the corresponding direction
696!! @param ks Index at which begins the loop for the corresponding direction
697!! @param ke Index at which finishes the loop for the corresponding direction
698!! @param DDENS Array to store density perturbation data at the time level n+1
699!! @param DDENS0 Array to store density perturbation data at the time level n
700!! @param DENS_hyd Array to hydrostatic part of density
701 subroutine timeint_rk_advance_trcvar3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke , &
702 DDENS, DDENS0, DENS_hyd )
703 implicit none
704 class(timeint_rk), intent(inout) :: this
705 integer, intent(in) :: nowstage
706 real(RP), intent(inout) :: q(:,:,:)
707 integer, intent(in) :: varID
708 integer, intent(in) :: is, ie ,js, je ,ks, ke
709 real(RP), intent(in) :: DDENS (:,:,:)
710 real(RP), intent(in) :: DDENS0(:,:,:)
711 real(RP), intent(in) :: DENS_hyd(:,:,:)
712 !----------------------------------------
713
714 if ( this%imex_flag ) then
715 if ( prc_ismaster ) write(*,*) "timeint_rk_advence_trcvar: IMEX is not supported. Check!"
716 call prc_abort
717 end if
718
719 if (this%low_storage_flag) then
720 call rk_advance_trcvar_low_storage3d( this, nowstage, q, varid, is, ie ,js, je ,ks, ke , this%size_each_var(1),this%size_each_var(2),this%size_each_var(3), this%var_num, &
721 ddens, ddens0, dens_hyd, this%var0_3D, this%varTmp_3d, this%tend_buf3D_ex )
722 else
723 call rk_advance_trcvar_general3d( this, nowstage, q, varid, is, ie ,js, je ,ks, ke , this%size_each_var(1),this%size_each_var(2),this%size_each_var(3), this%var_num, &
724 ddens, ddens0, dens_hyd, this%var0_3D, this%varTmp_3d, this%tend_buf3D_ex )
725 end if
726
727 return
728 end subroutine timeint_rk_advance_trcvar3d
729
730!> Store variable data at the time level n for the case of IMEX RK scheme
731!!
732!! @param q Array to store variable data at the time level n
733!! @param varID Index of the targeting variable
734!! @param is Index at which begins the loop for the corresponding direction
735!! @param ie Index at which finishes the loop for the corresponding direction
736!! @param js Index at which begins the loop for the corresponding direction
737!! @param je Index at which finishes the loop for the corresponding direction
738!! @param ks Index at which begins the loop for the corresponding direction
739!! @param ke Index at which finishes the loop for the corresponding direction
740 subroutine timeint_rk_store_var0_3d( this, q, varID, is, ie ,js, je ,ks, ke )
741 implicit none
742 class(timeint_rk), intent(inout) :: this
743 real(RP), intent(inout) :: q(:,:,:)
744 integer, intent(in) :: varID
745 integer, intent(in) :: is, ie ,js, je ,ks, ke
746
747 integer :: i,j,k
748 !----------------------------------------
749
750 !$omp parallel private(i,j,k)
751 !$omp do collapse(2)
752 do k=ks, ke
753 do j=js, je
754 do i=is, ie
755 this%var0_3D(i,j,k,varid) = q(i,j,k)
756 end do
757 end do
758 end do
759 !$omp end parallel
760
761 return
762 end subroutine timeint_rk_store_var0_3d
763
764!> Store variable data after the implicit part of current stage in IMEX RK scheme
765!!
766!! @param nowstage Current stage
767!! @param q Array to store variable data after the implicit part of current stage
768!! @param varID Index of the targeting variable
769!! @param is Index at which begins the loop for the corresponding direction
770!! @param ie Index at which finishes the loop for the corresponding direction
771!! @param js Index at which begins the loop for the corresponding direction
772!! @param je Index at which finishes the loop for the corresponding direction
773!! @param ks Index at which begins the loop for the corresponding direction
774!! @param ke Index at which finishes the loop for the corresponding direction
775 subroutine timeint_rk_storeimpl3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke )
776 implicit none
777 class(timeint_rk), intent(inout) :: this
778 integer, intent(in) :: nowstage
779 real(RP), intent(inout) :: q(:,:,:)
780 integer, intent(in) :: varID
781 integer, intent(in) :: is, ie ,js, je ,ks, ke
782
783 !----------------------------------------
784
785 call rk_storeimpl_general3d( this, nowstage, q, varid, is, ie ,js, je ,ks, ke , this%size_each_var(1),this%size_each_var(2),this%size_each_var(3), this%var_num, &
786 this%var0_3D, this%varTmp_3d, this%tend_buf3D_im )
787
788 return
789 end subroutine timeint_rk_storeimpl3d
790
791!-------------------
792
793
794!> Advance variable data at the current stage with explicit RK part (low-storage case)
795!!
796!! @param nowstage Current stage
797!! @param q Array to store variable data
798!! @param varID Index of the targeting variable
799!! @param is Index at which begins the loop for the corresponding direction
800!! @param ie Index at which finishes the loop for the corresponding direction
801!OCL SERIAL
802 subroutine rk_advance_low_storage1d( this, nowstage, q, varID, is, ie , IA, var_num, &
803 var0_1D, varTmp_1d, tend_buf1D_ex )
804 use scale_const, only: &
805 eps => const_eps
806 implicit none
807 class(timeint_rk), intent(inout) :: this
808 integer, intent(in) :: nowstage
809 integer, intent(in) :: IA
810 real(RP), intent(inout) :: q(IA)
811 integer, intent(in) :: varID
812 integer, intent(in) :: is, ie
813 integer, intent(in) :: var_num
814 real(RP), intent(inout) :: var0_1d(IA,var_num)
815 real(RP), intent(inout) :: varTmp_1d(IA,var_num)
816 real(RP), intent(inout) :: tend_buf1D_ex(IA,var_num,this%tend_buf_size)
817
818 integer :: i
819 real(RP) :: sig_ss
820 real(RP) :: gam_ss
821 real(RP) :: one_Minus_sig_ss
822 real(RP) :: sig_Ns
823 real(RP) :: gam_Ns
824
825 !----------------------------------------
826
827 call prof_rapstart( 'rk_advance_low_storage1D', 3)
828
829 sig_ss = this%coef_sig_ex(nowstage+1,nowstage)
830 sig_ns = this%coef_sig_ex(this%nstage+1,nowstage)
831 one_minus_sig_ss = 1.0_rp - sig_ss
832 gam_ss = this%dt * this%coef_gam_ex(nowstage+1,nowstage)
833 gam_ns = this%dt * this%coef_gam_ex(this%nstage+1,nowstage)
834
835 if ( nowstage == this%nstage ) then
836 !$omp parallel private(i)
837 !$omp do
838 do i=is, ie
839 q(i) = vartmp_1d(i,varid) &
840 + sig_ss * q(i) + gam_ss * tend_buf1d_ex(i,varid,1)
841 end do
842 !$omp end parallel
843 call prof_rapend( 'rk_advance_low_storage1D', 3)
844
845 return
846 end if
847
848 !$omp parallel private(i)
849 if (nowstage == 1) then
850 !$omp do
851 do i=is, ie
852 var0_1d(i,varid) = q(i)
853 vartmp_1d(i,varid) = 0.0_rp
854 end do
855 end if
856
857 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps ) then
858 !$omp do
859 do i=is, ie
860 vartmp_1d(i,varid) = vartmp_1d(i,varid) &
861 + sig_ns * q(i) + gam_ns * tend_buf1d_ex(i,varid,1)
862 end do
863 end if
864
865 !$omp do
866 do i=is, ie
867 q(i) = one_minus_sig_ss * var0_1d(i,varid) &
868 + sig_ss * q(i) + gam_ss * tend_buf1d_ex(i,varid,1)
869 end do
870
871 !$omp end parallel
872
873 call prof_rapend( 'rk_advance_low_storage1D', 3)
874
875 return
876 end subroutine rk_advance_low_storage1d
877
878!> Advance variable data at the current stage with explicit RK part (low-storage case)
879!!
880!! @param nowstage Current stage
881!! @param q1 Array to store variable data
882!! @param q2 Array to store variable data
883!! @param varID1 Index of the targeting variable
884!! @param varID2 Index of the targeting variable
885!! @param is Index at which begins the loop for the corresponding direction
886!! @param ie Index at which finishes the loop for the corresponding direction
887!OCL SERIAL
888 subroutine rk_advance_low_storage1d_var2( this, nowstage, q1, q2, varID1, varID2, is, ie , IA, var_num, &
889 var0_1D, varTmp_1d, tend_buf1D_ex )
890 use scale_const, only: &
891 eps => const_eps
892 implicit none
893 class(timeint_rk), intent(inout) :: this
894 integer, intent(in) :: nowstage
895 integer, intent(in) :: IA
896 real(RP), intent(inout) :: q1(IA)
897 real(RP), intent(inout) :: q2(IA)
898 integer, intent(in) :: varID1
899 integer, intent(in) :: varID2
900 integer, intent(in) :: is, ie
901 integer, intent(in) :: var_num
902 real(RP), intent(inout) :: var0_1d(IA,var_num)
903 real(RP), intent(inout) :: varTmp_1d(IA,var_num)
904 real(RP), intent(inout) :: tend_buf1D_ex(IA,var_num,this%tend_buf_size)
905
906 integer :: i
907 real(RP) :: sig_ss
908 real(RP) :: gam_ss
909 real(RP) :: one_Minus_sig_ss
910 real(RP) :: sig_Ns
911 real(RP) :: gam_Ns
912
913 !----------------------------------------
914
915 call prof_rapstart( 'rk_advance_low_storage1D', 3)
916
917 sig_ss = this%coef_sig_ex(nowstage+1,nowstage)
918 sig_ns = this%coef_sig_ex(this%nstage+1,nowstage)
919 one_minus_sig_ss = 1.0_rp - sig_ss
920 gam_ss = this%dt * this%coef_gam_ex(nowstage+1,nowstage)
921 gam_ns = this%dt * this%coef_gam_ex(this%nstage+1,nowstage)
922
923 if ( nowstage == this%nstage ) then
924 !$omp parallel private(i)
925 !$omp do
926 do i=is, ie
927 q1(i) = vartmp_1d(i,varid1) &
928 + sig_ss * q1(i) + gam_ss * tend_buf1d_ex(i,varid1,1)
929 q2(i) = vartmp_1d(i,varid2) &
930 + sig_ss * q2(i) + gam_ss * tend_buf1d_ex(i,varid2,1)
931 end do
932 !$omp end parallel
933 call prof_rapend( 'rk_advance_low_storage1D', 3)
934
935 return
936 end if
937
938 !$omp parallel private(i)
939 if (nowstage == 1) then
940 !$omp do
941 do i=is, ie
942 var0_1d(i,varid1) = q1(i)
943 var0_1d(i,varid2) = q2(i)
944
945 vartmp_1d(i,varid1) = 0.0_rp
946 vartmp_1d(i,varid2) = 0.0_rp
947 end do
948 end if
949
950 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps ) then
951 !$omp do
952 do i=is, ie
953 vartmp_1d(i,varid1) = vartmp_1d(i,varid1) &
954 + sig_ns * q1(i) + gam_ns * tend_buf1d_ex(i,varid1,1)
955 vartmp_1d(i,varid2) = vartmp_1d(i,varid2) &
956 + sig_ns * q2(i) + gam_ns * tend_buf1d_ex(i,varid2,1)
957 end do
958 end if
959
960 !$omp do
961 do i=is, ie
962 q1(i) = one_minus_sig_ss * var0_1d(i,varid1) &
963 + sig_ss * q1(i) + gam_ss * tend_buf1d_ex(i,varid1,1)
964 q2(i) = one_minus_sig_ss * var0_1d(i,varid2) &
965 + sig_ss * q2(i) + gam_ss * tend_buf1d_ex(i,varid2,1)
966 end do
967
968 !$omp end parallel
969
970 call prof_rapend( 'rk_advance_low_storage1D', 3)
971
972 return
973 end subroutine rk_advance_low_storage1d_var2
974
975!> Advance tracer variable data at the current stage with explicit RK part (low-storage case)
976!!
977!! @param nowstage Current stage
978!! @param q Array to store tracer variable data
979!! @param varID Index of the targeting variable
980!! @param is Index at which begins the loop for the corresponding direction
981!! @param ie Index at which finishes the loop for the corresponding direction
982!! @param DDENS Array to store density perturbation data at the time level n+1
983!! @param DDENS0 Array to store density perturbation data at the time level n
984!! @param DENS_hyd Array to hydrostatic part of density
985!OCL SERIAL
986 subroutine rk_advance_trcvar_low_storage1d( this, nowstage, q, varID, is, ie , IA, var_num, &
987 DDENS, DDENS0, DENS_hyd, &
988 var0_1D, varTmp_1d, tend_buf1D_ex )
989 use scale_const, only: &
990 eps => const_eps
991 implicit none
992 class(timeint_rk), intent(inout) :: this
993 integer, intent(in) :: nowstage
994 integer, intent(in) :: IA
995 real(RP), intent(inout) :: q(IA)
996 integer, intent(in) :: varID
997 integer, intent(in) :: is, ie
998 integer, intent(in) :: var_num
999 real(RP), intent(in) :: DDENS (IA)
1000 real(RP), intent(in) :: DDENS0(IA)
1001 real(RP), intent(in) :: DENS_hyd(IA)
1002 real(RP), intent(inout) :: var0_1d(IA,var_num)
1003 real(RP), intent(inout) :: varTmp_1d(IA,var_num)
1004 real(RP), intent(inout) :: tend_buf1D_ex(IA,var_num,this%tend_buf_size)
1005
1006 integer :: i
1007 real(RP) :: sig_ss
1008 real(RP) :: gam_ss
1009 real(RP) :: one_Minus_sig_ss
1010 real(RP) :: sig_Ns
1011 real(RP) :: gam_Ns
1012 real(RP) :: c_ssm1
1013 real(RP) :: c_ss
1014
1015 real(RP) :: dens_ssm1
1016 real(RP) :: dens_ss
1017
1018 !----------------------------------------
1019
1020 call prof_rapstart( 'rk_advance_trcvar_low_storage1D', 3)
1021
1022 sig_ss = this%coef_sig_ex(nowstage+1,nowstage)
1023 sig_ns = this%coef_sig_ex(this%nstage+1,nowstage)
1024 one_minus_sig_ss = 1.0_rp - sig_ss
1025 gam_ss = this%dt * this%coef_gam_ex(nowstage+1,nowstage)
1026 gam_ns = this%dt * this%coef_gam_ex(this%nstage+1,nowstage)
1027 c_ssm1 = this%coef_c_ex(nowstage)
1028
1029 if ( nowstage == this%nstage ) then
1030 !$omp parallel private(i)
1031 !$omp do
1032 do i=is, ie
1033 dens_ssm1 = dens_hyd(i) + ddens0(i) + c_ssm1 * ( ddens(i) - ddens0(i) )
1034 q(i) = ( vartmp_1d(i,varid) &
1035 + sig_ss * q(i) * dens_ssm1 + gam_ss * tend_buf1d_ex(i,varid,1) ) &
1036 / ( dens_hyd(i) + ddens(i) )
1037 end do
1038 !$omp end parallel
1039 call prof_rapend( 'rk_advance_trcvar_low_storage1D', 3)
1040
1041 return
1042 end if
1043
1044 c_ss = this%coef_c_ex(nowstage+1)
1045
1046 !$omp parallel private(i,dens_ss,dens_ssm1)
1047 if (nowstage == 1) then
1048 !$omp do
1049 do i=is, ie
1050 var0_1d(i,varid) = q(i) * ( dens_hyd(i) + ddens0(i) )
1051 vartmp_1d(i,varid) = 0.0_rp
1052 end do
1053 end if
1054
1055 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps ) then
1056 !$omp do
1057 do i=is, ie
1058 dens_ssm1 = dens_hyd(i) + ddens0(i) + c_ssm1 * ( ddens(i) - ddens0(i) )
1059 vartmp_1d(i,varid) = vartmp_1d(i,varid) &
1060 + sig_ns * q(i) * dens_ssm1 + gam_ns * tend_buf1d_ex(i,varid,1)
1061 end do
1062 end if
1063
1064 !$omp do
1065 do i=is, ie
1066 dens_ssm1 = dens_hyd(i) + ddens0(i) + c_ssm1 * ( ddens(i) - ddens0(i) )
1067 dens_ss = dens_hyd(i) + ddens0(i) + c_ss * ( ddens(i) - ddens0(i) )
1068
1069 q(i) = ( one_minus_sig_ss * var0_1d(i,varid) &
1070 + sig_ss * q(i) * dens_ssm1 + gam_ss * tend_buf1d_ex(i,varid,1) ) &
1071 / dens_ss
1072 end do
1073
1074 !$omp end parallel
1075
1076 call prof_rapend( 'rk_advance_trcvar_low_storage1D', 3)
1077
1078 return
1079 end subroutine rk_advance_trcvar_low_storage1d
1080
1081
1082!> Advance variable data at the current stage with explicit RK part (low-storage case)
1083!!
1084!! @param nowstage Current stage
1085!! @param q Array to store variable data
1086!! @param varID Index of the targeting variable
1087!! @param is Index at which begins the loop for the corresponding direction
1088!! @param ie Index at which finishes the loop for the corresponding direction
1089!! @param js Index at which begins the loop for the corresponding direction
1090!! @param je Index at which finishes the loop for the corresponding direction
1091!OCL SERIAL
1092 subroutine rk_advance_low_storage2d( this, nowstage, q, varID, is, ie ,js, je , IA,JA, var_num, &
1093 var0_2D, varTmp_2d, tend_buf2D_ex )
1094 use scale_const, only: &
1095 eps => const_eps
1096 implicit none
1097 class(timeint_rk), intent(inout) :: this
1098 integer, intent(in) :: nowstage
1099 integer, intent(in) :: IA,JA
1100 real(RP), intent(inout) :: q(IA,JA)
1101 integer, intent(in) :: varID
1102 integer, intent(in) :: is, ie ,js, je
1103 integer, intent(in) :: var_num
1104 real(RP), intent(inout) :: var0_2d(IA,JA,var_num)
1105 real(RP), intent(inout) :: varTmp_2d(IA,JA,var_num)
1106 real(RP), intent(inout) :: tend_buf2D_ex(IA,JA,var_num,this%tend_buf_size)
1107
1108 integer :: i,j
1109 real(RP) :: sig_ss
1110 real(RP) :: gam_ss
1111 real(RP) :: one_Minus_sig_ss
1112 real(RP) :: sig_Ns
1113 real(RP) :: gam_Ns
1114
1115 !----------------------------------------
1116
1117 call prof_rapstart( 'rk_advance_low_storage2D', 3)
1118
1119 sig_ss = this%coef_sig_ex(nowstage+1,nowstage)
1120 sig_ns = this%coef_sig_ex(this%nstage+1,nowstage)
1121 one_minus_sig_ss = 1.0_rp - sig_ss
1122 gam_ss = this%dt * this%coef_gam_ex(nowstage+1,nowstage)
1123 gam_ns = this%dt * this%coef_gam_ex(this%nstage+1,nowstage)
1124
1125 if ( nowstage == this%nstage ) then
1126 !$omp parallel private(i,j)
1127 !$omp do
1128 do j=js, je
1129 do i=is, ie
1130 q(i,j) = vartmp_2d(i,j,varid) &
1131 + sig_ss * q(i,j) + gam_ss * tend_buf2d_ex(i,j,varid,1)
1132 end do
1133 end do
1134 !$omp end parallel
1135 call prof_rapend( 'rk_advance_low_storage2D', 3)
1136
1137 return
1138 end if
1139
1140 !$omp parallel private(i,j)
1141 if (nowstage == 1) then
1142 !$omp do
1143 do j=js, je
1144 do i=is, ie
1145 var0_2d(i,j,varid) = q(i,j)
1146 vartmp_2d(i,j,varid) = 0.0_rp
1147 end do
1148 end do
1149 end if
1150
1151 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps ) then
1152 !$omp do
1153 do j=js, je
1154 do i=is, ie
1155 vartmp_2d(i,j,varid) = vartmp_2d(i,j,varid) &
1156 + sig_ns * q(i,j) + gam_ns * tend_buf2d_ex(i,j,varid,1)
1157 end do
1158 end do
1159 end if
1160
1161 !$omp do
1162 do j=js, je
1163 do i=is, ie
1164 q(i,j) = one_minus_sig_ss * var0_2d(i,j,varid) &
1165 + sig_ss * q(i,j) + gam_ss * tend_buf2d_ex(i,j,varid,1)
1166 end do
1167 end do
1168
1169 !$omp end parallel
1170
1171 call prof_rapend( 'rk_advance_low_storage2D', 3)
1172
1173 return
1174 end subroutine rk_advance_low_storage2d
1175
1176!> Advance variable data at the current stage with explicit RK part (low-storage case)
1177!!
1178!! @param nowstage Current stage
1179!! @param q1 Array to store variable data
1180!! @param q2 Array to store variable data
1181!! @param varID1 Index of the targeting variable
1182!! @param varID2 Index of the targeting variable
1183!! @param is Index at which begins the loop for the corresponding direction
1184!! @param ie Index at which finishes the loop for the corresponding direction
1185!! @param js Index at which begins the loop for the corresponding direction
1186!! @param je Index at which finishes the loop for the corresponding direction
1187!OCL SERIAL
1188 subroutine rk_advance_low_storage2d_var2( this, nowstage, q1, q2, varID1, varID2, is, ie ,js, je , IA,JA, var_num, &
1189 var0_2D, varTmp_2d, tend_buf2D_ex )
1190 use scale_const, only: &
1191 eps => const_eps
1192 implicit none
1193 class(timeint_rk), intent(inout) :: this
1194 integer, intent(in) :: nowstage
1195 integer, intent(in) :: IA,JA
1196 real(RP), intent(inout) :: q1(IA,JA)
1197 real(RP), intent(inout) :: q2(IA,JA)
1198 integer, intent(in) :: varID1
1199 integer, intent(in) :: varID2
1200 integer, intent(in) :: is, ie ,js, je
1201 integer, intent(in) :: var_num
1202 real(RP), intent(inout) :: var0_2d(IA,JA,var_num)
1203 real(RP), intent(inout) :: varTmp_2d(IA,JA,var_num)
1204 real(RP), intent(inout) :: tend_buf2D_ex(IA,JA,var_num,this%tend_buf_size)
1205
1206 integer :: i,j
1207 real(RP) :: sig_ss
1208 real(RP) :: gam_ss
1209 real(RP) :: one_Minus_sig_ss
1210 real(RP) :: sig_Ns
1211 real(RP) :: gam_Ns
1212
1213 !----------------------------------------
1214
1215 call prof_rapstart( 'rk_advance_low_storage2D', 3)
1216
1217 sig_ss = this%coef_sig_ex(nowstage+1,nowstage)
1218 sig_ns = this%coef_sig_ex(this%nstage+1,nowstage)
1219 one_minus_sig_ss = 1.0_rp - sig_ss
1220 gam_ss = this%dt * this%coef_gam_ex(nowstage+1,nowstage)
1221 gam_ns = this%dt * this%coef_gam_ex(this%nstage+1,nowstage)
1222
1223 if ( nowstage == this%nstage ) then
1224 !$omp parallel private(i,j)
1225 !$omp do
1226 do j=js, je
1227 do i=is, ie
1228 q1(i,j) = vartmp_2d(i,j,varid1) &
1229 + sig_ss * q1(i,j) + gam_ss * tend_buf2d_ex(i,j,varid1,1)
1230 q2(i,j) = vartmp_2d(i,j,varid2) &
1231 + sig_ss * q2(i,j) + gam_ss * tend_buf2d_ex(i,j,varid2,1)
1232 end do
1233 end do
1234 !$omp end parallel
1235 call prof_rapend( 'rk_advance_low_storage2D', 3)
1236
1237 return
1238 end if
1239
1240 !$omp parallel private(i,j)
1241 if (nowstage == 1) then
1242 !$omp do
1243 do j=js, je
1244 do i=is, ie
1245 var0_2d(i,j,varid1) = q1(i,j)
1246 var0_2d(i,j,varid2) = q2(i,j)
1247
1248 vartmp_2d(i,j,varid1) = 0.0_rp
1249 vartmp_2d(i,j,varid2) = 0.0_rp
1250 end do
1251 end do
1252 end if
1253
1254 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps ) then
1255 !$omp do
1256 do j=js, je
1257 do i=is, ie
1258 vartmp_2d(i,j,varid1) = vartmp_2d(i,j,varid1) &
1259 + sig_ns * q1(i,j) + gam_ns * tend_buf2d_ex(i,j,varid1,1)
1260 vartmp_2d(i,j,varid2) = vartmp_2d(i,j,varid2) &
1261 + sig_ns * q2(i,j) + gam_ns * tend_buf2d_ex(i,j,varid2,1)
1262 end do
1263 end do
1264 end if
1265
1266 !$omp do
1267 do j=js, je
1268 do i=is, ie
1269 q1(i,j) = one_minus_sig_ss * var0_2d(i,j,varid1) &
1270 + sig_ss * q1(i,j) + gam_ss * tend_buf2d_ex(i,j,varid1,1)
1271 q2(i,j) = one_minus_sig_ss * var0_2d(i,j,varid2) &
1272 + sig_ss * q2(i,j) + gam_ss * tend_buf2d_ex(i,j,varid2,1)
1273 end do
1274 end do
1275
1276 !$omp end parallel
1277
1278 call prof_rapend( 'rk_advance_low_storage2D', 3)
1279
1280 return
1281 end subroutine rk_advance_low_storage2d_var2
1282
1283!> Advance tracer variable data at the current stage with explicit RK part (low-storage case)
1284!!
1285!! @param nowstage Current stage
1286!! @param q Array to store tracer variable data
1287!! @param varID Index of the targeting variable
1288!! @param is Index at which begins the loop for the corresponding direction
1289!! @param ie Index at which finishes the loop for the corresponding direction
1290!! @param js Index at which begins the loop for the corresponding direction
1291!! @param je Index at which finishes the loop for the corresponding direction
1292!! @param DDENS Array to store density perturbation data at the time level n+1
1293!! @param DDENS0 Array to store density perturbation data at the time level n
1294!! @param DENS_hyd Array to hydrostatic part of density
1295!OCL SERIAL
1296 subroutine rk_advance_trcvar_low_storage2d( this, nowstage, q, varID, is, ie ,js, je , IA,JA, var_num, &
1297 DDENS, DDENS0, DENS_hyd, &
1298 var0_2D, varTmp_2d, tend_buf2D_ex )
1299 use scale_const, only: &
1300 eps => const_eps
1301 implicit none
1302 class(timeint_rk), intent(inout) :: this
1303 integer, intent(in) :: nowstage
1304 integer, intent(in) :: IA,JA
1305 real(RP), intent(inout) :: q(IA,JA)
1306 integer, intent(in) :: varID
1307 integer, intent(in) :: is, ie ,js, je
1308 integer, intent(in) :: var_num
1309 real(RP), intent(in) :: DDENS (IA,JA)
1310 real(RP), intent(in) :: DDENS0(IA,JA)
1311 real(RP), intent(in) :: DENS_hyd(IA,JA)
1312 real(RP), intent(inout) :: var0_2d(IA,JA,var_num)
1313 real(RP), intent(inout) :: varTmp_2d(IA,JA,var_num)
1314 real(RP), intent(inout) :: tend_buf2D_ex(IA,JA,var_num,this%tend_buf_size)
1315
1316 integer :: i,j
1317 real(RP) :: sig_ss
1318 real(RP) :: gam_ss
1319 real(RP) :: one_Minus_sig_ss
1320 real(RP) :: sig_Ns
1321 real(RP) :: gam_Ns
1322 real(RP) :: c_ssm1
1323 real(RP) :: c_ss
1324
1325 real(RP) :: dens_ssm1
1326 real(RP) :: dens_ss
1327
1328 !----------------------------------------
1329
1330 call prof_rapstart( 'rk_advance_trcvar_low_storage2D', 3)
1331
1332 sig_ss = this%coef_sig_ex(nowstage+1,nowstage)
1333 sig_ns = this%coef_sig_ex(this%nstage+1,nowstage)
1334 one_minus_sig_ss = 1.0_rp - sig_ss
1335 gam_ss = this%dt * this%coef_gam_ex(nowstage+1,nowstage)
1336 gam_ns = this%dt * this%coef_gam_ex(this%nstage+1,nowstage)
1337 c_ssm1 = this%coef_c_ex(nowstage)
1338
1339 if ( nowstage == this%nstage ) then
1340 !$omp parallel private(i,j)
1341 !$omp do
1342 do j=js, je
1343 do i=is, ie
1344 dens_ssm1 = dens_hyd(i,j) + ddens0(i,j) + c_ssm1 * ( ddens(i,j) - ddens0(i,j) )
1345 q(i,j) = ( vartmp_2d(i,j,varid) &
1346 + sig_ss * q(i,j) * dens_ssm1 + gam_ss * tend_buf2d_ex(i,j,varid,1) ) &
1347 / ( dens_hyd(i,j) + ddens(i,j) )
1348 end do
1349 end do
1350 !$omp end parallel
1351 call prof_rapend( 'rk_advance_trcvar_low_storage2D', 3)
1352
1353 return
1354 end if
1355
1356 c_ss = this%coef_c_ex(nowstage+1)
1357
1358 !$omp parallel private(i,j,dens_ss,dens_ssm1)
1359 if (nowstage == 1) then
1360 !$omp do
1361 do j=js, je
1362 do i=is, ie
1363 var0_2d(i,j,varid) = q(i,j) * ( dens_hyd(i,j) + ddens0(i,j) )
1364 vartmp_2d(i,j,varid) = 0.0_rp
1365 end do
1366 end do
1367 end if
1368
1369 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps ) then
1370 !$omp do
1371 do j=js, je
1372 do i=is, ie
1373 dens_ssm1 = dens_hyd(i,j) + ddens0(i,j) + c_ssm1 * ( ddens(i,j) - ddens0(i,j) )
1374 vartmp_2d(i,j,varid) = vartmp_2d(i,j,varid) &
1375 + sig_ns * q(i,j) * dens_ssm1 + gam_ns * tend_buf2d_ex(i,j,varid,1)
1376 end do
1377 end do
1378 end if
1379
1380 !$omp do
1381 do j=js, je
1382 do i=is, ie
1383 dens_ssm1 = dens_hyd(i,j) + ddens0(i,j) + c_ssm1 * ( ddens(i,j) - ddens0(i,j) )
1384 dens_ss = dens_hyd(i,j) + ddens0(i,j) + c_ss * ( ddens(i,j) - ddens0(i,j) )
1385
1386 q(i,j) = ( one_minus_sig_ss * var0_2d(i,j,varid) &
1387 + sig_ss * q(i,j) * dens_ssm1 + gam_ss * tend_buf2d_ex(i,j,varid,1) ) &
1388 / dens_ss
1389 end do
1390 end do
1391
1392 !$omp end parallel
1393
1394 call prof_rapend( 'rk_advance_trcvar_low_storage2D', 3)
1395
1396 return
1397 end subroutine rk_advance_trcvar_low_storage2d
1398
1399
1400!> Advance variable data at the current stage with explicit RK part (low-storage case)
1401!!
1402!! @param nowstage Current stage
1403!! @param q Array to store variable data
1404!! @param varID Index of the targeting variable
1405!! @param is Index at which begins the loop for the corresponding direction
1406!! @param ie Index at which finishes the loop for the corresponding direction
1407!! @param js Index at which begins the loop for the corresponding direction
1408!! @param je Index at which finishes the loop for the corresponding direction
1409!! @param ks Index at which begins the loop for the corresponding direction
1410!! @param ke Index at which finishes the loop for the corresponding direction
1411!OCL SERIAL
1412 subroutine rk_advance_low_storage3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke , IA,JA,KA, var_num, &
1413 var0_3D, varTmp_3d, tend_buf3D_ex )
1414 use scale_const, only: &
1415 eps => const_eps
1416 implicit none
1417 class(timeint_rk), intent(inout) :: this
1418 integer, intent(in) :: nowstage
1419 integer, intent(in) :: IA,JA,KA
1420 real(RP), intent(inout) :: q(IA,JA,KA)
1421 integer, intent(in) :: varID
1422 integer, intent(in) :: is, ie ,js, je ,ks, ke
1423 integer, intent(in) :: var_num
1424 real(RP), intent(inout) :: var0_3d(IA,JA,KA,var_num)
1425 real(RP), intent(inout) :: varTmp_3d(IA,JA,KA,var_num)
1426 real(RP), intent(inout) :: tend_buf3D_ex(IA,JA,KA,var_num,this%tend_buf_size)
1427
1428 integer :: i,j,k
1429 real(RP) :: sig_ss
1430 real(RP) :: gam_ss
1431 real(RP) :: one_Minus_sig_ss
1432 real(RP) :: sig_Ns
1433 real(RP) :: gam_Ns
1434
1435 !----------------------------------------
1436
1437 call prof_rapstart( 'rk_advance_low_storage3D', 3)
1438
1439 sig_ss = this%coef_sig_ex(nowstage+1,nowstage)
1440 sig_ns = this%coef_sig_ex(this%nstage+1,nowstage)
1441 one_minus_sig_ss = 1.0_rp - sig_ss
1442 gam_ss = this%dt * this%coef_gam_ex(nowstage+1,nowstage)
1443 gam_ns = this%dt * this%coef_gam_ex(this%nstage+1,nowstage)
1444
1445 if ( nowstage == this%nstage ) then
1446 !$omp parallel private(i,j,k)
1447 !$omp do collapse(2)
1448 do k=ks, ke
1449 do j=js, je
1450 do i=is, ie
1451 q(i,j,k) = vartmp_3d(i,j,k,varid) &
1452 + sig_ss * q(i,j,k) + gam_ss * tend_buf3d_ex(i,j,k,varid,1)
1453 end do
1454 end do
1455 end do
1456 !$omp end parallel
1457 call prof_rapend( 'rk_advance_low_storage3D', 3)
1458
1459 return
1460 end if
1461
1462 !$omp parallel private(i,j,k)
1463 if (nowstage == 1) then
1464 !$omp do collapse(2)
1465 do k=ks, ke
1466 do j=js, je
1467 do i=is, ie
1468 var0_3d(i,j,k,varid) = q(i,j,k)
1469 vartmp_3d(i,j,k,varid) = 0.0_rp
1470 end do
1471 end do
1472 end do
1473 end if
1474
1475 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps ) then
1476 !$omp do collapse(2)
1477 do k=ks, ke
1478 do j=js, je
1479 do i=is, ie
1480 vartmp_3d(i,j,k,varid) = vartmp_3d(i,j,k,varid) &
1481 + sig_ns * q(i,j,k) + gam_ns * tend_buf3d_ex(i,j,k,varid,1)
1482 end do
1483 end do
1484 end do
1485 end if
1486
1487 !$omp do collapse(2)
1488 do k=ks, ke
1489 do j=js, je
1490 do i=is, ie
1491 q(i,j,k) = one_minus_sig_ss * var0_3d(i,j,k,varid) &
1492 + sig_ss * q(i,j,k) + gam_ss * tend_buf3d_ex(i,j,k,varid,1)
1493 end do
1494 end do
1495 end do
1496
1497 !$omp end parallel
1498
1499 call prof_rapend( 'rk_advance_low_storage3D', 3)
1500
1501 return
1502 end subroutine rk_advance_low_storage3d
1503
1504!> Advance variable data at the current stage with explicit RK part (low-storage case)
1505!!
1506!! @param nowstage Current stage
1507!! @param q1 Array to store variable data
1508!! @param q2 Array to store variable data
1509!! @param varID1 Index of the targeting variable
1510!! @param varID2 Index of the targeting variable
1511!! @param is Index at which begins the loop for the corresponding direction
1512!! @param ie Index at which finishes the loop for the corresponding direction
1513!! @param js Index at which begins the loop for the corresponding direction
1514!! @param je Index at which finishes the loop for the corresponding direction
1515!! @param ks Index at which begins the loop for the corresponding direction
1516!! @param ke Index at which finishes the loop for the corresponding direction
1517!OCL SERIAL
1518 subroutine rk_advance_low_storage3d_var2( this, nowstage, q1, q2, varID1, varID2, is, ie ,js, je ,ks, ke , IA,JA,KA, var_num, &
1519 var0_3D, varTmp_3d, tend_buf3D_ex )
1520 use scale_const, only: &
1521 eps => const_eps
1522 implicit none
1523 class(timeint_rk), intent(inout) :: this
1524 integer, intent(in) :: nowstage
1525 integer, intent(in) :: IA,JA,KA
1526 real(RP), intent(inout) :: q1(IA,JA,KA)
1527 real(RP), intent(inout) :: q2(IA,JA,KA)
1528 integer, intent(in) :: varID1
1529 integer, intent(in) :: varID2
1530 integer, intent(in) :: is, ie ,js, je ,ks, ke
1531 integer, intent(in) :: var_num
1532 real(RP), intent(inout) :: var0_3d(IA,JA,KA,var_num)
1533 real(RP), intent(inout) :: varTmp_3d(IA,JA,KA,var_num)
1534 real(RP), intent(inout) :: tend_buf3D_ex(IA,JA,KA,var_num,this%tend_buf_size)
1535
1536 integer :: i,j,k
1537 real(RP) :: sig_ss
1538 real(RP) :: gam_ss
1539 real(RP) :: one_Minus_sig_ss
1540 real(RP) :: sig_Ns
1541 real(RP) :: gam_Ns
1542
1543 !----------------------------------------
1544
1545 call prof_rapstart( 'rk_advance_low_storage3D', 3)
1546
1547 sig_ss = this%coef_sig_ex(nowstage+1,nowstage)
1548 sig_ns = this%coef_sig_ex(this%nstage+1,nowstage)
1549 one_minus_sig_ss = 1.0_rp - sig_ss
1550 gam_ss = this%dt * this%coef_gam_ex(nowstage+1,nowstage)
1551 gam_ns = this%dt * this%coef_gam_ex(this%nstage+1,nowstage)
1552
1553 if ( nowstage == this%nstage ) then
1554 !$omp parallel private(i,j,k)
1555 !$omp do collapse(2)
1556 do k=ks, ke
1557 do j=js, je
1558 do i=is, ie
1559 q1(i,j,k) = vartmp_3d(i,j,k,varid1) &
1560 + sig_ss * q1(i,j,k) + gam_ss * tend_buf3d_ex(i,j,k,varid1,1)
1561 q2(i,j,k) = vartmp_3d(i,j,k,varid2) &
1562 + sig_ss * q2(i,j,k) + gam_ss * tend_buf3d_ex(i,j,k,varid2,1)
1563 end do
1564 end do
1565 end do
1566 !$omp end parallel
1567 call prof_rapend( 'rk_advance_low_storage3D', 3)
1568
1569 return
1570 end if
1571
1572 !$omp parallel private(i,j,k)
1573 if (nowstage == 1) then
1574 !$omp do collapse(2)
1575 do k=ks, ke
1576 do j=js, je
1577 do i=is, ie
1578 var0_3d(i,j,k,varid1) = q1(i,j,k)
1579 var0_3d(i,j,k,varid2) = q2(i,j,k)
1580
1581 vartmp_3d(i,j,k,varid1) = 0.0_rp
1582 vartmp_3d(i,j,k,varid2) = 0.0_rp
1583 end do
1584 end do
1585 end do
1586 end if
1587
1588 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps ) then
1589 !$omp do collapse(2)
1590 do k=ks, ke
1591 do j=js, je
1592 do i=is, ie
1593 vartmp_3d(i,j,k,varid1) = vartmp_3d(i,j,k,varid1) &
1594 + sig_ns * q1(i,j,k) + gam_ns * tend_buf3d_ex(i,j,k,varid1,1)
1595 vartmp_3d(i,j,k,varid2) = vartmp_3d(i,j,k,varid2) &
1596 + sig_ns * q2(i,j,k) + gam_ns * tend_buf3d_ex(i,j,k,varid2,1)
1597 end do
1598 end do
1599 end do
1600 end if
1601
1602 !$omp do collapse(2)
1603 do k=ks, ke
1604 do j=js, je
1605 do i=is, ie
1606 q1(i,j,k) = one_minus_sig_ss * var0_3d(i,j,k,varid1) &
1607 + sig_ss * q1(i,j,k) + gam_ss * tend_buf3d_ex(i,j,k,varid1,1)
1608 q2(i,j,k) = one_minus_sig_ss * var0_3d(i,j,k,varid2) &
1609 + sig_ss * q2(i,j,k) + gam_ss * tend_buf3d_ex(i,j,k,varid2,1)
1610 end do
1611 end do
1612 end do
1613
1614 !$omp end parallel
1615
1616 call prof_rapend( 'rk_advance_low_storage3D', 3)
1617
1618 return
1619 end subroutine rk_advance_low_storage3d_var2
1620
1621!> Advance tracer variable data at the current stage with explicit RK part (low-storage case)
1622!!
1623!! @param nowstage Current stage
1624!! @param q Array to store tracer variable data
1625!! @param varID Index of the targeting variable
1626!! @param is Index at which begins the loop for the corresponding direction
1627!! @param ie Index at which finishes the loop for the corresponding direction
1628!! @param js Index at which begins the loop for the corresponding direction
1629!! @param je Index at which finishes the loop for the corresponding direction
1630!! @param ks Index at which begins the loop for the corresponding direction
1631!! @param ke Index at which finishes the loop for the corresponding direction
1632!! @param DDENS Array to store density perturbation data at the time level n+1
1633!! @param DDENS0 Array to store density perturbation data at the time level n
1634!! @param DENS_hyd Array to hydrostatic part of density
1635!OCL SERIAL
1636 subroutine rk_advance_trcvar_low_storage3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke , IA,JA,KA, var_num, &
1637 DDENS, DDENS0, DENS_hyd, &
1638 var0_3D, varTmp_3d, tend_buf3D_ex )
1639 use scale_const, only: &
1640 eps => const_eps
1641 implicit none
1642 class(timeint_rk), intent(inout) :: this
1643 integer, intent(in) :: nowstage
1644 integer, intent(in) :: IA,JA,KA
1645 real(RP), intent(inout) :: q(IA,JA,KA)
1646 integer, intent(in) :: varID
1647 integer, intent(in) :: is, ie ,js, je ,ks, ke
1648 integer, intent(in) :: var_num
1649 real(RP), intent(in) :: DDENS (IA,JA,KA)
1650 real(RP), intent(in) :: DDENS0(IA,JA,KA)
1651 real(RP), intent(in) :: DENS_hyd(IA,JA,KA)
1652 real(RP), intent(inout) :: var0_3d(IA,JA,KA,var_num)
1653 real(RP), intent(inout) :: varTmp_3d(IA,JA,KA,var_num)
1654 real(RP), intent(inout) :: tend_buf3D_ex(IA,JA,KA,var_num,this%tend_buf_size)
1655
1656 integer :: i,j,k
1657 real(RP) :: sig_ss
1658 real(RP) :: gam_ss
1659 real(RP) :: one_Minus_sig_ss
1660 real(RP) :: sig_Ns
1661 real(RP) :: gam_Ns
1662 real(RP) :: c_ssm1
1663 real(RP) :: c_ss
1664
1665 real(RP) :: dens_ssm1
1666 real(RP) :: dens_ss
1667
1668 !----------------------------------------
1669
1670 call prof_rapstart( 'rk_advance_trcvar_low_storage3D', 3)
1671
1672 sig_ss = this%coef_sig_ex(nowstage+1,nowstage)
1673 sig_ns = this%coef_sig_ex(this%nstage+1,nowstage)
1674 one_minus_sig_ss = 1.0_rp - sig_ss
1675 gam_ss = this%dt * this%coef_gam_ex(nowstage+1,nowstage)
1676 gam_ns = this%dt * this%coef_gam_ex(this%nstage+1,nowstage)
1677 c_ssm1 = this%coef_c_ex(nowstage)
1678
1679 if ( nowstage == this%nstage ) then
1680 !$omp parallel private(i,j,k)
1681 !$omp do collapse(2)
1682 do k=ks, ke
1683 do j=js, je
1684 do i=is, ie
1685 dens_ssm1 = dens_hyd(i,j,k) + ddens0(i,j,k) + c_ssm1 * ( ddens(i,j,k) - ddens0(i,j,k) )
1686 q(i,j,k) = ( vartmp_3d(i,j,k,varid) &
1687 + sig_ss * q(i,j,k) * dens_ssm1 + gam_ss * tend_buf3d_ex(i,j,k,varid,1) ) &
1688 / ( dens_hyd(i,j,k) + ddens(i,j,k) )
1689 end do
1690 end do
1691 end do
1692 !$omp end parallel
1693 call prof_rapend( 'rk_advance_trcvar_low_storage3D', 3)
1694
1695 return
1696 end if
1697
1698 c_ss = this%coef_c_ex(nowstage+1)
1699
1700 !$omp parallel private(i,j,k,dens_ss,dens_ssm1)
1701 if (nowstage == 1) then
1702 !$omp do collapse(2)
1703 do k=ks, ke
1704 do j=js, je
1705 do i=is, ie
1706 var0_3d(i,j,k,varid) = q(i,j,k) * ( dens_hyd(i,j,k) + ddens0(i,j,k) )
1707 vartmp_3d(i,j,k,varid) = 0.0_rp
1708 end do
1709 end do
1710 end do
1711 end if
1712
1713 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps ) then
1714 !$omp do collapse(2)
1715 do k=ks, ke
1716 do j=js, je
1717 do i=is, ie
1718 dens_ssm1 = dens_hyd(i,j,k) + ddens0(i,j,k) + c_ssm1 * ( ddens(i,j,k) - ddens0(i,j,k) )
1719 vartmp_3d(i,j,k,varid) = vartmp_3d(i,j,k,varid) &
1720 + sig_ns * q(i,j,k) * dens_ssm1 + gam_ns * tend_buf3d_ex(i,j,k,varid,1)
1721 end do
1722 end do
1723 end do
1724 end if
1725
1726 !$omp do collapse(2)
1727 do k=ks, ke
1728 do j=js, je
1729 do i=is, ie
1730 dens_ssm1 = dens_hyd(i,j,k) + ddens0(i,j,k) + c_ssm1 * ( ddens(i,j,k) - ddens0(i,j,k) )
1731 dens_ss = dens_hyd(i,j,k) + ddens0(i,j,k) + c_ss * ( ddens(i,j,k) - ddens0(i,j,k) )
1732
1733 q(i,j,k) = ( one_minus_sig_ss * var0_3d(i,j,k,varid) &
1734 + sig_ss * q(i,j,k) * dens_ssm1 + gam_ss * tend_buf3d_ex(i,j,k,varid,1) ) &
1735 / dens_ss
1736 end do
1737 end do
1738 end do
1739
1740 !$omp end parallel
1741
1742 call prof_rapend( 'rk_advance_trcvar_low_storage3D', 3)
1743
1744 return
1745 end subroutine rk_advance_trcvar_low_storage3d
1746
1747
1748
1749!> Advance variable data at the current stage with explicit RK part (general case)
1750!!
1751!! @param nowstage Current stage
1752!! @param q Array to store variable data
1753!! @param varID Index of the targeting variable
1754!! @param is Index at which begins the loop for the corresponding direction
1755!! @param ie Index at which finishes the loop for the corresponding direction
1756!OCL SERIAL
1757 subroutine rk_advance_general1d( this, nowstage, q, varID, is, ie , IA, var_num, &
1758 var0_1D, varTmp_1d, tend_buf1D_ex, tend_buf1D_im )
1759 implicit none
1760 class(timeint_rk), intent(inout) :: this
1761 integer, intent(in) :: nowstage
1762 integer, intent(in) :: IA
1763 real(RP), intent(inout) :: q(IA)
1764 integer, intent(in) :: varID
1765 integer, intent(in) :: is, ie
1766 integer, intent(in) :: var_num
1767 real(RP), intent(inout) :: var0_1d(IA,var_num)
1768 real(RP), intent(inout) :: varTmp_1d(IA,var_num)
1769 real(RP), intent(inout) :: tend_buf1D_ex(IA,var_num,this%tend_buf_size)
1770 real(RP), intent(inout), optional :: tend_buf1D_im(IA,var_num,this%tend_buf_size)
1771
1772 integer :: i
1773 integer :: s
1774 integer :: tintbuf_ind
1775
1776 real(RP) :: coef_a_ex_dt
1777 real(RP) :: coef_a_im_dt
1778 real(RP) :: coef_b_ex_dt
1779 real(RP) :: coef_b_im_dt
1780 !----------------------------------------
1781
1782 call prof_rapstart( 'rk_advance_general1D', 3)
1783
1784 tintbuf_ind = this%tend_buf_indmap(nowstage)
1785
1786 if ( this%nstage == 1 ) then
1787 !$omp parallel do
1788 do i=is, ie
1789 vartmp_1d(i,varid) = q(i)
1790 end do
1791 end if
1792
1793 if ( nowstage == this%nstage ) then
1794 if ( this%imex_flag ) then
1795 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
1796 coef_b_im_dt = this%coef_b_im(nowstage) * this%dt
1797 !$omp parallel do
1798 do i=is, ie
1799 q(i) = vartmp_1d(i,varid) &
1800 + coef_b_ex_dt * tend_buf1d_ex(i,varid,tintbuf_ind) &
1801 + coef_b_im_dt * tend_buf1d_im(i,varid,tintbuf_ind)
1802 end do
1803 else
1804 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
1805
1806 !$omp parallel do
1807 do i=is, ie
1808 q(i) = vartmp_1d(i,varid) &
1809 + coef_b_ex_dt * tend_buf1d_ex(i,varid,tintbuf_ind)
1810 end do
1811 end if
1812 call prof_rapend( 'rk_advance_general1D', 3)
1813
1814 return
1815 end if
1816
1817 !$omp parallel private( s, coef_a_ex_dt, coef_a_im_dt, coef_b_ex_dt, coef_b_im_dt ) &
1818 !$omp private( i )
1819
1820 if ( nowstage == 1 .and. (.not. this%imex_flag) ) then
1821 !$omp do
1822 do i=is, ie
1823 var0_1d(i,varid) = q(i)
1824 vartmp_1d(i,varid) = q(i)
1825 end do
1826 end if
1827
1828 if ( this%imex_flag ) then
1829 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
1830 coef_b_im_dt = this%coef_b_im(nowstage) * this%dt
1831
1832 !$omp do
1833 do i=is, ie
1834 q(i) = var0_1d(i,varid)
1835 vartmp_1d(i,varid) = vartmp_1d(i,varid) &
1836 + coef_b_ex_dt * tend_buf1d_ex(i,varid,tintbuf_ind) &
1837 + coef_b_im_dt * tend_buf1d_im(i,varid,tintbuf_ind)
1838 end do
1839 else
1840 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
1841
1842 !$omp do
1843 do i=is, ie
1844 q(i) = var0_1d(i,varid)
1845 vartmp_1d(i,varid) = vartmp_1d(i,varid) &
1846 + coef_b_ex_dt * tend_buf1d_ex(i,varid,tintbuf_ind)
1847 end do
1848 end if
1849
1850 if ( this%tend_buf_size == 1 .and. (.not. this%imex_flag) ) then
1851 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,nowstage)
1852 !$omp do
1853 do i=is, ie
1854 q(i) = var0_1d(i,varid) &
1855 + coef_a_ex_dt * tend_buf1d_ex(i,varid,1)
1856 end do
1857 else if ( .not. this%imex_flag ) then
1858 do s=1, nowstage
1859 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,s)
1860 !$omp do
1861 do i=is, ie
1862 q(i) = q(i) &
1863 + coef_a_ex_dt * tend_buf1d_ex(i,varid,s)
1864 end do
1865 end do
1866 else ! IMEX
1867 do s=1, nowstage
1868 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,s)
1869 coef_a_im_dt = this%dt * this%coef_a_im(nowstage+1,s)
1870 !$omp do
1871 do i=is, ie
1872 q(i) = q(i) &
1873 + coef_a_ex_dt * tend_buf1d_ex(i,varid,s) &
1874 + coef_a_im_dt * tend_buf1d_im(i,varid,s)
1875 end do
1876 end do
1877 end if
1878
1879 !$omp end parallel
1880 call prof_rapend( 'rk_advance_general1D', 3)
1881
1882 return
1883 end subroutine rk_advance_general1d
1884
1885!> Advance tracer variable data at the current stage with explicit RK part (general case)
1886!!
1887!! @param nowstage Current stage
1888!! @param q Array to store tracer variable data
1889!! @param varID Index of the targeting variable
1890!! @param is Index at which begins the loop for the corresponding direction
1891!! @param ie Index at which finishes the loop for the corresponding direction
1892!! @param DDENS Array to store density perturbation data at the time level n+1
1893!! @param DDENS0 Array to store density perturbation data at the time level n
1894!! @param DENS_hyd Array to hydrostatic part of density
1895!OCL SERIAL
1896 subroutine rk_advance_trcvar_general1d( this, nowstage, q, varID, is, ie , IA, var_num, &
1897 DDENS, DDENS0, DENS_hyd, &
1898 var0_1D, varTmp_1d, tend_buf1D_ex )
1899 implicit none
1900 class(timeint_rk), intent(inout) :: this
1901 integer, intent(in) :: nowstage
1902 integer, intent(in) :: IA
1903 real(RP), intent(inout) :: q(IA)
1904 integer, intent(in) :: varID
1905 integer, intent(in) :: is, ie
1906 integer, intent(in) :: var_num
1907 real(RP), intent(in) :: DDENS (IA)
1908 real(RP), intent(in) :: DDENS0(IA)
1909 real(RP), intent(in) :: DENS_hyd(IA)
1910 real(RP), intent(inout) :: var0_1d(IA,var_num)
1911 real(RP), intent(inout) :: varTmp_1d(IA,var_num)
1912 real(RP), intent(inout) :: tend_buf1D_ex(IA,var_num,this%tend_buf_size)
1913
1914 integer :: i
1915 integer :: s
1916 integer :: tintbuf_ind
1917
1918 real(RP) :: dens_
1919
1920 real(RP) :: coef_a_ex
1921 real(RP) :: coef_a_ex_dt
1922 real(RP) :: coef_b_ex_dt
1923 real(RP) :: c_ss
1924 !----------------------------------------
1925
1926 call prof_rapstart( 'rk_advance_trcvar_general1D', 3)
1927
1928 tintbuf_ind = this%tend_buf_indmap(nowstage)
1929
1930 if ( this%nstage == 1 ) then
1931 !$omp parallel do
1932 do i=is, ie
1933 vartmp_1d(i,varid) = q(i) * ( dens_hyd(i) + ddens0(i) )
1934 end do
1935 end if
1936
1937 if ( nowstage == this%nstage ) then
1938 coef_b_ex_dt = this%dt * this%coef_b_ex(nowstage)
1939 !$omp parallel do
1940 do i=is, ie
1941 q(i) = ( vartmp_1d(i,varid) &
1942 + coef_b_ex_dt * tend_buf1d_ex(i,varid,tintbuf_ind) ) &
1943 / ( dens_hyd(i) + ddens(i) )
1944 end do
1945 call prof_rapend( 'rk_advance_trcvar_general1D', 3)
1946
1947 return
1948 end if
1949
1950 c_ss = this%coef_c_ex(nowstage+1)
1951
1952 !$omp parallel private( s, dens_, coef_a_ex, coef_a_ex_dt, coef_b_ex_dt ) &
1953 !$omp private( i )
1954
1955 if ( nowstage == 1 .and. (.not. this%imex_flag) ) then
1956 !$omp do
1957 do i=is, ie
1958 var0_1d(i,varid) = q(i) * ( dens_hyd(i) + ddens0(i) )
1959 vartmp_1d(i,varid) = var0_1d(i,varid)
1960 end do
1961 end if
1962
1963 coef_b_ex_dt = this%dt * this%coef_b_ex(nowstage)
1964 !$omp do
1965 do i=is, ie
1966 q(i) = this%var0_1d(i,varid)
1967 vartmp_1d(i,varid) = vartmp_1d(i,varid) &
1968 + coef_b_ex_dt * tend_buf1d_ex(i,varid,tintbuf_ind)
1969 end do
1970
1971 if ( this%tend_buf_size == 1 .and. (.not. this%imex_flag) ) then
1972 coef_a_ex = this%coef_a_ex(nowstage+1,nowstage)
1973 coef_b_ex_dt = this%dt * this%coef_a_ex(nowstage+1,nowstage)
1974 !$omp do
1975 do i=is, ie
1976 dens_ = ( dens_hyd(i) + ddens0(i) ) &
1977 + coef_a_ex * ( ddens(i) - ddens0(i) )
1978 q(i) = ( var0_1d(i,varid) &
1979 + coef_b_ex_dt * tend_buf1d_ex(i,varid,1) ) &
1980 / dens_
1981 end do
1982 else if ( .not. this%imex_flag ) then
1983 do s=1, nowstage
1984 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,s)
1985 !$omp do
1986 do i=is, ie
1987 q(i) = q(i) &
1988 + coef_a_ex_dt * tend_buf1d_ex(i,varid,s)
1989 end do
1990 end do
1991 !$omp do
1992 do i=is, ie
1993 dens_ = dens_hyd(i) + ddens0(i) + c_ss * ( ddens(i) - ddens0(i) )
1994 q(i) = q(i) / dens_
1995 end do
1996 end if
1997
1998 !$omp end parallel
1999 call prof_rapend( 'rk_advance_trcvar_generall1D', 3)
2000
2001 return
2002 end subroutine rk_advance_trcvar_general1d
2003
2004!> Store variable data after the implicit part of current stage in IMEX RK scheme (general case)
2005!!
2006!! @param nowstage Current stage
2007!! @param q Array to store variable data after the implicit part of current stage
2008!! @param varID Index of the targeting variable
2009!! @param is Index at which begins the loop for the corresponding direction
2010!! @param ie Index at which finishes the loop for the corresponding direction
2011!OCL SERIAL
2012 subroutine rk_storeimpl_general1d( this, nowstage, q, varID, is, ie , IA, var_num, &
2013 var0_1D, varTmp_1d, tend_buf1D_im )
2014 implicit none
2015 class(timeint_rk), intent(inout) :: this
2016 integer, intent(in) :: nowstage
2017 integer, intent(in) :: IA
2018 real(RP), intent(inout) :: q(IA)
2019 integer, intent(in) :: varID
2020 integer, intent(in) :: is, ie
2021 integer, intent(in) :: var_num
2022 real(RP), intent(inout) :: var0_1d(IA,var_num)
2023 real(RP), intent(inout) :: varTmp_1d(IA,var_num)
2024 real(RP), intent(inout) :: tend_buf1D_im(IA,var_num,this%tend_buf_size)
2025
2026 integer :: i
2027 integer :: tintbuf_ind
2028 real(RP) :: coef_a_im_dt
2029 !----------------------------------------
2030
2031 if (.not. this%imex_flag ) return
2032
2033 tintbuf_ind = this%tend_buf_indmap(nowstage)
2034 coef_a_im_dt = this%dt * this%coef_a_im(nowstage,nowstage)
2035
2036 !$omp parallel
2037
2038 if ( nowstage == 1 ) then
2039 !$omp do
2040 do i=is, ie
2041 var0_1d(i,varid) = q(i)
2042 vartmp_1d(i,varid) = q(i)
2043 end do
2044 end if
2045
2046 !$omp do
2047 do i=is, ie
2048 q(i) = q(i) &
2049 + coef_a_im_dt * tend_buf1d_im(i,varid,tintbuf_ind)
2050 end do
2051 !$omp end parallel
2052
2053 return
2054 end subroutine rk_storeimpl_general1d
2055
2056
2057!> Advance variable data at the current stage with explicit RK part (general case)
2058!!
2059!! @param nowstage Current stage
2060!! @param q Array to store variable data
2061!! @param varID Index of the targeting variable
2062!! @param is Index at which begins the loop for the corresponding direction
2063!! @param ie Index at which finishes the loop for the corresponding direction
2064!! @param js Index at which begins the loop for the corresponding direction
2065!! @param je Index at which finishes the loop for the corresponding direction
2066!OCL SERIAL
2067 subroutine rk_advance_general2d( this, nowstage, q, varID, is, ie ,js, je , IA,JA, var_num, &
2068 var0_2D, varTmp_2d, tend_buf2D_ex, tend_buf2D_im )
2069 implicit none
2070 class(timeint_rk), intent(inout) :: this
2071 integer, intent(in) :: nowstage
2072 integer, intent(in) :: IA,JA
2073 real(RP), intent(inout) :: q(IA,JA)
2074 integer, intent(in) :: varID
2075 integer, intent(in) :: is, ie ,js, je
2076 integer, intent(in) :: var_num
2077 real(RP), intent(inout) :: var0_2d(IA,JA,var_num)
2078 real(RP), intent(inout) :: varTmp_2d(IA,JA,var_num)
2079 real(RP), intent(inout) :: tend_buf2D_ex(IA,JA,var_num,this%tend_buf_size)
2080 real(RP), intent(inout), optional :: tend_buf2D_im(IA,JA,var_num,this%tend_buf_size)
2081
2082 integer :: i,j
2083 integer :: s
2084 integer :: tintbuf_ind
2085
2086 real(RP) :: coef_a_ex_dt
2087 real(RP) :: coef_a_im_dt
2088 real(RP) :: coef_b_ex_dt
2089 real(RP) :: coef_b_im_dt
2090 !----------------------------------------
2091
2092 call prof_rapstart( 'rk_advance_general2D', 3)
2093
2094 tintbuf_ind = this%tend_buf_indmap(nowstage)
2095
2096 if ( this%nstage == 1 ) then
2097 !$omp parallel do
2098 do j=js, je
2099 do i=is, ie
2100 vartmp_2d(i,j,varid) = q(i,j)
2101 end do
2102 end do
2103 end if
2104
2105 if ( nowstage == this%nstage ) then
2106 if ( this%imex_flag ) then
2107 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
2108 coef_b_im_dt = this%coef_b_im(nowstage) * this%dt
2109 !$omp parallel do
2110 do j=js, je
2111 do i=is, ie
2112 q(i,j) = vartmp_2d(i,j,varid) &
2113 + coef_b_ex_dt * tend_buf2d_ex(i,j,varid,tintbuf_ind) &
2114 + coef_b_im_dt * tend_buf2d_im(i,j,varid,tintbuf_ind)
2115 end do
2116 end do
2117 else
2118 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
2119
2120 !$omp parallel do
2121 do j=js, je
2122 do i=is, ie
2123 q(i,j) = vartmp_2d(i,j,varid) &
2124 + coef_b_ex_dt * tend_buf2d_ex(i,j,varid,tintbuf_ind)
2125 end do
2126 end do
2127 end if
2128 call prof_rapend( 'rk_advance_general2D', 3)
2129
2130 return
2131 end if
2132
2133 !$omp parallel private( s, coef_a_ex_dt, coef_a_im_dt, coef_b_ex_dt, coef_b_im_dt ) &
2134 !$omp private( i,j )
2135
2136 if ( nowstage == 1 .and. (.not. this%imex_flag) ) then
2137 !$omp do
2138 do j=js, je
2139 do i=is, ie
2140 var0_2d(i,j,varid) = q(i,j)
2141 vartmp_2d(i,j,varid) = q(i,j)
2142 end do
2143 end do
2144 end if
2145
2146 if ( this%imex_flag ) then
2147 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
2148 coef_b_im_dt = this%coef_b_im(nowstage) * this%dt
2149
2150 !$omp do
2151 do j=js, je
2152 do i=is, ie
2153 q(i,j) = var0_2d(i,j,varid)
2154 vartmp_2d(i,j,varid) = vartmp_2d(i,j,varid) &
2155 + coef_b_ex_dt * tend_buf2d_ex(i,j,varid,tintbuf_ind) &
2156 + coef_b_im_dt * tend_buf2d_im(i,j,varid,tintbuf_ind)
2157 end do
2158 end do
2159 else
2160 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
2161
2162 !$omp do
2163 do j=js, je
2164 do i=is, ie
2165 q(i,j) = var0_2d(i,j,varid)
2166 vartmp_2d(i,j,varid) = vartmp_2d(i,j,varid) &
2167 + coef_b_ex_dt * tend_buf2d_ex(i,j,varid,tintbuf_ind)
2168 end do
2169 end do
2170 end if
2171
2172 if ( this%tend_buf_size == 1 .and. (.not. this%imex_flag) ) then
2173 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,nowstage)
2174 !$omp do
2175 do j=js, je
2176 do i=is, ie
2177 q(i,j) = var0_2d(i,j,varid) &
2178 + coef_a_ex_dt * tend_buf2d_ex(i,j,varid,1)
2179 end do
2180 end do
2181 else if ( .not. this%imex_flag ) then
2182 do s=1, nowstage
2183 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,s)
2184 !$omp do
2185 do j=js, je
2186 do i=is, ie
2187 q(i,j) = q(i,j) &
2188 + coef_a_ex_dt * tend_buf2d_ex(i,j,varid,s)
2189 end do
2190 end do
2191 end do
2192 else ! IMEX
2193 do s=1, nowstage
2194 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,s)
2195 coef_a_im_dt = this%dt * this%coef_a_im(nowstage+1,s)
2196 !$omp do
2197 do j=js, je
2198 do i=is, ie
2199 q(i,j) = q(i,j) &
2200 + coef_a_ex_dt * tend_buf2d_ex(i,j,varid,s) &
2201 + coef_a_im_dt * tend_buf2d_im(i,j,varid,s)
2202 end do
2203 end do
2204 end do
2205 end if
2206
2207 !$omp end parallel
2208 call prof_rapend( 'rk_advance_general2D', 3)
2209
2210 return
2211 end subroutine rk_advance_general2d
2212
2213!> Advance tracer variable data at the current stage with explicit RK part (general case)
2214!!
2215!! @param nowstage Current stage
2216!! @param q Array to store tracer variable data
2217!! @param varID Index of the targeting variable
2218!! @param is Index at which begins the loop for the corresponding direction
2219!! @param ie Index at which finishes the loop for the corresponding direction
2220!! @param js Index at which begins the loop for the corresponding direction
2221!! @param je Index at which finishes the loop for the corresponding direction
2222!! @param DDENS Array to store density perturbation data at the time level n+1
2223!! @param DDENS0 Array to store density perturbation data at the time level n
2224!! @param DENS_hyd Array to hydrostatic part of density
2225!OCL SERIAL
2226 subroutine rk_advance_trcvar_general2d( this, nowstage, q, varID, is, ie ,js, je , IA,JA, var_num, &
2227 DDENS, DDENS0, DENS_hyd, &
2228 var0_2D, varTmp_2d, tend_buf2D_ex )
2229 implicit none
2230 class(timeint_rk), intent(inout) :: this
2231 integer, intent(in) :: nowstage
2232 integer, intent(in) :: IA,JA
2233 real(RP), intent(inout) :: q(IA,JA)
2234 integer, intent(in) :: varID
2235 integer, intent(in) :: is, ie ,js, je
2236 integer, intent(in) :: var_num
2237 real(RP), intent(in) :: DDENS (IA,JA)
2238 real(RP), intent(in) :: DDENS0(IA,JA)
2239 real(RP), intent(in) :: DENS_hyd(IA,JA)
2240 real(RP), intent(inout) :: var0_2d(IA,JA,var_num)
2241 real(RP), intent(inout) :: varTmp_2d(IA,JA,var_num)
2242 real(RP), intent(inout) :: tend_buf2D_ex(IA,JA,var_num,this%tend_buf_size)
2243
2244 integer :: i,j
2245 integer :: s
2246 integer :: tintbuf_ind
2247
2248 real(RP) :: dens_
2249
2250 real(RP) :: coef_a_ex
2251 real(RP) :: coef_a_ex_dt
2252 real(RP) :: coef_b_ex_dt
2253 real(RP) :: c_ss
2254 !----------------------------------------
2255
2256 call prof_rapstart( 'rk_advance_trcvar_general2D', 3)
2257
2258 tintbuf_ind = this%tend_buf_indmap(nowstage)
2259
2260 if ( this%nstage == 1 ) then
2261 !$omp parallel do
2262 do j=js, je
2263 do i=is, ie
2264 vartmp_2d(i,j,varid) = q(i,j) * ( dens_hyd(i,j) + ddens0(i,j) )
2265 end do
2266 end do
2267 end if
2268
2269 if ( nowstage == this%nstage ) then
2270 coef_b_ex_dt = this%dt * this%coef_b_ex(nowstage)
2271 !$omp parallel do
2272 do j=js, je
2273 do i=is, ie
2274 q(i,j) = ( vartmp_2d(i,j,varid) &
2275 + coef_b_ex_dt * tend_buf2d_ex(i,j,varid,tintbuf_ind) ) &
2276 / ( dens_hyd(i,j) + ddens(i,j) )
2277 end do
2278 end do
2279 call prof_rapend( 'rk_advance_trcvar_general2D', 3)
2280
2281 return
2282 end if
2283
2284 c_ss = this%coef_c_ex(nowstage+1)
2285
2286 !$omp parallel private( s, dens_, coef_a_ex, coef_a_ex_dt, coef_b_ex_dt ) &
2287 !$omp private( i,j )
2288
2289 if ( nowstage == 1 .and. (.not. this%imex_flag) ) then
2290 !$omp do
2291 do j=js, je
2292 do i=is, ie
2293 var0_2d(i,j,varid) = q(i,j) * ( dens_hyd(i,j) + ddens0(i,j) )
2294 vartmp_2d(i,j,varid) = var0_2d(i,j,varid)
2295 end do
2296 end do
2297 end if
2298
2299 coef_b_ex_dt = this%dt * this%coef_b_ex(nowstage)
2300 !$omp do
2301 do j=js, je
2302 do i=is, ie
2303 q(i,j) = this%var0_2d(i,j,varid)
2304 vartmp_2d(i,j,varid) = vartmp_2d(i,j,varid) &
2305 + coef_b_ex_dt * tend_buf2d_ex(i,j,varid,tintbuf_ind)
2306 end do
2307 end do
2308
2309 if ( this%tend_buf_size == 1 .and. (.not. this%imex_flag) ) then
2310 coef_a_ex = this%coef_a_ex(nowstage+1,nowstage)
2311 coef_b_ex_dt = this%dt * this%coef_a_ex(nowstage+1,nowstage)
2312 !$omp do
2313 do j=js, je
2314 do i=is, ie
2315 dens_ = ( dens_hyd(i,j) + ddens0(i,j) ) &
2316 + coef_a_ex * ( ddens(i,j) - ddens0(i,j) )
2317 q(i,j) = ( var0_2d(i,j,varid) &
2318 + coef_b_ex_dt * tend_buf2d_ex(i,j,varid,1) ) &
2319 / dens_
2320 end do
2321 end do
2322 else if ( .not. this%imex_flag ) then
2323 do s=1, nowstage
2324 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,s)
2325 !$omp do
2326 do j=js, je
2327 do i=is, ie
2328 q(i,j) = q(i,j) &
2329 + coef_a_ex_dt * tend_buf2d_ex(i,j,varid,s)
2330 end do
2331 end do
2332 end do
2333 !$omp do
2334 do j=js, je
2335 do i=is, ie
2336 dens_ = dens_hyd(i,j) + ddens0(i,j) + c_ss * ( ddens(i,j) - ddens0(i,j) )
2337 q(i,j) = q(i,j) / dens_
2338 end do
2339 end do
2340 end if
2341
2342 !$omp end parallel
2343 call prof_rapend( 'rk_advance_trcvar_generall2D', 3)
2344
2345 return
2346 end subroutine rk_advance_trcvar_general2d
2347
2348!> Store variable data after the implicit part of current stage in IMEX RK scheme (general case)
2349!!
2350!! @param nowstage Current stage
2351!! @param q Array to store variable data after the implicit part of current stage
2352!! @param varID Index of the targeting variable
2353!! @param is Index at which begins the loop for the corresponding direction
2354!! @param ie Index at which finishes the loop for the corresponding direction
2355!! @param js Index at which begins the loop for the corresponding direction
2356!! @param je Index at which finishes the loop for the corresponding direction
2357!OCL SERIAL
2358 subroutine rk_storeimpl_general2d( this, nowstage, q, varID, is, ie ,js, je , IA,JA, var_num, &
2359 var0_2D, varTmp_2d, tend_buf2D_im )
2360 implicit none
2361 class(timeint_rk), intent(inout) :: this
2362 integer, intent(in) :: nowstage
2363 integer, intent(in) :: IA,JA
2364 real(RP), intent(inout) :: q(IA,JA)
2365 integer, intent(in) :: varID
2366 integer, intent(in) :: is, ie ,js, je
2367 integer, intent(in) :: var_num
2368 real(RP), intent(inout) :: var0_2d(IA,JA,var_num)
2369 real(RP), intent(inout) :: varTmp_2d(IA,JA,var_num)
2370 real(RP), intent(inout) :: tend_buf2D_im(IA,JA,var_num,this%tend_buf_size)
2371
2372 integer :: i,j
2373 integer :: tintbuf_ind
2374 real(RP) :: coef_a_im_dt
2375 !----------------------------------------
2376
2377 if (.not. this%imex_flag ) return
2378
2379 tintbuf_ind = this%tend_buf_indmap(nowstage)
2380 coef_a_im_dt = this%dt * this%coef_a_im(nowstage,nowstage)
2381
2382 !$omp parallel
2383
2384 if ( nowstage == 1 ) then
2385 !$omp do
2386 do j=js, je
2387 do i=is, ie
2388 var0_2d(i,j,varid) = q(i,j)
2389 vartmp_2d(i,j,varid) = q(i,j)
2390 end do
2391 end do
2392 end if
2393
2394 !$omp do
2395 do j=js, je
2396 do i=is, ie
2397 q(i,j) = q(i,j) &
2398 + coef_a_im_dt * tend_buf2d_im(i,j,varid,tintbuf_ind)
2399 end do
2400 end do
2401 !$omp end parallel
2402
2403 return
2404 end subroutine rk_storeimpl_general2d
2405
2406
2407!> Advance variable data at the current stage with explicit RK part (general case)
2408!!
2409!! @param nowstage Current stage
2410!! @param q Array to store variable data
2411!! @param varID Index of the targeting variable
2412!! @param is Index at which begins the loop for the corresponding direction
2413!! @param ie Index at which finishes the loop for the corresponding direction
2414!! @param js Index at which begins the loop for the corresponding direction
2415!! @param je Index at which finishes the loop for the corresponding direction
2416!! @param ks Index at which begins the loop for the corresponding direction
2417!! @param ke Index at which finishes the loop for the corresponding direction
2418!OCL SERIAL
2419 subroutine rk_advance_general3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke , IA,JA,KA, var_num, &
2420 var0_3D, varTmp_3d, tend_buf3D_ex, tend_buf3D_im )
2421 implicit none
2422 class(timeint_rk), intent(inout) :: this
2423 integer, intent(in) :: nowstage
2424 integer, intent(in) :: IA,JA,KA
2425 real(RP), intent(inout) :: q(IA,JA,KA)
2426 integer, intent(in) :: varID
2427 integer, intent(in) :: is, ie ,js, je ,ks, ke
2428 integer, intent(in) :: var_num
2429 real(RP), intent(inout) :: var0_3d(IA,JA,KA,var_num)
2430 real(RP), intent(inout) :: varTmp_3d(IA,JA,KA,var_num)
2431 real(RP), intent(inout) :: tend_buf3D_ex(IA,JA,KA,var_num,this%tend_buf_size)
2432 real(RP), intent(inout), optional :: tend_buf3D_im(IA,JA,KA,var_num,this%tend_buf_size)
2433
2434 integer :: i,j,k
2435 integer :: s
2436 integer :: tintbuf_ind
2437
2438 real(RP) :: coef_a_ex_dt
2439 real(RP) :: coef_a_im_dt
2440 real(RP) :: coef_b_ex_dt
2441 real(RP) :: coef_b_im_dt
2442 !----------------------------------------
2443
2444 call prof_rapstart( 'rk_advance_general3D', 3)
2445
2446 tintbuf_ind = this%tend_buf_indmap(nowstage)
2447
2448 if ( this%nstage == 1 ) then
2449 !$omp parallel do collapse(2)
2450 do k=ks, ke
2451 do j=js, je
2452 do i=is, ie
2453 vartmp_3d(i,j,k,varid) = q(i,j,k)
2454 end do
2455 end do
2456 end do
2457 end if
2458
2459 if ( nowstage == this%nstage ) then
2460 if ( this%imex_flag ) then
2461 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
2462 coef_b_im_dt = this%coef_b_im(nowstage) * this%dt
2463 !$omp parallel do collapse(2)
2464 do k=ks, ke
2465 do j=js, je
2466 do i=is, ie
2467 q(i,j,k) = vartmp_3d(i,j,k,varid) &
2468 + coef_b_ex_dt * tend_buf3d_ex(i,j,k,varid,tintbuf_ind) &
2469 + coef_b_im_dt * tend_buf3d_im(i,j,k,varid,tintbuf_ind)
2470 end do
2471 end do
2472 end do
2473 else
2474 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
2475
2476 !$omp parallel do collapse(2)
2477 do k=ks, ke
2478 do j=js, je
2479 do i=is, ie
2480 q(i,j,k) = vartmp_3d(i,j,k,varid) &
2481 + coef_b_ex_dt * tend_buf3d_ex(i,j,k,varid,tintbuf_ind)
2482 end do
2483 end do
2484 end do
2485 end if
2486 call prof_rapend( 'rk_advance_general3D', 3)
2487
2488 return
2489 end if
2490
2491 !$omp parallel private( s, coef_a_ex_dt, coef_a_im_dt, coef_b_ex_dt, coef_b_im_dt ) &
2492 !$omp private( i,j,k )
2493
2494 if ( nowstage == 1 .and. (.not. this%imex_flag) ) then
2495 !$omp do collapse(2)
2496 do k=ks, ke
2497 do j=js, je
2498 do i=is, ie
2499 var0_3d(i,j,k,varid) = q(i,j,k)
2500 vartmp_3d(i,j,k,varid) = q(i,j,k)
2501 end do
2502 end do
2503 end do
2504 end if
2505
2506 if ( this%imex_flag ) then
2507 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
2508 coef_b_im_dt = this%coef_b_im(nowstage) * this%dt
2509
2510 !$omp do collapse(2)
2511 do k=ks, ke
2512 do j=js, je
2513 do i=is, ie
2514 q(i,j,k) = var0_3d(i,j,k,varid)
2515 vartmp_3d(i,j,k,varid) = vartmp_3d(i,j,k,varid) &
2516 + coef_b_ex_dt * tend_buf3d_ex(i,j,k,varid,tintbuf_ind) &
2517 + coef_b_im_dt * tend_buf3d_im(i,j,k,varid,tintbuf_ind)
2518 end do
2519 end do
2520 end do
2521 else
2522 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
2523
2524 !$omp do collapse(2)
2525 do k=ks, ke
2526 do j=js, je
2527 do i=is, ie
2528 q(i,j,k) = var0_3d(i,j,k,varid)
2529 vartmp_3d(i,j,k,varid) = vartmp_3d(i,j,k,varid) &
2530 + coef_b_ex_dt * tend_buf3d_ex(i,j,k,varid,tintbuf_ind)
2531 end do
2532 end do
2533 end do
2534 end if
2535
2536 if ( this%tend_buf_size == 1 .and. (.not. this%imex_flag) ) then
2537 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,nowstage)
2538 !$omp do collapse(2)
2539 do k=ks, ke
2540 do j=js, je
2541 do i=is, ie
2542 q(i,j,k) = var0_3d(i,j,k,varid) &
2543 + coef_a_ex_dt * tend_buf3d_ex(i,j,k,varid,1)
2544 end do
2545 end do
2546 end do
2547 else if ( .not. this%imex_flag ) then
2548 do s=1, nowstage
2549 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,s)
2550 !$omp do collapse(2)
2551 do k=ks, ke
2552 do j=js, je
2553 do i=is, ie
2554 q(i,j,k) = q(i,j,k) &
2555 + coef_a_ex_dt * tend_buf3d_ex(i,j,k,varid,s)
2556 end do
2557 end do
2558 end do
2559 end do
2560 else ! IMEX
2561 do s=1, nowstage
2562 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,s)
2563 coef_a_im_dt = this%dt * this%coef_a_im(nowstage+1,s)
2564 !$omp do collapse(2)
2565 do k=ks, ke
2566 do j=js, je
2567 do i=is, ie
2568 q(i,j,k) = q(i,j,k) &
2569 + coef_a_ex_dt * tend_buf3d_ex(i,j,k,varid,s) &
2570 + coef_a_im_dt * tend_buf3d_im(i,j,k,varid,s)
2571 end do
2572 end do
2573 end do
2574 end do
2575 end if
2576
2577 !$omp end parallel
2578 call prof_rapend( 'rk_advance_general3D', 3)
2579
2580 return
2581 end subroutine rk_advance_general3d
2582
2583!> Advance tracer variable data at the current stage with explicit RK part (general case)
2584!!
2585!! @param nowstage Current stage
2586!! @param q Array to store tracer variable data
2587!! @param varID Index of the targeting variable
2588!! @param is Index at which begins the loop for the corresponding direction
2589!! @param ie Index at which finishes the loop for the corresponding direction
2590!! @param js Index at which begins the loop for the corresponding direction
2591!! @param je Index at which finishes the loop for the corresponding direction
2592!! @param ks Index at which begins the loop for the corresponding direction
2593!! @param ke Index at which finishes the loop for the corresponding direction
2594!! @param DDENS Array to store density perturbation data at the time level n+1
2595!! @param DDENS0 Array to store density perturbation data at the time level n
2596!! @param DENS_hyd Array to hydrostatic part of density
2597!OCL SERIAL
2598 subroutine rk_advance_trcvar_general3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke , IA,JA,KA, var_num, &
2599 DDENS, DDENS0, DENS_hyd, &
2600 var0_3D, varTmp_3d, tend_buf3D_ex )
2601 implicit none
2602 class(timeint_rk), intent(inout) :: this
2603 integer, intent(in) :: nowstage
2604 integer, intent(in) :: IA,JA,KA
2605 real(RP), intent(inout) :: q(IA,JA,KA)
2606 integer, intent(in) :: varID
2607 integer, intent(in) :: is, ie ,js, je ,ks, ke
2608 integer, intent(in) :: var_num
2609 real(RP), intent(in) :: DDENS (IA,JA,KA)
2610 real(RP), intent(in) :: DDENS0(IA,JA,KA)
2611 real(RP), intent(in) :: DENS_hyd(IA,JA,KA)
2612 real(RP), intent(inout) :: var0_3d(IA,JA,KA,var_num)
2613 real(RP), intent(inout) :: varTmp_3d(IA,JA,KA,var_num)
2614 real(RP), intent(inout) :: tend_buf3D_ex(IA,JA,KA,var_num,this%tend_buf_size)
2615
2616 integer :: i,j,k
2617 integer :: s
2618 integer :: tintbuf_ind
2619
2620 real(RP) :: dens_
2621
2622 real(RP) :: coef_a_ex
2623 real(RP) :: coef_a_ex_dt
2624 real(RP) :: coef_b_ex_dt
2625 real(RP) :: c_ss
2626 !----------------------------------------
2627
2628 call prof_rapstart( 'rk_advance_trcvar_general3D', 3)
2629
2630 tintbuf_ind = this%tend_buf_indmap(nowstage)
2631
2632 if ( this%nstage == 1 ) then
2633 !$omp parallel do collapse(2)
2634 do k=ks, ke
2635 do j=js, je
2636 do i=is, ie
2637 vartmp_3d(i,j,k,varid) = q(i,j,k) * ( dens_hyd(i,j,k) + ddens0(i,j,k) )
2638 end do
2639 end do
2640 end do
2641 end if
2642
2643 if ( nowstage == this%nstage ) then
2644 coef_b_ex_dt = this%dt * this%coef_b_ex(nowstage)
2645 !$omp parallel do collapse(2)
2646 do k=ks, ke
2647 do j=js, je
2648 do i=is, ie
2649 q(i,j,k) = ( vartmp_3d(i,j,k,varid) &
2650 + coef_b_ex_dt * tend_buf3d_ex(i,j,k,varid,tintbuf_ind) ) &
2651 / ( dens_hyd(i,j,k) + ddens(i,j,k) )
2652 end do
2653 end do
2654 end do
2655 call prof_rapend( 'rk_advance_trcvar_general3D', 3)
2656
2657 return
2658 end if
2659
2660 c_ss = this%coef_c_ex(nowstage+1)
2661
2662 !$omp parallel private( s, dens_, coef_a_ex, coef_a_ex_dt, coef_b_ex_dt ) &
2663 !$omp private( i,j,k )
2664
2665 if ( nowstage == 1 .and. (.not. this%imex_flag) ) then
2666 !$omp do collapse(2)
2667 do k=ks, ke
2668 do j=js, je
2669 do i=is, ie
2670 var0_3d(i,j,k,varid) = q(i,j,k) * ( dens_hyd(i,j,k) + ddens0(i,j,k) )
2671 vartmp_3d(i,j,k,varid) = var0_3d(i,j,k,varid)
2672 end do
2673 end do
2674 end do
2675 end if
2676
2677 coef_b_ex_dt = this%dt * this%coef_b_ex(nowstage)
2678 !$omp do collapse(2)
2679 do k=ks, ke
2680 do j=js, je
2681 do i=is, ie
2682 q(i,j,k) = this%var0_3d(i,j,k,varid)
2683 vartmp_3d(i,j,k,varid) = vartmp_3d(i,j,k,varid) &
2684 + coef_b_ex_dt * tend_buf3d_ex(i,j,k,varid,tintbuf_ind)
2685 end do
2686 end do
2687 end do
2688
2689 if ( this%tend_buf_size == 1 .and. (.not. this%imex_flag) ) then
2690 coef_a_ex = this%coef_a_ex(nowstage+1,nowstage)
2691 coef_b_ex_dt = this%dt * this%coef_a_ex(nowstage+1,nowstage)
2692 !$omp do collapse(2)
2693 do k=ks, ke
2694 do j=js, je
2695 do i=is, ie
2696 dens_ = ( dens_hyd(i,j,k) + ddens0(i,j,k) ) &
2697 + coef_a_ex * ( ddens(i,j,k) - ddens0(i,j,k) )
2698 q(i,j,k) = ( var0_3d(i,j,k,varid) &
2699 + coef_b_ex_dt * tend_buf3d_ex(i,j,k,varid,1) ) &
2700 / dens_
2701 end do
2702 end do
2703 end do
2704 else if ( .not. this%imex_flag ) then
2705 do s=1, nowstage
2706 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,s)
2707 !$omp do collapse(2)
2708 do k=ks, ke
2709 do j=js, je
2710 do i=is, ie
2711 q(i,j,k) = q(i,j,k) &
2712 + coef_a_ex_dt * tend_buf3d_ex(i,j,k,varid,s)
2713 end do
2714 end do
2715 end do
2716 end do
2717 !$omp do collapse(2)
2718 do k=ks, ke
2719 do j=js, je
2720 do i=is, ie
2721 dens_ = dens_hyd(i,j,k) + ddens0(i,j,k) + c_ss * ( ddens(i,j,k) - ddens0(i,j,k) )
2722 q(i,j,k) = q(i,j,k) / dens_
2723 end do
2724 end do
2725 end do
2726 end if
2727
2728 !$omp end parallel
2729 call prof_rapend( 'rk_advance_trcvar_generall3D', 3)
2730
2731 return
2732 end subroutine rk_advance_trcvar_general3d
2733
2734!> Store variable data after the implicit part of current stage in IMEX RK scheme (general case)
2735!!
2736!! @param nowstage Current stage
2737!! @param q Array to store variable data after the implicit part of current stage
2738!! @param varID Index of the targeting variable
2739!! @param is Index at which begins the loop for the corresponding direction
2740!! @param ie Index at which finishes the loop for the corresponding direction
2741!! @param js Index at which begins the loop for the corresponding direction
2742!! @param je Index at which finishes the loop for the corresponding direction
2743!! @param ks Index at which begins the loop for the corresponding direction
2744!! @param ke Index at which finishes the loop for the corresponding direction
2745!OCL SERIAL
2746 subroutine rk_storeimpl_general3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke , IA,JA,KA, var_num, &
2747 var0_3D, varTmp_3d, tend_buf3D_im )
2748 implicit none
2749 class(timeint_rk), intent(inout) :: this
2750 integer, intent(in) :: nowstage
2751 integer, intent(in) :: IA,JA,KA
2752 real(RP), intent(inout) :: q(IA,JA,KA)
2753 integer, intent(in) :: varID
2754 integer, intent(in) :: is, ie ,js, je ,ks, ke
2755 integer, intent(in) :: var_num
2756 real(RP), intent(inout) :: var0_3d(IA,JA,KA,var_num)
2757 real(RP), intent(inout) :: varTmp_3d(IA,JA,KA,var_num)
2758 real(RP), intent(inout) :: tend_buf3D_im(IA,JA,KA,var_num,this%tend_buf_size)
2759
2760 integer :: i,j,k
2761 integer :: tintbuf_ind
2762 real(RP) :: coef_a_im_dt
2763 !----------------------------------------
2764
2765 if (.not. this%imex_flag ) return
2766
2767 tintbuf_ind = this%tend_buf_indmap(nowstage)
2768 coef_a_im_dt = this%dt * this%coef_a_im(nowstage,nowstage)
2769
2770 !$omp parallel
2771
2772 if ( nowstage == 1 ) then
2773 !$omp do collapse(2)
2774 do k=ks, ke
2775 do j=js, je
2776 do i=is, ie
2777 var0_3d(i,j,k,varid) = q(i,j,k)
2778 vartmp_3d(i,j,k,varid) = q(i,j,k)
2779 end do
2780 end do
2781 end do
2782 end if
2783
2784 !$omp do collapse(2)
2785 do k=ks, ke
2786 do j=js, je
2787 do i=is, ie
2788 q(i,j,k) = q(i,j,k) &
2789 + coef_a_im_dt * tend_buf3d_im(i,j,k,varid,tintbuf_ind)
2790 end do
2791 end do
2792 end do
2793 !$omp end parallel
2794
2795 return
2796 end subroutine rk_storeimpl_general3d
2797
2798end module scale_timeint_rk
Module common / Runge-Kutta scheme.
subroutine, public timeint_rk_butcher_tab_get_info(rk_scheme_name, nstage, tend_buf_size, low_storage_flag, imex_flag)
subroutine, public timeint_rk_butcher_tab_get(rk_scheme_name, nstage, imex_flag, coef_a_ex, coef_b_ex, coef_c_ex, coef_sig_ex, coef_gam_ex, coef_a_im, coef_b_im, coef_c_im, tend_buf_indmap)
Module common / Runge-Kutta scheme.
subroutine timeint_rk_init(this, rk_scheme_name, dt, var_num, ndim, size_each_var)
Initialize a object to provide RK scheme.
Derived type to provide RK scheme.