Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALSimulation 6.2.0.1-8a6b96f
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
LALSimIMRPrecessingNRSur.h
Go to the documentation of this file.
1/* * Copyright (C) 2017 Jonathan Blackman, Vijay Varma
2 * NRSur7dq2 and NRSur7dq4 NR surrogate models.
3 * Papers: https://arxiv.org/abs/1705.07089, https://arxiv.org/abs/1905.09300.
4 * Based on the python implementation found at:
5 * https://www.black-holes.org/data/surrogates/index.html
6 * which uses the same hdf5 data file.
7 *
8 * This program is free software; you can redistribute it and/or modify
9 * it under the terms of the GNU General Public License as published by
10 * the Free Software Foundation; either version 2 of the License, or
11 * (at your option) any later version.
12 *
13 * This program is distributed in the hope that it will be useful,
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 * GNU General Public License for more details.
17 *
18 * You should have received a copy of the GNU General Public License
19 * along with with program; see the file COPYING. If not, write to the
20 * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
21 * MA 02110-1301 USA
22 */
23
25
26
27/****************************** Constants ***********************************/
28
29///// NRSur7dq2 constants
30
31// These are to rescale the mass ratio fit range from [0.99, 2.01] to [-1, 1].
32static const double NRSUR7DQ2_Q_FIT_OFFSET = -2.9411764705882355; // = (0.5*(2.01 - 0.99)) * (2 / (2.01 - 0.99))
33static const double NRSUR7DQ2_Q_FIT_SLOPE = 1.9607843137254901; // = 2 / (2.01 - 0.99)
34
35// Allow extrapolation in mass ratio and spin beyond these values but throw
36// a warning
37static const double NRSUR7DQ2_Q_MAX_WARN = 2.01;
38static const double NRSUR7DQ2_CHI_MAX_WARN = 0.81;
39
40// Do not allow extrapolation beyond these, unless unlimited_extrapolation=1
41static const double NRSUR7DQ2_Q_MAX = 3.01;
42static const double NRSUR7DQ2_CHI_MAX = 1;
43
44static const double NRSUR7DQ2_START_TIME = -4499.99999999; // The first time node is ever so slightly larger than -4500; this avoids out-of-range issues;
45
46
47///// NRSur7dq4 constants
48
49// These are to rescale log(q) fit range from [-0.01, np.log(4.01)] to [-1, 1].
50// slope = 2./(np.log(4.01) + 0.01)
51static const double NRSUR7DQ4_Q_FIT_SLOPE = 1.4298059216576398 ;
52// offset = np.log(4.01)*1.4298059216576398 - 1
53static const double NRSUR7DQ4_Q_FIT_OFFSET = -0.9857019407834238;
54
55// Allow extrapolation in mass ratio and spin beyond these values but throw
56// a warning
57static const double NRSUR7DQ4_Q_MAX_WARN = 4.01;
58static const double NRSUR7DQ4_CHI_MAX_WARN = 0.81;
59
60// Do not allow extrapolation beyond these, unless unlimited_extrapolation=1
61static const double NRSUR7DQ4_Q_MAX = 6.01;
62static const double NRSUR7DQ4_CHI_MAX = 1;
63
64static const double NRSUR7DQ4_START_TIME = -4299.99999999; // The first time node is ever so slightly larger than -4300; this avoids out-of-range issues;
65
66
67
68
69
70static const int NRSUR_LMAX = 4;
71
72// Surrogate model data, make sure the files are in your LAL_DATA_PATH.
73//
74// The binary data for NRSur7dq2 (NRSur7dq2.h5) is available in lalsuite-extra or at
75// https://www.black-holes.org/surrogates
76// The binary data for NRSur7dq4 (NRSur7dq4_v1.0.h5) is available at:
77// https://git.ligo.org/waveforms/software/lalsuite-waveform-data
78// and on Zenodo at: https://zenodo.org/records/14999310.
79// Get the lalsuite-waveform-data repo or put the data into a location in your
80// LAL_DATA_PATH.
81// The NRSur7dq4 data is also available on CIT at /home/lalsimulation_data and and via CVMFS
82// at /cvmfs/shared.storage.igwn.org/igwn/shared/auxiliary/obs_sci/cbc/waveform/lalsimulation_data
83static const char NRSUR7DQ2_DATAFILE[] = "NRSur7dq2.h5";
84static const char NRSUR7DQ4_DATAFILE[] = "NRSur7dq4_v1.0.h5";
85
86/***********************************************************************************/
87/****************************** Type definitions ***********************************/
88/***********************************************************************************/
89
90/**
91 * Data used in a single scalar fit
92 */
93typedef struct tagFitData {
94 gsl_matrix_long *basisFunctionOrders; /**< matrix of (n_coefs x 7) basis function orders
95 giving the polynomial order in f(q), chiA components, and chiB components. */
96 gsl_vector *coefs; /**< coefficient vector of length n_coefs */
97 int n_coefs; /**< Number of coefficients in the fit */
98} FitData;
99
100/**
101 * Data used in a single vector fit
102 * NOTE: basisFunctionOrders, coefs, componentIndices, and n_coefs are only
103 * used by NRSur7dq2. While fit_data is only used by NRSur7dq4.
104 */
105typedef struct tagVectorFitData {
106 gsl_matrix_long *basisFunctionOrders; /**< matrix of (n_coefs x 7) basis function orders
107 giving the polynomial order in f(q), chiA components, and chiB components. */
108 gsl_vector *coefs; /**< coefficient vector of length n_coefs */
109 gsl_vector_long *componentIndices; /**< Each fit coefficient applies to a single component
110 of the vector; this gives the component indices. */
111 int n_coefs; /**< Number of coefficients in the fit */
112 int vec_dim; /**< Dimension of the vector */
113 FitData **fit_data; /**< Vector of FitData */
115
116/**
117 * Data for a single dynamics node
118 */
119typedef struct tagDynamicsNodeFitData {
120 FitData *omega_data; /**< A fit to the orbital angular frequency */
121 VectorFitData *omega_copr_data; /**< A 2d vector fit for the x and y components of
122 Omega^{coorb}(t) in arxiv 1705.07089 */
123 VectorFitData *chiA_dot_data; /**< A 3d vector fit for the coorbital components of the
124 time derivative of chiA taken in the coprecessing frame */
125 VectorFitData *chiB_dot_data; /**< A 3d vector fit for the coorbital components of the
126 time derivative of chiB taken in the coprecessing frame */
128
129/**
130 * Data for a single waveform data piece.
131 * Waveform data pieces are real time-dependent components of the waveform in the coorbital frame.
132 * See the beginning of section IV of https://arxiv.org/abs/1705.07089
133 */
134typedef struct tagWaveformDataPiece {
135 int n_nodes; /**< Number of empirical nodes */
136 FitData **fit_data; /**< FitData at each empirical node */
137 gsl_matrix *empirical_interpolant_basis; /**< The empirical interpolation matrix */
138 gsl_vector_long *empirical_node_indices; /**< The empirical node indices */
140
141/**
142 * All WaveformDataPieces needed to evaluate all modes with a fixed value of ell.
143 * For m=0 modes we model the real and imaginary parts separately.
144 * For m>0, we model the real and imaginary parts of
145 * X_pm^{ell, m} = frac{1}{2}( h^{ell, m} pm h^{ell, -m} )
146 * where this h is in the coorbital frame.
147 * (see https://arxiv.org/abs/1705.07089)
148 */
149typedef struct tagWaveformFixedEllModeData {
150 int ell; /**< The fixed value of ell */
151 WaveformDataPiece *m0_real_data; /**< The real (ell, 0) mode data piece */
152 WaveformDataPiece *m0_imag_data; /**< The imag (ell, 0) mode data piece */
153 WaveformDataPiece **X_real_plus_data; /**< One Re[X_+] for each 1 <= m <= ell */
154 WaveformDataPiece **X_real_minus_data; /**< One Re[X_-] for each 1 <= m <= ell */
155 WaveformDataPiece **X_imag_plus_data; /**< One Im[X_+] for each 1 <= m <= ell */
156 WaveformDataPiece **X_imag_minus_data; /**< One Im[X_-] for each 1 <= m <= ell */
158
159/**
160 * All data needed by the full surrogate model
161 */
162typedef struct tagPrecessingNRSurData {
163 UINT4 setup; /**< Indicates if this has been initialized */
164 int LMax; /**< Maximum ell mode that will ever be evaluated */
165 gsl_vector *t_ds; /** Vector of the dynamics surrogate node times, not including half times */
166 gsl_vector *t_ds_half_times; /**< t_1/2, t_3/2, and t_5/2 used to start up integration*/
167 gsl_vector *t_coorb; /**< Vector of the coorbital surrogate output times. */
168 DynamicsNodeFitData **ds_node_data; /** A DynamicsNodeFitData for each time in t_ds.*/
169 DynamicsNodeFitData **ds_half_node_data; /** A DynamicsNodeFitData for each time in t_ds_half_times. */
170 WaveformFixedEllModeData **coorbital_mode_data; /** One for each 2 <= ell <= LMax */
171 UINT4 PrecessingNRSurVersion; /**< 0 for NRSur7dq2, 1 for NRSur7dq4 */
173
174/***********************************************************************************/
175/****************************** Function declarations*******************************/
176/***********************************************************************************/
177#ifdef LAL_HDF5_ENABLED
178UNUSED static void NRSur7dq2_Init_LALDATA(void);
179UNUSED static void NRSur7dq4_Init_LALDATA(void);
180
181static int PrecessingNRSur_Init(PrecessingNRSurData *data, LALH5File *file, UINT4 PrecessingNRSurVersion);
182static void PrecessingNRSur_LoadFitData(FitData **fit_data, LALH5File *sub, const char *name);
183static void NRSur7dq4_LoadVectorFitData(VectorFitData **vector_fit_data, LALH5File *sub, const char *name, const size_t size);
184static void PrecessingNRSur_LoadDynamicsNode(DynamicsNodeFitData **ds_node_data, LALH5File *sub, int i, UINT4 PrecessingNRSurVersion);
185static void PrecessingNRSur_LoadCoorbitalEllModes(WaveformFixedEllModeData **coorbital_mode_data, LALH5File *file, int i);
186static void PrecessingNRSur_LoadWaveformDataPiece(LALH5File *sub, WaveformDataPiece **data, bool invert_sign);
187#endif
188
189static bool NRSur7dq2_IsSetup(void);
190static bool NRSur7dq4_IsSetup(void);
191static double ipow(double base, int exponent); // integer powers
192
193static double NRSur7dq2_eval_fit(FitData *data, double *x);
194
196 double *res, // Result
197 VectorFitData *data, // Data for fit
198 double *x // size 7, giving mass ratio q, and dimensionless spin components
199);
200
201static int NRSur7dq4_effective_spins(REAL8 *chiHat, REAL8 *chi_a,
202 const double q, const double chi1z, const double chi2z);
203static double NRSur7dq4_eval_fit(FitData *data, double *x);
204
206 double *res, // Result
207 VectorFitData *data, // Data for fit
208 double *x // size 7, giving mass ratio q, and dimensionless spin components
209);
210
211double PrecessingNRSur_eval_fit(FitData *data, double *x, PrecessingNRSurData *__sur_data);
212
213static void PrecessingNRSur_eval_vector_fit(double *res, VectorFitData *data, double *x, PrecessingNRSurData *__sur_data);
214
216 double chiANorm,
217 double chiBNorm,
218 double *y
219); // Helper to keep spin magnitude constant
220
222 double normA,
223 double normB,
224 gsl_vector **quat,
225 gsl_vector **chiA,
226 gsl_vector **chiB
227);
228
229static void PrecessingNRSur_rotate_spins(gsl_vector **chiA, gsl_vector **chiB, gsl_vector *phi);
230
232 double *x, //Result, length 7
233 double q, // Mass ratio
234 double *y // [q0, qx, qy, qz, phi, chiAx, chiAy, chiAz, chiBx, chiBy, chiBz]
235);
236
238 double *dydt, // Result, length 11
239 double *y, // ODE solution at the current time step
240 double *Omega_coorb_xy, // a form of time derivative of the coprecessing frame
241 double omega, // orbital angular frequency in the coprecessing frame
242 double *chiA_dot, // chiA time derivative
243 double *chiB_dot // chiB time derivative
244);
245
246
247static double cubic_interp(double xout, double *x, double *y);
248static gsl_vector *spline_array_interp(gsl_vector *xout, gsl_vector *x, gsl_vector *y);
249
250static double PrecessingNRSur_get_omega(size_t node_index, double q, double *y0, PrecessingNRSurData *__sur_data);
251
253 REAL8 omega_ref,
254 REAL8 q,
255 REAL8 *chiA0,
256 REAL8 *chiB0,
257 REAL8 *init_quat,
258 REAL8 init_orbphase,
259 PrecessingNRSurData *__sur_data
260);
261
263 double *dydt, // Output: dy/dt evaluated at the ODE time node with index i0. Must have space for 11 entries.
264 int i0, // Time node index. i0=-1, -2, and -3 are used for time nodes 1/2, 3/2, and 5/2 respectively.
265 double q, // Mass ratio
266 double *y, // Current ODE state: [q0, qx, qy, qz, orbphase, chiAx, chiAy, chiAz, chiBx, chiBy, chiBz]
267 PrecessingNRSurData *__sur_data
268);
269
271 double *dtdy, // Output: dy/dt evaluated at time t. Must have space for 11 entries.
272 double t, // Time at which the ODE should be evaluated.
273 double q, // Mass ratio
274 double *y, // Current ODE state
275 PrecessingNRSurData *__sur_data
276);
277
279 double *dynamics_data, // ODE output
280 double t_ref, // reference time. t_ds[i0] will be close to t_ref.
281 double q, // mass ratio
282 double *chiA0, // chiA at t_ref.
283 double *chiB0, // chiB at t_ref.
284 double init_orbphase, // orbital phase at t_ref.
285 double *init_quat, // quaternion at t_ref.
286 double normA, // |chiA|
287 double normB, // |chiB|
288 PrecessingNRSurData *__sur_data
289);
290
292 double *dynamics_data, // ODE output
293 double *time_steps, // Output: first three time steps. Should have size 3.
294 double *dydt0, // Output: dydt at node 0. Should have size 11.
295 double *dydt1, // Output: dydt at node 1. Should have size 11.
296 double *dydt2, // Output: dydt at node 2. Should have size 11.
297 double *dydt3, // Output: dydt at node 3. Should have size 11.
298 double normA, // |chiA|
299 double normB, // |chiB|
300 double q, // mass ratio
301 PrecessingNRSurData *__sur_data
302);
303
305 double *dynamics_data, // ODE output
306 double *time_steps, // Output: first three time steps. Should have size 3.
307 double *dydt0, // Output: dydt at node i0 + 0. Should have size 11.
308 double *dydt1, // Output: dydt at node i0 + 1. Should have size 11.
309 double *dydt2, // Output: dydt at node i0 + 2. Should have size 11.
310 double *dydt3, // Output: dydt at node i0 + 3. Should have size 11.
311 double normA, // |chiA|
312 double normB, // |chiB|
313 double q, // mass ratio
314 int i0, // the node that is already initialized
315 PrecessingNRSurData *__sur_data
316);
317
319 double *dynamics_data, // ODE output
320 double *time_steps, // The first three time steps beginning at i_start.
321 double *dydt0, // dydt at node i_start
322 double *dydt1, // dydt at node i_start + 1
323 double *dydt2, // dydt at node i_start + 2
324 double *dydt3, // dydt at node i_start + 3
325 double normA, // |chiA|
326 double normB, // |chiB|
327 double q, // mass ratio
328 int i_start, // nodes i_start through i_start+3 are already initialized
329 PrecessingNRSurData *__sur_data
330);
331
333 double *dynamics_data,
334 double q,
335 double *chiA0,
336 double *chiB0,
337 double omega_ref,
338 double init_orbphase,
339 double *init_quat,
340 LALDict* LALparams,
341 UINT4 PrecessingNRSurVersion
342);
343
345 gsl_vector *result,
346 double q,
347 gsl_vector **chiA,
348 gsl_vector **chiB,
350 PrecessingNRSurData *__sur_data
351);
352
354
357 double q,
358 double *chiA0,
359 double *chiB0,
360 double omega_ref,
361 double init_orbphase,
362 double *init_quat,
363 LALValue* ModeArray,
364 LALDict* LALparams,
366);
367
369 REAL8 m1, /**< mass of companion 1 (kg) */
370 REAL8 m2, /**< mass of companion 2 (kg) */
371 REAL8 s1x, /**< initial value of S1x */
372 REAL8 s1y, /**< initial value of S1y */
373 REAL8 s1z, /**< initial value of S1z */
374 REAL8 s2x, /**< initial value of S2x */
375 REAL8 s2y, /**< initial value of S2y */
376 REAL8 s2z, /**< initial value of S2z */
377 PrecessingNRSurData *__sur_data /**< 0 for NRSur7dq2, 1 for NRSur7dq4 */
378);
379
380
382 REAL8 *m1, /**< mass of companion 1 (kg) */
383 REAL8 *m2, /**< mass of companion 2 (kg) */
384 REAL8 *s1x, /**< initial value of S1x */
385 REAL8 *s1y, /**< initial value of S1y */
386 REAL8 *s1z, /**< initial value of S1z */
387 REAL8 *s2x, /**< initial value of S2x */
388 REAL8 *s2y, /**< initial value of S2y */
389 REAL8 *s2z /**< initial value of S2z */
390);
struct tagLALH5File LALH5File
static void NRSur7dq2_eval_vector_fit(double *res, VectorFitData *data, double *x)
static const double NRSUR7DQ4_Q_MAX
static const double NRSUR7DQ4_Q_FIT_OFFSET
static const char NRSUR7DQ4_DATAFILE[]
static bool NRSur7dq4_IsSetup(void)
static double cubic_interp(double xout, double *x, double *y)
static void PrecessingNRSur_ds_fit_x(double *x, double q, double *y)
static void PrecessingNRSur_initialize_RK4_with_half_nodes(double *dynamics_data, double *time_steps, double *dydt0, double *dydt1, double *dydt2, double *dydt3, double normA, double normB, double q, PrecessingNRSurData *__sur_data)
static const double NRSUR7DQ4_START_TIME
static const double NRSUR7DQ4_Q_FIT_SLOPE
static const double NRSUR7DQ2_Q_FIT_SLOPE
static void PrecessingNRSur_normalize_results(double normA, double normB, gsl_vector **quat, gsl_vector **chiA, gsl_vector **chiB)
static void PrecessingNRSur_rotate_spins(gsl_vector **chiA, gsl_vector **chiB, gsl_vector *phi)
static bool NRSur7dq2_IsSetup(void)
static void PrecessingNRSur_normalize_y(double chiANorm, double chiBNorm, double *y)
static int PrecessingNRSur_initialize_at_dynamics_node(double *dynamics_data, double t_ref, double q, double *chiA0, double *chiB0, double init_orbphase, double *init_quat, double normA, double normB, PrecessingNRSurData *__sur_data)
static const double NRSUR7DQ2_Q_MAX
static void PrecessingNRSur_get_time_deriv(double *dtdy, double t, double q, double *y, PrecessingNRSurData *__sur_data)
static void PrecessingNRSur_get_time_deriv_from_index(double *dydt, int i0, double q, double *y, PrecessingNRSurData *__sur_data)
static int PrecessingNRSur_IntegrateDynamics(double *dynamics_data, double q, double *chiA0, double *chiB0, double omega_ref, double init_orbphase, double *init_quat, LALDict *LALparams, UINT4 PrecessingNRSurVersion)
static void NRSur7dq4_eval_vector_fit(double *res, VectorFitData *data, double *x)
static const double NRSUR7DQ4_Q_MAX_WARN
static const double NRSUR7DQ2_CHI_MAX_WARN
static const int NRSUR_LMAX
static const double NRSUR7DQ2_CHI_MAX
static const double NRSUR7DQ2_Q_MAX_WARN
static const double NRSUR7DQ4_CHI_MAX_WARN
static const double NRSUR7DQ2_Q_FIT_OFFSET
static gsl_vector * spline_array_interp(gsl_vector *xout, gsl_vector *x, gsl_vector *y)
static PrecessingNRSurData * PrecessingNRSur_core(MultiModalWaveform **h, double q, double *chiA0, double *chiB0, double omega_ref, double init_orbphase, double *init_quat, LALValue *ModeArray, LALDict *LALparams, Approximant approximant)
static void PrecessingNRSur_integrate_AB4(double *dynamics_data, double *time_steps, double *dydt0, double *dydt1, double *dydt2, double *dydt3, double normA, double normB, double q, int i_start, PrecessingNRSurData *__sur_data)
static bool PrecessingNRSur_switch_labels_if_needed(REAL8 *m1, REAL8 *m2, REAL8 *s1x, REAL8 *s1y, REAL8 *s1z, REAL8 *s2x, REAL8 *s2y, REAL8 *s2z)
static const double NRSUR7DQ4_CHI_MAX
static void PrecessingNRSur_eval_data_piece(gsl_vector *result, double q, gsl_vector **chiA, gsl_vector **chiB, WaveformDataPiece *data, PrecessingNRSurData *__sur_data)
static int PrecessingNRSur_initialize_RK4(double *dynamics_data, double *time_steps, double *dydt0, double *dydt1, double *dydt2, double *dydt3, double normA, double normB, double q, int i0, PrecessingNRSurData *__sur_data)
static void PrecessingNRSur_eval_vector_fit(double *res, VectorFitData *data, double *x, PrecessingNRSurData *__sur_data)
static REAL8 PrecessingNRSur_get_t_ref(REAL8 omega_ref, REAL8 q, REAL8 *chiA0, REAL8 *chiB0, REAL8 *init_quat, REAL8 init_orbphase, PrecessingNRSurData *__sur_data)
static double PrecessingNRSur_get_omega(size_t node_index, double q, double *y0, PrecessingNRSurData *__sur_data)
static const double NRSUR7DQ2_START_TIME
static const char NRSUR7DQ2_DATAFILE[]
static PrecessingNRSurData * PrecessingNRSur_LoadData(Approximant approximant)
static void PrecessingNRSur_assemble_dydt(double *dydt, double *y, double *Omega_coorb_xy, double omega, double *chiA_dot, double *chiB_dot)
static int NRSur7dq4_effective_spins(REAL8 *chiHat, REAL8 *chi_a, const double q, const double chi1z, const double chi2z)
const char * name
double i
Definition: bh_ringdown.c:118
sigmaKerr data[0]
double REAL8
uint32_t UINT4
Approximant
Enum that specifies the PN approximant to be used in computing the waveform.
static double NRSur7dq4_eval_fit(FitData *data, double *x)
static double ipow(double base, int exponent)
Helper function for integer powers.
double PrecessingNRSur_eval_fit(FitData *data, double *x, PrecessingNRSurData *__sur_data)
static double PrecessingNRSur_StartFrequency(REAL8 m1, REAL8 m2, REAL8 s1x, REAL8 s1y, REAL8 s1z, REAL8 s2x, REAL8 s2y, REAL8 s2z, PrecessingNRSurData *__sur_data)
Computes the starting frequency of the NRSur7dq2 or NRSur7dq4 waveform approximant when evaluated usi...
static double NRSur7dq2_eval_fit(FitData *data, double *x)
static const INT4 q
Data for a single dynamics node.
FitData * omega_data
A fit to the orbital angular frequency.
VectorFitData * chiB_dot_data
A 3d vector fit for the coorbital components of the time derivative of chiB taken in the coprecessing...
VectorFitData * chiA_dot_data
A 3d vector fit for the coorbital components of the time derivative of chiA taken in the coprecessing...
VectorFitData * omega_copr_data
A 2d vector fit for the x and y components of Omega^{coorb}(t) in arxiv 1705.07089.
Data used in a single scalar fit.
int n_coefs
Number of coefficients in the fit.
gsl_matrix_long * basisFunctionOrders
matrix of (n_coefs x 7) basis function orders giving the polynomial order in f(q),...
gsl_vector * coefs
coefficient vector of length n_coefs
Structure for a multi-modal waveform.
All data needed by the full surrogate model.
gsl_vector * t_coorb
Vector of the coorbital surrogate output times.
UINT4 PrecessingNRSurVersion
One for each 2 <= ell <= LMax.
WaveformFixedEllModeData ** coorbital_mode_data
A DynamicsNodeFitData for each time in t_ds_half_times.
gsl_vector * t_ds_half_times
Vector of the dynamics surrogate node times, not including half times.
DynamicsNodeFitData ** ds_half_node_data
A DynamicsNodeFitData for each time in t_ds.
UINT4 setup
Indicates if this has been initialized.
DynamicsNodeFitData ** ds_node_data
int LMax
Maximum ell mode that will ever be evaluated.
Data used in a single vector fit NOTE: basisFunctionOrders, coefs, componentIndices,...
int n_coefs
Number of coefficients in the fit.
FitData ** fit_data
Vector of FitData.
gsl_matrix_long * basisFunctionOrders
matrix of (n_coefs x 7) basis function orders giving the polynomial order in f(q),...
int vec_dim
Dimension of the vector.
gsl_vector_long * componentIndices
Each fit coefficient applies to a single component of the vector; this gives the component indices.
gsl_vector * coefs
coefficient vector of length n_coefs
Data for a single waveform data piece.
gsl_matrix * empirical_interpolant_basis
The empirical interpolation matrix.
gsl_vector_long * empirical_node_indices
The empirical node indices.
int n_nodes
Number of empirical nodes.
FitData ** fit_data
FitData at each empirical node.
All WaveformDataPieces needed to evaluate all modes with a fixed value of ell.
WaveformDataPiece * m0_real_data
The real (ell, 0) mode data piece.
WaveformDataPiece ** X_real_minus_data
One Re[X_-] for each 1 <= m <= ell.
WaveformDataPiece ** X_imag_plus_data
One Im[X_+] for each 1 <= m <= ell.
WaveformDataPiece ** X_imag_minus_data
One Im[X_-] for each 1 <= m <= ell.
int ell
The fixed value of ell.
WaveformDataPiece * m0_imag_data
The imag (ell, 0) mode data piece.
WaveformDataPiece ** X_real_plus_data
One Re[X_+] for each 1 <= m <= ell.