aboutsummaryrefslogtreecommitdiff
path: root/src/ScalarBoundary.c
diff options
context:
space:
mode:
authorrideout <rideout@6a38eb6e-646e-4a02-a296-d141613ad6c4>2003-02-14 00:46:06 +0000
committerrideout <rideout@6a38eb6e-646e-4a02-a296-d141613ad6c4>2003-02-14 00:46:06 +0000
commitf93380441909a9c02fd7002a7071ee2f0b519f2f (patch)
treeeccc11dc7b7bc3eabd43b2a80d7f6f3ed2111438 /src/ScalarBoundary.c
parent7d93d9515961a769c22c32fd425d37439206f912 (diff)
Provide wrappers for all boundary conditions provided by thorn
Boundary, so that they can be 'selected for bc' under the new boundary interface. The details of default arguments and table keys will be documented in the thorn guide. git-svn-id: http://svn.cactuscode.org/arrangements/CactusBase/Boundary/trunk@198 6a38eb6e-646e-4a02-a296-d141613ad6c4
Diffstat (limited to 'src/ScalarBoundary.c')
-rw-r--r--src/ScalarBoundary.c146
1 files changed, 136 insertions, 10 deletions
diff --git a/src/ScalarBoundary.c b/src/ScalarBoundary.c
index 01d78ca..88fa8c9 100644
--- a/src/ScalarBoundary.c
+++ b/src/ScalarBoundary.c
@@ -17,6 +17,9 @@
#include <string.h>
#include "cctk.h"
+#include "cctk_Faces.h"
+#include "util_Table.h"
+#include "util_ErrorCodes.h"
#include "cctk_Parameters.h"
#include "cctk_FortranString.h"
@@ -25,12 +28,144 @@
/* the rcs ID and its dummy function to use it */
static const char *rcsid = "$Header$";
-CCTK_FILEVERSION(CactusBase_Boundary_ScalarBoundary_c)
+CCTK_FILEVERSION(CactusBase_Boundary_ScalarBoundary_c);
+static int ApplyBndScalar (const cGH *GH,
+ int stencil_dir,
+ const int *stencil_alldirs,
+ int dir,
+ CCTK_REAL scalar,
+ int first_var,
+ int num_vars);
/********************************************************************
******************** External Routines ************************
********************************************************************/
+
+/*@@
+ @routine BndScalar
+ @date 13 Feb 2003
+ @author David Rideout
+ @desc
+ Top level function which is registered as handling
+ the Scalar boundary condition
+ @enddesc
+ @calls ApplyBndScalar
+ @var GH
+ @vdesc Pointer to CCTK grid hierarchy
+ @vtype const cGH *
+ @vio in
+ @endvar
+ @var num_vars
+ @vdesc number of variables passed in through var_indices[]
+ @vtype int
+ @vio in
+ @endvar
+ @var var_indices
+ @vdesc array of variable indicies to which to apply this boundary
+ condition
+ @vtype int *
+ @vio in
+ @endvar
+ @var faces
+ @vdesc array of set of faces to which to apply the bc
+ @vtype int
+ @vio in
+ @endvar
+ @var table_handles
+ @vdesc array of table handles which hold extra arguments
+ @vtype int
+ @vio in
+ @endvar
+ @returntype int
+ @returndesc
+ return code of @seeroutine ApplyBndScalar
+ @endreturndesc
+@@*/
+
+int BndScalar(const cGH *GH, int num_vars, int *vars, int *faces, int *tables)
+{
+ int i, j, k, gdim, max_gdim, err, retval;
+
+ /* variables to pass to ApplyBndScalar */
+ /*int stencil_dir;*/ /* width of stencil in direction dir, not needed */
+ int *stencil_alldirs; /* width of stencil in all directions */
+ int dir; /* direction in which to apply bc */
+ CCTK_REAL scalar;
+
+ retval = 0; stencil_alldirs = NULL;
+
+ /* loop through variables, j at a time */
+ for (i=0; i<num_vars; i+=j) {
+
+ /* find other adjacent vars which are selected for identical bcs */
+ j=1;
+ while (i+j<num_vars && vars[i+j]==vars[i]+j && tables[i+j]==tables[i] &&
+ faces[i+j]==faces[i])
+ {
+ ++j;
+ }
+
+ /* Check to see if faces specification is valid */
+ if (faces[i] != CCTK_ALL_FACES)
+ {
+ CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
+ "Faces specification %d for Scalar boundary conditions on "
+ "%s is not implemented yet. "
+ "Applying radiative bcs to all (external) faces.", faces[i],
+ CCTK_VarName(vars[i]));
+ }
+ dir = 0; /* apply bc to all faces */
+
+ /* Set up default arguments for ApplyBndScalar */
+ /* Allocate memory and populate stencil width array */
+ gdim = CCTK_GroupDimFromVarI(vars[i]);
+ if (!stencil_alldirs)
+ {
+ stencil_alldirs = (int *) malloc(gdim*sizeof(int));
+ max_gdim = gdim;
+ } else if (gdim > max_gdim)
+ {
+ realloc(stencil_alldirs, gdim*sizeof(int));
+ max_gdim = gdim;
+ }
+ for (k=0; k<gdim; ++k)
+ {
+ stencil_alldirs[k] = 1;
+ }
+ /* Defaults for remainder of arguments */
+ scalar = 0.;
+
+ /* Look on table for possible non-default arguments
+ * (If any of these table look-ups fail, the value will be unchanged
+ * from its default value)
+ */
+ /* Scalar value */
+ err = Util_TableGetReal(tables[i], &scalar, "SCALAR");
+ if (err == UTIL_ERROR_BAD_HANDLE)
+ {
+ CCTK_VWarn(3, __LINE__, __FILE__, CCTK_THORNSTRING,
+ "Invalid table handle passed for Scalar boundary "
+ "conditions for %s. Using all default values.",
+ CCTK_VarName(vars[i]));
+ } else
+ {
+ /* Stencil width array */
+ Util_TableGetIntArray(tables[i], gdim, stencil_alldirs, "STENCIL WIDTH");
+ }
+
+ if ((retval = ApplyBndScalar(GH, 0, stencil_alldirs, dir, scalar, vars[i],
+ j)) < 0)
+ {
+ CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
+ "ApplyBndScalar() returned %d", retval);
+ }
+ }
+ free(stencil_alldirs);
+
+ return retval;
+}
+
/* prototypes for external C routines are declared in header Boundary.h
here only follow the fortran wrapper prototypes */
void CCTK_FCALL CCTK_FNAME (BndScalarDirVI)
@@ -89,15 +224,6 @@ void CCTK_FCALL CCTK_FNAME (BndScalarVN)
/********************************************************************
******************** Internal Routines ************************
********************************************************************/
-static int ApplyBndScalar (const cGH *GH,
- int stencil_dir,
- const int *stencil_alldirs,
- int dir,
- CCTK_REAL scalar,
- int first_var,
- int num_vars);
-
-
/*@@
@routine BndScalarDirVI
@date Tue Jan 16 2001