aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authordiener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a>2008-07-14 17:14:48 +0000
committerdiener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a>2008-07-14 17:14:48 +0000
commite3b9cd6b30215fb09a45910eafca7cd9121ceffa (patch)
tree3191911b23b09591699c4ace2bc1a2182a29ddc6 /src
parentb613d73f118c88bd116f0d8fe1c32cd3a7f855db (diff)
Make these routines aware of the offset applied to the underlying patch
structure so that integration done with these coefficients will leave out the overlap points. git-svn-id: https://svn.cct.lsu.edu/repos/numrel/LSUThorns/SummationByParts/trunk@115 f69c4107-0314-4c4f-9ad4-17e986b73f4a
Diffstat (limited to 'src')
-rw-r--r--src/GetScalProdDiag.F9027
-rw-r--r--src/get_offset.c10
2 files changed, 27 insertions, 10 deletions
diff --git a/src/GetScalProdDiag.F90 b/src/GetScalProdDiag.F90
index 00e3dff..c4ec3d1 100644
--- a/src/GetScalProdDiag.F90
+++ b/src/GetScalProdDiag.F90
@@ -20,6 +20,8 @@ subroutine SBP_GetScalProdDiag ( cctkGH, dir, nsize, sigmad )
CCTK_REAL, parameter :: zero = 0.0
integer, parameter :: wp = kind(zero)
integer, dimension(6) :: onesided
+ CCTK_INT, dimension(6) :: offset
+ CCTK_INT :: ol, or
CCTK_REAL, dimension(2), parameter :: bmask_2 = (/ 0.5_wp, 1.0_wp /)
CCTK_REAL, dimension(4), parameter :: bmask_4 = (/ 17.0_wp/48.0_wp, &
@@ -47,30 +49,37 @@ subroutine SBP_GetScalProdDiag ( cctkGH, dir, nsize, sigmad )
call SBP_determine_onesided_stencil (cctkGH, onesided)
- sigmad = 1.0_wp
+ call get_shiftout ( cctkGH, offset )
+
+ ol = offset(dir*2+1)
+ or = offset(dir*2+2)
+
+ sigmad(1:ol) = zero
+ sigmad(nsize+1-or:nsize) = zero
+ sigmad(1+ol:nsize-or) = 1.0_wp
if ( onesided(dir*2+1) == 1 ) then
select case (order)
case (2)
- sigmad(1:2) = bmask_2
+ sigmad(1+ol:2+ol) = bmask_2
case (4)
- sigmad(1:4) = bmask_4
+ sigmad(1+ol:4+ol) = bmask_4
case (6)
- sigmad(1:6) = bmask_6
+ sigmad(1+ol:6+ol) = bmask_6
case (8)
- sigmad(1:8) = bmask_8
+ sigmad(1+ol:8+ol) = bmask_8
end select
end if
if ( onesided(dir*2+2) == 1 ) then
select case (order)
case (2)
- sigmad(nsize:nsize-1:-1) = bmask_2
+ sigmad(nsize-or:nsize-1-or:-1) = bmask_2
case (4)
- sigmad(nsize:nsize-3:-1) = bmask_4
+ sigmad(nsize-or:nsize-3-or:-1) = bmask_4
case (6)
- sigmad(nsize:nsize-5:-1) = bmask_6
+ sigmad(nsize-or:nsize-5-or:-1) = bmask_6
case (8)
- sigmad(nsize:nsize-7:-1) = bmask_8
+ sigmad(nsize-or:nsize-7-or:-1) = bmask_8
end select
end if
diff --git a/src/get_offset.c b/src/get_offset.c
index c955681..ed2b024 100644
--- a/src/get_offset.c
+++ b/src/get_offset.c
@@ -14,7 +14,6 @@ void get_shiftout ( const CCTK_POINTER_TO_CONST cctkGH_,
cGH const * restrict const cctkGH = cctkGH_;
DECLARE_CCTK_PARAMETERS
- DECLARE_CCTK_ARGUMENTS
int i;
CCTK_INT nboundaryzones[6];
@@ -45,6 +44,15 @@ void get_shiftout ( const CCTK_POINTER_TO_CONST cctkGH_,
}
/***
+ * And it's fortran callable version...
+ */
+void CCTK_FCALL CCTK_FNAME(get_shiftout)( const CCTK_POINTER_TO_CONST cctkGH_,
+ CCTK_INT *shiftout )
+{
+ get_shiftout ( cctkGH_, shiftout );
+}
+
+/***
* This function returns the effective values of local index ranges,
* taking into account ghost zones and boundary_shiftout_* values.
* If the use_shiftout=no, boundary offsets are set to zero.