#include "cctk.h" #include "cctk_Functions.h" #include "cctk_Parameters.h" subroutine set_coeff_4_2 ( nsize, loc_order, bb, gsize, imin, imax, dd ) implicit none DECLARE_CCTK_FUNCTIONS DECLARE_CCTK_PARAMETERS CCTK_REAL, parameter :: zero = 0.0 integer, parameter :: wp = kind(zero) CCTK_INT, intent(IN) :: nsize, loc_order CCTK_INT, dimension(2), intent(IN) :: bb CCTK_INT, intent(IN) :: gsize CCTK_INT, dimension(nsize), intent(OUT) :: imin, imax CCTK_REAL, dimension(nsize,nsize), intent(OUT) :: dd CCTK_REAL, dimension(2), save :: a CCTK_REAL, dimension(6,4), save :: q CCTK_INT :: i, il, ir logical, save :: first = .true. if ( first ) then a(1) = 2.0_wp/3.0_wp; a(2) = -1.0_wp/12.0_wp q(1,1) = -24.0_wp/17.0_wp; q(2,1) = 59.0_wp/34.0_wp q(3,1) = -4.0_wp/17.0_wp; q(4,1) = -3.0_wp/34.0_wp q(5,1) = zero; q(6,1) = zero q(1,2) = -1.0_wp/2.0_wp; q(2,2) = zero q(3,2) = 1.0_wp/2.0_wp; q(4,2) = zero q(5,2) = zero; q(6,2) = zero q(1,3) = 4.0_wp/43.0_wp; q(2,3) = -59.0_wp/86.0_wp q(3,3) = zero; q(4,3) = 59.0_wp/86.0_wp q(5,3) = -4.0_wp/43.0_wp; q(6,3) = zero q(1,4) = 3.0_wp/98.0_wp; q(2,4) = zero q(3,4) = -59.0_wp/98.0_wp; q(4,4) = zero q(5,4) = 32.0_wp/49.0_wp; q(6,4) = -4.0_wp/49.0_wp first = .false. end if dd = zero imin = 0 imax = -1 if ( bb(1) == 0 ) then il = 1 + gsize else dd(1:4,1) = q(1:4,1) imin(1) = 1; imax(1) = 4 dd(1:3,2) = q(1:3,2) imin(2) = 1; imax(2) = 3 dd(1:5,3) = q(1:5,3) imin(3) = 1; imax(3) = 5 dd(1:6,4) = q(1:6,4) imin(4) = 1; imax(4) = 6 il = 5 end if if ( bb(2) == 0 ) then ir = nsize - gsize else dd(nsize-5:nsize,nsize-3) = -q(6:1:-1,4) imin(nsize-3) = nsize-5; imax(nsize-3) = nsize dd(nsize-4:nsize,nsize-2) = -q(5:1:-1,3) imin(nsize-2) = nsize-4; imax(nsize-2) = nsize dd(nsize-2:nsize,nsize-1) = -q(3:1:-1,2) imin(nsize-1) = nsize-2; imax(nsize-1) = nsize dd(nsize-3:nsize,nsize) = -q(4:1:-1,1) imin(nsize) = nsize-3; imax(nsize) = nsize ir = nsize - 4 end if do i = il, ir dd(i-2:i-1,i) = -a(2:1:-1); dd(i+1:i+2,i) = a(1:2) imin(i) = i-2; imax(i) = i+2 end do end subroutine set_coeff_4_2