aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authordiener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a>2010-03-19 22:24:27 +0000
committerdiener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a>2010-03-19 22:24:27 +0000
commit96a104f0d0703cbe827752189c5d7895f9ca4519 (patch)
tree7ce5e8835e7263393ac802432bff407dc17ac804
parentb041a5952a08c089d52a4c483fb27bec71b5061b (diff)
Add a routine to return the size of the boundary region, i.e. the region
where centered finite differences are not performed. This is necessary in order to mix SBP and hard coded stencils in CTGamma. git-svn-id: https://svn.cct.lsu.edu/repos/numrel/LSUThorns/SummationByParts/trunk@120 f69c4107-0314-4c4f-9ad4-17e986b73f4a
-rw-r--r--interface.ccl6
-rw-r--r--src/get_boundary_width.c104
-rw-r--r--src/make.code.defn3
3 files changed, 112 insertions, 1 deletions
diff --git a/interface.ccl b/interface.ccl
index 15fc072..58593b0 100644
--- a/interface.ccl
+++ b/interface.ccl
@@ -77,6 +77,12 @@ SUBROUTINE GetLshIndexRanges ( \
CCTK_INT OUT ARRAY imax )
PROVIDES FUNCTION GetLshIndexRanges WITH get_lsh_iranges LANGUAGE C
+SUBROUTINE GetBoundWidth (\
+ CCTK_POINTER_TO_CONST IN cctkGH, \
+ CCTK_INT OUT ARRAY bsize, \
+ CCTK_INT IN table_handle )
+PROVIDES FUNCTION GetBoundWidth WITH Get_Bound_Width LANGUAGE C
+
CCTK_INT FUNCTION GetDomainSpecification \
(CCTK_INT IN size, \
CCTK_REAL OUT ARRAY physical_min, \
diff --git a/src/get_boundary_width.c b/src/get_boundary_width.c
new file mode 100644
index 0000000..3c2686f
--- /dev/null
+++ b/src/get_boundary_width.c
@@ -0,0 +1,104 @@
+#include <assert.h>
+#include <stdio.h>
+
+#include "cctk.h"
+#include "cctk_Parameters.h"
+#include "cctk_Arguments.h"
+
+#include "util_ErrorCodes.h"
+#include "util_Table.h"
+
+#include "stencil.h"
+
+
+void Get_Bound_Width ( const CCTK_POINTER_TO_CONST cctkGH_, CCTK_INT *bsize,
+ const CCTK_INT table_handle )
+{
+ cGH const * restrict const cctkGH = cctkGH_;
+ DECLARE_CCTK_PARAMETERS
+ DECLARE_CCTK_ARGUMENTS
+
+ int onesided[6];
+ CCTK_INT loc_order;
+ CCTK_INT gsize[6];
+ int nelements;
+
+ SBP_determine_onesided_stencil (cctkGH, onesided);
+
+ if ( table_handle >=0 ) {
+ nelements = Util_TableGetInt ( table_handle, &loc_order, "order" );
+ if ( nelements == UTIL_ERROR_TABLE_NO_SUCH_KEY ) {
+ loc_order = order;
+ } else if ( nelements != 1) {
+ CCTK_WARN (0, "The options table has an entry \"order\", but it does not have the right properties");
+ }
+ } else {
+ loc_order = order;
+ }
+
+ gsize[0] = cctk_nghostzones[0];
+ gsize[1] = cctk_nghostzones[0];
+ gsize[2] = cctk_nghostzones[1];
+ gsize[3] = cctk_nghostzones[1];
+ gsize[4] = cctk_nghostzones[2];
+ gsize[5] = cctk_nghostzones[2];
+
+ if ( CCTK_Equals(norm_type,"Diagonal") ) {
+ if ( sbp_1st_deriv ) {
+ switch(loc_order) {
+ case 2: {
+ for (int d = 0; d<6; d++) bsize[d] = (onesided[d])?1:gsize[d];
+ break;
+ }
+ case 4: {
+ for (int d = 0; d<6; d++) bsize[d] = (onesided[d])?4:gsize[d];
+ break;
+ }
+ case 6: {
+ for (int d = 0; d<6; d++) bsize[d] = (onesided[d])?6:gsize[d];
+ break;
+ }
+ case 8: {
+ for (int d = 0; d<6; d++) bsize[d] = (onesided[d])?8:gsize[d];
+ break;
+ }
+ default:
+ CCTK_WARN (0, "Unknown stencil specified");
+ }
+ } else {
+ switch(loc_order) {
+ case 2: {
+ for (int d = 0; d<6; d++) bsize[d] = (onesided[d])?1:gsize[d];
+ break;
+ }
+ case 4: {
+ for (int d = 0; d<6; d++) bsize[d] = (onesided[d])?2:gsize[d];
+ break;
+ }
+ case 6: {
+ for (int d = 0; d<6; d++) bsize[d] = (onesided[d])?3:gsize[d];
+ break;
+ }
+ case 8: {
+ for (int d = 0; d<6; d++) bsize[d] = (onesided[d])?4:gsize[d];
+ break;
+ }
+ default:
+ CCTK_WARN (0, "Unknown 1st derivative stencil specified");
+ }
+ }
+ } else {
+ switch(loc_order) {
+ case 4: {
+ for (int d = 0; d<6; d++) bsize[d] = (onesided[d])?5:gsize[d];
+ break;
+ }
+ case 6: {
+ for (int d = 0; d<6; d++) bsize[d] = (onesided[d])?7:gsize[d];
+ break;
+ }
+ default:
+ CCTK_WARN (0, "Unknown stencil specified");
+ }
+ }
+}
diff --git a/src/make.code.defn b/src/make.code.defn
index 7fa5829..9536e0f 100644
--- a/src/make.code.defn
+++ b/src/make.code.defn
@@ -70,7 +70,8 @@ SRCS = call_derivs.c \
Coefficients_6_5_min_err_coeff.F90 \
Poisoning.F90 \
get_offset.c \
- stencil.c
+ stencil.c \
+ get_boundary_width.c
# Subdirectories containing source files
SUBDIRS =