From f05d032fcedb3026030217864a3ddcf13d3f27c8 Mon Sep 17 00:00:00 2001 From: rideout Date: Fri, 14 Feb 2003 14:29:23 +0000 Subject: Modified to use new boundary interface. git-svn-id: http://svn.cactuscode.org/arrangements/CactusWave/WaveToyF77/trunk@89 4451c3c6-1034-4891-99ea-21147727ccdf --- interface.ccl | 6 +++++- schedule.ccl | 3 +++ src/WaveToy.F77 | 41 +++++++++++++++++++++-------------------- 3 files changed, 29 insertions(+), 21 deletions(-) diff --git a/interface.ccl b/interface.ccl index ac58f9f..cdfebdc 100644 --- a/interface.ccl +++ b/interface.ccl @@ -4,7 +4,6 @@ implements: wavetoy inherits: Grid Boundary -USES INCLUDE: Boundary.h USES INCLUDE: Symmetry.h public: @@ -13,3 +12,8 @@ cctk_real scalarevolve type = GF Timelevels=3 { phi } "The evolved scalar field" + +# Function aliasing from Fortran does not work at the moment... +#CCTK_INT FUNCTION Boundary_SelectVarForBC CCTK_POINTER GH, CCTK_INT faces, \ +# CCTK_INT table_handle, CCTK_STRING var_name, CCTK_STRING bc_name +#USES FUNCTION Boundary_SelectVarForBC diff --git a/schedule.ccl b/schedule.ccl index 166388d..50ca3e0 100644 --- a/schedule.ccl +++ b/schedule.ccl @@ -24,3 +24,6 @@ schedule WaveToyF77_Boundaries as WaveToy_Boundaries at EVOL AFTER WaveToy_Evolu LANG: Fortran } "Boundaries of 3D wave equation" +schedule GROUP ApplyBCs at EVOL after WaveToy_Evolution +{ +} "Apply boundary conditions" diff --git a/src/WaveToy.F77 b/src/WaveToy.F77 index 4b287bc..91b985a 100644 --- a/src/WaveToy.F77 +++ b/src/WaveToy.F77 @@ -2,12 +2,13 @@ @file WaveToy.F77 @date @author Tom Goodale - @desc + @desc Evolution routines for the wave equation solver @enddesc @@*/ -#include "cctk.h" +#include "cctk.h" +#include "cctk_Faces.h" #include "cctk_Parameters.h" #include "cctk_Arguments.h" @@ -135,33 +136,33 @@ c ------------------------------------------------------------- c Apply the outer boundary conditions c ----------------------------------- +c Note: In each of the following calls to Boundary_SelectVarForBC, +c default arguments are used, so an invalid table handle of -1 can +c be passed if (CCTK_EQUALS(bound,"flat")) then - call BndFlatVN(ierr,cctkGH,sw,"wavetoy::phi") + call Boundary_SelectVarForBC(ierr, cctkGH, CCTK_ALL_FACES, -1, + $ "wavetoy::phi", "Flat"); else if (CCTK_EQUALS(bound,"static")) then - call BndStaticVN(ierr,cctkGH,sw,"wavetoy::phi") + call Boundary_SelectVarForBC(ierr, cctkGH, CCTK_ALL_FACES, -1, + $ "wavetoy::phi", "Static"); else if (CCTK_EQUALS(bound,"radiation")) then - call BndRadiativeVN(ierr,cctkGH,sw,zero,one,"wavetoy::phi", - & "wavetoy::phi") + call Boundary_SelectVarForBC(ierr, cctkGH, CCTK_ALL_FACES, -1, + $ "wavetoy::phi", "Radiative"); else if (CCTK_EQUALS(bound,"robin")) then - call BndRobinVN(ierr,cctkGH, sw, finf, npow,"wavetoy::phi") + call Boundary_SelectVarForBC(ierr, cctkGH, CCTK_ALL_FACES, -1, + $ "wavetoy::phi", "Robin"); else if (CCTK_EQUALS(bound,"zero")) then - call BndScalarVN(ierr,cctkGH,zero,sw,"wavetoy::phi") - else if (.NOT. CCTK_EQUALS(bound,"none")) then + call Boundary_SelectVarForBC(ierr, cctkGH, CCTK_ALL_FACES, -1, + $ "wavetoy::phi", "Scalar"); + else if (CCTK_EQUALS(bound,"none")) then + call Boundary_SelectVarForBC(ierr, cctkGH, CCTK_ALL_FACES, -1, + $ "wavetoy::phi", "None"); + else call CCTK_WARN(0,"Unrecognized boundary condition") 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 return end - - - - - - - - - -- cgit v1.2.3