2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-12-08 20:53:29 -06:00
|
|
|
!!
|
2012-07-22 12:30:20 -03:00
|
|
|
!! This file is part of the AMUN source code, a program to perform
|
|
|
|
!! Newtonian or relativistic magnetohydrodynamical simulations on uniform or
|
|
|
|
!! adaptive mesh.
|
2008-12-08 20:53:29 -06:00
|
|
|
!!
|
2013-12-10 15:23:28 -02:00
|
|
|
!! Copyright (C) 2008-2013 Grzegorz Kowal <grzegorz@amuncode.org>
|
2008-12-08 20:53:29 -06:00
|
|
|
!!
|
2012-07-22 12:30:20 -03:00
|
|
|
!! This program is free software: you can redistribute it and/or modify
|
|
|
|
!! it under the terms of the GNU General Public License as published by
|
|
|
|
!! the Free Software Foundation, either version 3 of the License, or
|
|
|
|
!! (at your option) any later version.
|
2008-12-08 20:53:29 -06:00
|
|
|
!!
|
2011-04-29 11:21:30 -03:00
|
|
|
!! This program is distributed in the hope that it will be useful,
|
2008-12-08 20:53:29 -06:00
|
|
|
!! but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
!! GNU General Public License for more details.
|
|
|
|
!!
|
|
|
|
!! You should have received a copy of the GNU General Public License
|
2012-07-22 12:30:20 -03:00
|
|
|
!! along with this program. If not, see <http://www.gnu.org/licenses/>.
|
2008-12-08 20:53:29 -06:00
|
|
|
!!
|
2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-12-08 20:53:29 -06:00
|
|
|
!!
|
2012-07-27 16:36:51 -03:00
|
|
|
!! module: INTERPOLATIONS
|
|
|
|
!!
|
2012-08-05 18:00:10 -03:00
|
|
|
!! This module provides subroutines to interpolate variables and reconstruct
|
|
|
|
!! the Riemann states.
|
|
|
|
!!
|
2012-07-22 12:30:20 -03:00
|
|
|
!!
|
|
|
|
!!******************************************************************************
|
2008-12-08 20:53:29 -06:00
|
|
|
!
|
2012-07-27 16:36:51 -03:00
|
|
|
module interpolations
|
2008-12-08 20:53:29 -06:00
|
|
|
|
2012-07-27 16:46:36 -03:00
|
|
|
! module variables are not implicit by default
|
|
|
|
!
|
2008-12-08 20:53:29 -06:00
|
|
|
implicit none
|
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
! pointers to the reconstruction and limiter procedures
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
procedure(reconstruct) , pointer, save :: reconstruct_states => null()
|
|
|
|
procedure(limiter_zero), pointer, save :: limiter => null()
|
2013-12-11 22:34:29 -02:00
|
|
|
|
2012-07-27 16:46:36 -03:00
|
|
|
! module parameters
|
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
real(kind=8), save :: eps = epsilon(1.0d+00)
|
|
|
|
|
|
|
|
! flags for reconstruction corrections
|
|
|
|
!
|
|
|
|
logical , save :: positivity = .false.
|
2012-07-27 16:46:36 -03:00
|
|
|
|
2012-08-05 19:42:06 -03:00
|
|
|
! by default everything is private
|
2012-07-27 16:46:36 -03:00
|
|
|
!
|
2012-08-05 19:42:06 -03:00
|
|
|
private
|
|
|
|
|
|
|
|
! declare public subroutines
|
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
public :: initialize_interpolations, finalize_interpolations
|
2013-12-12 12:34:05 -02:00
|
|
|
public :: reconstruct, limiter
|
2013-12-11 10:16:06 -02:00
|
|
|
public :: fix_positivity
|
2012-07-27 16:46:36 -03:00
|
|
|
|
|
|
|
!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
|
|
!
|
2008-12-08 20:53:29 -06:00
|
|
|
contains
|
|
|
|
!
|
2012-07-27 16:46:36 -03:00
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine INITIALIZE_INTERPOLATIONS:
|
|
|
|
! ------------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine initializes the interpolation module by reading the module
|
|
|
|
! parameters.
|
|
|
|
!
|
2012-08-05 18:00:10 -03:00
|
|
|
!
|
2012-07-27 16:46:36 -03:00
|
|
|
!===============================================================================
|
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
subroutine initialize_interpolations(verbose, iret)
|
2012-07-27 16:46:36 -03:00
|
|
|
|
|
|
|
! include external procedures
|
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
use parameters, only : get_parameter_string, get_parameter_real
|
2012-07-27 16:46:36 -03:00
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
logical, intent(in) :: verbose
|
|
|
|
integer, intent(inout) :: iret
|
|
|
|
|
2012-07-27 16:46:36 -03:00
|
|
|
! local variables
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
character(len=255) :: sreconstruction = "tvd"
|
|
|
|
character(len=255) :: slimiter = "mm"
|
|
|
|
character(len=255) :: positivity_fix = "off"
|
|
|
|
character(len=255) :: name_rec = ""
|
|
|
|
character(len=255) :: name_lim = ""
|
2012-07-27 16:46:36 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
! obtain the user defined interpolation methods and coefficients
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
call get_parameter_string("reconstruction", sreconstruction)
|
|
|
|
call get_parameter_string("limiter" , slimiter )
|
|
|
|
call get_parameter_string("fix_positivity", positivity_fix )
|
|
|
|
call get_parameter_real ("eps" , eps )
|
2013-12-11 22:34:29 -02:00
|
|
|
|
|
|
|
! select the reconstruction method
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
select case(trim(sreconstruction))
|
2013-12-11 22:34:29 -02:00
|
|
|
case ("tvd", "TVD")
|
|
|
|
name_rec = "2nd order TVD"
|
|
|
|
reconstruct_states => reconstruct_tvd
|
2013-12-12 14:32:27 -02:00
|
|
|
case ("weno3", "WENO3")
|
|
|
|
name_rec = "3rd order WENO"
|
|
|
|
reconstruct_states => reconstruct_weno3
|
2013-12-11 22:34:29 -02:00
|
|
|
case default
|
|
|
|
if (verbose) then
|
|
|
|
write (*,"(1x,a)") "The selected reconstruction method is not " // &
|
2013-12-12 12:34:05 -02:00
|
|
|
"implemented: " // trim(sreconstruction)
|
2013-12-11 22:34:29 -02:00
|
|
|
stop
|
|
|
|
end if
|
|
|
|
end select
|
|
|
|
|
|
|
|
! select the limiter
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
select case(trim(slimiter))
|
|
|
|
case ("mm", "minmod")
|
|
|
|
name_lim = "minmod"
|
|
|
|
limiter => limiter_minmod
|
|
|
|
case ("mc", "monotonized_central")
|
|
|
|
name_lim = "monotonized central"
|
|
|
|
limiter => limiter_monotonized_central
|
|
|
|
case ("sb", "superbee")
|
|
|
|
name_lim = "superbee"
|
|
|
|
limiter => limiter_superbee
|
|
|
|
case ("vl", "vanleer")
|
|
|
|
name_lim = "van Leer"
|
|
|
|
limiter => limiter_vanleer
|
|
|
|
case ("va", "vanalbada")
|
|
|
|
name_lim = "van Albada"
|
|
|
|
limiter => limiter_vanalbada
|
2013-12-11 22:34:29 -02:00
|
|
|
case default
|
2013-12-12 12:34:05 -02:00
|
|
|
name_lim = "zero derivative"
|
|
|
|
limiter => limiter_zero
|
2013-12-11 22:34:29 -02:00
|
|
|
end select
|
|
|
|
|
|
|
|
! check additional reconstruction limiting
|
2012-07-27 16:46:36 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
select case(trim(positivity_fix))
|
|
|
|
case ("on", "ON", "t", "T", "y", "Y", "true", "TRUE", "yes", "YES")
|
|
|
|
positivity = .true.
|
|
|
|
case default
|
|
|
|
positivity = .false.
|
|
|
|
end select
|
2012-07-27 16:46:36 -03:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
! print informations about the reconstruction methods and parameters
|
2012-07-27 16:46:36 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
if (verbose) then
|
|
|
|
|
2013-12-12 13:52:05 -02:00
|
|
|
write (*,"(4x,a14,9x,'=',1x,a)") "reconstruction", trim(name_rec)
|
|
|
|
write (*,"(4x,a14,9x,'=',1x,a)") "limiter ", trim(name_lim)
|
|
|
|
write (*,"(4x,a14,9x,'=',1x,a)") "fix positivity", trim(positivity_fix)
|
2013-12-11 22:34:29 -02:00
|
|
|
|
|
|
|
end if
|
2012-07-27 16:46:36 -03:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine initialize_interpolations
|
2013-12-11 22:34:29 -02:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine FINALIZE_INTERPOLATIONS:
|
|
|
|
! ----------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine finalizes the interpolation module by releasing all memory used
|
|
|
|
! by its module variables.
|
|
|
|
!
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine finalize_interpolations(iret)
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
integer, intent(inout) :: iret
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! release the procedure pointers
|
|
|
|
!
|
|
|
|
nullify(reconstruct_states)
|
2013-12-12 12:34:05 -02:00
|
|
|
nullify(limiter)
|
2013-12-11 22:34:29 -02:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine finalize_interpolations
|
2012-08-05 19:42:06 -03:00
|
|
|
!
|
2008-12-08 20:53:29 -06:00
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-08-05 19:42:06 -03:00
|
|
|
! subroutine RECONSTRUCT:
|
|
|
|
! ----------------------
|
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
! Subroutine calls a reconstruction procedure, depending on the compilation
|
|
|
|
! flag SPACE, in order to interpolate the left and right states from their
|
|
|
|
! cell integrals. These states are required by any approximate Riemann
|
|
|
|
! solver.
|
2012-08-05 19:42:06 -03:00
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! n - the length of the input vector;
|
|
|
|
! h - the spatial step; this is required for some reconstruction methods;
|
|
|
|
! f - the input vector of cell averaged values;
|
|
|
|
! fl - the left side state reconstructed for location (i+1/2);
|
|
|
|
! fr - the right side state reconstructed for location (i+1/2);
|
2008-12-08 20:53:29 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2010-12-06 17:05:44 -02:00
|
|
|
subroutine reconstruct(n, h, f, fl, fr)
|
2008-12-08 20:53:29 -06:00
|
|
|
|
2012-07-27 16:46:36 -03:00
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
2008-12-08 20:53:29 -06:00
|
|
|
implicit none
|
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
! subroutine arguments
|
2008-12-08 20:53:29 -06:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
integer , intent(in) :: n
|
|
|
|
real(kind=8) , intent(in) :: h
|
|
|
|
real(kind=8), dimension(n), intent(in) :: f
|
|
|
|
real(kind=8), dimension(n), intent(out) :: fl, fr
|
2008-12-08 20:53:29 -06:00
|
|
|
!
|
2010-10-13 03:32:10 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-08 20:53:29 -06:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
! reconstruct the states using the selected subroutine
|
2010-12-06 17:05:44 -02:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
call reconstruct_states(n, h, f(:), fl(:), fr(:))
|
2012-08-05 19:42:06 -03:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
!-------------------------------------------------------------------------------
|
2012-08-05 19:42:06 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
end subroutine reconstruct
|
2012-08-05 19:42:06 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
!===============================================================================
|
2012-08-05 19:42:06 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
! subroutine RECONSTRUCT_TVD:
|
|
|
|
! --------------------------
|
2012-08-05 19:42:06 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
! Subroutine reconstructs the interface states using the second order TVD
|
|
|
|
! method with a selected limiter.
|
2012-08-05 22:03:13 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
! Arguments are described in subroutine reconstruct().
|
2012-08-05 19:42:06 -03:00
|
|
|
!
|
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine reconstruct_tvd(n, h, f, fl, fr)
|
2011-05-28 09:49:35 -03:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
! local variables are not implicit by default
|
2012-08-05 19:42:06 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
implicit none
|
2012-08-05 19:42:06 -03:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
! subroutine arguments
|
2012-08-05 19:42:06 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
integer , intent(in) :: n
|
|
|
|
real(kind=8) , intent(in) :: h
|
|
|
|
real(kind=8), dimension(n), intent(in) :: f
|
|
|
|
real(kind=8), dimension(n), intent(out) :: fl, fr
|
2008-12-08 20:53:29 -06:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
! local variables
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
integer :: i, im1, ip1
|
|
|
|
real(kind=8) :: df, dfl, dfr
|
2013-12-11 22:34:29 -02:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! calculate the left- and right-side interface interpolations
|
2013-12-11 22:34:29 -02:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
do i = 1, n
|
2012-08-05 19:42:06 -03:00
|
|
|
|
2013-12-12 12:34:05 -02:00
|
|
|
! calculate left and right indices
|
2012-08-05 19:42:06 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
im1 = max(1, i - 1)
|
|
|
|
ip1 = min(n, i + 1)
|
|
|
|
|
|
|
|
! calculate left and right side derivatives
|
|
|
|
!
|
|
|
|
dfl = f(i ) - f(im1)
|
|
|
|
dfr = f(ip1) - f(i )
|
|
|
|
|
|
|
|
! obtain the TVD limited derivative
|
|
|
|
!
|
|
|
|
df = limiter(0.5d+00, dfl, dfr)
|
2012-08-05 19:42:06 -03:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
! update the left and right-side interpolation states
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
fl(i ) = f(i) + df
|
|
|
|
fr(im1) = f(i) - df
|
2012-08-05 19:42:06 -03:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
end do ! i = 1, n
|
2008-12-08 20:53:29 -06:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
! update the interpolation of the first and last points
|
2010-12-06 17:05:44 -02:00
|
|
|
!
|
2012-08-05 19:42:06 -03:00
|
|
|
fl(1) = f(1)
|
2010-12-06 17:05:44 -02:00
|
|
|
fr(n) = f(n)
|
2011-05-29 12:22:50 -03:00
|
|
|
|
2011-05-28 16:45:49 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
end subroutine reconstruct_tvd
|
2012-08-05 19:42:06 -03:00
|
|
|
!
|
2011-05-28 17:42:13 -03:00
|
|
|
!===============================================================================
|
|
|
|
!
|
2013-12-12 14:32:27 -02:00
|
|
|
! subroutine RECONSTRUCT_WENO3:
|
|
|
|
! ----------------------------
|
|
|
|
!
|
|
|
|
! Subroutine reconstructs the interface states using the third order
|
|
|
|
! Weighted Essentially Non-Oscillatory (WENO) method.
|
|
|
|
!
|
|
|
|
! Arguments are described in subroutine reconstruct().
|
|
|
|
!
|
|
|
|
! References:
|
|
|
|
!
|
|
|
|
! [1] Yamaleev & Carpenter, 2009, J. Comput. Phys., 228, 3025
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine reconstruct_weno3(n, h, f, fl, fr)
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
integer , intent(in) :: n
|
|
|
|
real(kind=8) , intent(in) :: h
|
|
|
|
real(kind=8), dimension(n), intent(in) :: f
|
|
|
|
real(kind=8), dimension(n), intent(out) :: fl, fr
|
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
|
|
|
integer :: i, im1, ip1
|
|
|
|
real(kind=8) :: bp, bm, ap, am, wp, wm, ww
|
|
|
|
real(kind=8) :: dfl, dfr, df, fp, fm, fc, h2
|
|
|
|
|
|
|
|
! selection weights
|
|
|
|
!
|
|
|
|
real, parameter :: dp = 2.0d+00 / 3.0d+00, dm = 1.0d+00 / 3.0d+00
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! prepare common parameters
|
|
|
|
!
|
|
|
|
h2 = h * h
|
|
|
|
|
|
|
|
! iterate along the vector
|
|
|
|
!
|
|
|
|
do i = 1, n
|
|
|
|
|
|
|
|
! prepare neighbour indices
|
|
|
|
!
|
|
|
|
im1 = max(1, i - 1)
|
|
|
|
ip1 = min(n, i + 1)
|
|
|
|
|
|
|
|
! calculate the left and right derivatives
|
|
|
|
!
|
|
|
|
dfl = f(i ) - f(im1)
|
|
|
|
dfr = f(ip1) - f(i )
|
|
|
|
|
|
|
|
! calculate coefficient omega
|
|
|
|
!
|
|
|
|
ww = (dfr - dfl)**2
|
|
|
|
|
|
|
|
! calculate corresponding betas
|
|
|
|
!
|
|
|
|
bp = dfr * dfr
|
|
|
|
bm = dfl * dfl
|
|
|
|
|
|
|
|
! calculate improved alphas
|
|
|
|
!
|
|
|
|
ap = 1.0d+00 + ww / (bp + h2)
|
|
|
|
am = 1.0d+00 + ww / (bm + h2)
|
|
|
|
|
|
|
|
! calculate weights
|
|
|
|
!
|
|
|
|
wp = dp * ap
|
|
|
|
wm = dm * am
|
|
|
|
ww = 2.0d+00 * (wp + wm)
|
|
|
|
|
|
|
|
! calculate central interpolation
|
|
|
|
!
|
|
|
|
fp = f(i ) + f(ip1)
|
|
|
|
|
|
|
|
! calculate left side interpolation
|
|
|
|
!
|
|
|
|
fm = - f(im1) + 3.0d+00 * f(i )
|
|
|
|
|
|
|
|
! calculate the left state
|
|
|
|
!
|
|
|
|
fl( i ) = (wp * fp + wm * fm) / ww
|
|
|
|
|
|
|
|
! calculate weights
|
|
|
|
!
|
|
|
|
wp = dp * am
|
|
|
|
wm = dm * ap
|
|
|
|
ww = 2.0d+00 * (wp + wm)
|
|
|
|
|
|
|
|
! calculate central interpolation
|
|
|
|
!
|
|
|
|
fp = f(i ) + f(im1)
|
|
|
|
|
|
|
|
! calculate right side interpolation
|
|
|
|
!
|
|
|
|
fm = - f(ip1) + 3.0d+00 * f(i )
|
|
|
|
|
|
|
|
! calculate the right state
|
|
|
|
!
|
|
|
|
fr(im1) = (wp * fp + wm * fm) / ww
|
|
|
|
|
|
|
|
end do ! i = 1, n
|
|
|
|
|
|
|
|
! update the interpolation of the first and last points
|
|
|
|
!
|
|
|
|
fl(1) = f (1)
|
|
|
|
fr(n) = fl(n)
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine reconstruct_weno3
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! function LIMITER_ZERO:
|
|
|
|
! ---------------------
|
|
|
|
!
|
|
|
|
! Function returns zero.
|
|
|
|
!
|
|
|
|
! Arguments:
|
2013-12-11 22:34:29 -02:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! x - scaling factor;
|
|
|
|
! a, b - the input values;
|
2011-05-28 17:42:13 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
function limiter_zero(x, a, b) result(c)
|
2011-05-28 17:42:13 -03:00
|
|
|
|
2013-12-12 12:34:05 -02:00
|
|
|
! local variables are not implicit by default
|
2011-05-28 17:42:13 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
implicit none
|
2011-05-28 17:42:13 -03:00
|
|
|
|
2013-12-12 12:34:05 -02:00
|
|
|
! input arguments
|
2011-05-28 17:42:13 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
real(kind=8), intent(in) :: x, a, b
|
|
|
|
real(kind=8) :: c
|
2011-05-28 17:42:13 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
!-------------------------------------------------------------------------------
|
2011-05-28 17:42:13 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
c = 0.0d+00
|
2011-05-28 17:42:13 -03:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
!-------------------------------------------------------------------------------
|
2011-05-28 17:42:13 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
end function limiter_zero
|
2011-05-28 17:42:13 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
!===============================================================================
|
2011-05-28 17:42:13 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! function LIMITER_MINMOD:
|
|
|
|
! -----------------------
|
|
|
|
!
|
|
|
|
! Function returns the minimum module value among two arguments using
|
|
|
|
! minmod limiter.
|
2011-05-28 17:42:13 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! x - scaling factor;
|
|
|
|
! a, b - the input values;
|
2011-05-28 17:42:13 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
!===============================================================================
|
2011-05-28 17:42:13 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
function limiter_minmod(x, a, b) result(c)
|
2013-12-11 22:34:29 -02:00
|
|
|
|
2013-12-12 12:34:05 -02:00
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
implicit none
|
2011-05-28 17:42:13 -03:00
|
|
|
|
2013-12-12 12:34:05 -02:00
|
|
|
! input arguments
|
2011-05-28 17:42:13 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
real(kind=8), intent(in) :: x, a, b
|
|
|
|
real(kind=8) :: c
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
c = 0.5d+00 * (sign(x, a) + sign(x, b)) * min(abs(a), abs(b))
|
2011-05-28 17:42:13 -03:00
|
|
|
|
2013-12-12 12:34:05 -02:00
|
|
|
!-------------------------------------------------------------------------------
|
2011-05-28 17:42:13 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
end function limiter_minmod
|
2013-12-11 22:34:29 -02:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
!===============================================================================
|
2013-12-11 22:34:29 -02:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! function LIMITER_MONOTONIZED_CENTRAL:
|
|
|
|
! ------------------------------------
|
2013-12-11 22:34:29 -02:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! Function returns the minimum module value among two arguments using
|
|
|
|
! the monotonized central TVD limiter.
|
2013-12-11 22:34:29 -02:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! x - scaling factor;
|
|
|
|
! a, b - the input values;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
function limiter_monotonized_central(x, a, b) result(c)
|
2013-12-11 22:34:29 -02:00
|
|
|
|
2013-12-12 12:34:05 -02:00
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
2013-12-11 22:34:29 -02:00
|
|
|
|
2013-12-12 12:34:05 -02:00
|
|
|
! input arguments
|
|
|
|
!
|
|
|
|
real(kind=8), intent(in) :: x, a, b
|
|
|
|
real(kind=8) :: c
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
c = (sign(x, a) + sign(x, b)) * min(abs(a), abs(b), 2.5d-01 * abs(a + b))
|
2011-05-29 12:22:50 -03:00
|
|
|
|
2011-05-28 17:42:13 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
end function limiter_monotonized_central
|
2011-05-28 17:42:13 -03:00
|
|
|
!
|
2011-05-28 16:45:49 -03:00
|
|
|
!===============================================================================
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! function LIMITER_SUPERBEE:
|
|
|
|
! -------------------------
|
2013-12-11 22:34:29 -02:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! Function returns the minimum module value among two arguments using
|
|
|
|
! superbee limiter.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! x - scaling factor;
|
|
|
|
! a, b - the input values;
|
2011-05-28 16:45:49 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
function limiter_superbee(x, a, b) result(c)
|
2011-05-28 16:45:49 -03:00
|
|
|
|
2013-12-12 12:34:05 -02:00
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
2011-05-28 16:45:49 -03:00
|
|
|
implicit none
|
|
|
|
|
2013-12-12 12:34:05 -02:00
|
|
|
! input arguments
|
2011-05-28 16:45:49 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
real(kind=8), intent(in) :: x, a, b
|
|
|
|
real(kind=8) :: c
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
c = 0.5d+00 * (sign(x, a) + sign(x, b)) &
|
|
|
|
* max(min(2.0d+00 * abs(a), abs(b)), min(abs(a), 2.0d+00 * abs(b)))
|
2011-05-28 16:45:49 -03:00
|
|
|
|
2013-12-12 12:34:05 -02:00
|
|
|
!-------------------------------------------------------------------------------
|
2011-05-28 16:45:49 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
end function limiter_superbee
|
2011-05-28 16:45:49 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
!===============================================================================
|
2011-05-28 16:45:49 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! function LIMITER_VANLEER:
|
|
|
|
! ------------------------
|
2010-12-06 20:54:42 -02:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! Function returns the minimum module value among two arguments using
|
|
|
|
! van Leer's limiter.
|
2010-12-06 20:54:42 -02:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! Arguments:
|
2011-03-25 01:05:58 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! x - scaling factor;
|
|
|
|
! a, b - the input values;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
2010-12-06 20:54:42 -02:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
function limiter_vanleer(x, a, b) result(c)
|
2010-12-06 20:54:42 -02:00
|
|
|
|
2013-12-12 12:34:05 -02:00
|
|
|
! local variables are not implicit by default
|
2010-12-06 20:54:42 -02:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
implicit none
|
2010-12-06 20:54:42 -02:00
|
|
|
|
2013-12-12 12:34:05 -02:00
|
|
|
! input arguments
|
|
|
|
!
|
|
|
|
real(kind=8), intent(in) :: x, a, b
|
|
|
|
real(kind=8) :: c
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
c = a * b
|
|
|
|
if (c > 0.0d+00) then
|
|
|
|
c = 2.0d+00 * x * c / (a + b)
|
|
|
|
else
|
|
|
|
c = 0.0d+00
|
|
|
|
end if
|
2011-05-29 12:22:50 -03:00
|
|
|
|
2011-05-28 16:45:49 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
end function limiter_vanleer
|
2008-12-12 16:39:03 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! function LIMITER_VANALBADA:
|
|
|
|
! --------------------------
|
2012-08-05 19:42:06 -03:00
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! Function returns the minimum module value among two arguments using
|
|
|
|
! van Albada's limiter.
|
2012-08-05 19:42:06 -03:00
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
! x - scaling factor;
|
2013-12-11 22:34:29 -02:00
|
|
|
! a, b - the input values;
|
2010-12-07 10:08:30 -02:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
function limiter_vanalbada(x, a, b) result(c)
|
2010-12-07 10:08:30 -02:00
|
|
|
|
2012-08-05 19:42:06 -03:00
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
2010-12-07 10:08:30 -02:00
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
real(kind=8), intent(in) :: x, a, b
|
|
|
|
real(kind=8) :: c
|
2010-12-07 10:08:30 -02:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
c = x * a * b * (a + b) / max(eps, a * a + b * b)
|
2011-06-09 14:12:33 -03:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2013-12-12 12:34:05 -02:00
|
|
|
end function limiter_vanalbada
|
2010-12-07 10:08:30 -02:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-08-05 18:00:10 -03:00
|
|
|
! subroutine FIX_POSITIVITY:
|
|
|
|
! -------------------------
|
2011-03-24 16:15:51 -03:00
|
|
|
!
|
2012-08-05 18:00:10 -03:00
|
|
|
! Subroutine scans the input arrays of the left and right states fl(:) and
|
|
|
|
! fr(:) for negative values. If it finds a negative value, it repeates the
|
|
|
|
! state reconstruction from f(:) using the zeroth order interpolation.
|
2011-03-24 16:15:51 -03:00
|
|
|
!
|
2012-08-05 18:00:10 -03:00
|
|
|
!===============================================================================
|
2011-03-24 16:15:51 -03:00
|
|
|
!
|
2012-08-05 18:00:10 -03:00
|
|
|
subroutine fix_positivity(n, f, fl, fr)
|
2011-06-09 14:47:59 -03:00
|
|
|
|
2012-08-05 18:00:10 -03:00
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
2011-06-09 14:47:59 -03:00
|
|
|
|
2012-08-05 18:00:10 -03:00
|
|
|
! input/output arguments
|
2011-06-09 14:47:59 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
integer , intent(in) :: n
|
|
|
|
real(kind=8), dimension(n), intent(in) :: f
|
|
|
|
real(kind=8), dimension(n), intent(inout) :: fl, fr
|
2012-08-05 18:00:10 -03:00
|
|
|
|
|
|
|
! local variables
|
2011-03-24 16:15:51 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
integer :: i, im1, ip1
|
|
|
|
real(kind=8) :: fmn, fmx
|
2011-03-24 16:15:51 -03:00
|
|
|
!
|
2012-08-05 18:00:10 -03:00
|
|
|
!------------------------------------------------------------------------------
|
2010-12-07 10:08:30 -02:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
! check positivity only if desired
|
2010-12-07 10:08:30 -02:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
if (positivity) then
|
2010-12-07 10:08:30 -02:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
! look for negative values in the states along the vector
|
2010-12-07 10:08:30 -02:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
do i = 1, n
|
2012-08-05 18:00:10 -03:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
! check if the left state has a negative value
|
2010-12-07 10:08:30 -02:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
if (fl(i) <= 0.0d+00) then
|
2010-12-07 10:08:30 -02:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
! calculate the left neighbour index
|
2010-12-07 10:08:30 -02:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
im1 = max(1, i - 1)
|
2010-07-03 23:37:02 -03:00
|
|
|
|
2012-08-05 18:00:10 -03:00
|
|
|
! limit the states using the zeroth-order reconstruction
|
2010-07-03 23:37:02 -03:00
|
|
|
!
|
2012-08-05 18:00:10 -03:00
|
|
|
fl(i ) = f(i)
|
|
|
|
fr(im1) = f(i)
|
2010-07-03 23:37:02 -03:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
end if ! fl ≤ 0
|
2010-02-11 23:30:46 -02:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
! check if the right state has a negative value
|
2010-03-30 17:44:48 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
if (fr(i) <= 0.0d+00) then
|
|
|
|
|
|
|
|
! calculate the right neighbour index
|
2010-03-30 17:44:48 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
ip1 = min(n, i + 1)
|
|
|
|
|
|
|
|
! limit the states using the zeroth-order reconstruction
|
2010-03-30 17:44:48 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
fl(ip1) = f(ip1)
|
|
|
|
fr(i ) = f(ip1)
|
2010-03-30 17:44:48 -03:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
end if ! fr ≤ 0
|
2010-03-30 17:44:48 -03:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
end do ! i = 1, n
|
2010-03-30 17:44:48 -03:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
end if ! positivity == .true.
|
2010-03-30 17:44:48 -03:00
|
|
|
|
2013-12-11 22:34:29 -02:00
|
|
|
!-------------------------------------------------------------------------------
|
2010-03-30 17:44:48 -03:00
|
|
|
!
|
2013-12-11 22:34:29 -02:00
|
|
|
end subroutine fix_positivity
|
2012-08-01 16:42:45 -03:00
|
|
|
|
2008-12-08 20:53:29 -06:00
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-27 16:36:51 -03:00
|
|
|
end module interpolations
|