aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorallen <allen@f80f6fb6-8356-4fd4-90bc-d84ad503c100>1999-09-21 11:21:20 +0000
committerallen <allen@f80f6fb6-8356-4fd4-90bc-d84ad503c100>1999-09-21 11:21:20 +0000
commitac9dedf87e0f1e99c3054eef4eb4b309671d3c51 (patch)
tree32a44e9e2c30837c01dde01f433754b69719f597 /src
parent22495226c82a63a15394df5b62b24243ad4770c4 (diff)
Mainly changing to new IO parameters and synchronising routines
git-svn-id: http://svn.cactuscode.org/arrangements/CactusWave/WaveToyF90/trunk@18 f80f6fb6-8356-4fd4-90bc-d84ad503c100
Diffstat (limited to 'src')
-rw-r--r--src/InitSymBound.F5
-rw-r--r--src/WaveToy.F39
2 files changed, 22 insertions, 22 deletions
diff --git a/src/InitSymBound.F b/src/InitSymBound.F
index 8f1dd64..a6c4685 100644
--- a/src/InitSymBound.F
+++ b/src/InitSymBound.F
@@ -3,7 +3,7 @@
@date
@author Gabrielle Allen
@desc
- Sets the symmetries for Wave Toy
+ Sets the symmetries across the coordinate axes
@enddesc
@@*/
@@ -34,8 +34,7 @@
INTEGER, DIMENSION(3) :: sym = 1
call SetCartSymmetry(cctkGH, sym, 'wavetoy::phi')
-
- return
+ call SetCartSymmetry(cctkGH, sym, 'wavetoyf90::phi_tmp')
end subroutine WaveToyF90_InitSymbound
diff --git a/src/WaveToy.F b/src/WaveToy.F
index eb9bc0f..1c9bc8f 100644
--- a/src/WaveToy.F
+++ b/src/WaveToy.F
@@ -32,20 +32,20 @@
@@*/
- subroutine wavetoyf90_evolution(CCTK_FARGUMENTS)
+ subroutine WaveToyf90_Evolution(CCTK_FARGUMENTS)
implicit none
DECLARE_CCTK_FARGUMENTS
DECLARE_CCTK_PARAMETERS
-c Declare local variables
+! Declare local variables
INTEGER :: i,j,k
INTEGER :: istart, jstart, kstart, iend, jend, kend
CCTK_REAL :: dx,dy,dz,dt
-c Set up shorthands
-c -----------------
+! Set up shorthands
+! -----------------
dx = CCTK_DELTA_SPACE(1)
dy = CCTK_DELTA_SPACE(2)
dz = CCTK_DELTA_SPACE(3)
@@ -59,8 +59,8 @@ c -----------------
jend = cctk_lsh(2)-1
kend = cctk_lsh(3)-1
-c Do the evolution
-c ----------------
+! Do the evolution
+! ----------------
do k = kstart, kend
do j = jstart, jend
do i = istart, iend
@@ -77,20 +77,21 @@ c ----------------
end do
end do
+! Synchronize
+! -----------
+ call CCTK_SyncGroup(cctkGH,"wavetoyf90::temps")
+
+! Apply boundary conditions
+! -------------------------
+ call WaveToyF90_Boundaries(CCTK_PASS_FTOF)
+
! Update timeslices
! -----------------
phi_old = phi
phi = phi_tmp
-! Apply boundary conditions
-! -------------------------
- call WaveToyF90_Boundaries(CCTK_FARGUMENTS)
-
-! Synchronize
-! -----------
- call CCTK_SyncGroup(cctkGH,"wavetoy::scalarevolve")
- end subroutine wavetoyf90_evolution
+ end subroutine WaveToyF90_Evolution
/*@@
@@ -121,20 +122,20 @@ c ----------------
CCTK_REAL,parameter :: zero = 0.0
CCTK_REAL,parameter :: one = 1.0
- call ApplySymmetry(cctkGH,"wavetoy::scalarevolve")
+ call ApplySymmetry(cctkGH,"wavetoyf90::temps")
if (CCTK_EQUALS(bound,"flat")) then
- call ApplyFlatBC(ierr,cctkGH,sw,"wavetoy::phi")
+ call ApplyFlatBC(ierr,cctkGH,sw,"wavetoyf90::phi_tmp")
else if (CCTK_EQUALS(bound,"radiation")) then
call ApplyRadiativeBC(ierr,cctkGH,zero,one,sw,
- & "wavetoy::phi","wavetoy::phi_old")
+ & "wavetoyf90::phi_tmp","wavetoy::phi")
end if
if (ierr < 0) then
- call CCTK_WARN(0,"Boundary conditions not applied - giving up!");
+ call CCTK_WARN(0,"Boundary conditions not applied - giving up!")
end if
- end subroutine wavetoyf90_boundaries
+ end subroutine WaveToyF90_Boundaries