aboutsummaryrefslogtreecommitdiff
path: root/CarpetExtra/FOWaveToyF77
diff options
context:
space:
mode:
authorschnetter <>2003-07-08 19:39:00 +0000
committerschnetter <>2003-07-08 19:39:00 +0000
commite012009aa12be53c6277cab9737812cd1d3070b9 (patch)
treedb4bd8bcc052ef106798b7c3d04911b03d0fd999 /CarpetExtra/FOWaveToyF77
parent062a5eb015829c65eaeff1b30d187ac9a09faf19 (diff)
Evolve up to the boundary.
darcs-hash:20030708193923-07bb3-4bf8f7536ac456bfbf92a3465f75a71841f4ff77.gz
Diffstat (limited to 'CarpetExtra/FOWaveToyF77')
-rw-r--r--CarpetExtra/FOWaveToyF77/schedule.ccl5
-rw-r--r--CarpetExtra/FOWaveToyF77/src/FOWaveToy.F7742
2 files changed, 31 insertions, 16 deletions
diff --git a/CarpetExtra/FOWaveToyF77/schedule.ccl b/CarpetExtra/FOWaveToyF77/schedule.ccl
index 01d8b5605..43b97e968 100644
--- a/CarpetExtra/FOWaveToyF77/schedule.ccl
+++ b/CarpetExtra/FOWaveToyF77/schedule.ccl
@@ -1,5 +1,5 @@
# Schedule definitions for thorn FOWaveToy77
-# $Header: /home/eschnett/C/carpet/Carpet/CarpetExtra/FOWaveToyF77/schedule.ccl,v 1.6 2003/06/27 16:44:03 schnetter Exp $
+# $Header: /home/eschnett/C/carpet/Carpet/CarpetExtra/FOWaveToyF77/schedule.ccl,v 1.7 2003/07/08 21:39:23 schnetter Exp $
STORAGE: scalarevolve[3] scalarevolve_derivs[3]
@@ -19,6 +19,7 @@ schedule FOWaveToyF77_InitSymBound at BASEGRID
schedule FOWaveToyF77_Boundaries as FOWaveToy_Boundaries at INITIAL after FOWaveToy_InitialData
{
LANG: Fortran
+ OPTIONS: level
SYNC: scalarevolve scalarevolve_derivs
} "Select boundary conditions of 3D wave equation"
@@ -36,6 +37,7 @@ schedule FOWaveToyF77_Evolution as FOWaveToy_Evolution at EVOL
schedule FOWaveToyF77_Boundaries as FOWaveToy_Boundaries at EVOL after FOWaveToy_Evolution
{
LANG: Fortran
+ OPTIONS: level
SYNC: scalarevolve scalarevolve_derivs
} "Select boundary conditions of 3D wave equation"
@@ -48,6 +50,7 @@ schedule group ApplyBCs as FOWaveToy_ApplyBCs at EVOL after FOWaveToy_Boundaries
schedule FOWaveToyF77_Boundaries as FOWaveToy_Boundaries in POSTRESTRICT
{
LANG: Fortran
+ OPTIONS: level
SYNC: scalarevolve scalarevolve_derivs
} "Select boundary conditions of 3D wave equation"
diff --git a/CarpetExtra/FOWaveToyF77/src/FOWaveToy.F77 b/CarpetExtra/FOWaveToyF77/src/FOWaveToy.F77
index 3d5a31469..3a0efb628 100644
--- a/CarpetExtra/FOWaveToyF77/src/FOWaveToy.F77
+++ b/CarpetExtra/FOWaveToyF77/src/FOWaveToy.F77
@@ -43,8 +43,12 @@ c Declare variables in argument list
DECLARE_CCTK_FUNCTIONS
DECLARE_CCTK_PARAMETERS
+ INTEGER bndwidth
+ PARAMETER (bndwidth=1)
+
INTEGER i,j,k
- INTEGER istart, jstart, kstart, iend, jend, kend
+ INTEGER istart(3), iend(3)
+ INTEGER d
CCTK_REAL dx,dy,dz,dt
CCTK_REAL dxi,dyi,dzi
@@ -61,19 +65,24 @@ c -----------------
dyi = 1/dy
dzi = 1/dz
- istart = 1+cctk_nghostzones(1)
- jstart = 1+cctk_nghostzones(2)
- kstart = 1+cctk_nghostzones(3)
-
- iend = cctk_lsh(1)-cctk_nghostzones(1)
- jend = cctk_lsh(2)-cctk_nghostzones(2)
- kend = cctk_lsh(3)-cctk_nghostzones(3)
+ do d = 1, 3
+ if (cctk_bbox(2*d-1).eq.0) then
+ istart(d) = 1+cctk_nghostzones(d)
+ else
+ istart(d) = 1+bndwidth
+ end if
+ if (cctk_bbox(2*d).eq.0) then
+ iend(d) = cctk_lsh(d)-cctk_nghostzones(d)
+ else
+ iend(d) = cctk_lsh(d)-bndwidth
+ end if
+ end do
c Do the evolution
c ----------------
- do k = kstart, kend
- do j = jstart, jend
- do i = istart, iend
+ do k = istart(3), iend(3)
+ do j = istart(2), iend(2)
+ do i = istart(1), iend(1)
pi(i,j,k) = pi_p_p(i,j,k) + dt * (
$ (phix_p(i+1,j,k) - phix_p(i-1,j,k))*dxi
@@ -128,6 +137,9 @@ c Local declarations
CCTK_REAL zero, one
parameter (zero=0, one=1)
+ INTEGER bndwidth
+ PARAMETER (bndwidth=1)
+
integer table
save table
data table /-1/
@@ -194,11 +206,11 @@ c -----------------------------------
call CCTK_FortranString (fboundlen, bound, fbound)
if (fboundlen.lt.0) call CCTK_WARN (0, "internal error")
- ierr = Boundary_SelectGroupForBC (cctkGH, CCTK_ALL_FACES, 1, table,
- $ "fowavetoy::scalarevolve", fbound)
+ ierr = Boundary_SelectGroupForBC (cctkGH, CCTK_ALL_FACES, bndwidth,
+ $ table, "fowavetoy::scalarevolve", fbound)
if (ierr.ne.0) call CCTK_WARN (0, "internal error")
- ierr = Boundary_SelectGroupForBC (cctkGH, CCTK_ALL_FACES, 1, table,
- $ "fowavetoy::scalarevolve_derivs", fbound)
+ ierr = Boundary_SelectGroupForBC (cctkGH, CCTK_ALL_FACES, bndwidth,
+ $ table, "fowavetoy::scalarevolve_derivs", fbound)
if (ierr.ne.0) call CCTK_WARN (0, "internal error")
end