aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorallen <allen@daab82bb-f315-4ad1-b6d0-9353ff8b6e27>2002-05-31 10:04:58 +0000
committerallen <allen@daab82bb-f315-4ad1-b6d0-9353ff8b6e27>2002-05-31 10:04:58 +0000
commit758b3946a61c260cfe5f7f988f361dfc5dffe7ac (patch)
tree70eab8feb0b080083e89fa767c1fb00e89437d87
parente0dd72623fb65ae1e5c4cdb305cd61593fbd8518 (diff)
Standardizing all boundary conditions across wavetoy thorns
git-svn-id: http://svn.cactuscode.org/arrangements/CactusWave/WaveToyFreeF90/trunk@43 daab82bb-f315-4ad1-b6d0-9353ff8b6e27
-rw-r--r--param.ccl4
-rw-r--r--src/WaveToy.F9023
2 files changed, 21 insertions, 6 deletions
diff --git a/param.ccl b/param.ccl
index 162d5b4..9e3d28f 100644
--- a/param.ccl
+++ b/param.ccl
@@ -7,7 +7,9 @@ KEYWORD bound "Type of boundary condition to use"
{
"none" :: "No boundary condition"
"flat" :: "Flat boundary condition"
+ "static" :: "Static boundary condition"
"radiation" :: "Radiation boundary condition"
+ "robin" :: "Robin boundary condition"
+ "zero" :: "Zero boundary condition"
} "none"
-
diff --git a/src/WaveToy.F90 b/src/WaveToy.F90
index 915df8a..f82aed7 100644
--- a/src/WaveToy.F90
+++ b/src/WaveToy.F90
@@ -116,18 +116,31 @@ subroutine WaveToyFreeF90_Boundaries(CCTK_ARGUMENTS)
integer, dimension(3):: sw(3)=1
CCTK_REAL, parameter :: zero = 0.0
CCTK_REAL, parameter :: one = 1.0
-
+ CCTK_REAL :: finf
+ integer :: npow
+
+ finf = 1.0d0
+ npow = 1
+
call CartSymGN(ierr,cctkGH,"wavetoy::scalarevolve")
if (CCTK_EQUALS(bound,"flat")) then
call BndFlatVN(ierr,cctkGH,sw,"wavetoy::phi")
+ else if (CCTK_EQUALS(bound,"static")) then
+ call BndStaticVN(ierr,cctkGH,sw,"wavetoy::phi")
else if (CCTK_EQUALS(bound,"radiation")) then
- call BndRadiativeVN(ierr,cctkGH,sw,zero,one, &
- "wavetoy::phi","wavetoy::phi")
+ call BndRadiativeVN(ierr,cctkGH,sw,zero,one,"wavetoy::phi", &
+ "wavetoy::phi")
+ else if (CCTK_EQUALS(bound,"robin")) then
+ call BndRobinVN(ierr,cctkGH, sw, finf, npow,"wavetoy::phi")
+ else if (CCTK_EQUALS(bound,"zero")) then
+ call BndScalarVN(ierr,cctkGH,zero,sw,"wavetoy::phi")
+ else if (.NOT. CCTK_EQUALS(bound,"none")) then
+ call CCTK_WARN(0,"Unrecognized boundary condition")
end if
-
+
if (ierr < 0) then
- call CCTK_WARN(0,"WaveToyFreeF90_Boundaries: Boundary conditions not applied - giving up!")
+ call CCTK_WARN(0,"WaveToyFreeF90_Boundaries: Boundary conditions not applied")
end if
end subroutine WaveToyFreeF90_Boundaries