aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorgoodale <goodale@4451c3c6-1034-4891-99ea-21147727ccdf>2000-03-29 07:42:21 +0000
committergoodale <goodale@4451c3c6-1034-4891-99ea-21147727ccdf>2000-03-29 07:42:21 +0000
commit168a3c0a7f89b9ec54504e47c18dd0ca7bcb6940 (patch)
tree695a56da4cdc0bd830c22c12f1270fd30ab48613 /src
parenta457121df89fc3569421ccf68534a74838877c1f (diff)
Changed to use timelevels - thanks for putting the rotation in PUGH Ed.
Note that the no-boundary condition tests assumed that phi_new had zero at its boundaries. This is, or course, no longer true with timelevel rotation. Tom git-svn-id: http://svn.cactuscode.org/arrangements/CactusWave/WaveToyF77/trunk@50 4451c3c6-1034-4891-99ea-21147727ccdf
Diffstat (limited to 'src')
-rw-r--r--src/InitSymBound.F771
-rw-r--r--src/WaveToy.F7733
2 files changed, 11 insertions, 23 deletions
diff --git a/src/InitSymBound.F77 b/src/InitSymBound.F77
index 423b35c..cf5c2b4 100644
--- a/src/InitSymBound.F77
+++ b/src/InitSymBound.F77
@@ -40,7 +40,6 @@
sym(3) = one
call SetCartSymmetry(cctkGH, sym,'wavetoy::phi')
- call SetCartSymmetry(cctkGH, sym,'wavetoyf77::phi_next')
return
end
diff --git a/src/WaveToy.F77 b/src/WaveToy.F77
index 6130447..98d21be 100644
--- a/src/WaveToy.F77
+++ b/src/WaveToy.F77
@@ -45,6 +45,8 @@ c Declare variables in argument list
CCTK_REAL dx2,dy2,dz2,dt2
CCTK_REAL dx2i,dy2i,dz2i
+ CCTK_REAL factor
+
c Set up shorthands
c -----------------
dx = CCTK_DELTA_SPACE(1)
@@ -69,16 +71,16 @@ c -----------------
jend = cctk_lsh(2)-1
kend = cctk_lsh(3)-1
+ factor = 2*(1 - (dt2)*(dx2i + dy2i + dz2i))
+
c Do the evolution
c ----------------
do k = kstart, kend
do j = jstart, jend
do i = istart, iend
- phi_next(i,j,k) =
- & 2.0*(1.0 - (dt2)*(dx2i +
- & dy2i +dz2i))*phi(i,j,k) -
- & phi_old(i,j,k) + (dt2) *
+ phi_n(i,j,k) = factor*phi(i,j,k) -
+ & phi_p(i,j,k) + (dt2) *
& ((phi(i+1,j,k)+phi(i-1,j,k))*dx2i
& +(phi(i,j+1,k)+phi(i,j-1,k))*dy2i
& +(phi(i,j,k+1)+phi(i,j,k-1))*dz2i)
@@ -89,25 +91,12 @@ c ----------------
c Synchronize
c -----------
- call CCTK_SyncGroup(cctkGH,"wavetoyf77::scalartmps")
+ call CCTK_SyncGroup(cctkGH,"wavetoy::scalarevolve")
c Apply boundary conditions
c -------------------------
call WaveToyF77_Boundaries(CCTK_PASS_FTOF)
-c Update timeslices
-c -----------------
- do k = 1, cctk_lsh(3)
- do j = 1, cctk_lsh(2)
- do i = 1, cctk_lsh(1)
-
- phi_old(i,j,k) = phi(i,j,k)
- phi(i,j,k) = phi_next(i,j,k)
-
- end do
- end do
- end do
-
return
end
@@ -155,16 +144,16 @@ c ---------------------
c Apply the symmetry boundary conditions on any coordinate axes
c -------------------------------------------------------------
- call CartSymBCGroup(ierr,cctkGH,"wavetoyf77::scalartmps")
+ call CartSymBCGroup(ierr,cctkGH,"wavetoy::scalarevolve")
c Apply the outer boundary conditions
c -----------------------------------
if (CCTK_EQUALS(bound,"flat")) then
- call FlatBCVar(ierr,cctkGH,sw,"wavetoyf77::phi_next")
+ call FlatBCVar(ierr,cctkGH,sw,"wavetoy::phi")
else if (CCTK_EQUALS(bound,"zero")) then
- call ConstantBCVar(ierr,cctkGH,zero,sw,"wavetoyf77::phi_next")
+ call ConstantBCVar(ierr,cctkGH,zero,sw,"wavetoy::phi")
else if (CCTK_Equals(bound,"radiation").eq.1) then
- call RadiativeBCVar(ierr,cctkGH,zero,one,sw,"wavetoyf77::phi_next",
+ call RadiativeBCVar(ierr,cctkGH,zero,one,sw,"wavetoy::phi",
& "wavetoy::phi")
end if