aboutsummaryrefslogtreecommitdiff
path: root/src/Schwarzschild.c
blob: 17b1c47b4c8f1a43e647569bf2a48fa340f72764 (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
 /*@@
   @file      Schwarzschild.c
   @date      Sun Oct 17 10:35:41 1999
   @author    Tom Goodale
   @desc 
   C version of Scwhwarzschild lapse routine
   @enddesc 
 @@*/

#include <math.h>

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


/* Need include file from Einstein */
#include "CactusEinstein/Einstein/src/Einstein.h"

void Schwarzschild(CCTK_CARGUMENTS)
{
  DECLARE_CCTK_CARGUMENTS
  DECLARE_CCTK_PARAMETERS

  CCTK_REAL zero,one,two,three;
  CCTK_REAL tmp;
  CCTK_REAL r_squared;

  int i,j,k;
  int nx,ny,nz;
  int index;

  zero = 0.0;
  one = 1.0;
  two = 2.0;
  three = 3.0;

  nx = cctk_lsh[0];
  ny = cctk_lsh[1];
  nz = cctk_lsh[2];
  
  /*     conformal metric flag */
  if(use_conformal == 1)
  {
  
    *conformal_state = CONFORMAL_METRIC;


    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);

          /*        Compute conformal factor */
          psi[index] = ( one + mass/two/r[index]);
         

          /*        derivatives of psi / psi */

          tmp = mass/two/pow(r[index],3) / psi[index];
          psix[index] = -x[index]*tmp;
          psiy[index] = -y[index]*tmp;
          psiz[index] = -z[index]*tmp;
         
          tmp = mass/two/pow(r[index],5)/psi[index];
          psixy[index] = three*x[index]*y[index]*tmp;
          psixz[index] = three*x[index]*z[index]*tmp;
          psiyz[index] = three*y[index]*z[index]*tmp;

          r_squared = r[index]*r[index];
          psixx[index]  = (three*x[index]*x[index] - r_squared)*tmp;
          psiyy[index]  = (three*y[index]*y[index] - r_squared)*tmp;
          psizz[index]  = (three*z[index]*z[index] - r_squared)*tmp;
         
          gxx[index] = one;
          gyy[index] = one;
          gzz[index] = one;
          gxy[index] = zero; 
          gxz[index] = zero;  
          gyz[index] = zero;  
        }
      }
    }
  }     
  else
  {
    *conformal_state = NOCONFORMAL_METRIC;

    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);

          gxx[index] = pow(( one + mass/two/r[index]), 4);
          gyy[index] = gxx[index];
          gzz[index] = gxx[index];
          gxy[index] = zero;
          gxz[index] = zero;
          gyz[index] = zero;
        }
      }
    }
  }

  /*     If the initial lapse is not one ... */
  if (CCTK_Equals(initial_lapse,"schwarz"))
  { 
    CCTK_INFO("Initialise with schwarzschild lapse");

    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);

          alp[index] = (2.*r[index] - mass)/(2.*r[index]+mass);
        }
      }
    }
  }
      
  /*     time symmetric initial slice */

  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);
        
        kxx[index] = zero;
        kxy[index] = zero;
        kxz[index] = zero;
        kyy[index] = zero;
        kyz[index] = zero;
        kzz[index] = zero;
      }
    }
  }
    
  return;

  
  

}