aboutsummaryrefslogtreecommitdiff
path: root/src/Dissipation_8_4.F90
diff options
context:
space:
mode:
authordiener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a>2005-12-06 17:20:11 +0000
committerdiener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a>2005-12-06 17:20:11 +0000
commit74de5846f7ddeb3bf181b3205f6539b4f0081705 (patch)
tree70bbe4f8269d32fadb6526eb58f810837d64a018 /src/Dissipation_8_4.F90
parent2f50e2f41868514ac23f5948caa4210607cf8917 (diff)
Added a new boolean parameter "scale_with_h" that allows the user to
choose to scale the dissipation with the grid spacing h. This results formally in the loss of one order of convergence. The default value is "no". git-svn-id: https://svn.cct.lsu.edu/repos/numrel/LSUThorns/SummationByParts/trunk@52 f69c4107-0314-4c4f-9ad4-17e986b73f4a
Diffstat (limited to 'src/Dissipation_8_4.F90')
-rw-r--r--src/Dissipation_8_4.F9025
1 files changed, 20 insertions, 5 deletions
diff --git a/src/Dissipation_8_4.F90 b/src/Dissipation_8_4.F90
index bc225af..ec98407 100644
--- a/src/Dissipation_8_4.F90
+++ b/src/Dissipation_8_4.F90
@@ -1,18 +1,21 @@
! $Header$
#include "cctk.h"
+#include "cctk_Parameters.h"
-subroutine dissipation_8_4 (var, ni, nj, nk, bb, gsize, delta, epsdis, rhs)
+subroutine dissipation_8_4 (var, ni, nj, nk, bb, gsize, delta, epsilon, rhs)
implicit none
+ DECLARE_CCTK_PARAMETERS
+
integer :: ni, nj, nk
CCTK_REAL, dimension(ni,nj,nk), intent(in) :: var
CCTK_REAL, dimension(ni,nj,nk), intent(inout) :: rhs
CCTK_INT, dimension(6), intent(in) :: bb
CCTK_INT, dimension(3), intent(in) :: gsize
CCTK_REAL, dimension(3), intent(in) :: delta
- CCTK_REAL, intent(in) :: epsdis
+ CCTK_REAL, intent(in) :: epsilon
CCTK_REAL :: zero = 0.0
integer, parameter :: wp = kind(zero)
@@ -23,7 +26,11 @@ subroutine dissipation_8_4 (var, ni, nj, nk, bb, gsize, delta, epsdis, rhs)
call set_coeff ( a )
- idel = epsdis / 256
+ if ( scale_with_h > 0 ) then
+ idel = epsilon / ( 256 * delta(1) )
+ else
+ idel = epsilon / 256
+ end if
if ( bb(1) == 0 ) then
il = 1 + gsize(1)
@@ -134,7 +141,11 @@ subroutine dissipation_8_4 (var, ni, nj, nk, bb, gsize, delta, epsdis, rhs)
call set_coeff ( a )
- idel = epsdis / 256
+ if ( scale_with_h > 0 ) then
+ idel = epsilon / ( 256 * delta(2) )
+ else
+ idel = epsilon / 256
+ end if
if ( bb(3) == 0 ) then
jl = 1 + gsize(2)
@@ -245,7 +256,11 @@ subroutine dissipation_8_4 (var, ni, nj, nk, bb, gsize, delta, epsdis, rhs)
call set_coeff ( a )
- idel = epsdis / 256
+ if ( scale_with_h > 0 ) then
+ idel = epsilon / ( 256 * delta(3) )
+ else
+ idel = epsilon / 256
+ end if
if ( bb(5) == 0 ) then
kl = 1 + gsize(3)