Merge branch 'master' into reconnection

This commit is contained in:
Grzegorz Kowal 2020-09-26 10:06:12 -03:00
commit 907e2afdab

View File

@ -361,6 +361,12 @@ module interpolations
reconstruct_states => reconstruct_crmp7
order = 7
nghosts = max(nghosts, 4)
case ("crmp7l", "crmp7ld", "CRMP7L", "CRMP7LD")
name_rec = "7th order Low-Dissipation Compact Monotonicity Preserving"
interfaces => interfaces_dir
reconstruct_states => reconstruct_crmp7ld
order = 7
nghosts = max(nghosts, 4)
case ("gp", "GP")
write(stmp, '(f16.1)') sgp
write(name_rec, '("Gaussian Process (",i1,"-point, δ=",a,")")') ngp &
@ -5129,6 +5135,189 @@ module interpolations
!
!===============================================================================
!
! subroutine RECONSTRUCT_CRMP7LD:
! ------------------------------
!
! Subroutine reconstructs the interface states using the seventh order
! low dissipation Compact Reconstruction Monotonicity Preserving (CRMP)
! method.
!
! Arguments are described in subroutine reconstruct().
!
! References:
!
! [1] Suresh, A. & Huynh, H. T.,
! "Accurate Monotonicity-Preserving Schemes with Runge-Kutta
! Time Stepping"
! Journal on Computational Physics,
! 1997, vol. 136, pp. 83-99,
! http://dx.doi.org/10.1006/jcph.1997.5745
! [2] He, ZhiWei, Li, XinLiang, Fu, DeXun, & Ma, YanWen,
! "A 5th order monotonicity-preserving upwind compact difference
! scheme",
! Science China Physics, Mechanics and Astronomy,
! Volume 54, Issue 3, pp. 511-522,
! http://dx.doi.org/10.1007/s11433-010-4220-x
!
!===============================================================================
!
subroutine reconstruct_crmp7ld(h, fc, fl, fr)
! include external procedures
!
use algebra , only : tridiag
! local variables are not implicit by default
!
implicit none
! subroutine arguments
!
real(kind=8) , intent(in) :: h
real(kind=8), dimension(:), intent(in) :: fc
real(kind=8), dimension(:), intent(out) :: fl, fr
! local variables
!
integer :: n, i, iret
! local arrays for derivatives
!
real(kind=8), dimension(size(fc)) :: fi
real(kind=8), dimension(size(fc)) :: a, b, c
real(kind=8), dimension(size(fc)) :: r
real(kind=8), dimension(size(fc)) :: u
! local parameters
!
real(kind=8), dimension(3), parameter :: &
di = [ 5.0d+00, 1.2d+01, 4.0d+00 ] / 2.1d+01
real(kind=8), dimension(6), parameter :: &
ci = [-2.00d+00, 4.20d+01, 6.37d+02, &
5.57d+02, 2.70d+01,-1.00d+00 ] / 1.26d+03
real(kind=8), dimension(7), parameter :: &
ce7 = [-3.0d+00, 2.5d+01,-1.01d+02, 3.19d+02, &
2.14d+02,-3.8d+01, 4.0d+00 ] / 4.2d+02
real(kind=8), dimension(5), parameter :: &
ce5 = [ 2.0d+00,-1.3d+01, 4.70d+01, &
2.70d+01,-3.0d+00 ] / 6.0d+01
real(kind=8), dimension(3), parameter :: &
ce3 = [-1.0d+00, 5.0d+00, 2.0d+00 ] / 6.0d+00
real(kind=8), dimension(2), parameter :: &
ce2 = [ 5.0d-01, 5.0d-01 ]
!
!-------------------------------------------------------------------------------
!
! get the input vector length
!
n = size(fc)
! prepare the diagonals of the tridiagonal matrix
!
do i = 1, ng
a(i) = 0.0d+00
b(i) = 1.0d+00
c(i) = 0.0d+00
end do
do i = ng + 1, n - ng - 1
a(i) = di(1)
b(i) = di(2)
c(i) = di(3)
end do
do i = n - ng, n
a(i) = 0.0d+00
b(i) = 1.0d+00
c(i) = 0.0d+00
end do
!! === left-side interpolation ===
!!
! prepare the tridiagonal system coefficients for the interior part
!
do i = ng, n - ng + 1
r(i) = sum( ci(:) * fc(i-2:i+3))
end do
! interpolate ghost zones using the explicit method
!
r( 1) = sum(ce2(:) * fc( 1: 2))
r( 2) = sum(ce3(:) * fc( 1: 3))
r( 3) = sum(ce5(:) * fc( 1: 5))
do i = 4, ng
r(i) = sum(ce7(:) * fc(i-3:i+3))
end do
do i = n - ng, n - 3
r(i) = sum(ce7(:) * fc(i-3:i+3))
end do
r(n-2) = sum(ce5(:) * fc(n-4: n))
r(n-1) = sum(ce3(:) * fc(n-2: n))
r(n ) = fc(n )
! solve the tridiagonal system of equations
!
call tridiag(n, a(1:n), b(1:n), c(1:n), r(1:n), u(1:n), iret)
! apply the monotonicity preserving limiting
!
call mp_limiting(fc(:), u(:))
! copy the interpolation to the respective vector
!
fl(1:n) = u(1:n)
!! === right-side interpolation ===
!!
! invert the cell-centered value vector
!
fi(1:n) = fc(n:1:-1)
! prepare the tridiagonal system coefficients for the interior part
!
do i = ng, n - ng + 1
r(i) = sum( ci(:) * fi(i-2:i+3))
end do ! i = ng, n - ng + 1
! interpolate ghost zones using the explicit method
!
r( 1) = sum(ce2(:) * fi( 1: 2))
r( 2) = sum(ce3(:) * fi( 1: 3))
r( 3) = sum(ce5(:) * fi( 1: 5))
do i = 4, ng
r(i) = sum(ce7(:) * fi(i-3:i+3))
end do
do i = n - ng, n - 3
r(i) = sum(ce7(:) * fi(i-3:i+3))
end do
r(n-2) = sum(ce5(:) * fi(n-4: n))
r(n-1) = sum(ce3(:) * fi(n-2: n))
r(n ) = fi(n )
! solve the tridiagonal system of equations
!
call tridiag(n, a(1:n), b(1:n), c(1:n), r(1:n), u(1:n), iret)
! apply the monotonicity preserving limiting
!
call mp_limiting(fi(:), u(:))
! copy the interpolation to the respective vector
!
fr(1:n-1) = u(n-1:1:-1)
! update the interpolation of the first and last points
!
i = n - 1
fl(1) = 0.5d+00 * (fc(1) + fc(2))
fr(i) = 0.5d+00 * (fc(i) + fc(n))
fl(n) = fc(n)
fr(n) = fc(n)
!-------------------------------------------------------------------------------
!
end subroutine reconstruct_crmp7ld
!
!===============================================================================
!
! subroutine PREPARE_GP:
! ---------------------
!