aboutsummaryrefslogtreecommitdiff
path: root/src/Bin_NS.cc
blob: 66dd9df4d96ce2ed5f0ad287518d8bdd785ff916 (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
#include <cstdio>
#include <cassert>
#include <vector>

#include <cctk.h>
#include <cctk_Arguments.h>
#include <cctk_Parameters.h>

#include <bin_ns.h>

using namespace std;


static void set_dt_from_domega (CCTK_ARGUMENTS,
                                CCTK_REAL const* const var,
                                CCTK_REAL      * const dtvar,
                                CCTK_REAL const& omega)
{
  DECLARE_CCTK_ARGUMENTS;
  
  int const npoints = cctk_lsh[0] * cctk_lsh[1] * cctk_lsh[2];
  vector<CCTK_REAL> dxvar(npoints), dyvar(npoints);
  
  Diff_gv (cctkGH, 0, var, &dxvar[0], -1);
  Diff_gv (cctkGH, 1, var, &dyvar[0], -1);
  
#pragma omp parallel for
  for (int i=0; i<npoints; ++i) {
    CCTK_REAL const ephix = +y[i];
    CCTK_REAL const ephiy = -x[i];
    CCTK_REAL const dphi_var = ephix * dxvar[i] + ephiy * dyvar[i];
    dtvar[i] = omega * dphi_var;
  }
}



extern "C"
void ID_Bin_NS_initialise (CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  CCTK_INFO ("Setting up LORENE Bin_NS initial data");
  
  // Meudon data are distributed in SI units (MKSA).  Here are some
  // conversion factors.
  
  // Defined constants
  CCTK_REAL const c_light = 299792458.0; // speed of light [m/s]
  
  // Constants of nature (IAU, CODATA):
  CCTK_REAL const G_grav = 6.67428e-11; // gravitational constant [m^3/kg/s^2]
  CCTK_REAL const M_sun  = 1.98892e+30; // solar mass [kg]
  
  // Cactus units in terms of SI units:
  // (These are derived from M = M_sun, c = G = 1, and using 1/M_sun
  // for the magnetic field)
  CCTK_REAL const cactusM = M_sun;
  CCTK_REAL const cactusL = cactusM * G_grav / pow(c_light,2);
  CCTK_REAL const cactusT = cactusL / c_light;
  
  // Other quantities in terms of Cactus units
  CCTK_REAL const coord_unit = cactusL / 1.0e+3;         // from km
  CCTK_REAL const rho_unit   = cactusM / pow(cactusL,3); // from kg/m^3
  CCTK_REAL const ener_unit  = pow(cactusL,2);           // from c^2
  CCTK_REAL const vel_unit   = cactusL / cactusT / c_light; // from c
  
  
  
  CCTK_INFO ("Setting up coordinates");
  
  int const npoints = cctk_lsh[0] * cctk_lsh[1] * cctk_lsh[2];
  vector<double> xx(npoints), yy(npoints), zz(npoints);
  
#pragma omp parallel for
  for (int i=0; i<npoints; ++i) {
    xx[i] = x[i] * coord_unit;
    yy[i] = y[i] * coord_unit;
    zz[i] = z[i] * coord_unit;
  }
  
  
  
  CCTK_VInfo (CCTK_THORNSTRING, "Reading from file \"%s\"", filename);
  
  Bin_NS bin_ns (npoints, &xx[0], &yy[0], &zz[0], filename);
  
  CCTK_VInfo (CCTK_THORNSTRING, "omega [rad/s]:       %g", bin_ns.omega);
  CCTK_VInfo (CCTK_THORNSTRING, "dist [km]:           %g", bin_ns.dist);
  CCTK_VInfo (CCTK_THORNSTRING, "dist_mass [km]:      %g", bin_ns.dist_mass);
  CCTK_VInfo (CCTK_THORNSTRING, "mass1_b [M_sun]:     %g", bin_ns.mass1_b);
  CCTK_VInfo (CCTK_THORNSTRING, "mass2_b [M_sun]:     %g", bin_ns.mass2_b);
  CCTK_VInfo (CCTK_THORNSTRING, "mass_ADM [M_sun]:    %g", bin_ns.mass_adm);
  CCTK_VInfo (CCTK_THORNSTRING, "L_tot [G M_sun^2/c]: %g", bin_ns.angu_mom);
  CCTK_VInfo (CCTK_THORNSTRING, "rad1_x_comp [km]:    %g", bin_ns.rad1_x_comp);
  CCTK_VInfo (CCTK_THORNSTRING, "rad1_y [km]:         %g", bin_ns.rad1_y);
  CCTK_VInfo (CCTK_THORNSTRING, "rad1_z [km]:         %g", bin_ns.rad1_z);
  CCTK_VInfo (CCTK_THORNSTRING, "rad1_x_opp [km]:     %g", bin_ns.rad1_x_opp);
  CCTK_VInfo (CCTK_THORNSTRING, "rad2_x_comp [km]:    %g", bin_ns.rad2_x_comp);
  CCTK_VInfo (CCTK_THORNSTRING, "rad2_y [km]:         %g", bin_ns.rad2_y);
  CCTK_VInfo (CCTK_THORNSTRING, "rad2_z [km]:         %g", bin_ns.rad2_z);
  CCTK_VInfo (CCTK_THORNSTRING, "rad2_x_opp [km]:     %g", bin_ns.rad2_x_opp);
  assert (bin_ns.np == npoints);
  
  
  
  CCTK_INFO ("Filling in Cactus grid points");
  
#pragma omp parallel for
  for (int i=0; i<npoints; ++i) {
    
    alp[i] = bin_ns.nnn[i];
    
    betax[i] = -bin_ns.beta_x[i];
    betay[i] = -bin_ns.beta_y[i];
    betaz[i] = -bin_ns.beta_z[i];
    
    CCTK_REAL g[3][3];
    g[0][0] = bin_ns.g_xx[i];
    g[0][1] = bin_ns.g_xy[i];
    g[0][2] = bin_ns.g_xz[i];
    g[1][1] = bin_ns.g_yy[i];
    g[1][2] = bin_ns.g_yz[i];
    g[2][2] = bin_ns.g_zz[i];
    g[1][0] = g[0][1];
    g[2][0] = g[0][2];
    g[2][1] = g[1][2];
    
    CCTK_REAL ku[3][3];
    ku[0][0] = bin_ns.k_xx[i];
    ku[0][1] = bin_ns.k_xy[i];
    ku[0][2] = bin_ns.k_xz[i];
    ku[1][1] = bin_ns.k_yy[i];
    ku[1][2] = bin_ns.k_yz[i];
    ku[2][2] = bin_ns.k_zz[i];
    ku[1][0] = ku[0][1];
    ku[2][0] = ku[0][2];
    ku[2][1] = ku[1][2];
    
    CCTK_REAL k[3][3];
    for (int a=0; a<3; ++a) {
      for (int b=0; b<3; ++b) {
        k[a][b] = 0.0;
        for (int c=0; c<3; ++c) {
          for (int d=0; d<3; ++d) {
            k[a][b] += g[a][c] * g[b][d] * ku[c][d];
          }
        }
      }
    }
    
    gxx[i] = g[0][0];
    gxy[i] = g[0][1];
    gxz[i] = g[0][2];
    gyy[i] = g[1][1];
    gyz[i] = g[1][2];
    gzz[i] = g[2][2];
    
    kxx[i] = ku[0][0] * coord_unit;
    kxy[i] = ku[0][1] * coord_unit;
    kxz[i] = ku[0][2] * coord_unit;
    kyy[i] = ku[1][1] * coord_unit;
    kyz[i] = ku[1][2] * coord_unit;
    kzz[i] = ku[2][2] * coord_unit;
    
    rho[i] = bin_ns.nbar[i] / rho_unit;

    eps[i] = rho[i] * bin_ns.ener_spec[i] / ener_unit;
    
    vel[i          ] = bin_ns.u_euler_x[i] / vel_unit;
    vel[i+  npoints] = bin_ns.u_euler_y[i] / vel_unit;
    vel[i+2*npoints] = bin_ns.u_euler_z[i] / vel_unit;

    if (rho[i] < 1.e-20) {
      rho[i          ] = 1.e-20;
      vel[i          ] = 0.0;
      vel[i+  npoints] = 0.0;
      vel[i+2*npoints] = 0.0;
    }

  } // for i
  
  CCTK_INFO ("Calculating time derivatives of lapse and shift");
  {
    // Angular velocity
    CCTK_REAL const omega = bin_ns.omega * cactusT;
    
    // These initial data assume a helical Killing vector field
    
    if (CCTK_EQUALS (initial_dtlapse, "ID_Bin_NS")) {
      set_dt_from_domega (CCTK_PASS_CTOC, alp, dtalp, omega);
    } else if (CCTK_EQUALS (initial_dtlapse, "none")) {
      // do nothing
    } else {
      CCTK_WARN (CCTK_WARN_ABORT, "internal error");
    }
    
    if (CCTK_EQUALS (initial_dtshift, "ID_Bin_NS")) {
      set_dt_from_domega (CCTK_PASS_CTOC, betax, dtbetax, omega);
      set_dt_from_domega (CCTK_PASS_CTOC, betay, dtbetay, omega);
      set_dt_from_domega (CCTK_PASS_CTOC, betaz, dtbetaz, omega);
    } else if (CCTK_EQUALS (initial_dtshift, "none")) {
      // do nothing
    } else {
      CCTK_WARN (CCTK_WARN_ABORT, "internal error");
    }
  }
  
  
  
  CCTK_INFO ("Done.");
}