aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_Weights.c
blob: 6e01c95c571326f283d312340e466695c3b5c0d1 (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
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
 /*@@
   @file      GRHydro_Weights.F90
   @date      Sat Jan 26 01:40:14 2002
   @author    Ian Hawke
   @desc 
   This routine calculates the "weights" of the cells and cell boundaries.
   This is only required if using FishEye.

   There is also a routine that calculates the "physical" velocity.

   @enddesc 
 @@*/
   
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"

#include "FishEye.h"

 /*@@
   @routine    GRHydro_Weights
   @date       Sat Jan 26 01:41:02 2002
   @author     Ian Hawke
   @desc 
   Calculates the weights from FishEye.
   @enddesc 
   @calls     
   @calledby   
   @history 

   @endhistory 

@@*/

void GRHydro_Weights(CCTK_ARGUMENTS)
{
  
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;

  CCTK_REAL  dx, dy, dz, xbound, ybound, zbound;
  CCTK_REAL  jdet[27], jacobian[9], cell_centre_slice, cell_surface_l;
  CCTK_INT  i, j, k, nx, ny, nz, l, m, n, index, totalsize;

  dx = cctk_delta_space[0] / cctk_levfac[0];
  dy = cctk_delta_space[1] / cctk_levfac[1];
  dz = cctk_delta_space[2] / cctk_levfac[2];
  nx = cctk_lsh[0];
  ny = cctk_lsh[1];
  nz = cctk_lsh[2];
  totalsize = nx*ny*nz;

  if (use_weighted_fluxes)
  {
    
    for (k = 0; k < nz; k++)
    {
      for (j = 0; j < ny; j++)
      {
        for (i = 0; i < nx; i++)
        {

          /*
            Calculate the jacobian of the transformation at the 
            27 points; 3 points in every direction
          */

          index = CCTK_GFINDEX3D(cctkGH, i, j, k);

          for (n = 0; n < 3; n++)
          {
            for (m = 0; m < 3; m++)
            {
              for (l = 0; l < 3; l++)
              {
                
                xbound = x[index] + 0.5 * (l-2) * dx;
                ybound = y[index] + 0.5 * (m-2) * dy;
                zbound = z[index] + 0.5 * (n-2) * dz;

#ifdef FISHEYE_ACTIVE
                activejacobian(xbound,ybound,zbound,
                               &jacobian[0],&jacobian[1],&jacobian[2],
                               &jacobian[3],&jacobian[4],&jacobian[5],
                               &jacobian[6],&jacobian[7],&jacobian[8]);
#else
                CCTK_WARN(0,"You must compile with FishEye for this to work!");
#endif
                jdet[l+3*m+9*n] = 
                  -jacobian[2]*jacobian[4]*jacobian[6] + 
                  jacobian[1]*jacobian[5]*jacobian[6] + 
                  jacobian[2]*jacobian[3]*jacobian[7] - 
                  jacobian[0]*jacobian[5]*jacobian[7] - 
                  jacobian[1]*jacobian[3]*jacobian[8] + 
                  jacobian[0]*jacobian[4]*jacobian[8];

              }
            }
          }

          cell_jdet[index] = jdet[0];
          
          /* 
             Integrate by Simpsons rule over both directions
          */

          /*
            x
          */

          cell_surface[index] = 
            1.0 / 36.0 * ( (jdet[2] + 4.0 * jdet[5] + jdet[8])+ 
                          4.0 *(jdet[11] + 4.0 * jdet[14] + jdet[17])+
                          (jdet[20] + 4.0 * jdet[23] + jdet[26]) );

          /*
            y
          */

          cell_surface[index + totalsize] = 
            1.0 / 36.0 * ( (jdet[6] + 4.0 * jdet[7] + jdet[8])+
                          4.0 *(jdet[15] + 4.0 * jdet[16] + jdet[17])+
                          (jdet[24] + 4.0 * jdet[25] + jdet[26]) );

          /*
            z
          */

          cell_surface[index + 2 * totalsize] = 
            1.0 / 36.0 * ( (jdet[18] + 4.0 * jdet[19] + jdet[20])+
                          4.0 *(jdet[21] + 4.0 * jdet[22] + jdet[23])+
                          (jdet[24] + 4.0 * jdet[25] + jdet[26]) );
          
          /* 
             scalar for the x direction at the left, needed for
             the volume calculations
          */

          cell_surface_l = 
            1.0 / 36.0 * ( (jdet[0] + 4.0 * jdet[3] + jdet[6])+ 
                          4.0 *(jdet[9] + 4.0 * jdet[12] + jdet[15])+
                          (jdet[18] + 4.0 * jdet[21] + jdet[24]));

          /*
            The slice through the centre holding x fixed.
          */

          cell_centre_slice = 
            1.0 / 36.0 * ( (jdet[1] + 4.0 * jdet[4] + jdet[7])+ 
                          4.0 *(jdet[10] + 4.0 * jdet[13] + jdet[16])+
                          (jdet[19] + 4.0 * jdet[22] + jdet[25]) );

          /*
            The cell volume.
            Again using Simpsons rule, only in the x direction.
          */

          cell_volume[index] = 
            1.0 / 6.0 * ( cell_surface_l + 4.0 * cell_centre_slice + 
                          cell_surface[index] );
          
          
        }
      }
    }
    
  }
  else
  {
    
    for (k = 0; k < nz; k++)
    {
      for (j = 0; j < ny; j++)
      {
        for (i = 0; i < nx; i++)
        {

          index = CCTK_GFINDEX3D(cctkGH, i, j, k);
    
          cell_volume[index] = 1.0;
          cell_surface[index] = 1.0;
          cell_surface[index + totalsize] = 1.0;
          cell_surface[index + 2 * totalsize] = 1.0;
          
        }
      }
    }
    
  }
  
  return;
}



 /*@@
   @routine    GRHydro_Fisheye_Analysis
   @date       
   @author     Christian Ott, Ian Hawke
   @desc 
   Prepare output of GRHydro variables in physical coordinates
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

int GRHydro_Fisheye_Analysis(CCTK_ARGUMENTS)
{

  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  double jacobian[9];

  int i,j,k, index;  

  /* let's first do the scalars: simple copying */

  for(i=0;
	i<cctkGH->cctk_lsh[0]*cctkGH->cctk_lsh[1]*cctkGH->cctk_lsh[2];
       i++) 
  {

    frho[i] = rho[i];
    fpress[i] = press[i];
    feps[i] = eps[i];

  }    

  /* now the vectors: more complicated... first we need to (numerically) 
     calculate the jacobian */


  for(i=0;i<cctkGH->cctk_lsh[0];i++)
  {
    for(j=0;j<cctkGH->cctk_lsh[1];j++)
    {
      for(k=0;k<cctkGH->cctk_lsh[2];k++)
      {
       
        index = CCTK_GFINDEX3D(cctkGH,i,j,k);
        
#ifdef FISHEYE_ACTIVE
        activejacobian(x[index], y[index], z[index],
                       &jacobian[0],&jacobian[1],&jacobian[2],
                       &jacobian[3],&jacobian[4],&jacobian[5],
                       &jacobian[6],&jacobian[7],&jacobian[8]);
#else
        CCTK_WARN(1,"The output of \"FishEye\" velocities is only meaningful if FishEye is active!");
        jacobian[0]=1.0;
        jacobian[1]=0.0;
        jacobian[2]=0.0;
        jacobian[3]=0.0;
        jacobian[4]=1.0;
        jacobian[5]=0.0;
        jacobian[6]=0.0;
        jacobian[7]=0.0;
        jacobian[8]=1.0;
#endif

        fvelx[index] =
          jacobian[0] * velx[index] +
          jacobian[1] * vely[index] +
          jacobian[2] * velz[index];
        
        fvely[index] =
          jacobian[3] * velx[index] +
          jacobian[4] * vely[index] +
          jacobian[5] * velz[index];
        
        fvelz[index] =
          jacobian[6] * velx[index] +
          jacobian[7] * vely[index] +
          jacobian[8] * velz[index];
        
      }
    }
  }
  
  return;

}