aboutsummaryrefslogtreecommitdiff
path: root/src/dissipation.c
diff options
context:
space:
mode:
authordiener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a>2006-02-13 19:22:04 +0000
committerdiener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a>2006-02-13 19:22:04 +0000
commitdf7c02eb3b05ab5360a62cf9da3c22df8c9c83a6 (patch)
tree1f0eea57f27a9f881fa69f56965bd93d47fc8028 /src/dissipation.c
parent6c9e944df63221316b029d28b7cee534819da525 (diff)
Only use one sided finite differences near multipatch boundaries. Should fix
the problems we have with periodic boundary conditions. Also added the missing declarations of fortran routines with Kreiss-Oliger type dissipation for certain cases. git-svn-id: https://svn.cct.lsu.edu/repos/numrel/LSUThorns/SummationByParts/trunk@58 f69c4107-0314-4c4f-9ad4-17e986b73f4a
Diffstat (limited to 'src/dissipation.c')
-rw-r--r--src/dissipation.c51
1 files changed, 50 insertions, 1 deletions
diff --git a/src/dissipation.c b/src/dissipation.c
index cc93b33..48f7ba4 100644
--- a/src/dissipation.c
+++ b/src/dissipation.c
@@ -7,6 +7,9 @@
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"
+#include "util_ErrorCodes.h"
+#include "util_Table.h"
+
void CCTK_FCALL CCTK_FNAME(dissipation_4_3_opt) (const CCTK_REAL *var,
const CCTK_INT *lsh,
const CCTK_INT *gsh,
@@ -17,6 +20,16 @@ void CCTK_FCALL CCTK_FNAME(dissipation_4_3_opt) (const CCTK_REAL *var,
const CCTK_REAL *epsdis,
const CCTK_REAL *diss_fraction,
CCTK_REAL *rhs);
+void CCTK_FCALL CCTK_FNAME(dissipation_4_3_alt) (const CCTK_REAL *var,
+ const CCTK_INT *lsh,
+ const CCTK_INT *gsh,
+ const CCTK_INT *lbnd,
+ const CCTK_INT *bbox,
+ const CCTK_INT *gsize,
+ const CCTK_REAL *dx,
+ const CCTK_REAL *epsdis,
+ const CCTK_REAL *diss_fraction,
+ CCTK_REAL *rhs);
void CCTK_FCALL CCTK_FNAME(dissipation_6_5_opt) (const CCTK_REAL *var,
const CCTK_INT *ni,
const CCTK_INT *nj,
@@ -64,6 +77,15 @@ void CCTK_FCALL CCTK_FNAME(dissipation_4_2) (const CCTK_REAL *var,
const CCTK_REAL *dx,
const CCTK_REAL *epsdis,
CCTK_REAL *rhs);
+void CCTK_FCALL CCTK_FNAME(dissipation_4_2_alt) (const CCTK_REAL *var,
+ const CCTK_INT *ni,
+ const CCTK_INT *nj,
+ const CCTK_INT *nk,
+ const CCTK_INT *bbox,
+ const CCTK_INT *gsize,
+ const CCTK_REAL *dx,
+ const CCTK_REAL *epsdis,
+ CCTK_REAL *rhs);
void CCTK_FCALL CCTK_FNAME(dissipation_6_3) (const CCTK_REAL *var,
const CCTK_INT *ni,
const CCTK_INT *nj,
@@ -73,6 +95,15 @@ void CCTK_FCALL CCTK_FNAME(dissipation_6_3) (const CCTK_REAL *var,
const CCTK_REAL *dx,
const CCTK_REAL *epsdis,
CCTK_REAL *rhs);
+void CCTK_FCALL CCTK_FNAME(dissipation_6_3_alt) (const CCTK_REAL *var,
+ const CCTK_INT *ni,
+ const CCTK_INT *nj,
+ const CCTK_INT *nk,
+ const CCTK_INT *bbox,
+ const CCTK_INT *gsize,
+ const CCTK_REAL *dx,
+ const CCTK_REAL *epsdis,
+ CCTK_REAL *rhs);
void CCTK_FCALL CCTK_FNAME(dissipation_8_4) (const CCTK_REAL *var,
const CCTK_INT *ni,
const CCTK_INT *nj,
@@ -120,6 +151,8 @@ apply (int const varindex, char const * const optstring, void * const arg)
int n;
int d;
int ierr;
+ CCTK_INT symtable, pen_sym_handle;
+ CCTK_INT symbnd[6];
assert (varindex >= 0);
@@ -128,8 +161,24 @@ apply (int const varindex, char const * const optstring, void * const arg)
gsize[d] = cctk_nghostzones[d];
}
+ symtable = SymmetryTableHandleForGrid (cctkGH);
+ if (symtable<0) {
+ CCTK_WARN(0,"symtable is out of bounds");
+ }
+
+ ierr = Util_TableGetIntArray (symtable, 6, symbnd, "symmetry_handle");
+ if (ierr!=6) {
+ CCTK_WARN(0,"Util_TableGetIntArray returned error");
+ }
+
+ pen_sym_handle = SymmetryHandleOfName ( "multipatch" );
+
for (d=0; d<6; ++d) {
- bbox[d] = cctk_bbox[d];
+ if ( symbnd[d] == pen_sym_handle && cctk_bbox[d] == 1 ) {
+ bbox[d] = 1;
+ } else {
+ bbox[d] = 0;
+ }
}
rhsindex = MoLQueryEvolvedRHS (varindex);