aboutsummaryrefslogtreecommitdiff
path: root/CarpetExtra/FOWaveToyF77
diff options
context:
space:
mode:
authorschnetter <>2003-06-27 14:26:00 +0000
committerschnetter <>2003-06-27 14:26:00 +0000
commitea582beb0d512a99f6e3658060488c107be4d8e5 (patch)
tree76c1cb3f242bae96785e6e9af211bd7373fe4783 /CarpetExtra/FOWaveToyF77
parent6b25e7c5e7cfb07a94102b172100598a95557e84 (diff)
Remove superfluous parameter checking.
Remove superfluous parameter checking. Provide correct asymptotic limits. darcs-hash:20030627142611-07bb3-106ad5807b1a2ab2feea367886e03cbe9ca2377d.gz
Diffstat (limited to 'CarpetExtra/FOWaveToyF77')
-rw-r--r--CarpetExtra/FOWaveToyF77/src/FOWaveToy.F7723
1 files changed, 18 insertions, 5 deletions
diff --git a/CarpetExtra/FOWaveToyF77/src/FOWaveToy.F77 b/CarpetExtra/FOWaveToyF77/src/FOWaveToy.F77
index b8cb5fe86..1450ff325 100644
--- a/CarpetExtra/FOWaveToyF77/src/FOWaveToy.F77
+++ b/CarpetExtra/FOWaveToyF77/src/FOWaveToy.F77
@@ -128,6 +128,8 @@ c Local declarations
CCTK_REAL zero, one
parameter (zero=0, one=1)
+ integer table, table_derivs
+
character fbound*100
CCTK_INT fboundlen
@@ -183,16 +185,27 @@ c -----------------------------------
call CCTK_FortranString (fboundlen, bound, fbound)
if (fboundlen.lt.0) call CCTK_WARN (0, "internal error")
- ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, -1,
+
+ call Util_TableCreateFromString (table, "LIMIT=1.0 SPEED=1.0")
+ if (table.lt.0) call CCTK_WARN (0, "internal error")
+ call Util_TableCreateFromString (table_derivs, "LIMIT=0.0 SPEED=1.0")
+ if (table_derivs.lt.0) call CCTK_WARN (0, "internal error")
+
+ ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, table,
$ "fowavetoy::phi", fbound)
- ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, -1,
+ ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, table_derivs,
$ "fowavetoy::pi", fbound)
- ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, -1,
+ ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, table_derivs,
$ "fowavetoy::phix", fbound)
- ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, -1,
+ ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, table_derivs,
$ "fowavetoy::phiy", fbound)
- ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, -1,
+ ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, 1, table_derivs,
$ "fowavetoy::phiz", fbound)
if (ierr.ne.0) call CCTK_WARN (0, "internal error")
+ call Util_TableDestroy (ierr, table)
+ if (ierr.ne.0) call CCTK_WARN (0, "internal error")
+ call Util_TableDestroy (ierr, table_derivs)
+ if (ierr.ne.0) call CCTK_WARN (0, "internal error")
+
end