aboutsummaryrefslogtreecommitdiff
path: root/src/Bin_NS.cc
blob: 621f41f37047731594bd18407427ce569b999591 (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
/* (c) 2009 Erik Schnetter
 * (c) 2010 Frank Loeffler */

#include <cstdio>
#include <cassert>
#include <vector>
#include <ios>

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

#include <bin_ns.h>
#include <unites.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 Meudon_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.
  // Be aware: these are the constants Lorene uses. They do differ from other
  // conventions, but they gave the best results in some tests.

  CCTK_REAL const c_light  = Unites::c_si;      // speed of light [m/s]
  CCTK_REAL const nuc_dens = Unites::rhonuc_si; // Nuclear density as used in Lorene units [kg/m^3]
  CCTK_REAL const G_grav   = Unites::g_si;      // gravitational constant [m^3/kg/s^2]
  CCTK_REAL const M_sun    = Unites::msol_si;   // 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 (~1.477)
  CCTK_REAL const rho_unit   = cactusM / pow(cactusL,3); // from kg/m^3


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

  try {
  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);
  double K = bin_ns.kappa_poly1 * pow(c_light, 6.0) /
             ( pow(G_grav, 3.0) * M_sun * M_sun *
               pow(nuc_dens, bin_ns.gamma_poly1-1.) );
  CCTK_VInfo (CCTK_THORNSTRING, "K [ET unit]:         %.15g", K);

  assert (bin_ns.np == npoints);

  CCTK_INFO ("Filling in Cactus grid points");

#pragma omp parallel for
  for (int i=0; i<npoints; ++i) {

    if (CCTK_EQUALS(initial_lapse, "Meudon_Bin_NS")) { 
      alp[i] = bin_ns.nnn[i];
    }

    if (CCTK_EQUALS(initial_shift, "Meudon_Bin_NS")) { 
      betax[i] = -bin_ns.beta_x[i];
      betay[i] = -bin_ns.beta_y[i];
      betaz[i] = -bin_ns.beta_z[i];
    }

    if (CCTK_EQUALS(initial_data, "Meudon_Bin_NS")) {
      gxx[i] = bin_ns.g_xx[i];
      gxy[i] = bin_ns.g_xy[i];
      gxz[i] = bin_ns.g_xz[i];
      gyy[i] = bin_ns.g_yy[i];
      gyz[i] = bin_ns.g_yz[i];
      gzz[i] = bin_ns.g_zz[i];

      kxx[i] = bin_ns.k_xx[i] * coord_unit;
      kxy[i] = bin_ns.k_xy[i] * coord_unit;
      kxz[i] = bin_ns.k_xz[i] * coord_unit;
      kyy[i] = bin_ns.k_yy[i] * coord_unit;
      kyz[i] = bin_ns.k_yz[i] * coord_unit;
      kzz[i] = bin_ns.k_zz[i] * coord_unit;
    }

    if (CCTK_EQUALS(initial_data, "Meudon_Bin_NS")) {
      rho[i] = bin_ns.nbar[i] / rho_unit;
      eps[i] = bin_ns.ener_spec[i];
      // The following would recompute eps using the polytopic EOS, but
      // setting eps directly seems to work as well
      // eps[i] = K * pow(rho[i], bin_ns.gamma_poly1-1.) / (bin_ns.gamma_poly1-1.);

      // TODO: we should really use some EOS calls for the pressure, but for the
      //       moment this works with polytropes at least
      press[i] = K * pow(rho[i], bin_ns.gamma_poly1);

      vel[i          ] = bin_ns.u_euler_x[i];
      vel[i+  npoints] = bin_ns.u_euler_y[i];
      vel[i+2*npoints] = bin_ns.u_euler_z[i];

      // Especially the velocity is set to strange values outside of the
      // matter region, so take care of this in the following way
      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;
        eps[i          ] = K * pow(rho[i], bin_ns.gamma_poly1-1.) / (bin_ns.gamma_poly1-1.);
        press[i        ] = K * pow(rho[i], bin_ns.gamma_poly1);
      }
    }

  } // for i

  {
    // Angular velocity
    CCTK_REAL const omega = bin_ns.omega * cactusT;

    // These initial data assume a helical Killing vector field

    if (CCTK_EQUALS(initial_lapse, "Meudon_Bin_NS")) { 
      if (CCTK_EQUALS (initial_dtlapse, "Meudon_Bin_NS")) {
        CCTK_INFO ("Calculating time derivatives of lapse");
        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_shift, "Meudon_Bin_NS")) { 
      if (CCTK_EQUALS (initial_dtshift, "Meudon_Bin_NS")) {
        CCTK_INFO ("Calculating time derivatives of shift");
        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.");
  } catch (ios::failure e) {
    CCTK_VWarn (CCTK_WARN_ABORT, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Could not read initial data from file '%s': %s", filename, e.what());
  }
}