diff options
author | diener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a> | 2005-12-15 03:54:19 +0000 |
---|---|---|
committer | diener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a> | 2005-12-15 03:54:19 +0000 |
commit | 4c7144d8e695475f86fec845e0730179a8959fc2 (patch) | |
tree | 5bc120016e74900734d0cd7967c1861f94a2fbad | |
parent | 74de5846f7ddeb3bf181b3205f6539b4f0081705 (diff) |
Introduce a new kind of dissipation operators. Kreiss-Oliger type dissipation
operator suitably modified near the boundary. Since the dissipation operator
in the interior uses a stencil that is one point wider than the derivative
operator, an additional check on the number of ghostzones in multi processor
runs have been added. Only implemented for the 8-4 and optimized 6-5
operators.
git-svn-id: https://svn.cct.lsu.edu/repos/numrel/LSUThorns/SummationByParts/trunk@53 f69c4107-0314-4c4f-9ad4-17e986b73f4a
-rw-r--r-- | param.ccl | 6 | ||||
-rw-r--r-- | src/CheckGridSizes.F90 | 10 | ||||
-rw-r--r-- | src/Dissipation_6_5_min_err_coeff_alt.F90 | 753 | ||||
-rw-r--r-- | src/Dissipation_8_4_alt.F90 | 513 | ||||
-rw-r--r-- | src/dissipation.c | 39 | ||||
-rw-r--r-- | src/make.code.defn | 2 |
6 files changed, 1319 insertions, 4 deletions
@@ -27,6 +27,12 @@ BOOLEAN scale_with_h "Should we scale the dissipation with the grid spacing h" { } "no" +KEYWORD dissipation_type "Type of dissipation operator" STEERABLE=always +{ + "Mattson-Svard-Nordstrom" :: "Mattson, Svaerd and Nordstroem type" + "Kreiss-Oliger" :: "Kreiss-Oliger modified near the boundaries" +} "Mattson-Svard-Nordstrom" + REAL epsdis "Dissipation strength" STEERABLE=always { *:* :: "Values typical between 0 and 1" diff --git a/src/CheckGridSizes.F90 b/src/CheckGridSizes.F90 index 3240f6a..b6af08e 100644 --- a/src/CheckGridSizes.F90 +++ b/src/CheckGridSizes.F90 @@ -25,6 +25,16 @@ subroutine SBP_CheckGridSizes(CCTK_ARGUMENTS) bc(2) = cctk_bbox(3)+cctk_bbox(4)+1 ! 1 processor boundary -> bc = 2 bc(3) = cctk_bbox(5)+cctk_bbox(6)+1 ! 2 processor boundaries -> bc = 1 + if ( ( use_dissipation > 0 ) .and. & + ( CCTK_EQUALS(dissipation_type,"Kreiss-Oliger" ) ) .and. & + ( any (cctk_nghostzones < order / 2 + 1) ) .and. & + ( CCTK_nProcs(cctkGH) > 1 ) ) then + write(info_message,'(a22,i1,a13,i2,a32)') 'You need ghostsize >= ', & + order/2+1, ' to run with ', order+2, & + ' order Kreiss-Oliger dissipation' + call CCTK_WARN(0,info_message) + end if + if ( any (cctk_nghostzones < order / 2) .and. CCTK_nProcs(cctkGH) > 1 ) then write(info_message,'(a22,i1,a13,i1,a25)') 'You need ghostsize >= ', & order/2, ' to run with ', order, & diff --git a/src/Dissipation_6_5_min_err_coeff_alt.F90 b/src/Dissipation_6_5_min_err_coeff_alt.F90 new file mode 100644 index 0000000..102ac5d --- /dev/null +++ b/src/Dissipation_6_5_min_err_coeff_alt.F90 @@ -0,0 +1,753 @@ +! $Header$ + +#include "cctk.h" +#include "cctk_Parameters.h" +#include "cctk_Functions.h" + + +subroutine dissipation_6_5_alt (var, lsh, gsh, lbnd, bb, gsize, & + delta, epsdisl, dfl, rhs) + + implicit none + + DECLARE_CCTK_PARAMETERS + DECLARE_CCTK_FUNCTIONS + + CCTK_INT, dimension(3), intent(in) :: lsh, gsh, lbnd + CCTK_INT :: ni, nj, nk + CCTK_REAL, dimension(lsh(1),lsh(2),lsh(3)), intent(in) :: var + CCTK_REAL, dimension(lsh(1),lsh(2),lsh(3)), 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) :: epsdisl, dfl + + CCTK_REAL :: zero = 0.0 + integer, parameter :: wp = kind(zero) + integer :: i, j, k, center + CCTK_REAL, dimension(:,:), allocatable :: a, d, b, h, tmp + CCTK_REAL :: idel + + CCTK_INT :: il, ir, jl, jr, kl, kr + + ni = lsh(1); nj = lsh(2); nk = lsh(3) + + allocate ( a(ni,ni), d(ni,ni), b(ni,ni), h(ni,ni), tmp(ni,ni) ) + a = zero; d = zero; b = zero; h = zero; tmp = zero + + call set_bmatrix ( b, bb(1:2), lsh(1), gsh(1), lbnd(1), delta(1), dfl ) + + call set_hmatrix ( h, bb(1:2) ) + + call set_dmatrix ( d, bb(1:2) ) + + a = -transpose ( matmul ( h, matmul ( transpose(d), matmul ( b, d ) ) ) ) + + idel = epsdisl / ( 256.0_wp * delta(1) ) + + if ( bb(1) == 0 ) then + + il = 1 + gsize(1) + + else + + rhs(1,:,:) = rhs(1,:,:) + & + ( a(1,1) * var(1,:,:) + a(2,1) * var(2,:,:) + & + a(3,1) * var(3,:,:) + a(4,1) * var(4,:,:) + & + a(5,1) * var(5,:,:) ) * idel + rhs(2,:,:) = rhs(2,:,:) + & + ( a(1,2) * var(1,:,:) + a(2,2) * var(2,:,:) + & + a(3,2) * var(3,:,:) + a(4,2) * var(4,:,:) + & + a(5,2) * var(5,:,:) + a(6,2) * var(6,:,:) + & + a(7,2) * var(7,:,:) + a(8,2) * var(8,:,:) + & + a(9,2) * var(9,:,:) + a(10,2) * var(10,:,:) + & + a(11,2) * var(11,:,:) ) * idel + rhs(3,:,:) = rhs(3,:,:) + & + ( a(1,3) * var(1,:,:) + a(2,3) * var(2,:,:) + & + a(3,3) * var(3,:,:) + a(4,3) * var(4,:,:) + & + a(5,3) * var(5,:,:) + a(6,3) * var(6,:,:) + & + a(7,3) * var(7,:,:) + a(8,3) * var(8,:,:) + & + a(9,3) * var(9,:,:) + a(10,3) * var(10,:,:) + & + a(11,3) * var(11,:,:) ) * idel + rhs(4,:,:) = rhs(4,:,:) + & + ( a(1,4) * var(1,:,:) + a(2,4) * var(2,:,:) + & + a(3,4) * var(3,:,:) + a(4,4) * var(4,:,:) + & + a(5,4) * var(5,:,:) + a(6,4) * var(6,:,:) + & + a(7,4) * var(7,:,:) + a(8,4) * var(8,:,:) + & + a(9,4) * var(9,:,:) + a(10,4) * var(10,:,:) + & + a(11,4) * var(11,:,:) ) * idel + rhs(5,:,:) = rhs(5,:,:) + & + ( a(1,5) * var(1,:,:) + a(2,5) * var(2,:,:) + & + a(3,5) * var(3,:,:) + a(4,5) * var(4,:,:) + & + a(5,5) * var(5,:,:) + a(6,5) * var(6,:,:) + & + a(7,5) * var(7,:,:) + a(8,5) * var(8,:,:) + & + a(9,5) * var(9,:,:) + a(10,5) * var(10,:,:) + & + a(11,5) * var(11,:,:) ) * idel + rhs(6,:,:) = rhs(6,:,:) + & + ( a(1,6) * var(1,:,:) + a(2,6) * var(2,:,:) + & + a(3,6) * var(3,:,:) + a(4,6) * var(4,:,:) + & + a(5,6) * var(5,:,:) + a(6,6) * var(6,:,:) + & + a(7,6) * var(7,:,:) + a(8,6) * var(8,:,:) + & + a(9,6) * var(9,:,:) + a(10,6) * var(10,:,:) + & + a(11,6) * var(11,:,:) ) * idel + rhs(7,:,:) = rhs(7,:,:) + & + ( a(1,7) * var(1,:,:) + a(2,7) * var(2,:,:) + & + a(3,7) * var(3,:,:) + a(4,7) * var(4,:,:) + & + a(5,7) * var(5,:,:) + a(6,7) * var(6,:,:) + & + a(7,7) * var(7,:,:) + a(8,7) * var(8,:,:) + & + a(9,7) * var(9,:,:) + a(10,7) * var(10,:,:) + & + a(11,7) * var(11,:,:) ) * idel + + il = 8 + + end if + + if ( bb(2) == 0 ) then + + ir = ni - gsize(1) + + else + + rhs(ni-6,:,:) = rhs(ni-6,:,:) + & + ( a(ni-10,ni-6) * var(ni-10,:,:) + & + a(ni-9,ni-6) * var(ni-9,:,:) + & + a(ni-8,ni-6) * var(ni-8,:,:) + & + a(ni-7,ni-6) * var(ni-7,:,:) + & + a(ni-6,ni-6) * var(ni-6,:,:) + & + a(ni-5,ni-6) * var(ni-5,:,:) + & + a(ni-4,ni-6) * var(ni-4,:,:) + & + a(ni-3,ni-6) * var(ni-3,:,:) + & + a(ni-2,ni-6) * var(ni-2,:,:) + & + a(ni-1,ni-6) * var(ni-1,:,:) + & + a(ni,ni-6) * var(ni,:,:) ) * idel + rhs(ni-5,:,:) = rhs(ni-5,:,:) + & + ( a(ni-10,ni-5) * var(ni-10,:,:) + & + a(ni-9,ni-5) * var(ni-9,:,:) + & + a(ni-8,ni-5) * var(ni-8,:,:) + & + a(ni-7,ni-5) * var(ni-7,:,:) + & + a(ni-6,ni-5) * var(ni-6,:,:) + & + a(ni-5,ni-5) * var(ni-5,:,:) + & + a(ni-4,ni-5) * var(ni-4,:,:) + & + a(ni-3,ni-5) * var(ni-3,:,:) + & + a(ni-2,ni-5) * var(ni-2,:,:) + & + a(ni-1,ni-5) * var(ni-1,:,:) + & + a(ni,ni-5) * var(ni,:,:) ) * idel + rhs(ni-4,:,:) = rhs(ni-4,:,:) + & + ( a(ni-10,ni-4) * var(ni-10,:,:) + & + a(ni-9,ni-4) * var(ni-9,:,:) + & + a(ni-8,ni-4) * var(ni-8,:,:) + & + a(ni-7,ni-4) * var(ni-7,:,:) + & + a(ni-6,ni-4) * var(ni-6,:,:) + & + a(ni-5,ni-4) * var(ni-5,:,:) + & + a(ni-4,ni-4) * var(ni-4,:,:) + & + a(ni-3,ni-4) * var(ni-3,:,:) + & + a(ni-2,ni-4) * var(ni-2,:,:) + & + a(ni-1,ni-4) * var(ni-1,:,:) + & + a(ni,ni-4) * var(ni,:,:) ) * idel + rhs(ni-3,:,:) = rhs(ni-3,:,:) + & + ( a(ni-10,ni-3) * var(ni-10,:,:) + & + a(ni-9,ni-3) * var(ni-9,:,:) + & + a(ni-8,ni-3) * var(ni-8,:,:) + & + a(ni-7,ni-3) * var(ni-7,:,:) + & + a(ni-6,ni-3) * var(ni-6,:,:) + & + a(ni-5,ni-3) * var(ni-5,:,:) + & + a(ni-4,ni-3) * var(ni-4,:,:) + & + a(ni-3,ni-3) * var(ni-3,:,:) + & + a(ni-2,ni-3) * var(ni-2,:,:) + & + a(ni-1,ni-3) * var(ni-1,:,:) + & + a(ni,ni-3) * var(ni,:,:) ) * idel + rhs(ni-2,:,:) = rhs(ni-2,:,:) + & + ( a(ni-10,ni-2) * var(ni-10,:,:) + & + a(ni-9,ni-2) * var(ni-9,:,:) + & + a(ni-8,ni-2) * var(ni-8,:,:) + & + a(ni-7,ni-2) * var(ni-7,:,:) + & + a(ni-6,ni-2) * var(ni-6,:,:) + & + a(ni-5,ni-2) * var(ni-5,:,:) + & + a(ni-4,ni-2) * var(ni-4,:,:) + & + a(ni-3,ni-2) * var(ni-3,:,:) + & + a(ni-2,ni-2) * var(ni-2,:,:) + & + a(ni-1,ni-2) * var(ni-1,:,:) + & + a(ni,ni-2) * var(ni,:,:) ) * idel + rhs(ni-1,:,:) = rhs(ni-1,:,:) + & + ( a(ni-10,ni-1) * var(ni-10,:,:) + & + a(ni-9,ni-1) * var(ni-9,:,:) + & + a(ni-8,ni-1) * var(ni-8,:,:) + & + a(ni-7,ni-1) * var(ni-7,:,:) + & + a(ni-6,ni-1) * var(ni-6,:,:) + & + a(ni-5,ni-1) * var(ni-5,:,:) + & + a(ni-4,ni-1) * var(ni-4,:,:) + & + a(ni-3,ni-1) * var(ni-3,:,:) + & + a(ni-2,ni-1) * var(ni-2,:,:) + & + a(ni-1,ni-1) * var(ni-1,:,:) + & + a(ni,ni-1) * var(ni,:,:) ) * idel + rhs(ni,:,:) = rhs(ni,:,:) + & + ( a(ni-4,ni) * var(ni-4,:,:) + & + a(ni-3,ni) * var(ni-3,:,:) + & + a(ni-2,ni) * var(ni-2,:,:) + & + a(ni-1,ni) * var(ni-1,:,:) + & + a(ni,ni) * var(ni,:,:) ) * idel + + ir = ni - 7 + + end if + + do i = il, ir + + rhs(i,:,:) = rhs(i,:,:) + & + ( a(i-4,i) * var(i-4,:,:) + & + a(i-3,i) * var(i-3,:,:) + & + a(i-2,i) * var(i-2,:,:) + & + a(i-1,i) * var(i-1,:,:) + & + a(i,i) * var(i,:,:) + & + a(i+1,i) * var(i+1,:,:) + & + a(i+2,i) * var(i+2,:,:) + & + a(i+3,i) * var(i+3,:,:) + & + a(i+4,i) * var(i+4,:,:) ) * idel + + end do + + deallocate ( a, d, b, h, tmp ) + + allocate ( a(nj,nj), d(nj,nj), b(nj,nj), h(nj,nj), tmp(nj,nj) ) + a = zero; d = zero; b = zero; h = zero; + + call set_bmatrix ( b, bb(3:4), lsh(2), gsh(2), lbnd(2), delta(2), dfl ) + + call set_hmatrix ( h, bb(3:4) ) + + call set_dmatrix ( d, bb(3:4) ) + + a = -transpose ( matmul ( h, matmul ( transpose(d), matmul ( b, d ) ) ) ) + + idel = epsdisl / ( 256.0_wp * delta(2) ) + + if ( bb(3) == 0 ) then + + jl = 1 + gsize(2) + + else + + rhs(:,1,:) = rhs(:,1,:) + & + ( a(1,1) * var(:,1,:) + a(2,1) * var(:,2,:) + & + a(3,1) * var(:,3,:) + a(4,1) * var(:,4,:) + & + a(5,1) * var(:,5,:) ) * idel + rhs(:,2,:) = rhs(:,2,:) + & + ( a(1,2) * var(:,1,:) + a(2,2) * var(:,2,:) + & + a(3,2) * var(:,3,:) + a(4,2) * var(:,4,:) + & + a(5,2) * var(:,5,:) + a(6,2) * var(:,6,:) + & + a(7,2) * var(:,7,:) + a(8,2) * var(:,8,:) + & + a(9,2) * var(:,9,:) + a(10,2) * var(:,10,:) + & + a(11,2) * var(:,11,:) ) * idel + rhs(:,3,:) = rhs(:,3,:) + & + ( a(1,3) * var(:,1,:) + a(2,3) * var(:,2,:) + & + a(3,3) * var(:,3,:) + a(4,3) * var(:,4,:) + & + a(5,3) * var(:,5,:) + a(6,3) * var(:,6,:) + & + a(7,3) * var(:,7,:) + a(8,3) * var(:,8,:) + & + a(9,3) * var(:,9,:) + a(10,3) * var(:,10,:) + & + a(11,3) * var(:,11,:) ) * idel + rhs(:,4,:) = rhs(:,4,:) + & + ( a(1,4) * var(:,1,:) + a(2,4) * var(:,2,:) + & + a(3,4) * var(:,3,:) + a(4,4) * var(:,4,:) + & + a(5,4) * var(:,5,:) + a(6,4) * var(:,6,:) + & + a(7,4) * var(:,7,:) + a(8,4) * var(:,8,:) + & + a(9,4) * var(:,9,:) + a(10,4) * var(:,10,:) + & + a(11,4) * var(:,11,:) ) * idel + rhs(:,5,:) = rhs(:,5,:) + & + ( a(1,5) * var(:,1,:) + a(2,5) * var(:,2,:) + & + a(3,5) * var(:,3,:) + a(4,5) * var(:,4,:) + & + a(5,5) * var(:,5,:) + a(6,5) * var(:,6,:) + & + a(7,5) * var(:,7,:) + a(8,5) * var(:,8,:) + & + a(9,5) * var(:,9,:) + a(10,5) * var(:,10,:) + & + a(11,5) * var(:,11,:) ) * idel + rhs(:,6,:) = rhs(:,6,:) + & + ( a(1,6) * var(:,1,:) + a(2,6) * var(:,2,:) + & + a(3,6) * var(:,3,:) + a(4,6) * var(:,4,:) + & + a(5,6) * var(:,5,:) + a(6,6) * var(:,6,:) + & + a(7,6) * var(:,7,:) + a(8,6) * var(:,8,:) + & + a(9,6) * var(:,9,:) + a(10,6) * var(:,10,:) + & + a(11,6) * var(:,11,:) ) * idel + rhs(:,7,:) = rhs(:,7,:) + & + ( a(1,7) * var(:,1,:) + a(2,7) * var(:,2,:) + & + a(3,7) * var(:,3,:) + a(4,7) * var(:,4,:) + & + a(5,7) * var(:,5,:) + a(6,7) * var(:,6,:) + & + a(7,7) * var(:,7,:) + a(8,7) * var(:,8,:) + & + a(9,7) * var(:,9,:) + a(10,7) * var(:,10,:) + & + a(11,7) * var(:,11,:) ) * idel + + jl = 8 + + end if + + if ( bb(4) == 0 ) then + + jr = nj - gsize(2) + + else + + rhs(:,nj-6,:) = rhs(:,nj-6,:) + & + ( a(nj-10,nj-6) * var(:,nj-10,:) + & + a(nj-9,nj-6) * var(:,nj-9,:) + & + a(nj-8,nj-6) * var(:,nj-8,:) + & + a(nj-7,nj-6) * var(:,nj-7,:) + & + a(nj-6,nj-6) * var(:,nj-6,:) + & + a(nj-5,nj-6) * var(:,nj-5,:) + & + a(nj-4,nj-6) * var(:,nj-4,:) + & + a(nj-3,nj-6) * var(:,nj-3,:) + & + a(nj-2,nj-6) * var(:,nj-2,:) + & + a(nj-1,nj-6) * var(:,nj-1,:) + & + a(nj,nj-6) * var(:,nj,:) ) * idel + rhs(:,nj-5,:) = rhs(:,nj-5,:) + & + ( a(nj-10,nj-5) * var(:,nj-10,:) + & + a(nj-9,nj-5) * var(:,nj-9,:) + & + a(nj-8,nj-5) * var(:,nj-8,:) + & + a(nj-7,nj-5) * var(:,nj-7,:) + & + a(nj-6,nj-5) * var(:,nj-6,:) + & + a(nj-5,nj-5) * var(:,nj-5,:) + & + a(nj-4,nj-5) * var(:,nj-4,:) + & + a(nj-3,nj-5) * var(:,nj-3,:) + & + a(nj-2,nj-5) * var(:,nj-2,:) + & + a(nj-1,nj-5) * var(:,nj-1,:) + & + a(nj,nj-5) * var(:,nj,:) ) * idel + rhs(:,nj-4,:) = rhs(:,nj-4,:) + & + ( a(nj-10,nj-4) * var(:,nj-10,:) + & + a(nj-9,nj-4) * var(:,nj-9,:) + & + a(nj-8,nj-4) * var(:,nj-8,:) + & + a(nj-7,nj-4) * var(:,nj-7,:) + & + a(nj-6,nj-4) * var(:,nj-6,:) + & + a(nj-5,nj-4) * var(:,nj-5,:) + & + a(nj-4,nj-4) * var(:,nj-4,:) + & + a(nj-3,nj-4) * var(:,nj-3,:) + & + a(nj-2,nj-4) * var(:,nj-2,:) + & + a(nj-1,nj-4) * var(:,nj-1,:) + & + a(nj,nj-4) * var(:,nj,:) ) * idel + rhs(:,nj-3,:) = rhs(:,nj-3,:) + & + ( a(nj-10,nj-3) * var(:,nj-10,:) + & + a(nj-9,nj-3) * var(:,nj-9,:) + & + a(nj-8,nj-3) * var(:,nj-8,:) + & + a(nj-7,nj-3) * var(:,nj-7,:) + & + a(nj-6,nj-3) * var(:,nj-6,:) + & + a(nj-5,nj-3) * var(:,nj-5,:) + & + a(nj-4,nj-3) * var(:,nj-4,:) + & + a(nj-3,nj-3) * var(:,nj-3,:) + & + a(nj-2,nj-3) * var(:,nj-2,:) + & + a(nj-1,nj-3) * var(:,nj-1,:) + & + a(nj,nj-3) * var(:,nj,:) ) * idel + rhs(:,nj-2,:) = rhs(:,nj-2,:) + & + ( a(nj-10,nj-2) * var(:,nj-10,:) + & + a(nj-9,nj-2) * var(:,nj-9,:) + & + a(nj-8,nj-2) * var(:,nj-8,:) + & + a(nj-7,nj-2) * var(:,nj-7,:) + & + a(nj-6,nj-2) * var(:,nj-6,:) + & + a(nj-5,nj-2) * var(:,nj-5,:) + & + a(nj-4,nj-2) * var(:,nj-4,:) + & + a(nj-3,nj-2) * var(:,nj-3,:) + & + a(nj-2,nj-2) * var(:,nj-2,:) + & + a(nj-1,nj-2) * var(:,nj-1,:) + & + a(nj,nj-2) * var(:,nj,:) ) * idel + rhs(:,nj-1,:) = rhs(:,nj-1,:) + & + ( a(nj-10,nj-1) * var(:,nj-10,:) + & + a(nj-9,nj-1) * var(:,nj-9,:) + & + a(nj-8,nj-1) * var(:,nj-8,:) + & + a(nj-7,nj-1) * var(:,nj-7,:) + & + a(nj-6,nj-1) * var(:,nj-6,:) + & + a(nj-5,nj-1) * var(:,nj-5,:) + & + a(nj-4,nj-1) * var(:,nj-4,:) + & + a(nj-3,nj-1) * var(:,nj-3,:) + & + a(nj-2,nj-1) * var(:,nj-2,:) + & + a(nj-1,nj-1) * var(:,nj-1,:) + & + a(nj,nj-1) * var(:,nj,:) ) * idel + rhs(:,nj,:) = rhs(:,nj,:) + & + ( a(nj-4,nj) * var(:,nj-4,:) + & + a(nj-3,nj) * var(:,nj-3,:) + & + a(nj-2,nj) * var(:,nj-2,:) + & + a(nj-1,nj) * var(:,nj-1,:) + & + a(nj,nj) * var(:,nj,:) ) * idel + + jr = nj - 7 + + end if + + do j = jl, jr + + rhs(:,j,:) = rhs(:,j,:) + & + ( a(j-4,j) * var(:,j-4,:) + & + a(j-3,j) * var(:,j-3,:) + & + a(j-2,j) * var(:,j-2,:) + & + a(j-1,j) * var(:,j-1,:) + & + a(j,j) * var(:,j,:) + & + a(j+1,j) * var(:,j+1,:) + & + a(j+2,j) * var(:,j+2,:) + & + a(j+3,j) * var(:,j+3,:) + & + a(j+4,j) * var(:,j+4,:) ) * idel + + end do + + deallocate ( a, d, b, h, tmp ) + + allocate ( a(nk,nk), d(nk,nk), b(nk,nk), h(nk,nk), tmp(nk,nk) ) + a = zero; d = zero; b = zero; h = zero; tmp = zero + + call set_bmatrix ( b, bb(5:6), lsh(3), gsh(3), lbnd(3), delta(3), dfl ) + + call set_hmatrix ( h, bb(5:6) ) + + call set_dmatrix ( d, bb(5:6) ) + + a = -transpose ( matmul ( h, matmul ( transpose(d), matmul ( b, d ) ) ) ) + + idel = epsdisl / ( 256.0_wp * delta(3) ) + + if ( bb(5) == 0 ) then + + kl = 1 + gsize(3) + + else + + rhs(:,:,1) = rhs(:,:,1) + & + ( a(1,1) * var(:,:,1) + a(2,1) * var(:,:,2) + & + a(3,1) * var(:,:,3) + a(4,1) * var(:,:,4) + & + a(5,1) * var(:,:,5) ) * idel + rhs(:,:,2) = rhs(:,:,2) + & + ( a(1,2) * var(:,:,1) + a(2,2) * var(:,:,2) + & + a(3,2) * var(:,:,3) + a(4,2) * var(:,:,4) + & + a(5,2) * var(:,:,5) + a(6,2) * var(:,:,6) + & + a(7,2) * var(:,:,7) + a(8,2) * var(:,:,8) + & + a(9,2) * var(:,:,9) + a(10,2) * var(:,:,10) + & + a(11,2) * var(:,:,11) ) * idel + rhs(:,:,3) = rhs(:,:,3) + & + ( a(1,3) * var(:,:,1) + a(2,3) * var(:,:,2) + & + a(3,3) * var(:,:,3) + a(4,3) * var(:,:,4) + & + a(5,3) * var(:,:,5) + a(6,3) * var(:,:,6) + & + a(7,3) * var(:,:,7) + a(8,3) * var(:,:,8) + & + a(9,3) * var(:,:,9) + a(10,3) * var(:,:,10) + & + a(11,3) * var(:,:,11) ) * idel + rhs(:,:,4) = rhs(:,:,4) + & + ( a(1,4) * var(:,:,1) + a(2,4) * var(:,:,2) + & + a(3,4) * var(:,:,3) + a(4,4) * var(:,:,4) + & + a(5,4) * var(:,:,5) + a(6,4) * var(:,:,6) + & + a(7,4) * var(:,:,7) + a(8,4) * var(:,:,8) + & + a(9,4) * var(:,:,9) + a(10,4) * var(:,:,10) + & + a(11,4) * var(:,:,11) ) * idel + rhs(:,:,5) = rhs(:,:,5) + & + ( a(1,5) * var(:,:,1) + a(2,5) * var(:,:,2) + & + a(3,5) * var(:,:,3) + a(4,5) * var(:,:,4) + & + a(5,5) * var(:,:,5) + a(6,5) * var(:,:,6) + & + a(7,5) * var(:,:,7) + a(8,5) * var(:,:,8) + & + a(9,5) * var(:,:,9) + a(10,5) * var(:,:,10) + & + a(11,5) * var(:,:,11) ) * idel + rhs(:,:,6) = rhs(:,:,6) + & + ( a(1,6) * var(:,:,1) + a(2,6) * var(:,:,2) + & + a(3,6) * var(:,:,3) + a(4,6) * var(:,:,4) + & + a(5,6) * var(:,:,5) + a(6,6) * var(:,:,6) + & + a(7,6) * var(:,:,7) + a(8,6) * var(:,:,8) + & + a(9,6) * var(:,:,9) + a(10,6) * var(:,:,10) + & + a(11,6) * var(:,:,11) ) * idel + rhs(:,:,7) = rhs(:,:,7) + & + ( a(1,7) * var(:,:,1) + a(2,7) * var(:,:,2) + & + a(3,7) * var(:,:,3) + a(4,7) * var(:,:,4) + & + a(5,7) * var(:,:,5) + a(6,7) * var(:,:,6) + & + a(7,7) * var(:,:,7) + a(8,7) * var(:,:,8) + & + a(9,7) * var(:,:,9) + a(10,7) * var(:,:,10) + & + a(11,7) * var(:,:,11) ) * idel + + kl = 8 + + end if + + if ( bb(6) == 0 ) then + + kr = nk - gsize(3) + + else + + rhs(:,:,nk-6) = rhs(:,:,nk-6) + & + ( a(nk-10,nk-6) * var(:,:,nk-10) + & + a(nk-9,nk-6) * var(:,:,nk-9) + & + a(nk-8,nk-6) * var(:,:,nk-8) + & + a(nk-7,nk-6) * var(:,:,nk-7) + & + a(nk-6,nk-6) * var(:,:,nk-6) + & + a(nk-5,nk-6) * var(:,:,nk-5) + & + a(nk-4,nk-6) * var(:,:,nk-4) + & + a(nk-3,nk-6) * var(:,:,nk-3) + & + a(nk-2,nk-6) * var(:,:,nk-2) + & + a(nk-1,nk-6) * var(:,:,nk-1) + & + a(nk,nk-6) * var(:,:,nk) ) * idel + rhs(:,:,nk-5) = rhs(:,:,nk-5) + & + ( a(nk-10,nk-5) * var(:,:,nk-10) + & + a(nk-9,nk-5) * var(:,:,nk-9) + & + a(nk-8,nk-5) * var(:,:,nk-8) + & + a(nk-7,nk-5) * var(:,:,nk-7) + & + a(nk-6,nk-5) * var(:,:,nk-6) + & + a(nk-5,nk-5) * var(:,:,nk-5) + & + a(nk-4,nk-5) * var(:,:,nk-4) + & + a(nk-3,nk-5) * var(:,:,nk-3) + & + a(nk-2,nk-5) * var(:,:,nk-2) + & + a(nk-1,nk-5) * var(:,:,nk-1) + & + a(nk,nk-5) * var(:,:,nk) ) * idel + rhs(:,:,nk-4) = rhs(:,:,nk-4) + & + ( a(nk-10,nk-4) * var(:,:,nk-10) + & + a(nk-9,nk-4) * var(:,:,nk-9) + & + a(nk-8,nk-4) * var(:,:,nk-8) + & + a(nk-7,nk-4) * var(:,:,nk-7) + & + a(nk-6,nk-4) * var(:,:,nk-6) + & + a(nk-5,nk-4) * var(:,:,nk-5) + & + a(nk-4,nk-4) * var(:,:,nk-4) + & + a(nk-3,nk-4) * var(:,:,nk-3) + & + a(nk-2,nk-4) * var(:,:,nk-2) + & + a(nk-1,nk-4) * var(:,:,nk-1) + & + a(nk,nk-4) * var(:,:,nk) ) * idel + rhs(:,:,nk-3) = rhs(:,:,nk-3) + & + ( a(nk-10,nk-3) * var(:,:,nk-10) + & + a(nk-9,nk-3) * var(:,:,nk-9) + & + a(nk-8,nk-3) * var(:,:,nk-8) + & + a(nk-7,nk-3) * var(:,:,nk-7) + & + a(nk-6,nk-3) * var(:,:,nk-6) + & + a(nk-5,nk-3) * var(:,:,nk-5) + & + a(nk-4,nk-3) * var(:,:,nk-4) + & + a(nk-3,nk-3) * var(:,:,nk-3) + & + a(nk-2,nk-3) * var(:,:,nk-2) + & + a(nk-1,nk-3) * var(:,:,nk-1) + & + a(nk,nk-3) * var(:,:,nk) ) * idel + rhs(:,:,nk-2) = rhs(:,:,nk-2) + & + ( a(nk-10,nk-2) * var(:,:,nk-10) + & + a(nk-9,nk-2) * var(:,:,nk-9) + & + a(nk-8,nk-2) * var(:,:,nk-8) + & + a(nk-7,nk-2) * var(:,:,nk-7) + & + a(nk-6,nk-2) * var(:,:,nk-6) + & + a(nk-5,nk-2) * var(:,:,nk-5) + & + a(nk-4,nk-2) * var(:,:,nk-4) + & + a(nk-3,nk-2) * var(:,:,nk-3) + & + a(nk-2,nk-2) * var(:,:,nk-2) + & + a(nk-1,nk-2) * var(:,:,nk-1) + & + a(nk,nk-2) * var(:,:,nk) ) * idel + rhs(:,:,nk-1) = rhs(:,:,nk-1) + & + ( a(nk-10,nk-1) * var(:,:,nk-10) + & + a(nk-9,nk-1) * var(:,:,nk-9) + & + a(nk-8,nk-1) * var(:,:,nk-8) + & + a(nk-7,nk-1) * var(:,:,nk-7) + & + a(nk-6,nk-1) * var(:,:,nk-6) + & + a(nk-5,nk-1) * var(:,:,nk-5) + & + a(nk-4,nk-1) * var(:,:,nk-4) + & + a(nk-3,nk-1) * var(:,:,nk-3) + & + a(nk-2,nk-1) * var(:,:,nk-2) + & + a(nk-1,nk-1) * var(:,:,nk-1) + & + a(nk,nk-1) * var(:,:,nk) ) * idel + rhs(:,:,nk) = rhs(:,:,nk) + & + ( a(nk-4,nk) * var(:,:,nk-4) + & + a(nk-3,nk) * var(:,:,nk-3) + & + a(nk-2,nk) * var(:,:,nk-2) + & + a(nk-1,nk) * var(:,:,nk-1) + & + a(nk,nk) * var(:,:,nk) ) * idel + + kr = nk - 7 + + end if + + do k = kl, kr + + rhs(:,:,k) = rhs(:,:,k) + & + ( a(k-4,k) * var(:,:,k-4) + & + a(k-3,k) * var(:,:,k-3) + & + a(k-2,k) * var(:,:,k-2) + & + a(k-1,k) * var(:,:,k-1) + & + a(k,k) * var(:,:,k) + & + a(k+1,k) * var(:,:,k+1) + & + a(k+2,k) * var(:,:,k+2) + & + a(k+3,k) * var(:,:,k+3) + & + a(k+4,k) * var(:,:,k+4) ) * idel + + end do + + deallocate ( a, d, b, h, tmp ) + +contains + + subroutine set_dmatrix ( d, bb ) + + implicit none + + CCTK_REAL, dimension(:,:), intent(out) :: d + CCTK_INT, dimension(2), intent(in) :: bb + CCTK_INT :: n + CCTK_REAL, dimension(5), save :: ac = (/ 1.0_wp, -4.0_wp, 6.0_wp, & + -4.0_wp, 1.0_wp /) + CCTK_INT :: i + + d = zero + + n = size(d,1) + + if ( bb(1) == 0 ) then + d(1,1:3) = ac(3:5) + d(2,1:4) = ac(2:5) + else + d(1,1:5) = ac + d(2,1:5) = ac + end if + if ( bb(2) == 0 ) then + d(n-1,n-3:n) = ac(1:4) + d(n,n-2:n) = ac(1:3) + else + d(n-1,n-4:n) = ac + d(n,n-4:n) = ac + end if + do i = 3, n-2 + d(i,i-2:i+2) = ac + end do + end subroutine set_dmatrix + + + subroutine set_bmatrix ( b, bb, lsh, gsh, lbnd, h, dfl ) + + implicit none + + CCTK_REAL, dimension(:,:), intent(out) :: b + CCTK_INT, dimension(2), intent(in) :: bb + CCTK_INT, intent(in) :: lsh, gsh, lbnd + CCTK_REAL, intent(in) :: h, dfl + CCTK_INT :: n, nwt + CCTK_INT :: i + CCTK_REAL :: xp, xmax + + CCTK_REAL :: zero = 0.0 + integer, parameter :: wp = kind(zero) + + n = size(b,1) + xp = ( gsh - 1 ) * dfl - 1 + + do i = 1, n + b(i,i) = bf ( real(i+lbnd-1,wp), xp, real(gsh-1,wp), h ) + end do + + end subroutine set_bmatrix + + + subroutine set_hmatrix ( h, bb ) + + implicit none + + CCTK_REAL, dimension(:,:), intent(out) :: h + CCTK_INT, dimension(2), intent(in) :: bb + CCTK_INT :: n, i + + CCTK_REAL :: zero = 0.0 + integer, parameter :: wp = kind(zero) + + n = size(h,1) + + h = zero + + do i = 1, n + h(i,i) = 1.0_wp + end do + + if ( bb(1) /= 0 ) then + h(1:7,1:7) = sigma ( ) + end if + + if ( bb(2) /= 0 ) then + h(n:n-6:-1,n:n-6:-1) = sigma ( ) + end if + + end subroutine set_hmatrix + + + function bf ( x, xp, xmax, h ) + + implicit none + + CCTK_REAL :: bf + CCTK_REAL, intent(in) :: x, xp, xmax, h + CCTK_REAL :: h2, h21, xp2, xp3 + + CCTK_REAL :: zero = 0.0 + integer, parameter :: wp = kind(zero) + + if ( x < xp ) then + + h2 = h**2 + h21 = 1.0_wp - h2 + xp2 = 1.0_wp / xp**2 + xp3 = xp2 / xp + bf = ( -2.0_wp * h21 * xp3 * x + 3.0_wp * h21 * xp2 ) * x**2 + h2 + + else if ( x > xmax - xp ) then + + h2 = h**2 + h21 = 1.0_wp - h2 + xp2 = 1.0_wp / xp**2 + xp3 = xp2 / xp + bf = ( 2.0_wp * h21 * xp3 * ( x - xmax ) + & + 3.0_wp * h21 * xp2 ) * ( x - xmax )**2 + h2 + + else + + bf = 1.0_wp + + end if + + end function bf + + + function sigma ( ) + + implicit none + + CCTK_REAL, dimension(7,7) :: sigma + + CCTK_REAL :: zero = 0.0 + integer, parameter :: wp = kind(zero) + + sigma(1,1) = 4.930709842221048047321555312222552006914_wp + sigma(2,1) = 0_wp + sigma(3,1) = 0_wp + sigma(4,1) = 0_wp + sigma(5,1) = 0_wp + sigma(6,1) = 0_wp + sigma(7,1) = 0_wp + sigma(1,2) = 0_wp + sigma(2,2) = 1.499128158967098993672808967791430710830_wp + sigma(3,2) = 0.8363320637711490458766332755492551638364_wp + sigma(4,2) = -0.2634194189384563273844251938374368097488_wp + sigma(5,2) = -0.1281142262970477369698649073120970227895_wp + sigma(6,2) = -0.1830121923043950334684235233648902467425_wp + sigma(7,2) = 0.01275879120879511857581469566228951245742_wp + sigma(1,3) = 0_wp + sigma(2,3) = 0.8363320637711490458766332755492551638364_wp + sigma(3,3) = 0.8785380112113613457715293301885586949116_wp + sigma(4,3) = 0.1454039706555408804118493444561963110723_wp + sigma(5,3) = -0.03346859376328961536059718572568704595897_wp + sigma(6,3) = -0.1759548942026651677201193466207997400242_wp + sigma(7,3) = 0.02020404709761823823384626194775259449799_wp + sigma(1,4) = 0_wp + sigma(2,4) = -0.2634194189384563273844251938374368097488_wp + sigma(3,4) = 0.1454039706555408804118493444561963110723_wp + sigma(4,4) = 0.5690707861055695191099807759535994201984_wp + sigma(5,4) = 0.3374169937858825026860519890688362591356_wp + sigma(6,4) = -0.03076801586199402080957339555604653302670_wp + sigma(7,4) = 0.004077760901760668929595481886645691229036_wp + sigma(1,5) = 0_wp + sigma(2,5) = -0.1281142262970477369698649073120970227895_wp + sigma(3,5) = -0.03346859376328961536059718572568704595897_wp + sigma(4,5) = 0.3374169937858825026860519890688362591356_wp + sigma(5,5) = 0.5001980842834785426891480539161939116986_wp + sigma(6,5) = 0.2904403943485980286998891564955850119600_wp + sigma(7,5) = -0.05655813410599246578174268632557910837116_wp + sigma(1,6) = 0_wp + sigma(2,6) = -0.1830121923043950334684235233648902467425_wp + sigma(3,6) = -0.1759548942026651677201193466207997400242_wp + sigma(4,6) = -0.03076801586199402080957339555604653302670_wp + sigma(5,6) = 0.2904403943485980286998891564955850119600_wp + sigma(6,6) = 0.8640467497799918355533313312208192675134_wp + sigma(7,6) = 0.03704582714017916999084939064121697271831_wp + sigma(1,7) = 0_wp + sigma(2,7) = 0.01275879120879511857581469566228951245742_wp + sigma(3,7) = 0.02020404709761823823384626194775259449799_wp + sigma(4,7) = 0.004077760901760668929595481886645691229036_wp + sigma(5,7) = -0.05655813410599246578174268632557910837116_wp + sigma(6,7) = 0.03704582714017916999084939064121697271831_wp + sigma(7,7) = 0.9909401061257062767072573096147576992963_wp + + end function sigma + +end subroutine dissipation_6_5_alt diff --git a/src/Dissipation_8_4_alt.F90 b/src/Dissipation_8_4_alt.F90 new file mode 100644 index 0000000..594b531 --- /dev/null +++ b/src/Dissipation_8_4_alt.F90 @@ -0,0 +1,513 @@ +! $Header$ + +#include "cctk.h" +#include "cctk_Parameters.h" + +subroutine dissipation_8_4_alt (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) :: epsilon + + CCTK_REAL :: zero = 0.0 + integer, parameter :: wp = kind(zero) + CCTK_REAL, dimension(13,8) :: a + CCTK_REAL :: idel + + CCTK_INT :: il, ir, jl, jr, kl, kr + + call set_coeff ( a ) + + idel = epsilon / ( 1024 * delta(1) ) + + if ( bb(1) == 0 ) then + il = 1 + gsize(1) + else + rhs(1,:,:) = rhs(1,:,:) + & + ( a(1,1) * var(1,:,:) + a(2,1) * var(2,:,:) + & + a(3,1) * var(3,:,:) + a(4,1) * var(4,:,:) + & + a(5,1) * var(5,:,:) + a(6,1) * var(6,:,:) ) * idel + rhs(2,:,:) = rhs(2,:,:) + & + ( a(1,2) * var(1,:,:) + a(2,2) * var(2,:,:) + & + a(3,2) * var(3,:,:) + a(4,2) * var(4,:,:) + & + a(5,2) * var(5,:,:) + a(6,2) * var(6,:,:) + & + a(7,2) * var(7,:,:) ) * idel + rhs(3,:,:) = rhs(3,:,:) + & + ( a(1,3) * var(1,:,:) + a(2,3) * var(2,:,:) + & + a(3,3) * var(3,:,:) + a(4,3) * var(4,:,:) + & + a(5,3) * var(5,:,:) + a(6,3) * var(6,:,:) + & + a(7,3) * var(7,:,:) + a(8,3) * var(8,:,:) ) * idel + rhs(4,:,:) = rhs(4,:,:) + & + ( a(1,4) * var(1,:,:) + a(2,4) * var(2,:,:) + & + a(3,4) * var(3,:,:) + a(4,4) * var(4,:,:) + & + a(5,4) * var(5,:,:) + a(6,4) * var(6,:,:) + & + a(7,4) * var(7,:,:) + a(8,4) * var(8,:,:) + & + a(9,4) * var(9,:,:) ) * idel + rhs(5,:,:) = rhs(5,:,:) + & + ( a(1,5) * var(1,:,:) + a(2,5) * var(2,:,:) + & + a(3,5) * var(3,:,:) + a(4,5) * var(4,:,:) + & + a(5,5) * var(5,:,:) + a(6,5) * var(6,:,:) + & + a(7,5) * var(7,:,:) + a(8,5) * var(8,:,:) + & + a(9,5) * var(9,:,:) + a(10,5) * var(10,:,:) ) * idel + rhs(6,:,:) = rhs(6,:,:) + & + ( a(1,6) * var(1,:,:) + a(2,6) * var(2,:,:) + & + a(3,6) * var(3,:,:) + a(4,6) * var(4,:,:) + & + a(5,6) * var(5,:,:) + a(6,6) * var(6,:,:) + & + a(7,6) * var(7,:,:) + a(8,6) * var(8,:,:) + & + a(9,6) * var(9,:,:) + a(10,6) * var(10,:,:) + & + a(11,6) * var(11,:,:) ) * idel + rhs(7,:,:) = rhs(7,:,:) + & + ( a(2,7) * var(2,:,:) + a(3,7) * var(3,:,:) + & + a(4,7) * var(4,:,:) + a(5,7) * var(5,:,:) + & + a(6,7) * var(6,:,:) + a(7,7) * var(7,:,:) + & + a(8,7) * var(8,:,:) + a(9,7) * var(9,:,:) + & + a(10,7) * var(10,:,:) + a(11,7) * var(11,:,:) + & + a(12,7) * var(12,:,:) ) * idel + rhs(8,:,:) = rhs(8,:,:) + & + ( a(3,8) * var(3,:,:) + a(4,8) * var(4,:,:) + & + a(5,8) * var(5,:,:) + a(6,8) * var(6,:,:) + & + a(7,8) * var(7,:,:) + a(8,8) * var(8,:,:) + & + a(9,8) * var(9,:,:) + a(10,8) * var(10,:,:) + & + a(11,8) * var(11,:,:) + a(12,8) * var(12,:,:) + & + a(13,8) * var(13,:,:) ) * idel + + il = 9 + end if + if ( bb(2) == 0 ) then + ir = ni - gsize(1) + else + rhs(ni-7,:,:) = rhs(ni-7,:,:) + & + ( a(3,8) * var(ni-2,:,:) + a(4,8) * var(ni-3,:,:) + & + a(5,8) * var(ni-4,:,:) + a(6,8) * var(ni-5,:,:) + & + a(7,8) * var(ni-6,:,:) + a(8,8) * var(ni-7,:,:) + & + a(9,8) * var(ni-8,:,:) + a(10,8) * var(ni-9,:,:) + & + a(11,8) * var(ni-10,:,:) + a(12,8) * var(ni-11,:,:) + & + a(13,8) * var(ni-12,:,:) ) * idel + rhs(ni-6,:,:) = rhs(ni-6,:,:) + & + ( a(2,7) * var(ni-1,:,:) + a(3,7) * var(ni-2,:,:) + & + a(4,7) * var(ni-3,:,:) + a(5,7) * var(ni-4,:,:) + & + a(6,7) * var(ni-5,:,:) + a(7,7) * var(ni-6,:,:) + & + a(8,7) * var(ni-7,:,:) + a(9,7) * var(ni-8,:,:) + & + a(10,7) * var(ni-9,:,:) + a(11,7) * var(ni-10,:,:) + & + a(12,7) * var(ni-11,:,:) ) * idel + rhs(ni-5,:,:) = rhs(ni-5,:,:) + & + ( a(1,6) * var(ni,:,:) + a(2,6) * var(ni-1,:,:) + & + a(3,6) * var(ni-2,:,:) + a(4,6) * var(ni-3,:,:) + & + a(5,6) * var(ni-4,:,:) + a(6,6) * var(ni-5,:,:) + & + a(7,6) * var(ni-6,:,:) + a(8,6) * var(ni-7,:,:) + & + a(9,6) * var(ni-8,:,:) + a(10,6) * var(ni-9,:,:) + & + a(11,6) * var(ni-10,:,:) ) * idel + rhs(ni-4,:,:) = rhs(ni-4,:,:) + & + ( a(1,5) * var(ni,:,:) + a(2,5) * var(ni-1,:,:) + & + a(3,5) * var(ni-2,:,:) + a(4,5) * var(ni-3,:,:) + & + a(5,5) * var(ni-4,:,:) + a(6,5) * var(ni-5,:,:) + & + a(7,5) * var(ni-6,:,:) + a(8,5) * var(ni-7,:,:) + & + a(9,5) * var(ni-8,:,:) + a(10,5) * var(ni-9,:,:) ) * idel + rhs(ni-3,:,:) = rhs(ni-3,:,:) + & + ( a(1,4) * var(ni,:,:) + a(2,4) * var(ni-1,:,:) + & + a(3,4) * var(ni-2,:,:) + a(4,4) * var(ni-3,:,:) + & + a(5,4) * var(ni-4,:,:) + a(6,4) * var(ni-5,:,:) + & + a(7,4) * var(ni-6,:,:) + a(8,4) * var(ni-7,:,:) + & + a(9,4) * var(ni-8,:,:) ) * idel + rhs(ni-2,:,:) = rhs(ni-2,:,:) + & + ( a(1,3) * var(ni,:,:) + a(2,3) * var(ni-1,:,:) + & + a(3,3) * var(ni-2,:,:) + a(4,3) * var(ni-3,:,:) + & + a(5,3) * var(ni-4,:,:) + a(6,3) * var(ni-5,:,:) + & + a(7,3) * var(ni-6,:,:) + a(8,3) * var(ni-7,:,:) ) * idel + rhs(ni-1,:,:) = rhs(ni-1,:,:) + & + ( a(1,2) * var(ni,:,:) + a(2,2) * var(ni-1,:,:) + & + a(3,2) * var(ni-2,:,:) + a(4,2) * var(ni-3,:,:) + & + a(5,2) * var(ni-4,:,:) + a(6,2) * var(ni-5,:,:) + & + a(7,2) * var(ni-6,:,:) ) * idel + rhs(ni,:,:) = rhs(ni,:,:) + & + ( a(1,1) * var(ni,:,:) + a(2,1) * var(ni-1,:,:) + & + a(3,1) * var(ni-2,:,:) + a(4,1) * var(ni-3,:,:) + & + a(5,1) * var(ni-4,:,:) + a(6,1) * var(ni-5,:,:) ) * idel + + ir = ni - 8 + end if + rhs(il:ir,:,:) = rhs(il:ir,:,:) + & + ( -252.0_wp * var(il:ir,:,:) + & + 210.0_wp * ( var(il-1:ir-1,:,:) + & + var(il+1:ir+1,:,:) ) - & + 120.0_wp * ( var(il-2:ir-2,:,:) + & + var(il+2:ir+2,:,:) ) + & + 45.0_wp * ( var(il-3:ir-3,:,:) + & + var(il+3:ir+3,:,:) ) - & + 10.0_wp * ( var(il-4:ir-4,:,:) + & + var(il+4:ir+4,:,:) ) + & + ( var(il-5:ir-5,:,:) + & + var(il+5:ir+5,:,:) ) ) * idel + + call set_coeff ( a ) + + idel = epsilon / ( 1024 * delta(2) ) + + if ( bb(3) == 0 ) then + jl = 1 + gsize(2) + else + rhs(:,1,:) = rhs(:,1,:) + & + ( a(1,1) * var(:,1,:) + a(2,1) * var(:,2,:) + & + a(3,1) * var(:,3,:) + a(4,1) * var(:,4,:) + & + a(5,1) * var(:,5,:) + a(6,1) * var(:,6,:) ) * idel + rhs(:,2,:) = rhs(:,2,:) + & + ( a(1,2) * var(:,1,:) + a(2,2) * var(:,2,:) + & + a(3,2) * var(:,3,:) + a(4,2) * var(:,4,:) + & + a(5,2) * var(:,5,:) + a(6,2) * var(:,6,:) + & + a(7,2) * var(:,7,:) ) * idel + rhs(:,3,:) = rhs(:,3,:) + & + ( a(1,3) * var(:,1,:) + a(2,3) * var(:,2,:) + & + a(3,3) * var(:,3,:) + a(4,3) * var(:,4,:) + & + a(5,3) * var(:,5,:) + a(6,3) * var(:,6,:) + & + a(7,3) * var(:,7,:) + a(8,3) * var(:,8,:) ) * idel + rhs(:,4,:) = rhs(:,4,:) + & + ( a(1,4) * var(:,1,:) + a(2,4) * var(:,2,:) + & + a(3,4) * var(:,3,:) + a(4,4) * var(:,4,:) + & + a(5,4) * var(:,5,:) + a(6,4) * var(:,6,:) + & + a(7,4) * var(:,7,:) + a(8,4) * var(:,8,:) + & + a(9,4) * var(:,9,:) ) * idel + rhs(:,5,:) = rhs(:,5,:) + & + ( a(1,5) * var(:,1,:) + a(2,5) * var(:,2,:) + & + a(3,5) * var(:,3,:) + a(4,5) * var(:,4,:) + & + a(5,5) * var(:,5,:) + a(6,5) * var(:,6,:) + & + a(7,5) * var(:,7,:) + a(8,5) * var(:,8,:) + & + a(9,5) * var(:,9,:) + a(10,5) * var(:,10,:) ) * idel + rhs(:,6,:) = rhs(:,6,:) + & + ( a(1,6) * var(:,1,:) + a(2,6) * var(:,2,:) + & + a(3,6) * var(:,3,:) + a(4,6) * var(:,4,:) + & + a(5,6) * var(:,5,:) + a(6,6) * var(:,6,:) + & + a(7,6) * var(:,7,:) + a(8,6) * var(:,8,:) + & + a(9,6) * var(:,9,:) + a(10,6) * var(:,10,:) + & + a(11,6) * var(:,11,:) ) * idel + rhs(:,7,:) = rhs(:,7,:) + & + ( a(2,7) * var(:,2,:) + a(3,7) * var(:,3,:) + & + a(4,7) * var(:,4,:) + a(5,7) * var(:,5,:) + & + a(6,7) * var(:,6,:) + a(7,7) * var(:,7,:) + & + a(8,7) * var(:,8,:) + a(9,7) * var(:,9,:) + & + a(10,7) * var(:,10,:) + a(11,7) * var(:,11,:) + & + a(12,7) * var(:,12,:) ) * idel + rhs(:,8,:) = rhs(:,8,:) + & + ( a(3,8) * var(:,3,:) + a(4,8) * var(:,4,:) + & + a(5,8) * var(:,5,:) + a(6,8) * var(:,6,:) + & + a(7,8) * var(:,7,:) + a(8,8) * var(:,8,:) + & + a(9,8) * var(:,9,:) + a(10,8) * var(:,10,:) + & + a(11,8) * var(:,11,:) + a(12,8) * var(:,12,:) + & + a(13,8) * var(:,13,:) ) * idel + + jl = 9 + end if + if ( bb(4) == 0 ) then + jr = nj - gsize(2) + else + rhs(:,nj-7,:) = rhs(:,nj-7,:) + & + ( a(3,8) * var(:,nj-2,:) + a(4,8) * var(:,nj-3,:) + & + a(5,8) * var(:,nj-4,:) + a(6,8) * var(:,nj-5,:) + & + a(7,8) * var(:,nj-6,:) + a(8,8) * var(:,nj-7,:) + & + a(9,8) * var(:,nj-8,:) + a(10,8) * var(:,nj-9,:) + & + a(11,8) * var(:,nj-10,:) + a(12,8) * var(:,nj-11,:) + & + a(13,8) * var(:,nj-12,:) ) * idel + rhs(:,nj-6,:) = rhs(:,nj-6,:) + & + ( a(2,7) * var(:,nj-1,:) + a(3,7) * var(:,nj-2,:) + & + a(4,7) * var(:,nj-3,:) + a(5,7) * var(:,nj-4,:) + & + a(6,7) * var(:,nj-5,:) + a(7,7) * var(:,nj-6,:) + & + a(8,7) * var(:,nj-7,:) + a(9,7) * var(:,nj-8,:) + & + a(10,7) * var(:,nj-9,:) + a(11,7) * var(:,nj-10,:) + & + a(12,7) * var(:,nj-11,:) ) * idel + rhs(:,nj-5,:) = rhs(:,nj-5,:) + & + ( a(1,6) * var(:,nj,:) + a(2,6) * var(:,nj-1,:) + & + a(3,6) * var(:,nj-2,:) + a(4,6) * var(:,nj-3,:) + & + a(5,6) * var(:,nj-4,:) + a(6,6) * var(:,nj-5,:) + & + a(7,6) * var(:,nj-6,:) + a(8,6) * var(:,nj-7,:) + & + a(9,6) * var(:,nj-8,:) + a(10,6) * var(:,nj-9,:) + & + a(11,6) * var(:,nj-10,:) ) * idel + rhs(:,nj-4,:) = rhs(:,nj-4,:) + & + ( a(1,5) * var(:,nj,:) + a(2,5) * var(:,nj-1,:) + & + a(3,5) * var(:,nj-2,:) + a(4,5) * var(:,nj-3,:) + & + a(5,5) * var(:,nj-4,:) + a(6,5) * var(:,nj-5,:) + & + a(7,5) * var(:,nj-6,:) + a(8,5) * var(:,nj-7,:) + & + a(9,5) * var(:,nj-8,:) + a(10,5) * var(:,nj-9,:) ) * idel + rhs(:,nj-3,:) = rhs(:,nj-3,:) + & + ( a(1,4) * var(:,nj,:) + a(2,4) * var(:,nj-1,:) + & + a(3,4) * var(:,nj-2,:) + a(4,4) * var(:,nj-3,:) + & + a(5,4) * var(:,nj-4,:) + a(6,4) * var(:,nj-5,:) + & + a(7,4) * var(:,nj-6,:) + a(8,4) * var(:,nj-7,:) + & + a(9,4) * var(:,nj-8,:) ) * idel + rhs(:,nj-2,:) = rhs(:,nj-2,:) + & + ( a(1,3) * var(:,nj,:) + a(2,3) * var(:,nj-1,:) + & + a(3,3) * var(:,nj-2,:) + a(4,3) * var(:,nj-3,:) + & + a(5,3) * var(:,nj-4,:) + a(6,3) * var(:,nj-5,:) + & + a(7,3) * var(:,nj-6,:) + a(8,3) * var(:,nj-7,:) ) * idel + rhs(:,nj-1,:) = rhs(:,nj-1,:) + & + ( a(1,2) * var(:,nj,:) + a(2,2) * var(:,nj-1,:) + & + a(3,2) * var(:,nj-2,:) + a(4,2) * var(:,nj-3,:) + & + a(5,2) * var(:,nj-4,:) + a(6,2) * var(:,nj-5,:) + & + a(7,2) * var(:,nj-6,:) ) * idel + rhs(:,nj,:) = rhs(:,nj,:) + & + ( a(1,1) * var(:,nj,:) + a(2,1) * var(:,nj-1,:) + & + a(3,1) * var(:,nj-2,:) + a(4,1) * var(:,nj-3,:) + & + a(5,1) * var(:,nj-4,:) + a(6,1) * var(:,nj-5,:) ) * idel + + jr = nj - 8 + end if + rhs(:,jl:jr,:) = rhs(:,jl:jr,:) + & + ( -252.0_wp * var(:,jl:jr,:) + & + 210.0_wp * ( var(:,jl-1:jr-1,:) + & + var(:,jl+1:jr+1,:) ) - & + 120.0_wp * ( var(:,jl-2:jr-2,:) + & + var(:,jl+2:jr+2,:) ) + & + 45.0_wp * ( var(:,jl-3:jr-3,:) + & + var(:,jl+3:jr+3,:) ) - & + 10.0_wp * ( var(:,jl-4:jr-4,:) + & + var(:,jl+4:jr+4,:) ) + & + ( var(:,jl-5:jr-5,:) + & + var(:,jl+5:jr+5,:) ) ) * idel + + call set_coeff ( a ) + + idel = epsilon / ( 1024 * delta(3) ) + + if ( bb(5) == 0 ) then + kl = 1 + gsize(3) + else + rhs(:,:,1) = rhs(:,:,1) + & + ( a(1,1) * var(:,:,1) + a(2,1) * var(:,:,2) + & + a(3,1) * var(:,:,3) + a(4,1) * var(:,:,4) + & + a(5,1) * var(:,:,5) + a(6,1) * var(:,:,6) ) * idel + rhs(:,:,2) = rhs(:,:,2) + & + ( a(1,2) * var(:,:,1) + a(2,2) * var(:,:,2) + & + a(3,2) * var(:,:,3) + a(4,2) * var(:,:,4) + & + a(5,2) * var(:,:,5) + a(6,2) * var(:,:,6) + & + a(7,2) * var(:,:,7) ) * idel + rhs(:,:,3) = rhs(:,:,3) + & + ( a(1,3) * var(:,:,1) + a(2,3) * var(:,:,2) + & + a(3,3) * var(:,:,3) + a(4,3) * var(:,:,4) + & + a(5,3) * var(:,:,5) + a(6,3) * var(:,:,6) + & + a(7,3) * var(:,:,7) + a(8,3) * var(:,:,8) ) * idel + rhs(:,:,4) = rhs(:,:,4) + & + ( a(1,4) * var(:,:,1) + a(2,4) * var(:,:,2) + & + a(3,4) * var(:,:,3) + a(4,4) * var(:,:,4) + & + a(5,4) * var(:,:,5) + a(6,4) * var(:,:,6) + & + a(7,4) * var(:,:,7) + a(8,4) * var(:,:,8) + & + a(9,4) * var(:,:,9) ) * idel + rhs(:,:,5) = rhs(:,:,5) + & + ( a(1,5) * var(:,:,1) + a(2,5) * var(:,:,2) + & + a(3,5) * var(:,:,3) + a(4,5) * var(:,:,4) + & + a(5,5) * var(:,:,5) + a(6,5) * var(:,:,6) + & + a(7,5) * var(:,:,7) + a(8,5) * var(:,:,8) + & + a(9,5) * var(:,:,9) + a(10,5) * var(:,:,10) ) * idel + rhs(:,:,6) = rhs(:,:,6) + & + ( a(1,6) * var(:,:,1) + a(2,6) * var(:,:,2) + & + a(3,6) * var(:,:,3) + a(4,6) * var(:,:,4) + & + a(5,6) * var(:,:,5) + a(6,6) * var(:,:,6) + & + a(7,6) * var(:,:,7) + a(8,6) * var(:,:,8) + & + a(9,6) * var(:,:,9) + a(10,6) * var(:,:,10) + & + a(11,6) * var(:,:,11) ) * idel + rhs(:,:,7) = rhs(:,:,7) + & + ( a(2,7) * var(:,:,2) + a(3,7) * var(:,:,3) + & + a(4,7) * var(:,:,4) + a(5,7) * var(:,:,5) + & + a(6,7) * var(:,:,6) + a(7,7) * var(:,:,7) + & + a(8,7) * var(:,:,8) + a(9,7) * var(:,:,9) + & + a(10,7) * var(:,:,10) + a(11,7) * var(:,:,11) + & + a(12,7) * var(:,:,12) ) * idel + rhs(:,:,8) = rhs(:,:,8) + & + ( a(3,8) * var(:,:,3) + a(4,8) * var(:,:,4) + & + a(5,8) * var(:,:,5) + a(6,8) * var(:,:,6) + & + a(7,8) * var(:,:,7) + a(8,8) * var(:,:,8) + & + a(9,8) * var(:,:,9) + a(10,8) * var(:,:,10) + & + a(11,8) * var(:,:,11) + a(12,8) * var(:,:,12) + & + a(13,8) * var(:,:,13) ) * idel + + kl = 9 + end if + if ( bb(6) == 0 ) then + kr = nk - gsize(3) + else + rhs(:,:,nk-7) = rhs(:,:,nk-7) + & + ( a(3,8) * var(:,:,nk-2) + a(4,8) * var(:,:,nk-3) + & + a(5,8) * var(:,:,nk-4) + a(6,8) * var(:,:,nk-5) + & + a(7,8) * var(:,:,nk-6) + a(8,8) * var(:,:,nk-7) + & + a(9,8) * var(:,:,nk-8) + a(10,8) * var(:,:,nk-9) + & + a(11,8) * var(:,:,nk-10) + a(12,8) * var(:,:,nk-11) + & + a(13,8) * var(:,:,nk-12) ) * idel + rhs(:,:,nk-6) = rhs(:,:,nk-6) + & + ( a(2,7) * var(:,:,nk-1) + a(3,7) * var(:,:,nk-2) + & + a(4,7) * var(:,:,nk-3) + a(5,7) * var(:,:,nk-4) + & + a(6,7) * var(:,:,nk-5) + a(7,7) * var(:,:,nk-6) + & + a(8,7) * var(:,:,nk-7) + a(9,7) * var(:,:,nk-8) + & + a(10,7) * var(:,:,nk-9) + a(11,7) * var(:,:,nk-10) + & + a(12,7) * var(:,:,nk-11) ) * idel + rhs(:,:,nk-5) = rhs(:,:,nk-5) + & + ( a(1,6) * var(:,:,nk) + a(2,6) * var(:,:,nk-1) + & + a(3,6) * var(:,:,nk-2) + a(4,6) * var(:,:,nk-3) + & + a(5,6) * var(:,:,nk-4) + a(6,6) * var(:,:,nk-5) + & + a(7,6) * var(:,:,nk-6) + a(8,6) * var(:,:,nk-7) + & + a(9,6) * var(:,:,nk-8) + a(10,6) * var(:,:,nk-9) + & + a(11,6) * var(:,:,nk-10) ) * idel + rhs(:,:,nk-4) = rhs(:,:,nk-4) + & + ( a(1,5) * var(:,:,nk) + a(2,5) * var(:,:,nk-1) + & + a(3,5) * var(:,:,nk-2) + a(4,5) * var(:,:,nk-3) + & + a(5,5) * var(:,:,nk-4) + a(6,5) * var(:,:,nk-5) + & + a(7,5) * var(:,:,nk-6) + a(8,5) * var(:,:,nk-7) + & + a(9,5) * var(:,:,nk-8) + a(10,5) * var(:,:,nk-9) ) * idel + rhs(:,:,nk-3) = rhs(:,:,nk-3) + & + ( a(1,4) * var(:,:,nk) + a(2,4) * var(:,:,nk-1) + & + a(3,4) * var(:,:,nk-2) + a(4,4) * var(:,:,nk-3) + & + a(5,4) * var(:,:,nk-4) + a(6,4) * var(:,:,nk-5) + & + a(7,4) * var(:,:,nk-6) + a(8,4) * var(:,:,nk-7) + & + a(9,4) * var(:,:,nk-8) ) * idel + rhs(:,:,nk-2) = rhs(:,:,nk-2) + & + ( a(1,3) * var(:,:,nk) + a(2,3) * var(:,:,nk-1) + & + a(3,3) * var(:,:,nk-2) + a(4,3) * var(:,:,nk-3) + & + a(5,3) * var(:,:,nk-4) + a(6,3) * var(:,:,nk-5) + & + a(7,3) * var(:,:,nk-6) + a(8,3) * var(:,:,nk-7) ) * idel + rhs(:,:,nk-1) = rhs(:,:,nk-1) + & + ( a(1,2) * var(:,:,nk) + a(2,2) * var(:,:,nk-1) + & + a(3,2) * var(:,:,nk-2) + a(4,2) * var(:,:,nk-3) + & + a(5,2) * var(:,:,nk-4) + a(6,2) * var(:,:,nk-5) + & + a(7,2) * var(:,:,nk-6) ) * idel + rhs(:,:,nk) = rhs(:,:,nk) + & + ( a(1,1) * var(:,:,nk) + a(2,1) * var(:,:,nk-1) + & + a(3,1) * var(:,:,nk-2) + a(4,1) * var(:,:,nk-3) + & + a(5,1) * var(:,:,nk-4) + a(6,1) * var(:,:,nk-5) ) * idel + + kr = nk - 8 + end if + rhs(:,:,kl:kr) = rhs(:,:,kl:kr) + & + ( -252.0_wp * var(:,:,kl:kr) + & + 210.0_wp * ( var(:,:,kl-1:kr-1) + & + var(:,:,kl+1:kr+1) ) - & + 120.0_wp * ( var(:,:,kl-2:kr-2) + & + var(:,:,kl+2:kr+2) ) + & + 45.0_wp * ( var(:,:,kl-3:kr-3) + & + var(:,:,kl+3:kr+3) ) - & + 10.0_wp * ( var(:,:,kl-4:kr-4) + & + var(:,:,kl+4:kr+4) ) + & + ( var(:,:,kl-5:kr-5) + & + var(:,:,kl+5:kr+5) ) ) * idel + +contains + + subroutine set_coeff ( a ) + + implicit none + + CCTK_REAL, dimension(13,8), intent(OUT) :: a + CCTK_REAL :: zero = 0.0 + integer, parameter :: wp = kind(zero) + + a(1,1) = -3.3910872088637970174997113084967416241083103770745_wp + a(2,1) = 16.955436044318985087498556542483708120541551885372_wp + a(3,1) = -33.910872088637970174997113084967416241083103770745_wp + a(4,1) = 33.910872088637970174997113084967416241083103770745_wp + a(5,1) = -16.955436044318985087498556542483708120541551885372_wp + a(6,1) = 3.3910872088637970174997113084967416241083103770745_wp + a(7,1) = 0.0_wp + a(8,1) = 0.0_wp + a(9,1) = 0.0_wp + a(10,1) = 0.0_wp + a(11,1) = 0.0_wp + a(12,1) = 0.0_wp + a(13,1) = 0.0_wp + a(1,2) = 3.2771399440263630592058029074141137010783820566473_wp + a(2,2) = -17.041127708937087907870175118553391245607586694566_wp + a(3,2) = 36.048539384289993651263831981555250711862202623121_wp + a(4,2) = -39.325679328316356710469634888969364412940584679768_wp + a(5,2) = 22.939979608184541414440620351898795907548674396531_wp + a(6,2) = -6.5542798880527261184116058148282274021567641132947_wp + a(7,2) = 0.65542798880527261184116058148282274021567641132947_wp + a(8,2) = 0.0_wp + a(9,2) = 0.0_wp + a(10,2) = 0.0_wp + a(11,2) = 0.0_wp + a(12,2) = 0.0_wp + a(13,2) = 0.0_wp + a(1,3) = -38.842059631038967294446317614758441308222147295410_wp + a(2,3) = 213.63132797071432011945474688117142719522181012475_wp + a(3,3) = -489.40995135109098791002360194595636048359905592216_wp + a(4,3) = 602.05192428110399306391792302875584027744328307885_wp + a(5,3) = -427.26265594142864023890949376234285439044362024951_wp + a(6,3) = 174.78926833967535282500842926641298588699966282934_wp + a(7,3) = -38.842059631038967294446317614758441308222147295410_wp + a(8,3) = 3.8842059631038967294446317614758441308222147295410_wp + a(9,3) = 0.0_wp + a(10,3) = 0.0_wp + a(11,3) = 0.0_wp + a(12,3) = 0.0_wp + a(13,3) = 0.0_wp + a(1,4) = 5.5613835719414344378807953109542612676331552744485_wp + a(2,4) = -33.368301431648606627284771865725567605798931646691_wp + a(3,4) = 86.201445365092233787152327319791049648313906753952_wp + a(4,4) = -125.68726872587641829610597402756630464850930920254_wp + a(5,4) = 114.00836322479940597655630387456235598647968312619_wp + a(6,4) = -66.736602863297213254569543731451135211597863293382_wp + a(7,4) = 25.026226073736454970463578899294175704349198735018_wp + a(8,4) = -5.5613835719414344378807953109542612676331552744485_wp + a(9,4) = 0.55613835719414344378807953109542612676331552744485_wp + a(10,4) = 0.0_wp + a(11,4) = 0.0_wp + a(12,4) = 0.0_wp + a(13,4) = 0.0_wp + a(1,5) = -12.115101476661536355654081268132755978592914829047_wp + a(2,5) = 84.805710336630754489578568876929291850150403803330_wp + a(3,5) = -266.53223248655379982438978789892063152904412623904_wp + a(4,5) = 496.71916054312299058181733199344299512230950799093_wp + a(5,5) = -608.17809412840912505383487966026435012536432441817_wp + a(6,5) = 508.83426201978452693747141326157575110090242281998_wp + a(7,5) = -290.76243543987687253569795043518614348622995589713_wp + a(8,5) = 109.03591328995382720088673141319480380733623346142_wp + a(9,5) = -24.230202953323072711308162536265511957185829658094_wp + a(10,5) = 2.4230202953323072711308162536265511957185829658094_wp + a(11,5) = 0.0_wp + a(12,5) = 0.0_wp + a(13,5) = 0.0_wp + a(1,6) = 0.78217600900123184961735065035840034142603567514089_wp + a(2,6) = -7.8217600900123184961735065035840034142603567514089_wp + a(3,6) = 35.197920405055433232780779266128015364171605381340_wp + a(4,6) = -93.861121080147821954082078043008040971124281016906_wp + a(5,6) = 164.25696189025868841964363657526407169946749177959_wp + a(6,6) = -197.10835426831042610357236389031688603936099013550_wp + a(7,6) = 164.25696189025868841964363657526407169946749177959_wp + a(8,6) = -93.861121080147821954082078043008040971124281016906_wp + a(9,6) = 35.197920405055433232780779266128015364171605381340_wp + a(10,6) = -7.8217600900123184961735065035840034142603567514089_wp + a(11,6) = 0.78217600900123184961735065035840034142603567514089_wp + a(12,6) = 0.0_wp + a(13,6) = 0.0_wp + a(1,7) = 0.0_wp + a(2,7) = 1.0830767761393601764536458481012280421614377748694_wp + a(3,7) = -10.830767761393601764536458481012280421614377748694_wp + a(4,7) = 48.738454926271207940414063164555261897264699869122_wp + a(5,7) = -129.96921313672322117443750177214736505937253298433_wp + a(6,7) = 227.44612298926563705526562810125788885390193272257_wp + a(7,7) = -272.93534758711876446631875372150946662468231926708_wp + a(8,7) = 227.44612298926563705526562810125788885390193272257_wp + a(9,7) = -129.96921313672322117443750177214736505937253298433_wp + a(10,7) = 48.738454926271207940414063164555261897264699869122_wp + a(11,7) = -10.830767761393601764536458481012280421614377748694_wp + a(12,7) = 1.0830767761393601764536458481012280421614377748694_wp + a(13,7) = 0.0_wp + a(1,8) = 0.0_wp + a(2,8) = 0.0_wp + a(3,8) = 0.99075245444434671889501396229410272246695863420506_wp + a(4,8) = -9.9075245444434671889501396229410272246695863420506_wp + a(5,8) = 44.583860449995602350275628303234622511013138539228_wp + a(6,8) = -118.89029453332160626740167547529232669603503610461_wp + a(7,8) = 208.05801543331281096795293208176157171806131318306_wp + a(8,8) = -249.66961851997537316154351849811388606167357581967_wp + a(9,8) = 208.05801543331281096795293208176157171806131318306_wp + a(10,8) = -118.89029453332160626740167547529232669603503610461_wp + a(11,8) = 44.583860449995602350275628303234622511013138539228_wp + a(12,8) = -9.9075245444434671889501396229410272246695863420506_wp + a(13,8) = 0.99075245444434671889501396229410272246695863420506_wp + + end subroutine set_coeff + +end subroutine dissipation_8_4_alt diff --git a/src/dissipation.c b/src/dissipation.c index 6426431..a1557cc 100644 --- a/src/dissipation.c +++ b/src/dissipation.c @@ -27,6 +27,16 @@ void CCTK_FCALL CCTK_FNAME(dissipation_6_5_opt) (const CCTK_REAL *var, const CCTK_REAL *epsdis, const CCTK_REAL *diss_fraction, CCTK_REAL *rhs); +void CCTK_FCALL CCTK_FNAME(dissipation_6_5_alt) (const CCTK_REAL *var, + const CCTK_INT *ni, + const CCTK_INT *nj, + const CCTK_INT *nk, + const CCTK_INT *bbox, + const CCTK_INT *gsize, + const CCTK_REAL *dx, + const CCTK_REAL *epsdis, + const CCTK_REAL *diss_fraction, + CCTK_REAL *rhs); void CCTK_FCALL CCTK_FNAME(dissipation_2_1) (const CCTK_REAL *var, const CCTK_INT *ni, const CCTK_INT *nj, @@ -63,6 +73,15 @@ void CCTK_FCALL CCTK_FNAME(dissipation_8_4) (const CCTK_REAL *var, const CCTK_REAL *dx, const CCTK_REAL *epsdis, CCTK_REAL *rhs); +void CCTK_FCALL CCTK_FNAME(dissipation_8_4_alt) (const CCTK_REAL *var, + const CCTK_INT *ni, + const CCTK_INT *nj, + const CCTK_INT *nk, + const CCTK_INT *bbox, + const CCTK_INT *gsize, + const CCTK_REAL *dx, + const CCTK_REAL *epsdis, + CCTK_REAL *rhs); static void apply (int const varindex, char const * const optstring, void * const arg); @@ -170,9 +189,15 @@ apply (int const varindex, char const * const optstring, void * const arg) break; } case 8: { - CCTK_FNAME(dissipation_8_4) - (varptr, &cctk_lsh[0], &cctk_lsh[1], &cctk_lsh[2], + if ( CCTK_Equals(dissipation_type,"Mattson-Svard-Nordstrom") ) { + CCTK_FNAME(dissipation_8_4) + (varptr, &cctk_lsh[0], &cctk_lsh[1], &cctk_lsh[2], + bbox, gsize, dx, &epsdis, rhsptr); + } else { + CCTK_FNAME(dissipation_8_4_alt) + (varptr, &cctk_lsh[0], &cctk_lsh[1], &cctk_lsh[2], bbox, gsize, dx, &epsdis, rhsptr); + } break; } default: @@ -187,9 +212,15 @@ apply (int const varindex, char const * const optstring, void * const arg) break; } case 6: { - CCTK_FNAME(dissipation_6_5_opt) - (varptr, cctk_lsh, cctk_gsh, cctk_lbnd, + if ( CCTK_Equals(dissipation_type,"Mattson-Svard-Nordstrom") ) { + CCTK_FNAME(dissipation_6_5_opt) + (varptr, cctk_lsh, cctk_gsh, cctk_lbnd, + bbox, gsize, dx, &epsdis, &diss_fraction, rhsptr); + } else { + CCTK_FNAME(dissipation_6_5_alt) + (varptr, cctk_lsh, cctk_gsh, cctk_lbnd, bbox, gsize, dx, &epsdis, &diss_fraction, rhsptr); + } break; } default: diff --git a/src/make.code.defn b/src/make.code.defn index 8767095..6ff753b 100644 --- a/src/make.code.defn +++ b/src/make.code.defn @@ -22,8 +22,10 @@ SRCS = call_derivs.c \ Dissipation_4_2.F90 \ Dissipation_6_3.F90 \ Dissipation_8_4.F90 \ + Dissipation_8_4_alt.F90 \ Dissipation_4_3_min_err_coeff.F90 \ Dissipation_6_5_min_err_coeff.F90 \ + Dissipation_6_5_min_err_coeff_alt.F90 \ get_coeffs.c \ Coefficients_2_1.F90 \ Coefficients_4_2.F90 \ |