! $Header$ #include "cctk.h" #include "cctk_Parameters.h" subroutine dissipation_8_4 (var, ni, nj, nk, bb, gsize, offset, 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_INT, dimension(6), intent(in) :: offset 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(12,8) :: a CCTK_REAL :: idel CCTK_INT :: il, ir, jl, jr, kl, kr, ol, or call set_coeff ( a ) 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) else ol = offset(1) rhs(1+ol,:,:) = rhs(1+ol,:,:) + & ( a(1,1) * var(1+ol,:,:) + a(2,1) * var(2+ol,:,:) + & a(3,1) * var(3+ol,:,:) + a(4,1) * var(4+ol,:,:) + & a(5,1) * var(5+ol,:,:) ) * idel rhs(2+ol,:,:) = rhs(2+ol,:,:) + & ( a(1,2) * var(1+ol,:,:) + a(2,2) * var(2+ol,:,:) + & a(3,2) * var(3+ol,:,:) + a(4,2) * var(4+ol,:,:) + & a(5,2) * var(5+ol,:,:) + a(6,2) * var(6+ol,:,:) ) * idel rhs(3+ol,:,:) = rhs(3+ol,:,:) + & ( a(1,3) * var(1+ol,:,:) + a(2,3) * var(2+ol,:,:) + & a(3,3) * var(3+ol,:,:) + a(4,3) * var(4+ol,:,:) + & a(5,3) * var(5+ol,:,:) + a(6,3) * var(6+ol,:,:) + & a(7,3) * var(7+ol,:,:) ) * idel rhs(4+ol,:,:) = rhs(4+ol,:,:) + & ( a(1,4) * var(1+ol,:,:) + a(2,4) * var(2+ol,:,:) + & a(3,4) * var(3+ol,:,:) + a(4,4) * var(4+ol,:,:) + & a(5,4) * var(5+ol,:,:) + a(6,4) * var(6+ol,:,:) + & a(7,4) * var(7+ol,:,:) + a(8,4) * var(8+ol,:,:) ) * idel rhs(5+ol,:,:) = rhs(5+ol,:,:) + & ( a(1,5) * var(1+ol,:,:) + a(2,5) * var(2+ol,:,:) + & a(3,5) * var(3+ol,:,:) + a(4,5) * var(4+ol,:,:) + & a(5,5) * var(5+ol,:,:) + a(6,5) * var(6+ol,:,:) + & a(7,5) * var(7+ol,:,:) + a(8,5) * var(8+ol,:,:) + & a(9,5) * var(9+ol,:,:) ) * idel rhs(6+ol,:,:) = rhs(6+ol,:,:) + & ( a(2,6) * var(2+ol,:,:) + a(3,6) * var(3+ol,:,:) + & a(4,6) * var(4+ol,:,:) + a(5,6) * var(5+ol,:,:) + & a(6,6) * var(6+ol,:,:) + a(7,6) * var(7+ol,:,:) + & a(8,6) * var(8+ol,:,:) + a(9,6) * var(9+ol,:,:) + & a(10,6) * var(10+ol,:,:) ) * idel rhs(7+ol,:,:) = rhs(7+ol,:,:) + & ( a(3,7) * var(3+ol,:,:) + a(4,7) * var(4+ol,:,:) + & a(5,7) * var(5+ol,:,:) + a(6,7) * var(6+ol,:,:) + & a(7,7) * var(7+ol,:,:) + a(8,7) * var(8+ol,:,:) + & a(9,7) * var(9+ol,:,:) + a(10,7) * var(10+ol,:,:) + & a(11,7) * var(11+ol,:,:) ) * idel rhs(8+ol,:,:) = rhs(8+ol,:,:) + & ( a(4,8) * var(4+ol,:,:) + a(5,8) * var(5+ol,:,:) + & a(6,8) * var(6+ol,:,:) + a(7,8) * var(7+ol,:,:) + & a(8,8) * var(8+ol,:,:) + a(9,8) * var(9+ol,:,:) + & a(10,8) * var(10+ol,:,:) + a(11,8) * var(11+ol,:,:) + & a(12,8) * var(12+ol,:,:) ) * idel il = 9 + ol end if if ( bb(2) == 0 ) then ir = ni - gsize(1) else or = ni - offset(2) rhs(or-7,:,:) = rhs(or-7,:,:) + & ( a(4,8) * var(or-3,:,:) + a(5,8) * var(or-4,:,:) + & a(6,8) * var(or-5,:,:) + a(7,8) * var(or-6,:,:) + & a(8,8) * var(or-7,:,:) + a(9,8) * var(or-8,:,:) + & a(10,8) * var(or-9,:,:) + a(11,8) * var(or-10,:,:) + & a(12,8) * var(or-11,:,:) ) * idel rhs(or-6,:,:) = rhs(or-6,:,:) + & ( a(3,7) * var(or-2,:,:) + a(4,7) * var(or-3,:,:) + & a(5,7) * var(or-4,:,:) + a(6,7) * var(or-5,:,:) + & a(7,7) * var(or-6,:,:) + a(8,7) * var(or-7,:,:) + & a(9,7) * var(or-8,:,:) + a(10,7) * var(or-9,:,:) + & a(11,7) * var(or-10,:,:) ) * idel rhs(or-5,:,:) = rhs(or-5,:,:) + & ( a(2,6) * var(or-1,:,:) + a(3,6) * var(or-2,:,:) + & a(4,6) * var(or-3,:,:) + a(5,6) * var(or-4,:,:) + & a(6,6) * var(or-5,:,:) + a(7,6) * var(or-6,:,:) + & a(8,6) * var(or-7,:,:) + a(9,6) * var(or-8,:,:) + & a(10,6) * var(or-9,:,:) ) * idel rhs(or-4,:,:) = rhs(or-4,:,:) + & ( a(1,5) * var(or,:,:) + a(2,5) * var(or-1,:,:) + & a(3,5) * var(or-2,:,:) + a(4,5) * var(or-3,:,:) + & a(5,5) * var(or-4,:,:) + a(6,5) * var(or-5,:,:) + & a(7,5) * var(or-6,:,:) + a(8,5) * var(or-7,:,:) + & a(9,5) * var(or-8,:,:) ) * idel rhs(or-3,:,:) = rhs(or-3,:,:) + & ( a(1,4) * var(or,:,:) + a(2,4) * var(or-1,:,:) + & a(3,4) * var(or-2,:,:) + a(4,4) * var(or-3,:,:) + & a(5,4) * var(or-4,:,:) + a(6,4) * var(or-5,:,:) + & a(7,4) * var(or-6,:,:) + a(8,4) * var(or-7,:,:) ) * idel rhs(or-2,:,:) = rhs(or-2,:,:) + & ( a(1,3) * var(or,:,:) + a(2,3) * var(or-1,:,:) + & a(3,3) * var(or-2,:,:) + a(4,3) * var(or-3,:,:) + & a(5,3) * var(or-4,:,:) + a(6,3) * var(or-5,:,:) + & a(7,3) * var(or-6,:,:) ) * idel rhs(or-1,:,:) = rhs(or-1,:,:) + & ( a(1,2) * var(or,:,:) + a(2,2) * var(or-1,:,:) + & a(3,2) * var(or-2,:,:) + a(4,2) * var(or-3,:,:) + & a(5,2) * var(or-4,:,:) + a(6,2) * var(or-5,:,:) ) * idel rhs(or,:,:) = rhs(or,:,:) + & ( a(1,1) * var(or,:,:) + a(2,1) * var(or-1,:,:) + & a(3,1) * var(or-2,:,:) + a(4,1) * var(or-3,:,:) + & a(5,1) * var(or-4,:,:) ) * idel ir = or - 8 end if rhs(il:ir,:,:) = rhs(il:ir,:,:) + & ( -70.0_wp * var(il:ir,:,:) + & 56.0_wp * ( var(il-1:ir-1,:,:) + & var(il+1:ir+1,:,:) ) - & 28.0_wp * ( var(il-2:ir-2,:,:) + & var(il+2:ir+2,:,:) ) + & 8.0_wp * ( var(il-3:ir-3,:,:) + & var(il+3:ir+3,:,:) ) - & ( var(il-4:ir-4,:,:) + & var(il+4:ir+4,:,:) ) ) * idel if ( zero_derivs_y == 0 ) then call set_coeff ( a ) 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) else ol = offset(3) rhs(:,1+ol,:) = rhs(:,1+ol,:) + & ( a(1,1) * var(:,1+ol,:) + a(2,1) * var(:,2+ol,:) + & a(3,1) * var(:,3+ol,:) + a(4,1) * var(:,4+ol,:) + & a(5,1) * var(:,5+ol,:) ) * idel rhs(:,2+ol,:) = rhs(:,2+ol,:) + & ( a(1,2) * var(:,1+ol,:) + a(2,2) * var(:,2+ol,:) + & a(3,2) * var(:,3+ol,:) + a(4,2) * var(:,4+ol,:) + & a(5,2) * var(:,5+ol,:) + a(6,2) * var(:,6+ol,:) ) * idel rhs(:,3+ol,:) = rhs(:,3+ol,:) + & ( a(1,3) * var(:,1+ol,:) + a(2,3) * var(:,2+ol,:) + & a(3,3) * var(:,3+ol,:) + a(4,3) * var(:,4+ol,:) + & a(5,3) * var(:,5+ol,:) + a(6,3) * var(:,6+ol,:) + & a(7,3) * var(:,7+ol,:) ) * idel rhs(:,4+ol,:) = rhs(:,4+ol,:) + & ( a(1,4) * var(:,1+ol,:) + a(2,4) * var(:,2+ol,:) + & a(3,4) * var(:,3+ol,:) + a(4,4) * var(:,4+ol,:) + & a(5,4) * var(:,5+ol,:) + a(6,4) * var(:,6+ol,:) + & a(7,4) * var(:,7+ol,:) + a(8,4) * var(:,8+ol,:) ) * idel rhs(:,5+ol,:) = rhs(:,5+ol,:) + & ( a(1,5) * var(:,1+ol,:) + a(2,5) * var(:,2+ol,:) + & a(3,5) * var(:,3+ol,:) + a(4,5) * var(:,4+ol,:) + & a(5,5) * var(:,5+ol,:) + a(6,5) * var(:,6+ol,:) + & a(7,5) * var(:,7+ol,:) + a(8,5) * var(:,8+ol,:) + & a(9,5) * var(:,9+ol,:) ) * idel rhs(:,6+ol,:) = rhs(:,6+ol,:) + & ( a(2,6) * var(:,2+ol,:) + a(3,6) * var(:,3+ol,:) + & a(4,6) * var(:,4+ol,:) + a(5,6) * var(:,5+ol,:) + & a(6,6) * var(:,6+ol,:) + a(7,6) * var(:,7+ol,:) + & a(8,6) * var(:,8+ol,:) + a(9,6) * var(:,9+ol,:) + & a(10,6) * var(:,10+ol,:) ) * idel rhs(:,7+ol,:) = rhs(:,7+ol,:) + & ( a(3,7) * var(:,3+ol,:) + a(4,7) * var(:,4+ol,:) + & a(5,7) * var(:,5+ol,:) + a(6,7) * var(:,6+ol,:) + & a(7,7) * var(:,7+ol,:) + a(8,7) * var(:,8+ol,:) + & a(9,7) * var(:,9+ol,:) + a(10,7) * var(:,10+ol,:) + & a(11,7) * var(:,11+ol,:) ) * idel rhs(:,8+ol,:) = rhs(:,8+ol,:) + & ( a(4,8) * var(:,4+ol,:) + a(5,8) * var(:,5+ol,:) + & a(6,8) * var(:,6+ol,:) + a(7,8) * var(:,7+ol,:) + & a(8,8) * var(:,8+ol,:) + a(9,8) * var(:,9+ol,:) + & a(10,8) * var(:,10+ol,:) + a(11,8) * var(:,11+ol,:) + & a(12,8) * var(:,12+ol,:) ) * idel jl = 9 + ol end if if ( bb(4) == 0 ) then jr = nj - gsize(2) else or = nj - offset(4) rhs(:,or-7,:) = rhs(:,or-7,:) + & ( a(4,8) * var(:,or-3,:) + a(5,8) * var(:,or-4,:) + & a(6,8) * var(:,or-5,:) + a(7,8) * var(:,or-6,:) + & a(8,8) * var(:,or-7,:) + a(9,8) * var(:,or-8,:) + & a(10,8) * var(:,or-9,:) + a(11,8) * var(:,or-10,:) + & a(12,8) * var(:,or-11,:) ) * idel rhs(:,or-6,:) = rhs(:,or-6,:) + & ( a(3,7) * var(:,or-2,:) + a(4,7) * var(:,or-3,:) + & a(5,7) * var(:,or-4,:) + a(6,7) * var(:,or-5,:) + & a(7,7) * var(:,or-6,:) + a(8,7) * var(:,or-7,:) + & a(9,7) * var(:,or-8,:) + a(10,7) * var(:,or-9,:) + & a(11,7) * var(:,or-10,:) ) * idel rhs(:,or-5,:) = rhs(:,or-5,:) + & ( a(2,6) * var(:,or-1,:) + a(3,6) * var(:,or-2,:) + & a(4,6) * var(:,or-3,:) + a(5,6) * var(:,or-4,:) + & a(6,6) * var(:,or-5,:) + a(7,6) * var(:,or-6,:) + & a(8,6) * var(:,or-7,:) + a(9,6) * var(:,or-8,:) + & a(10,6) * var(:,or-9,:) ) * idel rhs(:,or-4,:) = rhs(:,or-4,:) + & ( a(1,5) * var(:,or,:) + a(2,5) * var(:,or-1,:) + & a(3,5) * var(:,or-2,:) + a(4,5) * var(:,or-3,:) + & a(5,5) * var(:,or-4,:) + a(6,5) * var(:,or-5,:) + & a(7,5) * var(:,or-6,:) + a(8,5) * var(:,or-7,:) + & a(9,5) * var(:,or-8,:) ) * idel rhs(:,or-3,:) = rhs(:,or-3,:) + & ( a(1,4) * var(:,or,:) + a(2,4) * var(:,or-1,:) + & a(3,4) * var(:,or-2,:) + a(4,4) * var(:,or-3,:) + & a(5,4) * var(:,or-4,:) + a(6,4) * var(:,or-5,:) + & a(7,4) * var(:,or-6,:) + a(8,4) * var(:,or-7,:) ) * idel rhs(:,or-2,:) = rhs(:,or-2,:) + & ( a(1,3) * var(:,or,:) + a(2,3) * var(:,or-1,:) + & a(3,3) * var(:,or-2,:) + a(4,3) * var(:,or-3,:) + & a(5,3) * var(:,or-4,:) + a(6,3) * var(:,or-5,:) + & a(7,3) * var(:,or-6,:) ) * idel rhs(:,or-1,:) = rhs(:,or-1,:) + & ( a(1,2) * var(:,or,:) + a(2,2) * var(:,or-1,:) + & a(3,2) * var(:,or-2,:) + a(4,2) * var(:,or-3,:) + & a(5,2) * var(:,or-4,:) + a(6,2) * var(:,or-5,:) ) * idel rhs(:,or,:) = rhs(:,or,:) + & ( a(1,1) * var(:,or,:) + a(2,1) * var(:,or-1,:) + & a(3,1) * var(:,or-2,:) + a(4,1) * var(:,or-3,:) + & a(5,1) * var(:,or-4,:) ) * idel jr = or - 8 end if rhs(:,jl:jr,:) = rhs(:,jl:jr,:) + & ( -70.0_wp * var(:,jl:jr,:) + & 56.0_wp * ( var(:,jl-1:jr-1,:) + & var(:,jl+1:jr+1,:) ) - & 28.0_wp * ( var(:,jl-2:jr-2,:) + & var(:,jl+2:jr+2,:) ) + & 8.0_wp * ( var(:,jl-3:jr-3,:) + & var(:,jl+3:jr+3,:) ) - & ( var(:,jl-4:jr-4,:) + & var(:,jl+4:jr+4,:) ) ) * idel end if if ( zero_derivs_z == 0 ) then call set_coeff ( a ) 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) else ol = offset(5) rhs(:,:,1+ol) = rhs(:,:,1+ol) + & ( a(1,1) * var(:,:,1+ol) + a(2,1) * var(:,:,2+ol) + & a(3,1) * var(:,:,3+ol) + a(4,1) * var(:,:,4+ol) + & a(5,1) * var(:,:,5+ol) ) * idel rhs(:,:,2+ol) = rhs(:,:,2+ol) + & ( a(1,2) * var(:,:,1+ol) + a(2,2) * var(:,:,2+ol) + & a(3,2) * var(:,:,3+ol) + a(4,2) * var(:,:,4+ol) + & a(5,2) * var(:,:,5+ol) + a(6,2) * var(:,:,6+ol) ) * idel rhs(:,:,3+ol) = rhs(:,:,3+ol) + & ( a(1,3) * var(:,:,1+ol) + a(2,3) * var(:,:,2+ol) + & a(3,3) * var(:,:,3+ol) + a(4,3) * var(:,:,4+ol) + & a(5,3) * var(:,:,5+ol) + a(6,3) * var(:,:,6+ol) + & a(7,3) * var(:,:,7+ol) ) * idel rhs(:,:,4+ol) = rhs(:,:,4+ol) + & ( a(1,4) * var(:,:,1+ol) + a(2,4) * var(:,:,2+ol) + & a(3,4) * var(:,:,3+ol) + a(4,4) * var(:,:,4+ol) + & a(5,4) * var(:,:,5+ol) + a(6,4) * var(:,:,6+ol) + & a(7,4) * var(:,:,7+ol) + a(8,4) * var(:,:,8+ol) ) * idel rhs(:,:,5+ol) = rhs(:,:,5+ol) + & ( a(1,5) * var(:,:,1+ol) + a(2,5) * var(:,:,2+ol) + & a(3,5) * var(:,:,3+ol) + a(4,5) * var(:,:,4+ol) + & a(5,5) * var(:,:,5+ol) + a(6,5) * var(:,:,6+ol) + & a(7,5) * var(:,:,7+ol) + a(8,5) * var(:,:,8+ol) + & a(9,5) * var(:,:,9+ol) ) * idel rhs(:,:,6+ol) = rhs(:,:,6+ol) + & ( a(2,6) * var(:,:,2+ol) + a(3,6) * var(:,:,3+ol) + & a(4,6) * var(:,:,4+ol) + a(5,6) * var(:,:,5+ol) + & a(6,6) * var(:,:,6+ol) + a(7,6) * var(:,:,7+ol) + & a(8,6) * var(:,:,8+ol) + a(9,6) * var(:,:,9+ol) + & a(10,6) * var(:,:,10+ol) ) * idel rhs(:,:,7+ol) = rhs(:,:,7+ol) + & ( a(3,7) * var(:,:,3+ol) + a(4,7) * var(:,:,4+ol) + & a(5,7) * var(:,:,5+ol) + a(6,7) * var(:,:,6+ol) + & a(7,7) * var(:,:,7+ol) + a(8,7) * var(:,:,8+ol) + & a(9,7) * var(:,:,9+ol) + a(10,7) * var(:,:,10+ol) + & a(11,7) * var(:,:,11+ol) ) * idel rhs(:,:,8+ol) = rhs(:,:,8+ol) + & ( a(4,8) * var(:,:,4+ol) + a(5,8) * var(:,:,5+ol) + & a(6,8) * var(:,:,6+ol) + a(7,8) * var(:,:,7+ol) + & a(8,8) * var(:,:,8+ol) + a(9,8) * var(:,:,9+ol) + & a(10,8) * var(:,:,10+ol) + a(11,8) * var(:,:,11+ol) + & a(12,8) * var(:,:,12+ol) ) * idel kl = 9 + ol end if if ( bb(6) == 0 ) then kr = nk - gsize(3) else or = nk - offset(6) rhs(:,:,or-7) = rhs(:,:,or-7) + & ( a(4,8) * var(:,:,or-3) + a(5,8) * var(:,:,or-4) + & a(6,8) * var(:,:,or-5) + a(7,8) * var(:,:,or-6) + & a(8,8) * var(:,:,or-7) + a(9,8) * var(:,:,or-8) + & a(10,8) * var(:,:,or-9) + a(11,8) * var(:,:,or-10) + & a(12,8) * var(:,:,or-11) ) * idel rhs(:,:,or-6) = rhs(:,:,or-6) + & ( a(3,7) * var(:,:,or-2) + a(4,7) * var(:,:,or-3) + & a(5,7) * var(:,:,or-4) + a(6,7) * var(:,:,or-5) + & a(7,7) * var(:,:,or-6) + a(8,7) * var(:,:,or-7) + & a(9,7) * var(:,:,or-8) + a(10,7) * var(:,:,or-9) + & a(11,7) * var(:,:,or-10) ) * idel rhs(:,:,or-5) = rhs(:,:,or-5) + & ( a(2,6) * var(:,:,or-1) + a(3,6) * var(:,:,or-2) + & a(4,6) * var(:,:,or-3) + a(5,6) * var(:,:,or-4) + & a(6,6) * var(:,:,or-5) + a(7,6) * var(:,:,or-6) + & a(8,6) * var(:,:,or-7) + a(9,6) * var(:,:,or-8) + & a(10,6) * var(:,:,or-9) ) * idel rhs(:,:,or-4) = rhs(:,:,or-4) + & ( a(1,5) * var(:,:,or) + a(2,5) * var(:,:,or-1) + & a(3,5) * var(:,:,or-2) + a(4,5) * var(:,:,or-3) + & a(5,5) * var(:,:,or-4) + a(6,5) * var(:,:,or-5) + & a(7,5) * var(:,:,or-6) + a(8,5) * var(:,:,or-7) + & a(9,5) * var(:,:,or-8) ) * idel rhs(:,:,or-3) = rhs(:,:,or-3) + & ( a(1,4) * var(:,:,or) + a(2,4) * var(:,:,or-1) + & a(3,4) * var(:,:,or-2) + a(4,4) * var(:,:,or-3) + & a(5,4) * var(:,:,or-4) + a(6,4) * var(:,:,or-5) + & a(7,4) * var(:,:,or-6) + a(8,4) * var(:,:,or-7) ) * idel rhs(:,:,or-2) = rhs(:,:,or-2) + & ( a(1,3) * var(:,:,or) + a(2,3) * var(:,:,or-1) + & a(3,3) * var(:,:,or-2) + a(4,3) * var(:,:,or-3) + & a(5,3) * var(:,:,or-4) + a(6,3) * var(:,:,or-5) + & a(7,3) * var(:,:,or-6) ) * idel rhs(:,:,or-1) = rhs(:,:,or-1) + & ( a(1,2) * var(:,:,or) + a(2,2) * var(:,:,or-1) + & a(3,2) * var(:,:,or-2) + a(4,2) * var(:,:,or-3) + & a(5,2) * var(:,:,or-4) + a(6,2) * var(:,:,or-5) ) * idel rhs(:,:,or) = rhs(:,:,or) + & ( a(1,1) * var(:,:,or) + a(2,1) * var(:,:,or-1) + & a(3,1) * var(:,:,or-2) + a(4,1) * var(:,:,or-3) + & a(5,1) * var(:,:,or-4) ) * idel kr = or - 8 end if rhs(:,:,kl:kr) = rhs(:,:,kl:kr) + & ( -70.0_wp * var(:,:,kl:kr) + & 56.0_wp * ( var(:,:,kl-1:kr-1) + & var(:,:,kl+1:kr+1) ) - & 28.0_wp * ( var(:,:,kl-2:kr-2) + & var(:,:,kl+2:kr+2) ) + & 8.0_wp * ( var(:,:,kl-3:kr-3) + & var(:,:,kl+3:kr+3) ) - & ( var(:,:,kl-4:kr-4) + & var(:,:,kl+4:kr+4) ) ) * idel end if contains subroutine set_coeff ( a ) implicit none CCTK_REAL, dimension(12,8), intent(OUT) :: a CCTK_REAL :: zero = 0.0 integer, parameter :: wp = kind(zero) a(1,1) = -3.3910872088637970174997113084967416241083103770745_wp a(2,1) = 13.564348835455188069998845233986966496433241508298_wp a(3,1) = -20.346523253182782104998267850980449744649862262447_wp a(4,1) = 13.564348835455188069998845233986966496433241508298_wp a(5,1) = -3.3910872088637970174997113084967416241083103770745_wp a(6,1) = 0.0_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(1,2) = 2.6217119552210904473646423259312909608627056453179_wp a(2,2) = -11.142275809689634401299729885207986583666498992601_wp a(3,2) = 18.351983686547633131552496281519036726038939517225_wp a(4,2) = -14.419415753715997460505532792622100284744881049248_wp a(5,2) = 5.2434239104421808947292846518625819217254112906357_wp a(6,2) = -0.65542798880527261184116058148282274021567641132947_wp a(7,2) = 0.0_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(1,3) = -23.305235778623380376667790568855064784933288377246_wp a(2,3) = 108.75776696690910842444968932132363566302201242715_wp a(3,3) = -205.86291604450652666056548335821973893357738066567_wp a(4,3) = 201.97871008140262993112085159674389480275516593613_wp a(5,3) = -108.75776696690910842444968932132363566302201242715_wp a(6,3) = 31.073647704831173835557054091806753046577717836328_wp a(7,3) = -3.8842059631038967294446317614758441308222147295410_wp a(8,3) = 0.0_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(1,4) = 2.2245534287765737751523181243817045070532621097794_wp a(2,4) = -12.235043858271155763337749684099374788792941603787_wp a(3,4) = 28.919194574095459076980135616962158591692407427132_wp a(4,4) = -38.373546646395897621377487645584402746668771393695_wp a(5,4) = 31.143748002872032852132453741343863098745669536912_wp a(6,4) = -15.571874001436016426066226870671931549372834768456_wp a(7,4) = 4.4491068575531475503046362487634090141065242195588_wp a(8,4) = -0.55613835719414344378807953109542612676331552744485_wp a(9,4) = 0.0_wp a(10,4) = 0.0_wp a(11,4) = 0.0_wp a(12,4) = 0.0_wp a(1,5) = -2.4230202953323072711308162536265511957185829658094_wp a(2,5) = 19.384162362658458169046530029012409565748663726475_wp a(3,5) = -67.844568269304603591662855101543433480120323042664_wp a(4,5) = 135.68913653860920718332571020308686696024064608533_wp a(5,5) = -169.61142067326150897915713775385858370030080760666_wp a(6,5) = 135.68913653860920718332571020308686696024064608533_wp a(7,5) = -67.844568269304603591662855101543433480120323042664_wp a(8,5) = 19.384162362658458169046530029012409565748663726475_wp a(9,5) = -2.4230202953323072711308162536265511957185829658094_wp a(10,5) = 0.0_wp a(11,5) = 0.0_wp a(12,5) = 0.0_wp a(1,6) = 0.0_wp a(2,6) = -0.78217600900123184961735065035840034142603567514089_wp a(3,6) = 6.2574080720098547969388052028672027314082854011271_wp a(4,6) = -21.900928252034491789285818210035209559928998903945_wp a(5,6) = 43.801856504068983578571636420070419119857997807890_wp a(6,6) = -54.752320630086229473214545525088023899822497259862_wp a(7,6) = 43.801856504068983578571636420070419119857997807890_wp a(8,6) = -21.900928252034491789285818210035209559928998903945_wp a(9,6) = 6.2574080720098547969388052028672027314082854011271_wp a(10,6) = -0.78217600900123184961735065035840034142603567514089_wp a(11,6) = 0.0_wp a(12,6) = 0.0_wp a(1,7) = 0.0_wp a(2,7) = 0.0_wp a(3,7) = -1.0830767761393601764536458481012280421614377748694_wp a(4,7) = 8.6646142091148814116291667848098243372915021989551_wp a(5,7) = -30.326149731902084940702083746834385180520257696343_wp a(6,7) = 60.652299463804169881404167493668770361040515392685_wp a(7,7) = -75.815374329755212351755209367085962951300644240857_wp a(8,7) = 60.652299463804169881404167493668770361040515392685_wp a(9,7) = -30.326149731902084940702083746834385180520257696343_wp a(10,7) = 8.6646142091148814116291667848098243372915021989551_wp a(11,7) = -1.0830767761393601764536458481012280421614377748694_wp a(12,7) = 0.0_wp a(1,8) = 0.0_wp a(2,8) = 0.0_wp a(3,8) = 0.0_wp a(4,8) = -0.99075245444434671889501396229410272246695863420506_wp a(5,8) = 7.9260196355547737511601116983528217797356690736404_wp a(6,8) = -27.741068724441708129060390944234876229074841757742_wp a(7,8) = 55.482137448883416258120781888469752458149683515483_wp a(8,8) = -69.352671811104270322650977360587190572687104394354_wp a(9,8) = 55.482137448883416258120781888469752458149683515483_wp a(10,8) = -27.741068724441708129060390944234876229074841757742_wp a(11,8) = 7.9260196355547737511601116983528217797356690736404_wp a(12,8) = -0.99075245444434671889501396229410272246695863420506_wp end subroutine set_coeff end subroutine dissipation_8_4