13#include "scaleFElib.h"
36 real(RP),
private :: dt
38 integer,
private :: rk_scheme_id_ex
39 integer,
private :: rk_scheme_id_im
41 integer,
private,
allocatable :: size_each_var(:)
42 integer,
public :: nstage
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(:)
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(:,:,:,:)
77 logical,
private :: low_storage_flag
78 logical,
public :: imex_flag
79 integer,
public :: ndim
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
108 real(rp),
pointer :: var1d(:)
109 real(rp),
pointer :: var2d(:,:)
110 real(rp),
pointer :: var3d(:,:,:)
135 rk_scheme_name, dt, var_num, ndim, size_each_var )
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)
151 this%var_num = var_num
152 allocate( this%size_each_var(ndim) )
153 this%size_each_var(:) = size_each_var(:)
156 this%nstage, this%tend_buf_size, &
157 this%low_storage_flag, this%imex_flag )
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) )
162 allocate ( this%coef_a_im(this%nstage,this%nstage), this%coef_b_im(this%nstage), this%coef_c_im(this%nstage) )
165 select case(this%ndim)
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) )
171 allocate( this%var0_1D(size_each_var(1), var_num) )
172 allocate( this%varTmp_1D(size_each_var(1), var_num) )
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) )
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) )
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) )
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) )
188 allocate( this%tend_buf_indmap(this%nstage) )
191 rk_scheme_name, this%nstage, this%imex_flag, &
192 this%coef_a_ex, this%coef_b_ex, this%coef_c_ex, &
193 this%coef_sig_ex, this%coef_gam_ex, &
194 this%coef_a_im, this%coef_b_im, this%coef_c_im, &
195 this%tend_buf_indmap )
202 subroutine timeint_rk_final( this )
207 deallocate( this%coef_a_ex, this%coef_b_ex, this%coef_c_ex )
208 deallocate( this%coef_sig_ex, this%coef_gam_ex )
210 deallocate( this%coef_a_im, this%coef_b_im, this%coef_c_im )
213 deallocate( this%size_each_var, this%tend_buf_indmap )
215 select case(this%ndim)
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 )
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 )
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 )
234 end subroutine timeint_rk_final
238 elemental function timeint_rk_get_deltime( this )
result(dt)
247 end function timeint_rk_get_deltime
252 elemental function timeint_rk_get_implicit_diagfac( this, nowstage )
result(fac)
255 integer,
intent(in) :: nowstage
259 fac = this%coef_a_im(nowstage,nowstage) * this%dt
262 end function timeint_rk_get_implicit_diagfac
273 subroutine timeint_rk_advance1d( this, nowstage, q, varID, is, ie )
276 integer,
intent(in) :: nowstage
277 real(RP),
intent(inout) :: q(:)
278 integer,
intent(in) :: varID
279 integer,
intent(in) :: is, ie
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 )
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 )
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 )
295 end subroutine timeint_rk_advance1d
305 subroutine timeint_rk_advance1d_varlist( this, nowstage, var_list, varIDs, is, ie )
308 integer,
intent(in) :: nowstage
309 integer,
intent(in) :: varIDs(:)
311 integer,
intent(in) :: is, ie
319 if (this%low_storage_flag)
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 )
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 )
334 call this%Advance1D( nowstage, var_list(i)%var1D, varids(i), is, ie )
339 end subroutine timeint_rk_advance1d_varlist
351 subroutine timeint_rk_advance_trcvar1d( this, nowstage, q, varID, is, ie , &
352 DDENS, DDENS0, DENS_hyd )
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(:)
364 if ( this%imex_flag )
then
365 if ( prc_ismaster )
write(*,*)
"timeint_rk_advence_trcvar: IMEX is not supported. Check!"
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 )
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 )
378 end subroutine timeint_rk_advance_trcvar1d
386 subroutine timeint_rk_store_var0_1d( this, q, varID, is, ie )
389 real(RP),
intent(inout) :: q(:)
390 integer,
intent(in) :: varID
391 integer,
intent(in) :: is, ie
399 this%var0_1D(i,varid) = q(i)
404 end subroutine timeint_rk_store_var0_1d
413 subroutine timeint_rk_storeimpl1d( this, nowstage, q, varID, is, ie )
416 integer,
intent(in) :: nowstage
417 real(RP),
intent(inout) :: q(:)
418 integer,
intent(in) :: varID
419 integer,
intent(in) :: is, ie
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 )
427 end subroutine timeint_rk_storeimpl1d
438 subroutine timeint_rk_advance2d( this, nowstage, q, varID, is, ie ,js, je )
441 integer,
intent(in) :: nowstage
442 real(RP),
intent(inout) :: q(:,:)
443 integer,
intent(in) :: varID
444 integer,
intent(in) :: is, ie ,js, je
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 )
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 )
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 )
460 end subroutine timeint_rk_advance2d
472 subroutine timeint_rk_advance2d_varlist( this, nowstage, var_list, varIDs, is, ie ,js, je )
475 integer,
intent(in) :: nowstage
476 integer,
intent(in) :: varIDs(:)
478 integer,
intent(in) :: is, ie ,js, je
486 if (this%low_storage_flag)
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 )
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 )
501 call this%Advance2D( nowstage, var_list(i)%var2D, varids(i), is, ie ,js, je )
506 end subroutine timeint_rk_advance2d_varlist
520 subroutine timeint_rk_advance_trcvar2d( this, nowstage, q, varID, is, ie ,js, je , &
521 DDENS, DDENS0, DENS_hyd )
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(:,:)
533 if ( this%imex_flag )
then
534 if ( prc_ismaster )
write(*,*)
"timeint_rk_advence_trcvar: IMEX is not supported. Check!"
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 )
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 )
547 end subroutine timeint_rk_advance_trcvar2d
557 subroutine timeint_rk_store_var0_2d( this, q, varID, is, ie ,js, je )
560 real(RP),
intent(inout) :: q(:,:)
561 integer,
intent(in) :: varID
562 integer,
intent(in) :: is, ie ,js, je
571 this%var0_2D(i,j,varid) = q(i,j)
577 end subroutine timeint_rk_store_var0_2d
588 subroutine timeint_rk_storeimpl2d( this, nowstage, q, varID, is, ie ,js, je )
591 integer,
intent(in) :: nowstage
592 real(RP),
intent(inout) :: q(:,:)
593 integer,
intent(in) :: varID
594 integer,
intent(in) :: is, ie ,js, je
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 )
602 end subroutine timeint_rk_storeimpl2d
615 subroutine timeint_rk_advance3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke )
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
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 )
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 )
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 )
637 end subroutine timeint_rk_advance3d
651 subroutine timeint_rk_advance3d_varlist( this, nowstage, var_list, varIDs, is, ie ,js, je ,ks, ke )
654 integer,
intent(in) :: nowstage
655 integer,
intent(in) :: varIDs(:)
657 integer,
intent(in) :: is, ie ,js, je ,ks, ke
665 if (this%low_storage_flag)
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 )
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 )
680 call this%Advance3D( nowstage, var_list(i)%var3D, varids(i), is, ie ,js, je ,ks, ke )
685 end subroutine timeint_rk_advance3d_varlist
701 subroutine timeint_rk_advance_trcvar3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke , &
702 DDENS, DDENS0, DENS_hyd )
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(:,:,:)
714 if ( this%imex_flag )
then
715 if ( prc_ismaster )
write(*,*)
"timeint_rk_advence_trcvar: IMEX is not supported. Check!"
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 )
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 )
728 end subroutine timeint_rk_advance_trcvar3d
740 subroutine timeint_rk_store_var0_3d( this, q, varID, is, ie ,js, je ,ks, ke )
743 real(RP),
intent(inout) :: q(:,:,:)
744 integer,
intent(in) :: varID
745 integer,
intent(in) :: is, ie ,js, je ,ks, ke
755 this%var0_3D(i,j,k,varid) = q(i,j,k)
762 end subroutine timeint_rk_store_var0_3d
775 subroutine timeint_rk_storeimpl3d( this, nowstage, q, varID, is, ie ,js, je ,ks, ke )
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
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 )
789 end subroutine timeint_rk_storeimpl3d
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: &
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)
821 real(RP) :: one_Minus_sig_ss
827 call prof_rapstart(
'rk_advance_low_storage1D', 3)
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)
835 if ( nowstage == this%nstage )
then
839 q(i) = vartmp_1d(i,varid) &
840 + sig_ss * q(i) + gam_ss * tend_buf1d_ex(i,varid,1)
843 call prof_rapend(
'rk_advance_low_storage1D', 3)
849 if (nowstage == 1)
then
852 var0_1d(i,varid) = q(i)
853 vartmp_1d(i,varid) = 0.0_rp
857 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps )
then
860 vartmp_1d(i,varid) = vartmp_1d(i,varid) &
861 + sig_ns * q(i) + gam_ns * tend_buf1d_ex(i,varid,1)
867 q(i) = one_minus_sig_ss * var0_1d(i,varid) &
868 + sig_ss * q(i) + gam_ss * tend_buf1d_ex(i,varid,1)
873 call prof_rapend(
'rk_advance_low_storage1D', 3)
876 end subroutine rk_advance_low_storage1d
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: &
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)
909 real(RP) :: one_Minus_sig_ss
915 call prof_rapstart(
'rk_advance_low_storage1D', 3)
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)
923 if ( nowstage == this%nstage )
then
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)
933 call prof_rapend(
'rk_advance_low_storage1D', 3)
939 if (nowstage == 1)
then
942 var0_1d(i,varid1) = q1(i)
943 var0_1d(i,varid2) = q2(i)
945 vartmp_1d(i,varid1) = 0.0_rp
946 vartmp_1d(i,varid2) = 0.0_rp
950 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps )
then
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)
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)
970 call prof_rapend(
'rk_advance_low_storage1D', 3)
973 end subroutine rk_advance_low_storage1d_var2
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: &
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)
1009 real(RP) :: one_Minus_sig_ss
1015 real(RP) :: dens_ssm1
1020 call prof_rapstart(
'rk_advance_trcvar_low_storage1D', 3)
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)
1029 if ( nowstage == this%nstage )
then
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) )
1039 call prof_rapend(
'rk_advance_trcvar_low_storage1D', 3)
1044 c_ss = this%coef_c_ex(nowstage+1)
1047 if (nowstage == 1)
then
1050 var0_1d(i,varid) = q(i) * ( dens_hyd(i) + ddens0(i) )
1051 vartmp_1d(i,varid) = 0.0_rp
1055 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps )
then
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)
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) )
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) ) &
1076 call prof_rapend(
'rk_advance_trcvar_low_storage1D', 3)
1079 end subroutine rk_advance_trcvar_low_storage1d
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: &
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)
1111 real(RP) :: one_Minus_sig_ss
1117 call prof_rapstart(
'rk_advance_low_storage2D', 3)
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)
1125 if ( nowstage == this%nstage )
then
1130 q(i,j) = vartmp_2d(i,j,varid) &
1131 + sig_ss * q(i,j) + gam_ss * tend_buf2d_ex(i,j,varid,1)
1135 call prof_rapend(
'rk_advance_low_storage2D', 3)
1141 if (nowstage == 1)
then
1145 var0_2d(i,j,varid) = q(i,j)
1146 vartmp_2d(i,j,varid) = 0.0_rp
1151 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps )
then
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)
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)
1171 call prof_rapend(
'rk_advance_low_storage2D', 3)
1174 end subroutine rk_advance_low_storage2d
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: &
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)
1209 real(RP) :: one_Minus_sig_ss
1215 call prof_rapstart(
'rk_advance_low_storage2D', 3)
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)
1223 if ( nowstage == this%nstage )
then
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)
1235 call prof_rapend(
'rk_advance_low_storage2D', 3)
1241 if (nowstage == 1)
then
1245 var0_2d(i,j,varid1) = q1(i,j)
1246 var0_2d(i,j,varid2) = q2(i,j)
1248 vartmp_2d(i,j,varid1) = 0.0_rp
1249 vartmp_2d(i,j,varid2) = 0.0_rp
1254 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps )
then
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)
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)
1278 call prof_rapend(
'rk_advance_low_storage2D', 3)
1281 end subroutine rk_advance_low_storage2d_var2
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: &
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)
1319 real(RP) :: one_Minus_sig_ss
1325 real(RP) :: dens_ssm1
1330 call prof_rapstart(
'rk_advance_trcvar_low_storage2D', 3)
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)
1339 if ( nowstage == this%nstage )
then
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) )
1351 call prof_rapend(
'rk_advance_trcvar_low_storage2D', 3)
1356 c_ss = this%coef_c_ex(nowstage+1)
1359 if (nowstage == 1)
then
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
1369 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps )
then
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)
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) )
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) ) &
1394 call prof_rapend(
'rk_advance_trcvar_low_storage2D', 3)
1397 end subroutine rk_advance_trcvar_low_storage2d
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: &
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)
1431 real(RP) :: one_Minus_sig_ss
1437 call prof_rapstart(
'rk_advance_low_storage3D', 3)
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)
1445 if ( nowstage == this%nstage )
then
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)
1457 call prof_rapend(
'rk_advance_low_storage3D', 3)
1463 if (nowstage == 1)
then
1468 var0_3d(i,j,k,varid) = q(i,j,k)
1469 vartmp_3d(i,j,k,varid) = 0.0_rp
1475 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps )
then
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)
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)
1499 call prof_rapend(
'rk_advance_low_storage3D', 3)
1502 end subroutine rk_advance_low_storage3d
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: &
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)
1539 real(RP) :: one_Minus_sig_ss
1545 call prof_rapstart(
'rk_advance_low_storage3D', 3)
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)
1553 if ( nowstage == this%nstage )
then
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)
1567 call prof_rapend(
'rk_advance_low_storage3D', 3)
1573 if (nowstage == 1)
then
1578 var0_3d(i,j,k,varid1) = q1(i,j,k)
1579 var0_3d(i,j,k,varid2) = q2(i,j,k)
1581 vartmp_3d(i,j,k,varid1) = 0.0_rp
1582 vartmp_3d(i,j,k,varid2) = 0.0_rp
1588 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps )
then
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)
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)
1616 call prof_rapend(
'rk_advance_low_storage3D', 3)
1619 end subroutine rk_advance_low_storage3d_var2
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: &
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)
1659 real(RP) :: one_Minus_sig_ss
1665 real(RP) :: dens_ssm1
1670 call prof_rapstart(
'rk_advance_trcvar_low_storage3D', 3)
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)
1679 if ( nowstage == this%nstage )
then
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) )
1693 call prof_rapend(
'rk_advance_trcvar_low_storage3D', 3)
1698 c_ss = this%coef_c_ex(nowstage+1)
1701 if (nowstage == 1)
then
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
1713 if ( abs(sig_ns) > eps .or. abs(gam_ns) > eps )
then
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)
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) )
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) ) &
1742 call prof_rapend(
'rk_advance_trcvar_low_storage3D', 3)
1745 end subroutine rk_advance_trcvar_low_storage3d
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 )
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)
1774 integer :: tintbuf_ind
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
1782 call prof_rapstart(
'rk_advance_general1D', 3)
1784 tintbuf_ind = this%tend_buf_indmap(nowstage)
1786 if ( this%nstage == 1 )
then
1789 vartmp_1d(i,varid) = q(i)
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
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)
1804 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
1808 q(i) = vartmp_1d(i,varid) &
1809 + coef_b_ex_dt * tend_buf1d_ex(i,varid,tintbuf_ind)
1812 call prof_rapend(
'rk_advance_general1D', 3)
1820 if ( nowstage == 1 .and. (.not. this%imex_flag) )
then
1823 var0_1d(i,varid) = q(i)
1824 vartmp_1d(i,varid) = q(i)
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
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)
1840 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
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)
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)
1854 q(i) = var0_1d(i,varid) &
1855 + coef_a_ex_dt * tend_buf1d_ex(i,varid,1)
1857 else if ( .not. this%imex_flag )
then
1859 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,s)
1863 + coef_a_ex_dt * tend_buf1d_ex(i,varid,s)
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)
1873 + coef_a_ex_dt * tend_buf1d_ex(i,varid,s) &
1874 + coef_a_im_dt * tend_buf1d_im(i,varid,s)
1880 call prof_rapend(
'rk_advance_general1D', 3)
1883 end subroutine rk_advance_general1d
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 )
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)
1916 integer :: tintbuf_ind
1920 real(RP) :: coef_a_ex
1921 real(RP) :: coef_a_ex_dt
1922 real(RP) :: coef_b_ex_dt
1926 call prof_rapstart(
'rk_advance_trcvar_general1D', 3)
1928 tintbuf_ind = this%tend_buf_indmap(nowstage)
1930 if ( this%nstage == 1 )
then
1933 vartmp_1d(i,varid) = q(i) * ( dens_hyd(i) + ddens0(i) )
1937 if ( nowstage == this%nstage )
then
1938 coef_b_ex_dt = this%dt * this%coef_b_ex(nowstage)
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) )
1945 call prof_rapend(
'rk_advance_trcvar_general1D', 3)
1950 c_ss = this%coef_c_ex(nowstage+1)
1955 if ( nowstage == 1 .and. (.not. this%imex_flag) )
then
1958 var0_1d(i,varid) = q(i) * ( dens_hyd(i) + ddens0(i) )
1959 vartmp_1d(i,varid) = var0_1d(i,varid)
1963 coef_b_ex_dt = this%dt * this%coef_b_ex(nowstage)
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)
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)
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) ) &
1982 else if ( .not. this%imex_flag )
then
1984 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,s)
1988 + coef_a_ex_dt * tend_buf1d_ex(i,varid,s)
1993 dens_ = dens_hyd(i) + ddens0(i) + c_ss * ( ddens(i) - ddens0(i) )
1999 call prof_rapend(
'rk_advance_trcvar_generall1D', 3)
2002 end subroutine rk_advance_trcvar_general1d
2012 subroutine rk_storeimpl_general1d( this, nowstage, q, varID, is, ie , IA, var_num, &
2013 var0_1D, varTmp_1d, tend_buf1D_im )
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)
2027 integer :: tintbuf_ind
2028 real(RP) :: coef_a_im_dt
2031 if (.not. this%imex_flag )
return
2033 tintbuf_ind = this%tend_buf_indmap(nowstage)
2034 coef_a_im_dt = this%dt * this%coef_a_im(nowstage,nowstage)
2038 if ( nowstage == 1 )
then
2041 var0_1d(i,varid) = q(i)
2042 vartmp_1d(i,varid) = q(i)
2049 + coef_a_im_dt * tend_buf1d_im(i,varid,tintbuf_ind)
2054 end subroutine rk_storeimpl_general1d
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 )
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)
2084 integer :: tintbuf_ind
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
2092 call prof_rapstart(
'rk_advance_general2D', 3)
2094 tintbuf_ind = this%tend_buf_indmap(nowstage)
2096 if ( this%nstage == 1 )
then
2100 vartmp_2d(i,j,varid) = q(i,j)
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
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)
2118 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
2123 q(i,j) = vartmp_2d(i,j,varid) &
2124 + coef_b_ex_dt * tend_buf2d_ex(i,j,varid,tintbuf_ind)
2128 call prof_rapend(
'rk_advance_general2D', 3)
2136 if ( nowstage == 1 .and. (.not. this%imex_flag) )
then
2140 var0_2d(i,j,varid) = q(i,j)
2141 vartmp_2d(i,j,varid) = q(i,j)
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
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)
2160 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
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)
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)
2177 q(i,j) = var0_2d(i,j,varid) &
2178 + coef_a_ex_dt * tend_buf2d_ex(i,j,varid,1)
2181 else if ( .not. this%imex_flag )
then
2183 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,s)
2188 + coef_a_ex_dt * tend_buf2d_ex(i,j,varid,s)
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)
2200 + coef_a_ex_dt * tend_buf2d_ex(i,j,varid,s) &
2201 + coef_a_im_dt * tend_buf2d_im(i,j,varid,s)
2208 call prof_rapend(
'rk_advance_general2D', 3)
2211 end subroutine rk_advance_general2d
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 )
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)
2246 integer :: tintbuf_ind
2250 real(RP) :: coef_a_ex
2251 real(RP) :: coef_a_ex_dt
2252 real(RP) :: coef_b_ex_dt
2256 call prof_rapstart(
'rk_advance_trcvar_general2D', 3)
2258 tintbuf_ind = this%tend_buf_indmap(nowstage)
2260 if ( this%nstage == 1 )
then
2264 vartmp_2d(i,j,varid) = q(i,j) * ( dens_hyd(i,j) + ddens0(i,j) )
2269 if ( nowstage == this%nstage )
then
2270 coef_b_ex_dt = this%dt * this%coef_b_ex(nowstage)
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) )
2279 call prof_rapend(
'rk_advance_trcvar_general2D', 3)
2284 c_ss = this%coef_c_ex(nowstage+1)
2289 if ( nowstage == 1 .and. (.not. this%imex_flag) )
then
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)
2299 coef_b_ex_dt = this%dt * this%coef_b_ex(nowstage)
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)
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)
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) ) &
2322 else if ( .not. this%imex_flag )
then
2324 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,s)
2329 + coef_a_ex_dt * tend_buf2d_ex(i,j,varid,s)
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_
2343 call prof_rapend(
'rk_advance_trcvar_generall2D', 3)
2346 end subroutine rk_advance_trcvar_general2d
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 )
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)
2373 integer :: tintbuf_ind
2374 real(RP) :: coef_a_im_dt
2377 if (.not. this%imex_flag )
return
2379 tintbuf_ind = this%tend_buf_indmap(nowstage)
2380 coef_a_im_dt = this%dt * this%coef_a_im(nowstage,nowstage)
2384 if ( nowstage == 1 )
then
2388 var0_2d(i,j,varid) = q(i,j)
2389 vartmp_2d(i,j,varid) = q(i,j)
2398 + coef_a_im_dt * tend_buf2d_im(i,j,varid,tintbuf_ind)
2404 end subroutine rk_storeimpl_general2d
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 )
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)
2436 integer :: tintbuf_ind
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
2444 call prof_rapstart(
'rk_advance_general3D', 3)
2446 tintbuf_ind = this%tend_buf_indmap(nowstage)
2448 if ( this%nstage == 1 )
then
2453 vartmp_3d(i,j,k,varid) = q(i,j,k)
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
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)
2474 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
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)
2486 call prof_rapend(
'rk_advance_general3D', 3)
2494 if ( nowstage == 1 .and. (.not. this%imex_flag) )
then
2499 var0_3d(i,j,k,varid) = q(i,j,k)
2500 vartmp_3d(i,j,k,varid) = q(i,j,k)
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
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)
2522 coef_b_ex_dt = this%coef_b_ex(nowstage) * this%dt
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)
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)
2542 q(i,j,k) = var0_3d(i,j,k,varid) &
2543 + coef_a_ex_dt * tend_buf3d_ex(i,j,k,varid,1)
2547 else if ( .not. this%imex_flag )
then
2549 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,s)
2554 q(i,j,k) = q(i,j,k) &
2555 + coef_a_ex_dt * tend_buf3d_ex(i,j,k,varid,s)
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)
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)
2578 call prof_rapend(
'rk_advance_general3D', 3)
2581 end subroutine rk_advance_general3d
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 )
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)
2618 integer :: tintbuf_ind
2622 real(RP) :: coef_a_ex
2623 real(RP) :: coef_a_ex_dt
2624 real(RP) :: coef_b_ex_dt
2628 call prof_rapstart(
'rk_advance_trcvar_general3D', 3)
2630 tintbuf_ind = this%tend_buf_indmap(nowstage)
2632 if ( this%nstage == 1 )
then
2637 vartmp_3d(i,j,k,varid) = q(i,j,k) * ( dens_hyd(i,j,k) + ddens0(i,j,k) )
2643 if ( nowstage == this%nstage )
then
2644 coef_b_ex_dt = this%dt * this%coef_b_ex(nowstage)
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) )
2655 call prof_rapend(
'rk_advance_trcvar_general3D', 3)
2660 c_ss = this%coef_c_ex(nowstage+1)
2665 if ( nowstage == 1 .and. (.not. this%imex_flag) )
then
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)
2677 coef_b_ex_dt = this%dt * this%coef_b_ex(nowstage)
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)
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)
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) ) &
2704 else if ( .not. this%imex_flag )
then
2706 coef_a_ex_dt = this%dt * this%coef_a_ex(nowstage+1,s)
2711 q(i,j,k) = q(i,j,k) &
2712 + coef_a_ex_dt * tend_buf3d_ex(i,j,k,varid,s)
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_
2729 call prof_rapend(
'rk_advance_trcvar_generall3D', 3)
2732 end subroutine rk_advance_trcvar_general3d
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 )
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)
2761 integer :: tintbuf_ind
2762 real(RP) :: coef_a_im_dt
2765 if (.not. this%imex_flag )
return
2767 tintbuf_ind = this%tend_buf_indmap(nowstage)
2768 coef_a_im_dt = this%dt * this%coef_a_im(nowstage,nowstage)
2772 if ( nowstage == 1 )
then
2777 var0_3d(i,j,k,varid) = q(i,j,k)
2778 vartmp_3d(i,j,k,varid) = q(i,j,k)
2788 q(i,j,k) = q(i,j,k) &
2789 + coef_a_im_dt * tend_buf3d_im(i,j,k,varid,tintbuf_ind)
2796 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)
Initialize a object to provide RK scheme.
Derived type to provide RK scheme.