13#include "scaleFElib.h"
34 real(RP),
private :: dt
36 integer,
private :: rk_scheme_id_ex
37 integer,
private :: rk_scheme_id_im
39 integer,
public :: var_num
40 integer,
private,
allocatable :: size_each_var(:)
41 integer,
public :: nstage
42 integer,
public :: var_buf_size
43 integer,
public :: tend_buf_size
46 real(rp),
public,
allocatable :: coef_a_ex(:,:)
47 real(rp),
public,
allocatable :: coef_b_ex(:)
48 real(rp),
public,
allocatable :: coef_c_ex(:)
51 real(rp),
public,
allocatable :: coef_sig_ex(:,:)
52 real(rp),
public,
allocatable :: coef_gam_ex(:,:)
55 real(rp),
public,
allocatable :: coef_a_im(:,:)
56 real(rp),
public,
allocatable :: coef_b_im(:)
57 real(rp),
public,
allocatable :: coef_c_im(:)
59 integer,
allocatable :: tend_buf_indmap(:)
60 real(rp),
public,
allocatable :: var_buf1d_ex(:,:,:)
61 real(rp),
public,
allocatable :: tend_buf1d_ex(:,:,:)
62 real(rp),
public,
allocatable :: tend_buf1d_im(:,:,:)
63 real(rp),
public,
allocatable :: var0_1d(:,:)
64 real(rp),
private,
allocatable :: vartmp_1d(:,:)
65 real(rp),
public,
allocatable :: var_buf2d_ex(:,:,:,:)
66 real(rp),
public,
allocatable :: tend_buf2d_ex(:,:,:,:)
67 real(rp),
public,
allocatable :: tend_buf2d_im(:,:,:,:)
68 real(rp),
public,
allocatable :: var0_2d(:,:,:)
69 real(rp),
private,
allocatable :: vartmp_2d(:,:,:)
70 real(rp),
public,
allocatable :: var_buf3d_ex(:,:,:,:,:)
71 real(rp),
public,
allocatable :: tend_buf3d_ex(:,:,:,:,:)
72 real(rp),
public,
allocatable :: tend_buf3d_im(:,:,:,:,:)
73 real(rp),
public,
allocatable :: var0_3d(:,:,:,:)
74 real(rp),
private,
allocatable :: vartmp_3d(:,:,:,:)
76 logical,
private :: low_storage_flag
77 logical,
public :: imex_flag
78 integer,
public :: ndim
81 procedure,
public :: final => timeint_rk_final
82 procedure,
public :: get_implicit_diagfac => timeint_rk_get_implicit_diagfac
83 procedure,
public :: get_deltime => timeint_rk_get_deltime
84 procedure,
public :: advance1d => timeint_rk_advance1d
85 procedure,
public :: advance_trcvar_1d => timeint_rk_advance_trcvar1d
86 procedure,
public :: storevar0_1d => timeint_rk_store_var0_1d
87 procedure,
public :: storeimplicit1d => timeint_rk_storeimpl1d
88 procedure,
public :: advance2d => timeint_rk_advance2d
89 procedure,
public :: advance_trcvar_2d => timeint_rk_advance_trcvar2d
90 procedure,
public :: storevar0_2d => timeint_rk_store_var0_2d
91 procedure,
public :: storeimplicit2d => timeint_rk_storeimpl2d
92 procedure,
public :: advance3d => timeint_rk_advance3d
93 procedure,
public :: advance_trcvar_3d => timeint_rk_advance_trcvar3d
94 procedure,
public :: storevar0_3d => timeint_rk_store_var0_3d
95 procedure,
public :: storeimplicit3d => timeint_rk_storeimpl3d
96 generic,
public :: advance => advance1d, advance2d, advance3d
97 generic,
public :: advance_trcvar => advance_trcvar_1d, advance_trcvar_2d, advance_trcvar_3d
98 generic,
public :: storevar0 => storevar0_1d, storevar0_2d, storevar0_3d
99 generic,
public :: storeimplicit => storeimplicit1d, storeimplicit2d, storeimplicit3d
117 rk_scheme_name, dt, var_num, ndim, size_each_var )
124 character(*),
intent(in) :: rk_scheme_name
125 real(RP),
intent(in) :: dt
126 integer,
intent(in) :: var_num
127 integer,
intent(in) :: ndim
128 integer,
intent(in) :: size_each_var(ndim)
134 allocate( this%size_each_var(ndim) )
135 this%size_each_var(:) = size_each_var(:)
138 this%nstage, this%tend_buf_size, &
139 this%low_storage_flag, this%imex_flag )
141 allocate ( this%coef_a_ex(this%nstage,this%nstage), this%coef_b_ex(this%nstage), this%coef_c_ex(this%nstage) )
142 allocate ( this%coef_sig_ex(this%nstage+1,this%nstage), this%coef_gam_ex(this%nstage+1,this%nstage) )
144 allocate ( this%coef_a_im(this%nstage,this%nstage), this%coef_b_im(this%nstage), this%coef_c_im(this%nstage) )
147 select case(this%ndim)
149 allocate( this%tend_buf1D_ex(size_each_var(1), var_num, this%tend_buf_size) )
150 if ( this%imex_flag )
then
151 allocate( this%tend_buf1D_im(size_each_var(1), var_num, this%tend_buf_size) )
153 allocate( this%var0_1D(size_each_var(1), var_num) )
154 allocate( this%varTmp_1D(size_each_var(1), var_num) )
156 allocate( this%tend_buf2D_ex(size_each_var(1),size_each_var(2), var_num, this%tend_buf_size) )
157 if ( this%imex_flag )
then
158 allocate( this%tend_buf2D_im(size_each_var(1),size_each_var(2), var_num, this%tend_buf_size) )
160 allocate( this%var0_2D(size_each_var(1),size_each_var(2), var_num) )
161 allocate( this%varTmp_2D(size_each_var(1),size_each_var(2), var_num) )
163 allocate( this%tend_buf3D_ex(size_each_var(1),size_each_var(2),size_each_var(3), var_num, this%tend_buf_size) )
164 if ( this%imex_flag )
then
165 allocate( this%tend_buf3D_im(size_each_var(1),size_each_var(2),size_each_var(3), var_num, this%tend_buf_size) )
167 allocate( this%var0_3D(size_each_var(1),size_each_var(2),size_each_var(3), var_num) )
168 allocate( this%varTmp_3D(size_each_var(1),size_each_var(2),size_each_var(3), var_num) )
170 allocate( this%tend_buf_indmap(this%nstage) )
173 rk_scheme_name, this%nstage, this%imex_flag, &
174 this%coef_a_ex, this%coef_b_ex, this%coef_c_ex, &
175 this%coef_sig_ex, this%coef_gam_ex, &
176 this%coef_a_im, this%coef_b_im, this%coef_c_im, &
177 this%tend_buf_indmap )
182 subroutine timeint_rk_final( this )
187 deallocate( this%coef_a_ex, this%coef_b_ex, this%coef_c_ex )
188 deallocate( this%coef_sig_ex, this%coef_gam_ex )
190 deallocate( this%coef_a_im, this%coef_b_im, this%coef_c_im )
193 deallocate( this%size_each_var, this%tend_buf_indmap )
195 select case(this%ndim)
197 deallocate( this%var0_1d )
198 deallocate( this%tend_buf1D_ex )
199 if (this%imex_flag)
deallocate( this%tend_buf1D_im )
200 deallocate( this%varTmp_1d )
202 deallocate( this%var0_2d )
203 deallocate( this%tend_buf2D_ex )
204 if (this%imex_flag)
deallocate( this%tend_buf2D_im )
205 deallocate( this%varTmp_2d )
207 deallocate( this%var0_3d )
208 deallocate( this%tend_buf3D_ex )
209 if (this%imex_flag)
deallocate( this%tend_buf3D_im )
210 deallocate( this%varTmp_3d )
214 end subroutine timeint_rk_final
216 elemental function timeint_rk_get_deltime( this )
result(dt)
225 end function timeint_rk_get_deltime
227 elemental function timeint_rk_get_implicit_diagfac( this, nowstage )
result(fac)
230 integer,
intent(in) :: nowstage
234 fac = this%coef_a_im(nowstage,nowstage) * this%dt
237 end function timeint_rk_get_implicit_diagfac
241 subroutine timeint_rk_advance1d( this, nowstage, q, varID, is, ie )
244 integer,
intent(in) :: nowstage
245 real(RP),
intent(inout) :: q(:)
246 integer,
intent(in) :: varID
247 integer,
intent(in) :: is, ie
250 if (this%low_storage_flag)
then
251 call rk_advance_low_storage1d( this, nowstage, q, varid, is, ie )
253 call rk_advance_general1d( this, nowstage, q, varid, is, ie )
257 end subroutine timeint_rk_advance1d
259 subroutine timeint_rk_advance_trcvar1d( this, nowstage, q, varID, is, ie , &
260 DENS, DENS0, DENS_hyd )
263 integer,
intent(in) :: nowstage
264 real(RP),
intent(inout) :: q(:)
265 integer,
intent(in) :: varID
266 integer,
intent(in) :: is, ie
267 real(RP),
intent(in) :: DENS (:)
268 real(RP),
intent(in) :: DENS0(:)
269 real(RP),
intent(in) :: DENS_hyd(:)
272 if ( this%imex_flag )
then
273 if ( prc_ismaster )
write(*,*)
"timeint_rk_advence_trcvar: IMEX is not supported. Check!"
277 if (this%low_storage_flag)
then
278 call rk_advance_trcvar_low_storage1d( this, nowstage, q, varid, is, ie , &
279 dens, dens0, dens_hyd )
281 call rk_advance_trcvar_general1d( this, nowstage, q, varid, is, ie , &
282 dens, dens0, dens_hyd )
286 end subroutine timeint_rk_advance_trcvar1d
288 subroutine timeint_rk_store_var0_1d( this, q, varID, is, ie )
291 real(RP),
intent(inout) :: q(:)
292 integer,
intent(in) :: varID
293 integer,
intent(in) :: is, ie
301 this%var0_1D(i,varid) = q(i)
306 end subroutine timeint_rk_store_var0_1d
308 subroutine timeint_rk_storeimpl1d( this, nowstage, q, varID, is, ie )
311 integer,
intent(in) :: nowstage
312 real(RP),
intent(inout) :: q(:)
313 integer,
intent(in) :: varID
314 integer,
intent(in) :: is, ie
318 call rk_storeimpl_general1d( this, nowstage, q, varid, is, ie )
321 end subroutine timeint_rk_storeimpl1d
323 subroutine timeint_rk_advance2d( this, nowstage, q, varID, is, ie ,js, je )
326 integer,
intent(in) :: nowstage
327 real(RP),
intent(inout) :: q(:,:)
328 integer,
intent(in) :: varID
329 integer,
intent(in) :: is, ie ,js, je
332 if (this%low_storage_flag)
then
333 call rk_advance_low_storage2d( this, nowstage, q, varid, is, ie ,js, je )
335 call rk_advance_general2d( this, nowstage, q, varid, is, ie ,js, je )
339 end subroutine timeint_rk_advance2d
341 subroutine timeint_rk_advance_trcvar2d( this, nowstage, q, varID, is, ie ,js, je , &
342 DENS, DENS0, DENS_hyd )
345 integer,
intent(in) :: nowstage
346 real(RP),
intent(inout) :: q(:,:)
347 integer,
intent(in) :: varID
348 integer,
intent(in) :: is, ie ,js, je
349 real(RP),
intent(in) :: DENS (:,:)
350 real(RP),
intent(in) :: DENS0(:,:)
351 real(RP),
intent(in) :: DENS_hyd(:,:)
354 if ( this%imex_flag )
then
355 if ( prc_ismaster )
write(*,*)
"timeint_rk_advence_trcvar: IMEX is not supported. Check!"
359 if (this%low_storage_flag)
then
360 call rk_advance_trcvar_low_storage2d( this, nowstage, q, varid, is, ie ,js, je , &
361 dens, dens0, dens_hyd )
363 call rk_advance_trcvar_general2d( this, nowstage, q, varid, is, ie ,js, je , &
364 dens, dens0, dens_hyd )
368 end subroutine timeint_rk_advance_trcvar2d
370 subroutine timeint_rk_store_var0_2d( this, q, varID, is, ie ,js, je )
373 real(RP),
intent(inout) :: q(:,:)
374 integer,
intent(in) :: varID
375 integer,
intent(in) :: is, ie ,js, je
384 this%var0_2D(i,j,varid) = q(i,j)
390 end subroutine timeint_rk_store_var0_2d
392 subroutine timeint_rk_storeimpl2d( this, nowstage, q, varID, is, ie ,js, je )
395 integer,
intent(in) :: nowstage
396 real(RP),
intent(inout) :: q(:,:)
397 integer,
intent(in) :: varID
398 integer,
intent(in) :: is, ie ,js, je
402 call rk_storeimpl_general2d( this, nowstage, q, varid, is, ie ,js, je )
405 end subroutine timeint_rk_storeimpl2d
407 subroutine timeint_rk_advance3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke )
410 integer,
intent(in) :: nowstage
411 real(RP),
intent(inout) :: q(:,:,:)
412 integer,
intent(in) :: varID
413 integer,
intent(in) :: is, ie ,js, je ,ks, ke
416 if (this%low_storage_flag)
then
417 call rk_advance_low_storage3d( this, nowstage, q, varid, is, ie ,js, je ,ks, ke )
419 call rk_advance_general3d( this, nowstage, q, varid, is, ie ,js, je ,ks, ke )
423 end subroutine timeint_rk_advance3d
425 subroutine timeint_rk_advance_trcvar3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke , &
426 DENS, DENS0, DENS_hyd )
429 integer,
intent(in) :: nowstage
430 real(RP),
intent(inout) :: q(:,:,:)
431 integer,
intent(in) :: varID
432 integer,
intent(in) :: is, ie ,js, je ,ks, ke
433 real(RP),
intent(in) :: DENS (:,:,:)
434 real(RP),
intent(in) :: DENS0(:,:,:)
435 real(RP),
intent(in) :: DENS_hyd(:,:,:)
438 if ( this%imex_flag )
then
439 if ( prc_ismaster )
write(*,*)
"timeint_rk_advence_trcvar: IMEX is not supported. Check!"
443 if (this%low_storage_flag)
then
444 call rk_advance_trcvar_low_storage3d( this, nowstage, q, varid, is, ie ,js, je ,ks, ke , &
445 dens, dens0, dens_hyd )
447 call rk_advance_trcvar_general3d( this, nowstage, q, varid, is, ie ,js, je ,ks, ke , &
448 dens, dens0, dens_hyd )
452 end subroutine timeint_rk_advance_trcvar3d
454 subroutine timeint_rk_store_var0_3d( this, q, varID, is, ie ,js, je ,ks, ke )
457 real(RP),
intent(inout) :: q(:,:,:)
458 integer,
intent(in) :: varID
459 integer,
intent(in) :: is, ie ,js, je ,ks, ke
469 this%var0_3D(i,j,k,varid) = q(i,j,k)
476 end subroutine timeint_rk_store_var0_3d
478 subroutine timeint_rk_storeimpl3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke )
481 integer,
intent(in) :: nowstage
482 real(RP),
intent(inout) :: q(:,:,:)
483 integer,
intent(in) :: varID
484 integer,
intent(in) :: is, ie ,js, je ,ks, ke
488 call rk_storeimpl_general3d( this, nowstage, q, varid, is, ie ,js, je ,ks, ke )
491 end subroutine timeint_rk_storeimpl3d
498 subroutine rk_advance_low_storage1d( this, nowstage, q, varID, is, ie )
499 use scale_const,
only: &
503 integer,
intent(in) :: nowstage
504 real(RP),
intent(inout) :: q(:)
505 integer,
intent(in) :: varID
506 integer,
intent(in) :: is, ie
511 real(RP) :: one_Minus_sig_ss
517 call prof_rapstart(
'rk_advance_low_storage1D', 3)
519 sig_ss = this%coef_sig_ex(nowstage+1,nowstage)
520 sig_ns = this%coef_sig_ex(this%nstage+1,nowstage)
521 one_minus_sig_ss = 1.0_rp - sig_ss
522 gam_ss = this%dt * this%coef_gam_ex(nowstage+1,nowstage)
523 gam_ns = this%dt * this%coef_gam_ex(this%nstage+1,nowstage)
525 if ( nowstage == this%nstage )
then
529 q(i) = this%varTmp_1d(i,varid) &
530 + sig_ss * q(i) + gam_ss * this%tend_buf1D_ex(i,varid,1)
533 call prof_rapend(
'rk_advance_low_storage1D', 3)
539 if (nowstage == 1)
then
542 this%var0_1D(i,varid) = q(i)
543 this%varTmp_1D(i,varid) = 0.0_rp
547 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps )
then
550 this%varTmp_1d(i,varid) = this%varTmp_1d(i,varid) &
551 + sig_ns * q(i) + gam_ns * this%tend_buf1D_ex(i,varid,1)
557 q(i) = one_minus_sig_ss * this%var0_1d(i,varid) &
558 + sig_ss * q(i) + gam_ss * this%tend_buf1D_ex(i,varid,1)
563 call prof_rapend(
'rk_advance_low_storage1D', 3)
566 end subroutine rk_advance_low_storage1d
569 subroutine rk_advance_trcvar_low_storage1d( this, nowstage, q, varID, is, ie , &
570 DDENS, DDENS0, DENS_hyd )
571 use scale_const,
only: &
575 integer,
intent(in) :: nowstage
576 real(RP),
intent(inout) :: q(:)
577 integer,
intent(in) :: varID
578 integer,
intent(in) :: is, ie
579 real(RP),
intent(in) :: DDENS (:)
580 real(RP),
intent(in) :: DDENS0(:)
581 real(RP),
intent(in) :: DENS_hyd(:)
586 real(RP) :: one_Minus_sig_ss
592 real(RP) :: dens_ssm1
597 call prof_rapstart(
'rk_advance_trcvar_low_storage1D', 3)
599 sig_ss = this%coef_sig_ex(nowstage+1,nowstage)
600 sig_ns = this%coef_sig_ex(this%nstage+1,nowstage)
601 one_minus_sig_ss = 1.0_rp - sig_ss
602 gam_ss = this%dt * this%coef_gam_ex(nowstage+1,nowstage)
603 gam_ns = this%dt * this%coef_gam_ex(this%nstage+1,nowstage)
604 c_ssm1 = this%coef_c_ex(nowstage)
606 if ( nowstage == this%nstage )
then
610 dens_ssm1 = dens_hyd(i) + ddens0(i) + c_ssm1 * ( ddens(i) - ddens0(i) )
611 q(i) = ( this%varTmp_1d(i,varid) &
612 + sig_ss * q(i) * dens_ssm1 + gam_ss * this%tend_buf1D_ex(i,varid,1) ) &
613 / ( dens_hyd(i) + ddens(i) )
616 call prof_rapend(
'rk_advance_trcvar_low_storage1D', 3)
621 c_ss = this%coef_c_ex(nowstage+1)
624 if (nowstage == 1)
then
627 this%var0_1D(i,varid) = q(i) * ( dens_hyd(i) + ddens0(i) )
628 this%varTmp_1D(i,varid) = 0.0_rp
632 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps )
then
635 dens_ssm1 = dens_hyd(i) + ddens0(i) + c_ssm1 * ( ddens(i) - ddens0(i) )
636 this%varTmp_1d(i,varid) = this%varTmp_1d(i,varid) &
637 + sig_ns * q(i) * dens_ssm1 + gam_ns * this%tend_buf1D_ex(i,varid,1)
643 dens_ssm1 = dens_hyd(i) + ddens0(i) + c_ssm1 * ( ddens(i) - ddens0(i) )
644 dens_ss = dens_hyd(i) + ddens0(i) + c_ss * ( ddens(i) - ddens0(i) )
646 q(i) = ( one_minus_sig_ss * this%var0_1d(i,varid) &
647 + sig_ss * q(i) * dens_ssm1 + gam_ss * this%tend_buf1D_ex(i,varid,1) ) &
653 call prof_rapend(
'rk_advance_trcvar_low_storage1D', 3)
656 end subroutine rk_advance_trcvar_low_storage1d
661 subroutine rk_advance_low_storage2d( this, nowstage, q, varID, is, ie ,js, je )
662 use scale_const,
only: &
666 integer,
intent(in) :: nowstage
667 real(RP),
intent(inout) :: q(:,:)
668 integer,
intent(in) :: varID
669 integer,
intent(in) :: is, ie ,js, je
674 real(RP) :: one_Minus_sig_ss
680 call prof_rapstart(
'rk_advance_low_storage2D', 3)
682 sig_ss = this%coef_sig_ex(nowstage+1,nowstage)
683 sig_ns = this%coef_sig_ex(this%nstage+1,nowstage)
684 one_minus_sig_ss = 1.0_rp - sig_ss
685 gam_ss = this%dt * this%coef_gam_ex(nowstage+1,nowstage)
686 gam_ns = this%dt * this%coef_gam_ex(this%nstage+1,nowstage)
688 if ( nowstage == this%nstage )
then
693 q(i,j) = this%varTmp_2d(i,j,varid) &
694 + sig_ss * q(i,j) + gam_ss * this%tend_buf2D_ex(i,j,varid,1)
698 call prof_rapend(
'rk_advance_low_storage2D', 3)
704 if (nowstage == 1)
then
708 this%var0_2D(i,j,varid) = q(i,j)
709 this%varTmp_2D(i,j,varid) = 0.0_rp
714 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps )
then
718 this%varTmp_2d(i,j,varid) = this%varTmp_2d(i,j,varid) &
719 + sig_ns * q(i,j) + gam_ns * this%tend_buf2D_ex(i,j,varid,1)
727 q(i,j) = one_minus_sig_ss * this%var0_2d(i,j,varid) &
728 + sig_ss * q(i,j) + gam_ss * this%tend_buf2D_ex(i,j,varid,1)
734 call prof_rapend(
'rk_advance_low_storage2D', 3)
737 end subroutine rk_advance_low_storage2d
740 subroutine rk_advance_trcvar_low_storage2d( this, nowstage, q, varID, is, ie ,js, je , &
741 DDENS, DDENS0, DENS_hyd )
742 use scale_const,
only: &
746 integer,
intent(in) :: nowstage
747 real(RP),
intent(inout) :: q(:,:)
748 integer,
intent(in) :: varID
749 integer,
intent(in) :: is, ie ,js, je
750 real(RP),
intent(in) :: DDENS (:,:)
751 real(RP),
intent(in) :: DDENS0(:,:)
752 real(RP),
intent(in) :: DENS_hyd(:,:)
757 real(RP) :: one_Minus_sig_ss
763 real(RP) :: dens_ssm1
768 call prof_rapstart(
'rk_advance_trcvar_low_storage2D', 3)
770 sig_ss = this%coef_sig_ex(nowstage+1,nowstage)
771 sig_ns = this%coef_sig_ex(this%nstage+1,nowstage)
772 one_minus_sig_ss = 1.0_rp - sig_ss
773 gam_ss = this%dt * this%coef_gam_ex(nowstage+1,nowstage)
774 gam_ns = this%dt * this%coef_gam_ex(this%nstage+1,nowstage)
775 c_ssm1 = this%coef_c_ex(nowstage)
777 if ( nowstage == this%nstage )
then
782 dens_ssm1 = dens_hyd(i,j) + ddens0(i,j) + c_ssm1 * ( ddens(i,j) - ddens0(i,j) )
783 q(i,j) = ( this%varTmp_2d(i,j,varid) &
784 + sig_ss * q(i,j) * dens_ssm1 + gam_ss * this%tend_buf2D_ex(i,j,varid,1) ) &
785 / ( dens_hyd(i,j) + ddens(i,j) )
789 call prof_rapend(
'rk_advance_trcvar_low_storage2D', 3)
794 c_ss = this%coef_c_ex(nowstage+1)
797 if (nowstage == 1)
then
801 this%var0_2D(i,j,varid) = q(i,j) * ( dens_hyd(i,j) + ddens0(i,j) )
802 this%varTmp_2D(i,j,varid) = 0.0_rp
807 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps )
then
811 dens_ssm1 = dens_hyd(i,j) + ddens0(i,j) + c_ssm1 * ( ddens(i,j) - ddens0(i,j) )
812 this%varTmp_2d(i,j,varid) = this%varTmp_2d(i,j,varid) &
813 + sig_ns * q(i,j) * dens_ssm1 + gam_ns * this%tend_buf2D_ex(i,j,varid,1)
821 dens_ssm1 = dens_hyd(i,j) + ddens0(i,j) + c_ssm1 * ( ddens(i,j) - ddens0(i,j) )
822 dens_ss = dens_hyd(i,j) + ddens0(i,j) + c_ss * ( ddens(i,j) - ddens0(i,j) )
824 q(i,j) = ( one_minus_sig_ss * this%var0_2d(i,j,varid) &
825 + sig_ss * q(i,j) * dens_ssm1 + gam_ss * this%tend_buf2D_ex(i,j,varid,1) ) &
832 call prof_rapend(
'rk_advance_trcvar_low_storage2D', 3)
835 end subroutine rk_advance_trcvar_low_storage2d
840 subroutine rk_advance_low_storage3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke )
841 use scale_const,
only: &
845 integer,
intent(in) :: nowstage
846 real(RP),
intent(inout) :: q(:,:,:)
847 integer,
intent(in) :: varID
848 integer,
intent(in) :: is, ie ,js, je ,ks, ke
853 real(RP) :: one_Minus_sig_ss
859 call prof_rapstart(
'rk_advance_low_storage3D', 3)
861 sig_ss = this%coef_sig_ex(nowstage+1,nowstage)
862 sig_ns = this%coef_sig_ex(this%nstage+1,nowstage)
863 one_minus_sig_ss = 1.0_rp - sig_ss
864 gam_ss = this%dt * this%coef_gam_ex(nowstage+1,nowstage)
865 gam_ns = this%dt * this%coef_gam_ex(this%nstage+1,nowstage)
867 if ( nowstage == this%nstage )
then
873 q(i,j,k) = this%varTmp_3d(i,j,k,varid) &
874 + sig_ss * q(i,j,k) + gam_ss * this%tend_buf3D_ex(i,j,k,varid,1)
879 call prof_rapend(
'rk_advance_low_storage3D', 3)
885 if (nowstage == 1)
then
890 this%var0_3D(i,j,k,varid) = q(i,j,k)
891 this%varTmp_3D(i,j,k,varid) = 0.0_rp
897 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps )
then
902 this%varTmp_3d(i,j,k,varid) = this%varTmp_3d(i,j,k,varid) &
903 + sig_ns * q(i,j,k) + gam_ns * this%tend_buf3D_ex(i,j,k,varid,1)
913 q(i,j,k) = one_minus_sig_ss * this%var0_3d(i,j,k,varid) &
914 + sig_ss * q(i,j,k) + gam_ss * this%tend_buf3D_ex(i,j,k,varid,1)
921 call prof_rapend(
'rk_advance_low_storage3D', 3)
924 end subroutine rk_advance_low_storage3d
927 subroutine rk_advance_trcvar_low_storage3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke , &
928 DDENS, DDENS0, DENS_hyd )
929 use scale_const,
only: &
933 integer,
intent(in) :: nowstage
934 real(RP),
intent(inout) :: q(:,:,:)
935 integer,
intent(in) :: varID
936 integer,
intent(in) :: is, ie ,js, je ,ks, ke
937 real(RP),
intent(in) :: DDENS (:,:,:)
938 real(RP),
intent(in) :: DDENS0(:,:,:)
939 real(RP),
intent(in) :: DENS_hyd(:,:,:)
944 real(RP) :: one_Minus_sig_ss
950 real(RP) :: dens_ssm1
955 call prof_rapstart(
'rk_advance_trcvar_low_storage3D', 3)
957 sig_ss = this%coef_sig_ex(nowstage+1,nowstage)
958 sig_ns = this%coef_sig_ex(this%nstage+1,nowstage)
959 one_minus_sig_ss = 1.0_rp - sig_ss
960 gam_ss = this%dt * this%coef_gam_ex(nowstage+1,nowstage)
961 gam_ns = this%dt * this%coef_gam_ex(this%nstage+1,nowstage)
962 c_ssm1 = this%coef_c_ex(nowstage)
964 if ( nowstage == this%nstage )
then
970 dens_ssm1 = dens_hyd(i,j,k) + ddens0(i,j,k) + c_ssm1 * ( ddens(i,j,k) - ddens0(i,j,k) )
971 q(i,j,k) = ( this%varTmp_3d(i,j,k,varid) &
972 + sig_ss * q(i,j,k) * dens_ssm1 + gam_ss * this%tend_buf3D_ex(i,j,k,varid,1) ) &
973 / ( dens_hyd(i,j,k) + ddens(i,j,k) )
978 call prof_rapend(
'rk_advance_trcvar_low_storage3D', 3)
983 c_ss = this%coef_c_ex(nowstage+1)
986 if (nowstage == 1)
then
991 this%var0_3D(i,j,k,varid) = q(i,j,k) * ( dens_hyd(i,j,k) + ddens0(i,j,k) )
992 this%varTmp_3D(i,j,k,varid) = 0.0_rp
998 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps )
then
1003 dens_ssm1 = dens_hyd(i,j,k) + ddens0(i,j,k) + c_ssm1 * ( ddens(i,j,k) - ddens0(i,j,k) )
1004 this%varTmp_3d(i,j,k,varid) = this%varTmp_3d(i,j,k,varid) &
1005 + sig_ns * q(i,j,k) * dens_ssm1 + gam_ns * this%tend_buf3D_ex(i,j,k,varid,1)
1015 dens_ssm1 = dens_hyd(i,j,k) + ddens0(i,j,k) + c_ssm1 * ( ddens(i,j,k) - ddens0(i,j,k) )
1016 dens_ss = dens_hyd(i,j,k) + ddens0(i,j,k) + c_ss * ( ddens(i,j,k) - ddens0(i,j,k) )
1018 q(i,j,k) = ( one_minus_sig_ss * this%var0_3d(i,j,k,varid) &
1019 + sig_ss * q(i,j,k) * dens_ssm1 + gam_ss * this%tend_buf3D_ex(i,j,k,varid,1) ) &
1027 call prof_rapend(
'rk_advance_trcvar_low_storage3D', 3)
1030 end subroutine rk_advance_trcvar_low_storage3d
1035 subroutine rk_advance_general1d( this, nowstage, q, varID, is, ie )
1038 integer,
intent(in) :: nowstage
1039 real(RP),
intent(inout) :: q(:)
1040 integer,
intent(in) :: varID
1041 integer,
intent(in) :: is, ie
1045 integer :: tintbuf_ind
1048 call prof_rapstart(
'rk_advance_general1D', 3)
1050 tintbuf_ind = this%tend_buf_indmap(nowstage)
1052 if ( this%nstage == 1 )
then
1055 this%varTmp_1d(i,varid) = q(i)
1059 if ( nowstage == this%nstage )
then
1060 if ( this%imex_flag )
then
1064 q(i) = this%varTmp_1d(i,varid) + this%dt * ( &
1065 + this%coef_b_ex(nowstage) * this%tend_buf1D_ex(i,varid,tintbuf_ind) &
1066 + this%coef_b_im(nowstage) * this%tend_buf1D_im(i,varid,tintbuf_ind) )
1072 q(i) = this%varTmp_1d(i,varid) &
1073 + this%dt * this%coef_b_ex(nowstage) * this%tend_buf1D_ex(i,varid,tintbuf_ind)
1076 call prof_rapend(
'rk_advance_general1D', 3)
1084 if ( nowstage == 1 .and. (.not. this%imex_flag) )
then
1087 this%var0_1D(i,varid) = q(i)
1088 this%varTmp_1D(i,varid) = q(i)
1092 if ( this%imex_flag )
then
1095 q(i) = this%var0_1d(i,varid)
1096 this%varTmp_1d(i,varid) = this%varTmp_1d(i,varid) + this%dt * ( &
1097 + this%coef_b_ex(nowstage) * this%tend_buf1D_ex(i,varid,tintbuf_ind) &
1098 + this%coef_b_im(nowstage) * this%tend_buf1D_im(i,varid,tintbuf_ind) )
1103 q(i) = this%var0_1d(i,varid)
1104 this%varTmp_1d(i,varid) = this%varTmp_1d(i,varid) &
1105 + this%dt * this%coef_b_ex(nowstage) * this%tend_buf1D_ex(i,varid,tintbuf_ind)
1109 if ( this%tend_buf_size == 1 .and. (.not. this%imex_flag) )
then
1112 q(i) = this%var0_1d(i,varid) &
1113 + this%dt * this%coef_a_ex(nowstage+1,nowstage)*this%tend_buf1D_ex(i,varid,1)
1115 else if ( .not. this%imex_flag )
then
1120 + this%dt * this%coef_a_ex(nowstage+1,s)*this%tend_buf1D_ex(i,varid,s)
1128 + this%dt * ( this%coef_a_ex(nowstage+1,s)*this%tend_buf1D_ex(i,varid,s) &
1129 + this%coef_a_im(nowstage+1,s)*this%tend_buf1D_im(i,varid,s) )
1135 call prof_rapend(
'rk_advance_general1D', 3)
1138 end subroutine rk_advance_general1d
1141 subroutine rk_advance_trcvar_general1d( this, nowstage, q, varID, is, ie , &
1142 DDENS, DDENS0, DENS_hyd )
1145 integer,
intent(in) :: nowstage
1146 real(RP),
intent(inout) :: q(:)
1147 integer,
intent(in) :: varID
1148 integer,
intent(in) :: is, ie
1149 real(RP),
intent(in) :: DDENS (:)
1150 real(RP),
intent(in) :: DDENS0(:)
1151 real(RP),
intent(in) :: DENS_hyd(:)
1155 integer :: tintbuf_ind
1161 call prof_rapstart(
'rk_advance_trcvar_general1D', 3)
1163 tintbuf_ind = this%tend_buf_indmap(nowstage)
1165 if ( this%nstage == 1 )
then
1168 this%varTmp_1d(i,varid) = q(i) * ( dens_hyd(i) + ddens0(i) )
1172 if ( nowstage == this%nstage )
then
1175 q(i) = ( this%varTmp_1d(i,varid) &
1176 + this%dt * this%coef_b_ex(nowstage) * this%tend_buf1D_ex(i,varid,tintbuf_ind) ) &
1177 / ( dens_hyd(i) + ddens(i) )
1179 call prof_rapend(
'rk_advance_trcvar_general1D', 3)
1184 c_ss = this%coef_c_ex(nowstage+1)
1189 if ( nowstage == 1 .and. (.not. this%imex_flag) )
then
1192 this%var0_1D(i,varid) = q(i) * ( dens_hyd(i) + ddens0(i) )
1193 this%varTmp_1D(i,varid) = this%var0_1D(i,varid)
1199 q(i) = this%var0_1d(i,varid)
1200 this%varTmp_1d(i,varid) = this%varTmp_1d(i,varid) &
1201 + this%dt * this%coef_b_ex(nowstage) * this%tend_buf1D_ex(i,varid,tintbuf_ind)
1204 if ( this%tend_buf_size == 1 .and. (.not. this%imex_flag) )
then
1207 dens_ = ( dens_hyd(i) + ddens0(i) ) &
1208 + this%coef_a_ex(nowstage+1,nowstage)*( ddens(i) - ddens0(i) )
1209 q(i) = ( this%var0_1d(i,varid) &
1210 + this%dt * this%coef_a_ex(nowstage+1,nowstage)*this%tend_buf1D_ex(i,varid,1) ) &
1213 else if ( .not. this%imex_flag )
then
1218 + this%dt * this%coef_a_ex(nowstage+1,s)*this%tend_buf1D_ex(i,varid,s)
1223 dens_ = dens_hyd(i) + ddens0(i) + c_ss * ( ddens(i) - ddens0(i) )
1229 call prof_rapend(
'rk_advance_trcvar_generall1D', 3)
1232 end subroutine rk_advance_trcvar_general1d
1235 subroutine rk_storeimpl_general1d( this, nowstage, q, varID, is, ie )
1238 integer,
intent(in) :: nowstage
1239 real(RP),
intent(inout) :: q(:)
1240 integer,
intent(in) :: varID
1241 integer,
intent(in) :: is, ie
1244 integer :: tintbuf_ind
1247 if (.not. this%imex_flag )
return
1249 tintbuf_ind = this%tend_buf_indmap(nowstage)
1253 if ( nowstage == 1 )
then
1256 this%var0_1D(i,varid) = q(i)
1257 this%varTmp_1D(i,varid) = q(i)
1264 + this%dt * this%coef_a_im(nowstage,nowstage)*this%tend_buf1D_im(i,varid,tintbuf_ind)
1269 end subroutine rk_storeimpl_general1d
1273 subroutine rk_advance_general2d( this, nowstage, q, varID, is, ie ,js, je )
1276 integer,
intent(in) :: nowstage
1277 real(RP),
intent(inout) :: q(:,:)
1278 integer,
intent(in) :: varID
1279 integer,
intent(in) :: is, ie ,js, je
1283 integer :: tintbuf_ind
1286 call prof_rapstart(
'rk_advance_general2D', 3)
1288 tintbuf_ind = this%tend_buf_indmap(nowstage)
1290 if ( this%nstage == 1 )
then
1294 this%varTmp_2d(i,j,varid) = q(i,j)
1299 if ( nowstage == this%nstage )
then
1300 if ( this%imex_flag )
then
1305 q(i,j) = this%varTmp_2d(i,j,varid) + this%dt * ( &
1306 + this%coef_b_ex(nowstage) * this%tend_buf2D_ex(i,j,varid,tintbuf_ind) &
1307 + this%coef_b_im(nowstage) * this%tend_buf2D_im(i,j,varid,tintbuf_ind) )
1315 q(i,j) = this%varTmp_2d(i,j,varid) &
1316 + this%dt * this%coef_b_ex(nowstage) * this%tend_buf2D_ex(i,j,varid,tintbuf_ind)
1320 call prof_rapend(
'rk_advance_general2D', 3)
1328 if ( nowstage == 1 .and. (.not. this%imex_flag) )
then
1332 this%var0_2D(i,j,varid) = q(i,j)
1333 this%varTmp_2D(i,j,varid) = q(i,j)
1338 if ( this%imex_flag )
then
1342 q(i,j) = this%var0_2d(i,j,varid)
1343 this%varTmp_2d(i,j,varid) = this%varTmp_2d(i,j,varid) + this%dt * ( &
1344 + this%coef_b_ex(nowstage) * this%tend_buf2D_ex(i,j,varid,tintbuf_ind) &
1345 + this%coef_b_im(nowstage) * this%tend_buf2D_im(i,j,varid,tintbuf_ind) )
1352 q(i,j) = this%var0_2d(i,j,varid)
1353 this%varTmp_2d(i,j,varid) = this%varTmp_2d(i,j,varid) &
1354 + this%dt * this%coef_b_ex(nowstage) * this%tend_buf2D_ex(i,j,varid,tintbuf_ind)
1359 if ( this%tend_buf_size == 1 .and. (.not. this%imex_flag) )
then
1363 q(i,j) = this%var0_2d(i,j,varid) &
1364 + this%dt * this%coef_a_ex(nowstage+1,nowstage)*this%tend_buf2D_ex(i,j,varid,1)
1367 else if ( .not. this%imex_flag )
then
1373 + this%dt * this%coef_a_ex(nowstage+1,s)*this%tend_buf2D_ex(i,j,varid,s)
1383 + this%dt * ( this%coef_a_ex(nowstage+1,s)*this%tend_buf2D_ex(i,j,varid,s) &
1384 + this%coef_a_im(nowstage+1,s)*this%tend_buf2D_im(i,j,varid,s) )
1391 call prof_rapend(
'rk_advance_general2D', 3)
1394 end subroutine rk_advance_general2d
1397 subroutine rk_advance_trcvar_general2d( this, nowstage, q, varID, is, ie ,js, je , &
1398 DDENS, DDENS0, DENS_hyd )
1401 integer,
intent(in) :: nowstage
1402 real(RP),
intent(inout) :: q(:,:)
1403 integer,
intent(in) :: varID
1404 integer,
intent(in) :: is, ie ,js, je
1405 real(RP),
intent(in) :: DDENS (:,:)
1406 real(RP),
intent(in) :: DDENS0(:,:)
1407 real(RP),
intent(in) :: DENS_hyd(:,:)
1411 integer :: tintbuf_ind
1417 call prof_rapstart(
'rk_advance_trcvar_general2D', 3)
1419 tintbuf_ind = this%tend_buf_indmap(nowstage)
1421 if ( this%nstage == 1 )
then
1425 this%varTmp_2d(i,j,varid) = q(i,j) * ( dens_hyd(i,j) + ddens0(i,j) )
1430 if ( nowstage == this%nstage )
then
1434 q(i,j) = ( this%varTmp_2d(i,j,varid) &
1435 + this%dt * this%coef_b_ex(nowstage) * this%tend_buf2D_ex(i,j,varid,tintbuf_ind) ) &
1436 / ( dens_hyd(i,j) + ddens(i,j) )
1439 call prof_rapend(
'rk_advance_trcvar_general2D', 3)
1444 c_ss = this%coef_c_ex(nowstage+1)
1449 if ( nowstage == 1 .and. (.not. this%imex_flag) )
then
1453 this%var0_2D(i,j,varid) = q(i,j) * ( dens_hyd(i,j) + ddens0(i,j) )
1454 this%varTmp_2D(i,j,varid) = this%var0_2D(i,j,varid)
1462 q(i,j) = this%var0_2d(i,j,varid)
1463 this%varTmp_2d(i,j,varid) = this%varTmp_2d(i,j,varid) &
1464 + this%dt * this%coef_b_ex(nowstage) * this%tend_buf2D_ex(i,j,varid,tintbuf_ind)
1468 if ( this%tend_buf_size == 1 .and. (.not. this%imex_flag) )
then
1472 dens_ = ( dens_hyd(i,j) + ddens0(i,j) ) &
1473 + this%coef_a_ex(nowstage+1,nowstage)*( ddens(i,j) - ddens0(i,j) )
1474 q(i,j) = ( this%var0_2d(i,j,varid) &
1475 + this%dt * this%coef_a_ex(nowstage+1,nowstage)*this%tend_buf2D_ex(i,j,varid,1) ) &
1479 else if ( .not. this%imex_flag )
then
1485 + this%dt * this%coef_a_ex(nowstage+1,s)*this%tend_buf2D_ex(i,j,varid,s)
1492 dens_ = dens_hyd(i,j) + ddens0(i,j) + c_ss * ( ddens(i,j) - ddens0(i,j) )
1493 q(i,j) = q(i,j) / dens_
1499 call prof_rapend(
'rk_advance_trcvar_generall2D', 3)
1502 end subroutine rk_advance_trcvar_general2d
1505 subroutine rk_storeimpl_general2d( this, nowstage, q, varID, is, ie ,js, je )
1508 integer,
intent(in) :: nowstage
1509 real(RP),
intent(inout) :: q(:,:)
1510 integer,
intent(in) :: varID
1511 integer,
intent(in) :: is, ie ,js, je
1514 integer :: tintbuf_ind
1517 if (.not. this%imex_flag )
return
1519 tintbuf_ind = this%tend_buf_indmap(nowstage)
1523 if ( nowstage == 1 )
then
1527 this%var0_2D(i,j,varid) = q(i,j)
1528 this%varTmp_2D(i,j,varid) = q(i,j)
1537 + this%dt * this%coef_a_im(nowstage,nowstage)*this%tend_buf2D_im(i,j,varid,tintbuf_ind)
1543 end subroutine rk_storeimpl_general2d
1547 subroutine rk_advance_general3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke )
1550 integer,
intent(in) :: nowstage
1551 real(RP),
intent(inout) :: q(:,:,:)
1552 integer,
intent(in) :: varID
1553 integer,
intent(in) :: is, ie ,js, je ,ks, ke
1557 integer :: tintbuf_ind
1560 call prof_rapstart(
'rk_advance_general3D', 3)
1562 tintbuf_ind = this%tend_buf_indmap(nowstage)
1564 if ( this%nstage == 1 )
then
1569 this%varTmp_3d(i,j,k,varid) = q(i,j,k)
1575 if ( nowstage == this%nstage )
then
1576 if ( this%imex_flag )
then
1582 q(i,j,k) = this%varTmp_3d(i,j,k,varid) + this%dt * ( &
1583 + this%coef_b_ex(nowstage) * this%tend_buf3D_ex(i,j,k,varid,tintbuf_ind) &
1584 + this%coef_b_im(nowstage) * this%tend_buf3D_im(i,j,k,varid,tintbuf_ind) )
1594 q(i,j,k) = this%varTmp_3d(i,j,k,varid) &
1595 + this%dt * this%coef_b_ex(nowstage) * this%tend_buf3D_ex(i,j,k,varid,tintbuf_ind)
1600 call prof_rapend(
'rk_advance_general3D', 3)
1608 if ( nowstage == 1 .and. (.not. this%imex_flag) )
then
1613 this%var0_3D(i,j,k,varid) = q(i,j,k)
1614 this%varTmp_3D(i,j,k,varid) = q(i,j,k)
1620 if ( this%imex_flag )
then
1625 q(i,j,k) = this%var0_3d(i,j,k,varid)
1626 this%varTmp_3d(i,j,k,varid) = this%varTmp_3d(i,j,k,varid) + this%dt * ( &
1627 + this%coef_b_ex(nowstage) * this%tend_buf3D_ex(i,j,k,varid,tintbuf_ind) &
1628 + this%coef_b_im(nowstage) * this%tend_buf3D_im(i,j,k,varid,tintbuf_ind) )
1637 q(i,j,k) = this%var0_3d(i,j,k,varid)
1638 this%varTmp_3d(i,j,k,varid) = this%varTmp_3d(i,j,k,varid) &
1639 + this%dt * this%coef_b_ex(nowstage) * this%tend_buf3D_ex(i,j,k,varid,tintbuf_ind)
1645 if ( this%tend_buf_size == 1 .and. (.not. this%imex_flag) )
then
1650 q(i,j,k) = this%var0_3d(i,j,k,varid) &
1651 + this%dt * this%coef_a_ex(nowstage+1,nowstage)*this%tend_buf3D_ex(i,j,k,varid,1)
1655 else if ( .not. this%imex_flag )
then
1661 q(i,j,k) = q(i,j,k) &
1662 + this%dt * this%coef_a_ex(nowstage+1,s)*this%tend_buf3D_ex(i,j,k,varid,s)
1673 q(i,j,k) = q(i,j,k) &
1674 + this%dt * ( this%coef_a_ex(nowstage+1,s)*this%tend_buf3D_ex(i,j,k,varid,s) &
1675 + this%coef_a_im(nowstage+1,s)*this%tend_buf3D_im(i,j,k,varid,s) )
1683 call prof_rapend(
'rk_advance_general3D', 3)
1686 end subroutine rk_advance_general3d
1689 subroutine rk_advance_trcvar_general3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke , &
1690 DDENS, DDENS0, DENS_hyd )
1693 integer,
intent(in) :: nowstage
1694 real(RP),
intent(inout) :: q(:,:,:)
1695 integer,
intent(in) :: varID
1696 integer,
intent(in) :: is, ie ,js, je ,ks, ke
1697 real(RP),
intent(in) :: DDENS (:,:,:)
1698 real(RP),
intent(in) :: DDENS0(:,:,:)
1699 real(RP),
intent(in) :: DENS_hyd(:,:,:)
1703 integer :: tintbuf_ind
1709 call prof_rapstart(
'rk_advance_trcvar_general3D', 3)
1711 tintbuf_ind = this%tend_buf_indmap(nowstage)
1713 if ( this%nstage == 1 )
then
1718 this%varTmp_3d(i,j,k,varid) = q(i,j,k) * ( dens_hyd(i,j,k) + ddens0(i,j,k) )
1724 if ( nowstage == this%nstage )
then
1729 q(i,j,k) = ( this%varTmp_3d(i,j,k,varid) &
1730 + this%dt * this%coef_b_ex(nowstage) * this%tend_buf3D_ex(i,j,k,varid,tintbuf_ind) ) &
1731 / ( dens_hyd(i,j,k) + ddens(i,j,k) )
1735 call prof_rapend(
'rk_advance_trcvar_general3D', 3)
1740 c_ss = this%coef_c_ex(nowstage+1)
1745 if ( nowstage == 1 .and. (.not. this%imex_flag) )
then
1750 this%var0_3D(i,j,k,varid) = q(i,j,k) * ( dens_hyd(i,j,k) + ddens0(i,j,k) )
1751 this%varTmp_3D(i,j,k,varid) = this%var0_3D(i,j,k,varid)
1761 q(i,j,k) = this%var0_3d(i,j,k,varid)
1762 this%varTmp_3d(i,j,k,varid) = this%varTmp_3d(i,j,k,varid) &
1763 + this%dt * this%coef_b_ex(nowstage) * this%tend_buf3D_ex(i,j,k,varid,tintbuf_ind)
1768 if ( this%tend_buf_size == 1 .and. (.not. this%imex_flag) )
then
1773 dens_ = ( dens_hyd(i,j,k) + ddens0(i,j,k) ) &
1774 + this%coef_a_ex(nowstage+1,nowstage)*( ddens(i,j,k) - ddens0(i,j,k) )
1775 q(i,j,k) = ( this%var0_3d(i,j,k,varid) &
1776 + this%dt * this%coef_a_ex(nowstage+1,nowstage)*this%tend_buf3D_ex(i,j,k,varid,1) ) &
1781 else if ( .not. this%imex_flag )
then
1787 q(i,j,k) = q(i,j,k) &
1788 + this%dt * this%coef_a_ex(nowstage+1,s)*this%tend_buf3D_ex(i,j,k,varid,s)
1797 dens_ = dens_hyd(i,j,k) + ddens0(i,j,k) + c_ss * ( ddens(i,j,k) - ddens0(i,j,k) )
1798 q(i,j,k) = q(i,j,k) / dens_
1805 call prof_rapend(
'rk_advance_trcvar_generall3D', 3)
1808 end subroutine rk_advance_trcvar_general3d
1811 subroutine rk_storeimpl_general3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke )
1814 integer,
intent(in) :: nowstage
1815 real(RP),
intent(inout) :: q(:,:,:)
1816 integer,
intent(in) :: varID
1817 integer,
intent(in) :: is, ie ,js, je ,ks, ke
1820 integer :: tintbuf_ind
1823 if (.not. this%imex_flag )
return
1825 tintbuf_ind = this%tend_buf_indmap(nowstage)
1829 if ( nowstage == 1 )
then
1834 this%var0_3D(i,j,k,varid) = q(i,j,k)
1835 this%varTmp_3D(i,j,k,varid) = q(i,j,k)
1845 q(i,j,k) = q(i,j,k) &
1846 + this%dt * this%coef_a_im(nowstage,nowstage)*this%tend_buf3D_im(i,j,k,varid,tintbuf_ind)
1853 end subroutine rk_storeimpl_general3d
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)