aboutsummaryrefslogtreecommitdiff
path: root/src/Dissipation_2_1.F90
diff options
context:
space:
mode:
Diffstat (limited to 'src/Dissipation_2_1.F90')
-rw-r--r--src/Dissipation_2_1.F9025
1 files changed, 20 insertions, 5 deletions
diff --git a/src/Dissipation_2_1.F90 b/src/Dissipation_2_1.F90
index 62ba030..0904b33 100644
--- a/src/Dissipation_2_1.F90
+++ b/src/Dissipation_2_1.F90
@@ -1,18 +1,21 @@
! $Header$
#include "cctk.h"
+#include "cctk_Parameters.h"
-subroutine dissipation_2_1 (var, ni, nj, nk, bb, gsize, delta, epsdis, rhs)
+subroutine dissipation_2_1 (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)
@@ -24,7 +27,11 @@ subroutine dissipation_2_1 (var, ni, nj, nk, bb, gsize, delta, epsdis, rhs)
call set_coeff ( a )
- idel = epsdis / 4
+ if ( scale_with_h > 0 ) then
+ idel = epsilon / ( 4 * delta(1) )
+ else
+ idel = epsilon / 4
+ end if
if ( bb(1) == 0 ) then
il = 1 + gsize(1)
@@ -55,7 +62,11 @@ subroutine dissipation_2_1 (var, ni, nj, nk, bb, gsize, delta, epsdis, rhs)
call set_coeff ( a )
- idel = epsdis / 4
+ if ( scale_with_h > 0 ) then
+ idel = epsilon / ( 4 * delta(2) )
+ else
+ idel = epsilon / 4
+ end if
if ( bb(3) == 0 ) then
jl = 1 + gsize(2)
@@ -86,7 +97,11 @@ subroutine dissipation_2_1 (var, ni, nj, nk, bb, gsize, delta, epsdis, rhs)
call set_coeff ( a )
- idel = epsdis / 4
+ if ( scale_with_h > 0 ) then
+ idel = epsilon / ( 4 * delta(3) )
+ else
+ idel = epsilon / 4
+ end if
if ( bb(5) == 0 ) then
kl = 1 + gsize(3)