! $Header$ #include "cctk.h" #include "cctk_Parameters.h" 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) :: epsilon CCTK_REAL :: zero = 0.0 integer, parameter :: wp = kind(zero) CCTK_REAL, dimension(3,2) :: a CCTK_REAL :: idel CCTK_INT :: il, ir, jl, jr, kl, kr call set_coeff ( a ) 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) else rhs(1,:,:) = rhs(1,:,:) + & ( a(1,1) * var(1,:,:) + a(2,1) * var(2,:,:) ) * idel rhs(2,:,:) = rhs(2,:,:) + & ( a(1,2) * var(1,:,:) + a(2,2) * var(2,:,:) + & a(3,2) * var(3,:,:) ) * idel il = 3 end if if ( bb(2) == 0 ) then ir = ni - gsize(1) else rhs(ni-1,:,:) = rhs(ni-1,:,:) + & ( a(1,2) * var(ni,:,:) + a(2,2) * var(ni-1,:,:) + & a(3,2) * var(ni-2,:,:) ) * idel rhs(ni,:,:) = rhs(ni,:,:) + & ( a(1,1) * var(ni,:,:) + a(2,1) * var(ni-1,:,:) ) * idel ir = ni - 2 end if !$omp parallel workshare rhs(il:ir,:,:) = rhs(il:ir,:,:) + & ( -2.0_wp * var(il:ir,:,:) + & ( var(il-1:ir-1,:,:) + & var(il+1:ir+1,:,:) ) ) * idel !$omp end parallel workshare if ( zero_derivs_y == 0 ) then call set_coeff ( a ) 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) else rhs(:,1,:) = rhs(:,1,:) + & ( a(1,1) * var(:,1,:) + a(2,1) * var(:,2,:) ) * idel rhs(:,2,:) = rhs(:,2,:) + & ( a(1,2) * var(:,1,:) + a(2,2) * var(:,2,:) + & a(3,2) * var(:,3,:) ) * idel jl = 3 end if if ( bb(4) == 0 ) then jr = nj - gsize(2) else rhs(:,nj-1,:) = rhs(:,nj-1,:) + & ( a(1,2) * var(:,nj,:) + a(2,2) * var(:,nj-1,:) + & a(3,2) * var(:,nj-2,:) ) * idel rhs(:,nj,:) = rhs(:,nj,:) + & ( a(1,1) * var(:,nj,:) + a(2,1) * var(:,nj-1,:) ) * idel jr = nj - 2 end if !$omp parallel workshare rhs(:,jl:jr,:) = rhs(:,jl:jr,:) + & ( -2.0_wp * var(:,jl:jr,:) + & ( var(:,jl-1:jr-1,:) + & var(:,jl+1:jr+1,:) ) ) * idel !$omp end parallel workshare end if if ( zero_derivs_z == 0 ) then call set_coeff ( a ) 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) else rhs(:,:,1) = rhs(:,:,1) + & ( a(1,1) * var(:,:,1) + a(2,1) * var(:,:,2) ) * idel rhs(:,:,2) = rhs(:,:,2) + & ( a(1,2) * var(:,:,1) + a(2,2) * var(:,:,2) + & a(3,2) * var(:,:,3) ) * idel kl = 3 end if if ( bb(6) == 0 ) then kr = nk - gsize(3) else rhs(:,:,nk-1) = rhs(:,:,nk-1) + & ( a(1,2) * var(:,:,nk) + a(2,2) * var(:,:,nk-1) + & a(3,2) * var(:,:,nk-2) ) * idel rhs(:,:,nk) = rhs(:,:,nk) + & ( a(1,1) * var(:,:,nk) + a(2,1) * var(:,:,nk-1) ) * idel kr = nk - 2 end if !$omp parallel workshare rhs(:,:,kl:kr) = rhs(:,:,kl:kr) + & ( -2.0_wp * var(:,:,kl:kr) + & ( var(:,:,kl-1:kr-1) + & var(:,:,kl+1:kr+1) ) ) * idel !$omp end parallel workshare end if contains subroutine set_coeff ( a ) implicit none CCTK_REAL, dimension(3,2), intent(OUT) :: a CCTK_REAL :: zero = 0.0 integer, parameter :: wp = kind(zero) a(1,1) = -2.0_wp a(2,1) = 2.0_wp a(3,1) = 0.0_wp a(1,2) = 1.0_wp a(2,2) = -2.0_wp a(3,2) = 1.0_wp end subroutine set_coeff end subroutine dissipation_2_1