aboutsummaryrefslogtreecommitdiff
path: root/src/Misner_multiple.c
blob: a42414bb931efdb57400909dd2b90855c915b3af (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
 /*@@
   @file      Misner_multiple.F
   @date      
   @author    Carsten Gundlach
   @desc 
              Set up initial data for multiple Misner black holes
   @enddesc 
   @version   $Header$
 @@*/

#include <string.h>

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

#include "IDAnalyticBH.h"

static const char *rcsid = "$Header$";

CCTK_FILEVERSION(CactusEinstein_IDAnalyticBH_Misner_multiple_c)

 /*@@
   @routine    Misner_multiple
   @date       
   @author     Carsten Gundlach
   @desc 
               Set up initial data for multiple Misner black holes
   @enddesc 
   @calls      MisnerEvalPsi
   @history
   @hdate Fri Apr 26 10:04:05 2002 @hauthor Tom Goodale
   @hdesc Changed to use new StaticConformal stuff
   @hdate 6.May.2003   @hauthor Jonathan Thornburg
   @hdesc code cleanups
   @endhistory 
@@*/
void Misner_multiple(CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  int i, npoints;
  CCTK_REAL xval, yval, zval;
  CCTK_REAL inv_psi, tmp0, tmp1, tmp2, tmp3, tmp4;
  const CCTK_REAL nm_eps             = 1e-6;   /* finite differencing step*/
  const CCTK_REAL halved_inv_nm_eps  = 0.5 / nm_eps;
  const CCTK_REAL inv_nm_eps_squared = 1.0 / SQR(nm_eps);
  int make_conformal_derivs = 0;

  CCTK_VInfo(CCTK_THORNSTRING,
             "setting up Misner initial data for %d black holes",
             (int)misner_nbh);

  /* Check if we should create and store conformal factor stuff */
  if(CCTK_EQUALS(metric_type, "static conformal"))
  {
    if      (CCTK_EQUALS(conformal_storage,"factor"))
    {
      *conformal_state = 1;
      make_conformal_derivs = 0;
    }
    else if (CCTK_EQUALS(conformal_storage,"factor+derivs"))
    {
      *conformal_state = 2;
      make_conformal_derivs = 1;
    }
    else if (CCTK_EQUALS(conformal_storage,"factor+derivs+2nd derivs"))
    {
      *conformal_state = 3;
      make_conformal_derivs = 1;
    }
    else
    {
      CCTK_VWarn(0, __LINE__, __FILE__, CCTK_THORNSTRING,
"Misner_multiple(): impossible value for conformal_storage=\"%s\"!",
                 conformal_storage);                            /*NOTREACHED*/
    }
  }


  npoints = cctk_lsh[0] * cctk_lsh[1] * cctk_lsh[2];

  /*     Initialize C global variables
   *     -----------------------------
   */

  Misner_init(misner_nbh, mu, nmax);

  /*     Get value of psi pointwise from a C function
   *     --------------------------------------------
   */

  for(i = 0; i < npoints; i++)
  {
    xval = x[i];
    yval = y[i];
    zval = z[i];

    tmp0 = MisnerEvalPsi(xval, yval, zval);
    psi[i] = tmp0;
           
    /*              Only calculate derivatives of psi if required
     *              ---------------------------------------------
     */
    if (make_conformal_derivs)
    {
      tmp1 = MisnerEvalPsi(xval+nm_eps,yval,zval);
      tmp2 = MisnerEvalPsi(xval-nm_eps,yval,zval);
      psix[i] = (tmp1-tmp2) * halved_inv_nm_eps;

      if(*conformal_state > 2)
      {
        psixx[i] = (tmp1+tmp2-2.0*tmp0) * inv_nm_eps_squared;
      }

      tmp1 = MisnerEvalPsi(xval,yval+nm_eps,zval);
      tmp2 = MisnerEvalPsi(xval,yval-nm_eps,zval);
      psiy[i]  = (tmp1-tmp2) * halved_inv_nm_eps;

      if(*conformal_state > 2)
      {
        psiyy[i] = (tmp1+tmp2-2.0*tmp0) * inv_nm_eps_squared;
      }

      tmp1 = MisnerEvalPsi(xval,yval,zval+nm_eps);
      tmp2 = MisnerEvalPsi(xval,yval,zval-nm_eps);
      psiz[i] = (tmp1-tmp2) * halved_inv_nm_eps;

      if(*conformal_state > 2)
      {
        psizz[i] = (tmp1+tmp2-2.0*tmp0) * inv_nm_eps_squared;
      }

      if(*conformal_state > 2)
      {
        tmp1 = MisnerEvalPsi(xval+nm_eps,yval+nm_eps,zval);
        tmp2 = MisnerEvalPsi(xval+nm_eps,yval-nm_eps,zval);
        tmp3 = MisnerEvalPsi(xval-nm_eps,yval+nm_eps,zval);
        tmp4 = MisnerEvalPsi(xval-nm_eps,yval-nm_eps,zval);
        psixy[i] = 0.25*(tmp1-tmp2-tmp3+tmp4) * inv_nm_eps_squared;

        tmp1 = MisnerEvalPsi(xval,yval+nm_eps,zval+nm_eps);
        tmp2 = MisnerEvalPsi(xval,yval-nm_eps,zval+nm_eps);
        tmp3 = MisnerEvalPsi(xval,yval+nm_eps,zval-nm_eps);
        tmp4 = MisnerEvalPsi(xval,yval-nm_eps,zval-nm_eps);
        psiyz[i] = 0.25*(tmp1-tmp2-tmp3+tmp4) * inv_nm_eps_squared;

        tmp1 = MisnerEvalPsi(xval+nm_eps,yval,zval+nm_eps);
        tmp2 = MisnerEvalPsi(xval+nm_eps,yval,zval-nm_eps);
        tmp3 = MisnerEvalPsi(xval-nm_eps,yval,zval+nm_eps);
        tmp4 = MisnerEvalPsi(xval-nm_eps,yval,zval-nm_eps);
        psixz[i] = 0.25*(tmp1-tmp2-tmp3+tmp4) * inv_nm_eps_squared;
      }
    }
  }

  /*     Cactus conventions
   *     ------------------
   */

  if (make_conformal_derivs)
  {
    for(i = 0; i < npoints; i++)
    {
      inv_psi = 1.0 / psi[i];

      psix[i]  *= inv_psi;
      psiy[i]  *= inv_psi;
      psiz[i]  *= inv_psi;
      if(*conformal_state > 2)
      {
        psixx[i] *= inv_psi;
        psixy[i] *= inv_psi;
        psixz[i] *= inv_psi;
        psiyy[i] *= inv_psi;
        psiyz[i] *= inv_psi;
        psizz[i] *= inv_psi;
      }
    }
  }

  /*     Metric depends on conformal state
   *     ---------------------------------
   */

  if (CCTK_EQUALS(metric_type, "static conformal"))
  {
    for(i = 0; i < npoints; i++)
    {
      gxx[i] = 1.0;
      gyy[i] = 1.0;
      gzz[i] = 1.0;
    }
  }
  else 
  {
    for(i = 0; i < npoints; i++)
    {
      gxx[i] = psi[i] * psi[i] * psi[i] * psi[i];
      gyy[i] = gxx[i];
      gzz[i] = gxx[i];
    }
  }
  IDAnalyticBH_zero_CCTK_REAL_array(npoints, gxy);
  IDAnalyticBH_zero_CCTK_REAL_array(npoints, gxz);
  IDAnalyticBH_zero_CCTK_REAL_array(npoints, gyz);

  /*     Time-symmetric data
   *     -------------------
   */

  IDAnalyticBH_zero_CCTK_REAL_array(npoints, kxx);
  IDAnalyticBH_zero_CCTK_REAL_array(npoints, kxy);
  IDAnalyticBH_zero_CCTK_REAL_array(npoints, kxz);
  IDAnalyticBH_zero_CCTK_REAL_array(npoints, kyy);
  IDAnalyticBH_zero_CCTK_REAL_array(npoints, kyz);
  IDAnalyticBH_zero_CCTK_REAL_array(npoints, kzz);

}