aboutsummaryrefslogtreecommitdiff
path: root/src/dissipation.c
blob: 50d99a8f99dea2f9e1e52c4a001e0c8776eb27f0 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
/* $Header$ */

#include <assert.h>
#include <stdlib.h>

#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"

void CCTK_FCALL CCTK_FNAME(dissipation_6_5) (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);

static void
apply (int const varindex, char const * const optstring, void * const arg);

void
SBP_DissipationAdd (CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  CCTK_TraverseString (vars, apply, cctkGH, CCTK_GROUP_OR_VAR);
}

void
apply (int const varindex, char const * const optstring, void * const arg)
{
  cGH const * const cctkGH = (cGH const *) arg;
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  int rhsindex;
  int vargroup, rhsgroup;
  cGroup vardata, rhsdata;
  CCTK_REAL const * varptr;
  CCTK_REAL       * rhsptr;
  CCTK_INT ni, nj, nk, bbox[6], gsize[3];
  CCTK_REAL dx[3];
  int n;
  int d;
  int ierr;
  
  assert (varindex >= 0);
  
  for (d=0; d<3; ++d) {
    dx[d] = CCTK_DELTA_SPACE(d);
    gsize[d] = cctk_nghostzones[d];
  }
  
  for (d=0; d<6; ++d) {
    bbox[d] = cctk_bbox[d];
  }

  rhsindex = MoLQueryEvolvedRHS (varindex);
  if (rhsindex < 0) {
    char * const fullvarname = CCTK_FullName (varindex);
    assert (fullvarname);
    CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
                "There is no RHS variable registered with MoL for the evolved variable \"%s\"",
                fullvarname);
    free (fullvarname);
  }
  assert (rhsindex >= 0);
  
/*  if (verbose) {
    char * const fullvarname = CCTK_FullName (varindex);
    char * const fullrhsname = CCTK_FullName (rhsindex);
    assert (fullvarname);
    assert (fullrhsname);
    CCTK_VInfo (CCTK_THORNSTRING,
                "Applying dissipation to \"%s\" (RHS \"%s\")",
                fullvarname, fullrhsname);
    free (fullvarname);
    free (fullrhsname);
  } */
  
  vargroup = CCTK_GroupIndexFromVarI (varindex);
  assert (vargroup >= 0);
  rhsgroup = CCTK_GroupIndexFromVarI (rhsindex);
  assert (rhsgroup >= 0);
  
  ierr = CCTK_GroupData (vargroup, &vardata);
  assert (!ierr);
  ierr = CCTK_GroupData (rhsgroup, &rhsdata);
  assert (!ierr);
  
  assert (vardata.grouptype == CCTK_GF);
  assert (vardata.vartype == CCTK_VARIABLE_REAL);
  assert (vardata.dim == cctk_dim);
  assert (rhsdata.grouptype == CCTK_GF);
  assert (rhsdata.vartype == CCTK_VARIABLE_REAL);
  assert (rhsdata.dim == cctk_dim);
  
  varptr = CCTK_VarDataPtrI (cctkGH, 0, varindex);
  assert (varptr);
  rhsptr = CCTK_VarDataPtrI (cctkGH, 0, rhsindex);
  assert (rhsptr);
  
  if ( CCTK_Equals(norm_type,"Diagonal") ) {
    assert(0);
  } else {
    switch(order) {
    case 6: {
      CCTK_FNAME(dissipation_6_5)
        (varptr, &cctk_lsh[0], &cctk_lsh[1], &cctk_lsh[2],
                            bbox, gsize, dx, &epsdis, rhsptr);
      break;
    }
    default:
      assert(0);
    }
  }
}