aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorgoodale <goodale@daab82bb-f315-4ad1-b6d0-9353ff8b6e27>2004-05-06 02:10:18 +0000
committergoodale <goodale@daab82bb-f315-4ad1-b6d0-9353ff8b6e27>2004-05-06 02:10:18 +0000
commit63a9307937610ee7cbe380199f906fc67ff95765 (patch)
tree44c2c5af3413dbb5231d34817e9309459794c772
parent904d43e4b15fa9a6fb2d302424e43737ac5d6582 (diff)
Making WaveToy versions more consistent. This change also fixes the bug
in some of the implementations whereby "zero" rather than "scalar" was passed to the boundary condition routines - PR 1676. git-svn-id: http://svn.cactuscode.org/arrangements/CactusWave/WaveToyFreeF90/trunk@73 daab82bb-f315-4ad1-b6d0-9353ff8b6e27
-rw-r--r--src/WaveToy.F9034
1 files changed, 16 insertions, 18 deletions
diff --git a/src/WaveToy.F90 b/src/WaveToy.F90
index da4080a..41a650c 100644
--- a/src/WaveToy.F90
+++ b/src/WaveToy.F90
@@ -5,6 +5,7 @@
@desc
Evolution routines for the wave equation solver
@enddesc
+ @version $Header$
@@*/
#include "cctk.h"
@@ -112,35 +113,32 @@ subroutine WaveToyFreeF90_Boundaries(CCTK_ARGUMENTS)
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_FUNCTIONS
-
+
+! Local declarations
CCTK_INT :: ierr
+ CHARACTER (len=100) :: boundary
+ INTEGER :: length
+
ierr = 0
+! The "bound" parameter needs to be converted into a Fortran string.
+ call CCTK_FortranString(length,bound,boundary)
+
! Apply the outer boundary conditions
! -----------------------------------
! Note: In each of the following calls to Boundary_SelectVarForBC,
! default arguments are used, so an invalid table handle of -1 can
! be passed
- if (CCTK_EQUALS(bound,"flat")) then
- ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, -1, &
- "wavetoy::phi", "Flat")
- else if (CCTK_EQUALS(bound,"static")) then
- ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, -1, &
- "wavetoy::phi", "Static")
- else if (CCTK_EQUALS(bound,"radiation")) then
- ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, -1, &
- "wavetoy::phi", "Radiation")
- else if (CCTK_EQUALS(bound,"robin")) then
- ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, -1, &
- "wavetoy::phi", "Robin")
+ if (CCTK_EQUALS(bound,"flat") .or. CCTK_EQUALS(bound,"static") .or. &
+ CCTK_EQUALS(bound,"radiation") .or. CCTK_EQUALS(bound,"robin") .or. &
+ CCTK_EQUALS(bound,"none") ) then
+ ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, -1, &
+ "wavetoy::phi", boundary)
else if (CCTK_EQUALS(bound,"zero")) then
- ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, -1, &
- "wavetoy::phi", "Scalar")
- else if (CCTK_EQUALS(bound,"none")) then
- ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, -1, &
- "wavetoy::phi", "None")
+ ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, -1, &
+ "wavetoy::phi", "Scalar")
end if
if (ierr < 0) then