! $Header$ #include "cctk.h" #include "cctk_Parameters.h" #include "cctk_Functions.h" subroutine dissipation_6_5_alt (var, lsh, gsh, lbnd, bb, gsize, & delta, epsilon, dfl, npatches, patch, rhs) use dissipation_coeff 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) :: epsilon CCTK_REAL, dimension(3), intent(in) :: dfl CCTK_INT, intent(in) :: npatches, patch CCTK_REAL :: zero = 0.0 integer, parameter :: wp = kind(zero) integer :: i, j, k, center CCTK_REAL, dimension(:,:), allocatable :: atmp, d, b, h CCTK_REAL, dimension(:,:), pointer :: a CCTK_REAL :: idel CCTK_INT :: il, ir, jl, jr, kl, kr ni = lsh(1); nj = lsh(2); nk = lsh(3) if ( first ) then allocate ( savedx(npatches), savedy(npatches), savedz(npatches) ) savedx = .false.; savedy = .false.; savedz = .false. allocate ( xcoeff(npatches), ycoeff(npatches), zcoeff(npatches) ) first = .false. end if if ( .not. savedx(patch) ) then allocate ( atmp(ni,ni), d(ni,ni), b(ni,ni), h(ni,ni) ) atmp = zero; d = zero; b = zero; h = zero call set_bmatrix ( b, bb(1:2), lsh(1), gsh(1), lbnd(1), delta(1), dfl(1) ) call set_hmatrix ( h, bb(1:2) ) call set_dmatrix ( d, bb(1:2) ) atmp = -transpose ( matmul ( h, matmul ( transpose(d), matmul ( b, d ) ) ) ) allocate ( xcoeff(patch)%coeff(ni,ni) ) xcoeff(patch)%coeff = atmp savedx(patch) = .true. deallocate ( atmp, d, b, h ) end if a => xcoeff(patch)%coeff idel = epsilon / ( 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 if ( zero_derivs_y == 0 ) then if ( .not. savedy(patch) ) then allocate ( atmp(nj,nj), d(nj,nj), b(nj,nj), h(nj,nj) ) atmp = zero; d = zero; b = zero; h = zero call set_bmatrix ( b, bb(3:4), lsh(2), gsh(2), lbnd(2), delta(2), dfl(2) ) call set_hmatrix ( h, bb(3:4) ) call set_dmatrix ( d, bb(3:4) ) atmp = -transpose ( matmul ( h, matmul ( transpose(d), matmul ( b, d ) ) ) ) allocate ( ycoeff(patch)%coeff(nj,nj) ) ycoeff(patch)%coeff = atmp savedy(patch) = .true. deallocate ( atmp, d, b, h ) end if a => ycoeff(patch)%coeff idel = epsilon / ( 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 end if if ( zero_derivs_z == 0 ) then if ( .not. savedz(patch) ) then allocate ( atmp(nk,nk), d(nk,nk), b(nk,nk), h(nk,nk) ) atmp = zero; d = zero; b = zero; h = zero call set_bmatrix ( b, bb(5:6), lsh(3), gsh(3), lbnd(3), delta(3), dfl(3) ) call set_hmatrix ( h, bb(5:6) ) call set_dmatrix ( d, bb(5:6) ) atmp = -transpose ( matmul ( h, matmul ( transpose(d), matmul ( b, d ) ) ) ) allocate ( zcoeff(patch)%coeff(nk,nk) ) zcoeff(patch)%coeff = atmp savedz(patch) = .true. deallocate ( atmp, d, b, h ) end if a => zcoeff(patch)%coeff idel = epsilon / ( 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 end if 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(inout) :: 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.0_wp sigma(3,1) = 0.0_wp sigma(4,1) = 0.0_wp sigma(5,1) = 0.0_wp sigma(6,1) = 0.0_wp sigma(7,1) = 0.0_wp sigma(1,2) = 0.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.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.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.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.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.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