LALSimulation  5.4.0.1-fe68b98
LALSimIMRPrecessingNRSur.c
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 files.
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 
24 #ifdef __GNUC__
25 #define UNUSED __attribute__ ((unused))
26 #else
27 #define UNUSED
28 #endif
29 
30 #include <stdio.h>
31 #include <stdlib.h>
32 #include <math.h>
33 #include <complex.h>
34 #include <stdbool.h>
35 #include <alloca.h>
36 #include <string.h>
37 #include <libgen.h>
38 
39 #include <gsl/gsl_errno.h>
40 #include <gsl/gsl_bspline.h>
41 #include <gsl/gsl_blas.h>
42 #include <gsl/gsl_min.h>
43 #include <gsl/gsl_spline.h>
44 #include <gsl/gsl_complex_math.h>
45 #include <lal/Units.h>
46 #include <lal/SeqFactories.h>
47 #include <lal/LALConstants.h>
48 #include <lal/XLALError.h>
49 #include <lal/FrequencySeries.h>
50 #include <lal/Date.h>
51 #include <lal/StringInput.h>
52 #include <lal/Sequence.h>
53 #include <lal/LALStdio.h>
54 #include <lal/FileIO.h>
55 #include <lal/SphericalHarmonics.h>
56 #include <lal/LALSimInspiral.h>
57 #include <lal/LALSimSphHarmMode.h>
58 #include <lal/LALSimIMR.h>
59 #include <lal/H5FileIO.h>
60 
63 
64 #include <lal/LALConfig.h>
65 #ifdef LAL_PTHREAD_LOCK
66 #include <pthread.h>
67 #endif
68 
69 
70 #ifdef LAL_PTHREAD_LOCK
71 static pthread_once_t NRSur7dq2_is_initialized = PTHREAD_ONCE_INIT;
72 #endif
73 
74 #ifdef LAL_PTHREAD_LOCK
75 static pthread_once_t NRSur7dq4_is_initialized = PTHREAD_ONCE_INIT;
76 #endif
77 
78 
79 /**
80  * Global surrogate data.
81  * This data will be loaded at most once. Any executable which calls
82  * NRSur7dq2_Init_LALDATA or NRSur7dq4_Init_LALDATA directly or by calling any
83  * XLAL function will have a memory leak according to valgrind, because we
84  * never free this memory.
85  */
88 
89 /**
90  * @addtogroup LALSimNRSur7dq2_c
91  * @{
92  *
93  * @name Routines for NR surrogate models "NRSur7dq2" and "NRSur7dq4"
94  * @{
95  *
96  * @author Jonathan Blackman, Vijay Varma
97  *
98  * @brief C code for NRSur7dq2 and NRSur7dq4 NR surrogate waveform models.
99  *
100  *
101  *
102  * NRSur7dq2:
103  * This is a fully precessing time domain model including all subdominant modes up to ell=4.
104  * See Blackman et al \cite Blackman:q27d for details.
105  * Any studies that use this waveform model should include a reference to that paper.
106  * Using this model requires the file lalsuite-extra/data/lalsimulation/NRSur7dq2.h5
107  * Make sure your $LAL_DATA_PATH points to lalsuite-extra/data/lalsimulation/.
108  * The lalsuite-extra commit hash at the time of review was 77613e7f5f01d5ea11829ded5677783cafc0d298
109  *
110  * @note The range of validity of the model is:
111  * * Mass ratios 1 <= q <= 2
112  * * Spin magnitudes |chi_i| <= 0.8
113  * * Total time before merger <= 4500M, which in practice leads to a parameter-dependent lower bound for fmin.
114  *
115  * @note Additional notes:
116  * * Passing in a non-trivial ModeArray controls which co-orbital frame modes are evaluated.
117  * * A python version of this model can be installed with "pip install NRSur7dq2".
118  * * This lalsimulation implementation has been verified to agree with version 1.0.3 up to
119  * very small differences due to slightly differing time interpolation methods.
120  * * Note that for conventions to agree with ChooseTDWaveform (and XLALSimInspiralNRSur7dq2Polarizations),
121  * you must pass use_lalsimulation_conventions=True when calling the pip version of NRSur7dq2.
122  *
123  * @review NRSur7dq2 model and routines reviewed by Sebastian Khan, Harald Pfeiffer, Geraint Pratten, and Michael PĆ¼rrer.
124  * Reviewees were Jonathan Blackman, Scott Field, and Vijay Varma.
125  * The review page can be found at https://git.ligo.org/waveforms/reviews/nrsur/wikis/home
126  *
127  *
128  *
129  * NRSur7dq4:
130  * This is a q=4 extension of NRSur7dq2.
131  * See Varma et al. (arxiv:1905.09300) for details.
132  * Any studies that use this waveform model should include a reference to that
133  * paper.
134  * Using this model requires the file lalsuite-extra/data/lalsimulation/NRSur7dq4.h5
135  * Make sure your $LAL_DATA_PATH points to lalsuite-extra/data/lalsimulation/.
136  */
137 
138 
139 /***********************************************************************************/
140 /****************************** Function Definitions *******************************/
141 /***********************************************************************************/
142 
143 
144 
145 
146 /**
147  * This needs to be called once, before __lalsim_NRSur7dq2_data is used.
148  * It finds the hdf5 data file with the NRSur7dq2 data and calls PrecessingNRSur_Init.
149  */
150 static void NRSur7dq2_Init_LALDATA(void) {
151  if (NRSur7dq2_IsSetup()) return;
152 
154  if (path==NULL)
155  XLAL_ERROR_VOID(XLAL_EIO, "Unable to resolve data file %s in $LAL_DATA_PATH\n", NRSUR7DQ2_DATAFILE);
156  char *dir = dirname(path);
157  size_t size = strlen(dir) + strlen(NRSUR7DQ2_DATAFILE) + 2;
158  char *file_path = XLALMalloc(size);
159  snprintf(file_path, size, "%s/%s", dir, NRSUR7DQ2_DATAFILE);
160 
161  LALH5File *file = XLALH5FileOpen(file_path, "r");
162 
163  // 0 is for NRSur7dq2
165 
166  if (ret != XLAL_SUCCESS)
167  XLAL_ERROR_VOID(XLAL_FAILURE, "Failure loading data from %s\n", file_path);
168 
169  XLALFree(path);
170  XLALFree(file_path);
171 }
172 
173 /**
174  * This needs to be called once, before __lalsim_NRSur7dq4_data is used.
175  * It finds the hdf5 data file with the NRSur7dq4 data and calls PrecessingNRSur_Init.
176  */
177 static void NRSur7dq4_Init_LALDATA(void) {
178 
179  if (NRSur7dq4_IsSetup()) return;
180 
182  if (path==NULL)
183  XLAL_ERROR_VOID(XLAL_EIO, "Unable to resolve data file %s in $LAL_DATA_PATH\n", NRSUR7DQ4_DATAFILE);
184  char *dir = dirname(path);
185  size_t size = strlen(dir) + strlen(NRSUR7DQ4_DATAFILE) + 2;
186  char *file_path = XLALMalloc(size);
187  snprintf(file_path, size, "%s/%s", dir, NRSUR7DQ4_DATAFILE);
188 
189  LALH5File *file = XLALH5FileOpen(file_path, "r");
190 
191  // 1 is for NRSur7dq4
193 
194  if (ret != XLAL_SUCCESS)
195  XLAL_ERROR_VOID(XLAL_FAILURE, "Failure loading data from %s\n", file_path);
196 
197  XLALFree(path);
198  XLALFree(file_path);
199 }
200 
201 
202 /**
203  * Initialize a PrecessingNRSurData structure from an open hdf5 file.
204  * This will typically only be called once, from NRSur7dq2_Init_LALDATA or
205  * NRSur7dq4_Init_LALDATA.
206  */
208  PrecessingNRSurData *data, /**< Output: Surrogate data structure. */
209  LALH5File *file, /**< hdf5 file with surrogate data. */
210  UINT4 PrecessingNRSurVersion /**< 0 for NRSur7dq2, 1 for NRSur7dq4 */
211 ) {
212 
213  size_t i;
214 
215  if (data->setup) {
216  XLALPrintError("You tried to setup a NRSurrogate model that was already initialized. Ignoring.\n");
217  return XLAL_FAILURE;
218  }
219 
220  // Get the dynamics time nodes
221  gsl_vector *t_ds_with_halves = NULL;
222  ReadHDF5RealVectorDataset(file, "t_ds", &t_ds_with_halves);
223  gsl_vector *t_ds = gsl_vector_alloc(t_ds_with_halves->size - 3);
224  gsl_vector *t_ds_half_times = gsl_vector_alloc(3);
225  for (i=0; i < 3; i++) {
226  gsl_vector_set(t_ds, i, gsl_vector_get(t_ds_with_halves, 2*i));
227  gsl_vector_set(t_ds_half_times, i, gsl_vector_get(t_ds_with_halves, 2*i + 1));
228  }
229  for (i=3; i < t_ds->size; i++) {
230  gsl_vector_set(t_ds, i, gsl_vector_get(t_ds_with_halves, i+3));
231  }
232  gsl_vector_free(t_ds_with_halves);
233  data->t_ds = t_ds;
234  data->t_ds_half_times = t_ds_half_times;
235 
236  // Load the dynamics time node data. The nodes are stored in HDF5 groups "ds_node_%d"%(i)
237  // where i=0, 1, 2, ...
238  // These INCLUDE the half-nodes at t_1/2, t_3/2, and t_5/2, so the (node time, index) pairs are
239  // (t_0, 0), (t_1/2, 1), (t_1, 2), (t_3/2, 3), (t_2, 4), (t_5/2, 5), (t_3, 6), (t_4, 7), ...
240  // (t_n, n+3) for n >= 3
241  DynamicsNodeFitData **ds_node_data = XLALMalloc( (t_ds->size) * sizeof(*ds_node_data));
242  DynamicsNodeFitData **ds_half_node_data = XLALMalloc( 3 * sizeof(*ds_node_data) );
243  for (i=0; i < (t_ds->size); i++) ds_node_data[i] = NULL;
244  for (i=0; i < 3; i++) ds_half_node_data[i] = NULL;
245  LALH5File *sub;
246  char *sub_name = XLALMalloc(20);
247  int j;
248  for (i=0; i < (t_ds->size); i++) {
249  if (i < 3) {j = 2*i;} else {j = i+3;}
250  snprintf(sub_name, 20, "ds_node_%d", j);
251  sub = XLALH5GroupOpen(file, sub_name);
252  PrecessingNRSur_LoadDynamicsNode(ds_node_data, sub, i, PrecessingNRSurVersion);
253 
254  if (i < 3) {
255  snprintf(sub_name, 20, "ds_node_%d", j+1);
256  sub = XLALH5GroupOpen(file, sub_name);
257  PrecessingNRSur_LoadDynamicsNode(ds_half_node_data, sub, i, PrecessingNRSurVersion);
258  }
259  }
260  XLALFree(sub_name);
261  data->ds_node_data = ds_node_data;
262  data->ds_half_node_data = ds_half_node_data;
263 
264  // Get the coorbital time array
265  gsl_vector *t_coorb = NULL;
266  ReadHDF5RealVectorDataset(file, "t_coorb", &t_coorb);
267  data->t_coorb = t_coorb;
268 
269  // Load coorbital waveform surrogate data
270  WaveformFixedEllModeData **coorbital_mode_data = XLALMalloc( (NRSUR_LMAX - 1) * sizeof(*coorbital_mode_data) );
271  for (int ell_idx=0; ell_idx < NRSUR_LMAX-1; ell_idx++) {
272  PrecessingNRSur_LoadCoorbitalEllModes(coorbital_mode_data, file, ell_idx);
273  }
274  data->coorbital_mode_data = coorbital_mode_data;
275 
276  data->LMax = NRSUR_LMAX;
277  data->setup = 1;
278 
279  return XLAL_SUCCESS;
280 }
281 
282 /**
283  * Loads a single fit for NRSur7dq2 or NRSur7dq4.
284  */
286  FitData **fit_data, /**< Output: Data struct for fit data. Should be NULL; Will malloc space and load data into it. */
287  LALH5File *sub, /**< Subgroup containing fit data. */
288  const char *name /**< fit name. */
289 ) {
290  *fit_data = XLALMalloc(sizeof(FitData));
291 
292  const size_t str_size = 30;
293  char *tmp_name = XLALMalloc(str_size);
294  UNUSED size_t nwritten;
295 
296  nwritten = snprintf(tmp_name, str_size, "%s_coefs", name);
297  XLAL_CHECK_ABORT(nwritten < str_size);
298  (*fit_data)->coefs = NULL;
299  ReadHDF5RealVectorDataset(sub, tmp_name, &((*fit_data)->coefs));
300 
301  nwritten = snprintf(tmp_name, str_size, "%s_bfOrders", name);
302  XLAL_CHECK_ABORT(nwritten < str_size);
303  (*fit_data)->basisFunctionOrders = NULL;
304  ReadHDF5LongMatrixDataset(sub, tmp_name,
305  &((*fit_data)->basisFunctionOrders));
306 
307  (*fit_data)->n_coefs = (*fit_data)->coefs->size;
308 }
309 
310 /**
311  * Loads a vector fit for NRSur7dq4.
312  *
313  * For this model, vector fits are constructed as a vector of scalar fits.
314  */
316  VectorFitData **vector_fit_data, /**< Output: Data struct for vector fit data. Should be NULL; Will malloc space and load data into it. */
317  LALH5File *sub, /**< Subgroup containing fit data. */
318  const char *name, /**< fit name. */
319  const size_t size /**< size of vector. */
320 ) {
321  const size_t str_size = 20;
322  char *tmp_name = XLALMalloc(str_size);
323  UNUSED size_t nwritten;
324 
325  *vector_fit_data = XLALMalloc(sizeof(VectorFitData));
326  (*vector_fit_data)->vec_dim = size;
327  (*vector_fit_data)->fit_data = XLALMalloc(size * sizeof(FitData *) );
328 
329  for (size_t i=0; i<size; i++) {
330  nwritten = snprintf(tmp_name, str_size, "%s_%zu", name, i);
331  XLAL_CHECK_ABORT(nwritten < str_size);
332  FitData *fit_data = NULL;
333  PrecessingNRSur_LoadFitData(&fit_data, sub, tmp_name);
334  (*vector_fit_data)->fit_data[i] = fit_data;
335  }
336 }
337 
338 /**
339  * Loads the data for a single dynamics node into a DynamicsNodeFitData struct.
340  * This is only called during the initialization of the surrogate data through PrecessingNRSur_Init.
341  */
343  DynamicsNodeFitData **ds_node_data, /**< Entry i should be NULL; Will malloc space and load data into it. */
344  LALH5File *sub, /**< Subgroup containing data for dynamics node i. */
345  int i, /**< Dynamics node index. */
346  UINT4 PrecessingNRSurVersion /**< 0 for NRSur7dq2, 1 for NRSur7dq4 */
347 ) {
348  ds_node_data[i] = XLALMalloc( sizeof(*ds_node_data[i]) );
349 
350 
351  // omega
352  FitData *omega_data = NULL;
353  PrecessingNRSur_LoadFitData(&omega_data, sub, "omega");
354  ds_node_data[i]->omega_data = omega_data;
355 
356  // For NRSur7dq2 vector fits are done for omega_orb, chiA_dot and chiB_dot
357  if (PrecessingNRSurVersion == 0){
358  // omega_copr
359  VectorFitData *omega_copr_data = XLALMalloc(sizeof(VectorFitData));
360  omega_copr_data->coefs = NULL;
361  omega_copr_data->basisFunctionOrders = NULL;
362  omega_copr_data->componentIndices = NULL;
363  ReadHDF5RealVectorDataset(sub, "omega_orb_coefs", &(omega_copr_data->coefs));
364  ReadHDF5LongMatrixDataset(sub, "omega_orb_bfOrders", &(omega_copr_data->basisFunctionOrders));
365  ReadHDF5LongVectorDataset(sub, "omega_orb_bVecIndices", &(omega_copr_data->componentIndices));
366  omega_copr_data->n_coefs = omega_copr_data->coefs->size;
367  omega_copr_data->vec_dim = 2;
368  ds_node_data[i]->omega_copr_data = omega_copr_data;
369 
370  // chiA_dot
371  VectorFitData *chiA_dot_data = XLALMalloc(sizeof(VectorFitData));
372  chiA_dot_data->coefs = NULL;
373  chiA_dot_data->basisFunctionOrders = NULL;
374  chiA_dot_data->componentIndices = NULL;
375  ReadHDF5RealVectorDataset(sub, "chiA_coefs", &(chiA_dot_data->coefs));
376  ReadHDF5LongMatrixDataset(sub, "chiA_bfOrders", &(chiA_dot_data->basisFunctionOrders));
377  ReadHDF5LongVectorDataset(sub, "chiA_bVecIndices", &(chiA_dot_data->componentIndices));
378  chiA_dot_data->n_coefs = chiA_dot_data->coefs->size;
379  chiA_dot_data->vec_dim = 3;
380  ds_node_data[i]->chiA_dot_data = chiA_dot_data;
381 
382  // chiB_dot
383  // One chiB_dot node has 0 coefficients, and ReadHDF5RealVectorDataset fails.
384  VectorFitData *chiB_dot_data = XLALMalloc(sizeof(VectorFitData));
385  chiB_dot_data->coefs = NULL;
386  chiB_dot_data->basisFunctionOrders = NULL;
387  chiB_dot_data->componentIndices = NULL;
388 
389  UINT4Vector *dimLength;
390  size_t n;
391  LALH5Dataset *dset;
392  dset = XLALH5DatasetRead(sub, "chiB_coefs");
393  dimLength = XLALH5DatasetQueryDims(dset);
394  n = dimLength->data[0];
395  if (n==0) {
396  chiB_dot_data->n_coefs = 0;
397  } else {
398  ReadHDF5RealVectorDataset(sub, "chiB_coefs", &(chiB_dot_data->coefs));
399  ReadHDF5LongMatrixDataset(sub, "chiB_bfOrders", &(chiB_dot_data->basisFunctionOrders));
400  ReadHDF5LongVectorDataset(sub, "chiB_bVecIndices", &(chiB_dot_data->componentIndices));
401  chiB_dot_data->n_coefs = chiB_dot_data->coefs->size;
402  }
403  chiB_dot_data->vec_dim = 3;
404  ds_node_data[i]->chiB_dot_data = chiB_dot_data;
405 
406  // For NRSur7dq4 the vector fits are done simply as a vector of scalar
407  // fits, so we just need to loop over the indices.
408  } else if (PrecessingNRSurVersion == 1) {
409  // omega_copr
410  VectorFitData *omega_copr_data = NULL;
411  NRSur7dq4_LoadVectorFitData(&omega_copr_data, sub, "omega_orb", 2);
412  ds_node_data[i]->omega_copr_data = omega_copr_data;
413 
414  // chiA_dot
415  VectorFitData *chiA_dot_data = NULL;
416  NRSur7dq4_LoadVectorFitData(&chiA_dot_data, sub, "chiA", 3);
417  ds_node_data[i]->chiA_dot_data = chiA_dot_data;
418 
419  // chiB_dot
420  VectorFitData *chiB_dot_data = NULL;
421  NRSur7dq4_LoadVectorFitData(&chiB_dot_data, sub, "chiB", 3);
422  ds_node_data[i]->chiB_dot_data = chiB_dot_data;
423  }
424 }
425 
426 
427 /**
428  * Load the WaveformFixedEllModeData from file for a single value of ell.
429  * This is only called during the initialization of the surrogate data through PrecessingNRSur_Init.
430  */
432  WaveformFixedEllModeData **coorbital_mode_data, /**< Entry i should be NULL; will malloc space and load data into it.*/
433  LALH5File *file, /**< The open hdf5 file */
434  int i /**< The index of coorbital_mode_data. Equivalently, ell-2. */
435 ) {
436  WaveformFixedEllModeData *mode_data = XLALMalloc( sizeof(*coorbital_mode_data[i]) );
437  mode_data->ell = i+2;
438 
439  LALH5File *sub;
440  int str_size = 30; // Enough for L with 15 digits...
441  char *sub_name = XLALMalloc(str_size);
442 
443  // Real part of m=0 mode
444  snprintf(sub_name, str_size, "hCoorb_%d_0_real", i+2);
445  sub = XLALH5GroupOpen(file, sub_name);
446  PrecessingNRSur_LoadWaveformDataPiece(sub, &(mode_data->m0_real_data), false);
447 
448  // Imag part of m=0 mode
449  snprintf(sub_name, str_size, "hCoorb_%d_0_imag", i+2);
450  sub = XLALH5GroupOpen(file, sub_name);
451  PrecessingNRSur_LoadWaveformDataPiece(sub, &(mode_data->m0_imag_data), false);
452 
453  // NOTE:
454  // In the paper https://arxiv.org/abs/1705.07089, Eq. 16 uses
455  // h_\pm^{\ell, m} = \frac{1}{2} \left( h_\mathrm{coorb}^{\ell, m} \pm h_\mathrm{coorb}^{\ell, -m\, *} \right)
456  // and we often denote h_\pm^{\ell, m} == X_\pm^{\ell, m} in this file and other documentation,
457  // but when actually building the surrogate the (\ell, -m) mode was taken as the reference mode, so we store
458  // Y_\pm^{\ell, m} = \frac{1}{2} \left( h_\mathrm{coorb}^{\ell, -m} \pm h_\mathrm{coorb}^{\ell, m\, *} \right)
459  // Note that we still have everything we want:
460  // X_\pm^{\ell, m} = \pm Y_\pm^{\ell, m\, *}
461  // or:
462  // Re[X_+] = Re[Y_+],
463  // Im[X_+] = -Im[Y_+],
464  // Re[X_-] = -Re[Y_-],
465  // Im[X_-] = Im[Y_-]
466  // We work with X rather than Y in this file, so we need the minus signs when loading from Y.
467  mode_data->X_real_plus_data = XLALMalloc( (i+2) * sizeof(WaveformDataPiece *) );
468  mode_data->X_real_minus_data = XLALMalloc( (i+2) * sizeof(WaveformDataPiece *) );
469  mode_data->X_imag_plus_data = XLALMalloc( (i+2) * sizeof(WaveformDataPiece *) );
470  mode_data->X_imag_minus_data = XLALMalloc( (i+2) * sizeof(WaveformDataPiece *) );
471  for (int m=1; m<=(i+2); m++) {
472  snprintf(sub_name, str_size, "hCoorb_%d_%d_Re+", i+2, m);
473  sub = XLALH5GroupOpen(file, sub_name);
474  PrecessingNRSur_LoadWaveformDataPiece(sub, &(mode_data->X_real_plus_data[m-1]), false);
475  snprintf(sub_name, str_size, "hCoorb_%d_%d_Re-", i+2, m);
476  sub = XLALH5GroupOpen(file, sub_name);
477  PrecessingNRSur_LoadWaveformDataPiece(sub, &(mode_data->X_real_minus_data[m-1]), true);
478  snprintf(sub_name, str_size, "hCoorb_%d_%d_Im+", i+2, m);
479  sub = XLALH5GroupOpen(file, sub_name);
480  PrecessingNRSur_LoadWaveformDataPiece(sub, &(mode_data->X_imag_plus_data[m-1]), true);
481  snprintf(sub_name, str_size, "hCoorb_%d_%d_Im-", i+2, m);
482  sub = XLALH5GroupOpen(file, sub_name);
483  PrecessingNRSur_LoadWaveformDataPiece(sub, &(mode_data->X_imag_minus_data[m-1]), false);
484  }
485  XLALFree(sub_name);
486  coorbital_mode_data[i] = mode_data;
487 }
488 
489 /**
490  * Loads a single NRSur coorbital waveform data piece from file into a WaveformDataPiece.
491  * This is only called during the initialization of the surrogate data through PrecessingNRSur_Init.
492  */
494  LALH5File *sub, /**< HDF5 group containing data for this waveform data piece */
495  WaveformDataPiece **data, /**< Output - *data should be NULL. Space will be allocated. */
496  bool invert_sign /**< If true, multiply the empirical interpolation matrix by -1. */
497 ) {
498  *data = XLALMalloc(sizeof(WaveformDataPiece));
499 
500  gsl_matrix *EI_basis = NULL;
501  ReadHDF5RealMatrixDataset(sub, "EIBasis", &EI_basis);
502  if (invert_sign) {
503  gsl_matrix_scale(EI_basis, -1);
504  }
505  (*data)->empirical_interpolant_basis = EI_basis;
506 
507  gsl_vector_long *node_indices = NULL;
508  ReadHDF5LongVectorDataset(sub, "nodeIndices", &node_indices);
509  (*data)->empirical_node_indices = node_indices;
510 
511  int n_nodes = (*data)->empirical_node_indices->size;
512  (*data)->n_nodes = n_nodes;
513  (*data)->fit_data = XLALMalloc( n_nodes * sizeof(FitData *) );
514 
515  LALH5File *nodeModelers = XLALH5GroupOpen(sub, "nodeModelers");
516  int str_size = 20; // Enough for L with 11 digits...
517  char *sub_name = XLALMalloc(str_size);
518  for (int i=0; i<n_nodes; i++) {
519  FitData *node_data = XLALMalloc(sizeof(FitData));
520  node_data->coefs = NULL;
521  node_data->basisFunctionOrders = NULL;
522  snprintf(sub_name, str_size, "coefs_%d", i);
523  ReadHDF5RealVectorDataset(nodeModelers, sub_name, &(node_data->coefs));
524  snprintf(sub_name, str_size, "bfOrders_%d", i);
525  ReadHDF5LongMatrixDataset(nodeModelers, sub_name, &(node_data->basisFunctionOrders));
526  node_data->n_coefs = node_data->coefs->size;
527  (*data)->fit_data[i] = node_data;
528  }
529 }
530 
531 
532 /**
533  * Helper function which returns whether or not the global NRSur7dq2 surrogate data has been initialized.
534  */
535 static bool NRSur7dq2_IsSetup(void) {
536  if(__lalsim_NRSur7dq2_data.setup) return true;
537  else return false;
538 }
539 
540 /**
541  * Helper function which returns whether or not the global NRSur7dq4 surrogate data has been initialized.
542  */
543 static bool NRSur7dq4_IsSetup(void) {
544  if(__lalsim_NRSur7dq4_data.setup) return true;
545  else return false;
546 }
547 
548 /**
549  * Helper function for integer powers
550  */
551 REAL8 ipow(REAL8 base, int exponent) {
552  if (exponent == 0) return 1.0;
553  REAL8 res = base;
554  while (exponent > 1) {
555  res = res*base;
556  exponent -= 1;
557  }
558  return res;
559 }
560 
561 
562 /*
563  * Evaluate a NRSur7dq2 scalar fit.
564  * The fit result is given by
565  * \sum_{i=1}^{n} c_i * \prod_{j=1}^7 B_j(k_{i, j}; x_j)
566  * where i runs over fit coefficients, j runs over the 7 dimensional parameter
567  * space, and B_j is a basis function, taking an integer order k_{i, j} and
568  * the parameter component x_j. For this surrogate, B_j are monomials in the spin
569  * components, and monomials in an affine transformation of the mass ratio.
570  */
572  FitData *data, /**< Data for fit */
573  REAL8 *x /**< size 7, giving mass ratio q, and dimensionless spin components */
574 ) {
575  REAL8 x_powers[22]; // 3 per spin component, 4 for mass ratio
576  REAL8 res = 0.0;
577  REAL8 prod;
578  int i, j;
579 
580  // The fits were constructed using this rather than using q directly
582 
583  // Compute powers of components of x
584  for (i=0; i<22; i++) {
585  if (i%7==0) {
586  x_powers[i] = ipow(q_fit, i/7);
587  } else {
588  x_powers[i] = ipow(x[i%7], i/7);
589  }
590  }
591 
592  // Sum up fit terms
593  for (i=0; i < data->n_coefs; i++) {
594  // Initialize with q basis function:
595  prod = x_powers[7 * gsl_matrix_long_get(data->basisFunctionOrders, i, 0)];
596  // Multiply with spin basis functions:
597  for (j=1; j<7; j++) {
598  prod *= x_powers[7 * gsl_matrix_long_get(data->basisFunctionOrders, i, j) + j];
599  }
600  res += gsl_vector_get(data->coefs, i) * prod;
601  }
602 
603  return res;
604 }
605 
606 /*
607  * This is very similar to NRSur7dq2_eval_fit except that the result is a
608  * 2d or 3d vector instead of a scalar. Each fit coefficient now applies to
609  * just a single component of the result.
610  */
612  REAL8 *res, /**< Result */
613  VectorFitData *data, /**< Data for fit */
614  REAL8 *x /**< size 7, giving mass ratio q, and dimensionless spin components */
615 ) {
616  REAL8 x_powers[22]; // 3 per spin component, 4 for mass ratio
617  REAL8 prod;
618  int i, j;
619 
620  // Initialize the result
621  for (i=0; i < data->vec_dim; i++) {
622  res[i] = 0.0;
623  }
624 
625  // The fits were constructed using this rather than using q directly
627 
628  // Compute powers of components of x
629  for (i=0; i<22; i++) {
630  if (i%7==0) {
631  x_powers[i] = ipow(q_fit, i/7);
632  } else {
633  x_powers[i] = ipow(x[i%7], i/7);
634  }
635  }
636 
637  // Sum up fit terms
638  for (i=0; i < data->n_coefs; i++) {
639  // Initialize with q basis function:
640  prod = x_powers[7 * gsl_matrix_long_get(data->basisFunctionOrders, i, 0)];
641  // Multiply with spin basis functions:
642  for (j=1; j<7; j++) {
643  prod *= x_powers[7 * gsl_matrix_long_get(data->basisFunctionOrders, i, j) + j];
644  }
645  res[gsl_vector_long_get(data->componentIndices, i)] += gsl_vector_get(data->coefs, i) * prod;
646  }
647 }
648 
649 /*
650  * Computes effective spins chiHat and chi_a.
651  * chiHat is defined in Eq.(3) of 1508.07253.
652  * and chi_a = (chi1z - chi2z)/2.
653  */
655  REAL8 *chiHat, /**< Output: chiHat */
656  REAL8 *chi_a, /**< Output: chi_a */
657  const REAL8 q, /**< Mass ratio >= 1 */
658  const REAL8 chi1z, /**< Dimensionless z-spin of heavier BH */
659  const REAL8 chi2z /**< Dimensionless z-spin of lighter BH */
660 ) {
661  const REAL8 eta = q/(1.+q)/(1.+q);
662  const REAL8 chi_wtAvg = (q*chi1z+chi2z)/(1+q);
663  REAL8 chiHat_val = (chi_wtAvg - 38.*eta/113.*(chi1z
664  + chi2z))/(1. - 76.*eta/113.);
665  REAL8 chi_a_val = (chi1z - chi2z)/2.;
666  *chiHat = chiHat_val;
667  *chi_a = chi_a_val;
668  return XLAL_SUCCESS;
669 }
670 
671 /*
672  * Evaluate a NRSur7dq4 scalar fit.
673  */
675  FitData *data, /**< Data for fit */
676  REAL8 *x /**< size 7, giving mass ratio q, and dimensionless spin components */
677 ) {
678  REAL8 x_powers[22]; // 3 per spin component, 4 for mass ratio
679  REAL8 res = 0.0;
680  REAL8 prod;
681  int i, j;
682 
683  // get effective spins chiHat and chi_a
684  // chiHat is defined in Eq.(3) of 1508.07253.
685  // and chi_a = (chi1z - chi2z)/2.
686  REAL8 chiHat, chi_a;
687  NRSur7dq4_effective_spins(&chiHat, &chi_a, x[0], x[3], x[6]);
688 
689  // Convert from [q, chi1x, chi1y, chi1z, chi2x, chi2y, chi2z]
690  // to [log(q), chi1x, chi1y, chiHat, chi2x, chi2y, chi_a]
691  REAL8 fit_params[7];
692  fit_params[0] = log(x[0]);
693  fit_params[1] = x[1];
694  fit_params[2] = x[2];
695  fit_params[3] = chiHat;
696  fit_params[4] = x[4];
697  fit_params[5] = x[5];
698  fit_params[6] = chi_a;
699 
700  // The fits were constructed using this rather than using q directly
702  + NRSUR7DQ4_Q_FIT_SLOPE*fit_params[0];
703 
704  // Compute powers of components of fit_params
705  for (i=0; i<22; i++) {
706  if (i%7==0) {
707  x_powers[i] = ipow(q_fit, i/7);
708  } else {
709  x_powers[i] = ipow(fit_params[i%7], i/7);
710  }
711  }
712 
713  // Sum up fit terms
714  for (i=0; i < data->n_coefs; i++) {
715  // Initialize with q basis function:
716  prod = x_powers[7 * gsl_matrix_long_get(data->basisFunctionOrders, i, 0)];
717  // Multiply with spin basis functions:
718  for (j=1; j<7; j++) {
719  prod *= x_powers[7 * gsl_matrix_long_get(data->basisFunctionOrders, i, j) + j];
720  }
721  res += gsl_vector_get(data->coefs, i) * prod;
722  }
723 
724  return res;
725 }
726 
727 /*
728  * This is very similar to NRSur7dq4_eval_fit except that the result is a
729  * 2d or 3d vector instead of a scalar. Each fit coefficient now applies to
730  * just a single component of the result.
731  */
733  REAL8 *res, /**< Result */
734  VectorFitData *data, /**< Data for fit */
735  REAL8 *x /**< size 7, giving mass ratio q, and dimensionless spin components */
736 ) {
737  // loop over vector indices
738  for (int i=0; i < data->vec_dim; i++) {
739  res[i] = NRSur7dq4_eval_fit((*data).fit_data[i], x);
740  }
741 }
742 
743 
744 /*
745  * Wrapper for NRSur7dq2_eval_fit and NRSur7dq4_eval_fit
746  */
748  FitData *data, /**< Data for fit */
749  REAL8 *x, /**< size 7, giving mass ratio q, and dimensionless spin components */
750  PrecessingNRSurData *__sur_data /**< Loaded surrogate data */
751 ) {
752  if (__sur_data->PrecessingNRSurVersion == 0) {
753  return NRSur7dq2_eval_fit(data, x);
754  } else if (__sur_data->PrecessingNRSurVersion == 1) {
755  return NRSur7dq4_eval_fit(data, x);
756  } else {
757  XLAL_ERROR_REAL8(XLAL_FAILURE, "Only 0 or 1 are currently allowed for PrecessingNRSurVersion\n");
758  }
759 }
760 
761 /*
762  * Wrapper for NRSur7dq2_eval_vector_fit and NRSur7dq4_eval_vector_fit
763  */
765  REAL8 *res, /**< Result */
766  VectorFitData *data, /**< Data for fit */
767  REAL8 *x, /**< size 7, giving mass ratio q, and dimensionless spin components */
768  PrecessingNRSurData *__sur_data /**< Loaded surrogate data */
769 ) {
770  if (__sur_data->PrecessingNRSurVersion == 0) {
771  return NRSur7dq2_eval_vector_fit(res, data, x);
772  } else if (__sur_data->PrecessingNRSurVersion == 1) {
773  return NRSur7dq4_eval_vector_fit(res, data, x);
774  } else {
775  XLAL_ERROR_VOID(XLAL_FAILURE, "Only 0 or 1 are currently allowed for PrecessingNRSurVersion\n");
776  }
777 }
778 
779 
780 /* During the ODE integration, the norm of the spins will change due to
781  * integration errors and fit modeling errors. Keep them normalized.
782  * Normalizes in-place
783  */
785  REAL8 chiANorm, /**< ||vec{chi}_A|| */
786  REAL8 chiBNorm, /**< ||vec{chi}_B|| */
787  REAL8 *y /**< [q0, qx, qy, qz, phi, chiAx, chiAy, chiAz, chiBx, chiBy, chiBz] */
788 ) {
789 
790  REAL8 nQ, nA, nB, sum;
791  int i;
792 
793  // Compute current norms
794  sum = 0.0;
795  for (i=0; i<4; i++) {
796  sum += y[i]*y[i];
797  }
798  nQ = sqrt(sum);
799 
800  sum = 0.0;
801  for (i=5; i<8; i++) {
802  sum += y[i]*y[i];
803  }
804  nA = sqrt(sum);
805 
806  sum = 0.0;
807  for (i=8; i<11; i++) {
808  sum += y[i]*y[i];
809  }
810  nB = sqrt(sum);
811 
812  // Compute normalized output
813  for (i=0; i<4; i++) {
814  y[i] = y[i] / nQ;
815  }
816  y[4] = y[4];
817  for (i=5; i<8; i++) {
818  y[i] = y[i] * chiANorm / nA;
819  }
820  for (i=8; i<11; i++) {
821  y[i] = y[i] * chiBNorm / nB;
822  }
823 }
824 
825 /*
826  * Similar to normalize_y, but components given individually as arrays with n samples
827  */
829  REAL8 normA, /**< ||vec{chi}_A|| */
830  REAL8 normB, /**< ||vec{chi}_B|| */
831  gsl_vector **quat, /**< The four quaternion time-dependent components */
832  gsl_vector **chiA, /**< Time-dependent components of chiA */
833  gsl_vector **chiB /**< Time-dependent components of chiB */
834 ) {
835  REAL8 nA, nB, nQ;
836  int i;
837  int n = quat[0]->size;
838  REAL8 *chiAx = chiA[0]->data;
839  REAL8 *chiAy = chiA[1]->data;
840  REAL8 *chiAz = chiA[2]->data;
841  REAL8 *chiBx = chiB[0]->data;
842  REAL8 *chiBy = chiB[1]->data;
843  REAL8 *chiBz = chiB[2]->data;
844  REAL8 *q0 = quat[0]->data;
845  REAL8 *qx = quat[1]->data;
846  REAL8 *qy = quat[2]->data;
847  REAL8 *qz = quat[3]->data;
848 
849  if (normA > 0.0) {
850  for (i=0; i<n; i++) {
851  nA = sqrt(chiAx[i]*chiAx[i] + chiAy[i]*chiAy[i] + chiAz[i]*chiAz[i]);
852  chiAx[i] *= normA / nA;
853  chiAy[i] *= normA / nA;
854  chiAz[i] *= normA / nA;
855  }
856  }
857 
858  if (normB > 0.0) {
859  for (i=0; i<n; i++) {
860  nB = sqrt(chiBx[i]*chiBx[i] + chiBy[i]*chiBy[i] + chiBz[i]*chiBz[i]);
861  chiBx[i] *= normB / nB;
862  chiBy[i] *= normB / nB;
863  chiBz[i] *= normB / nB;
864  }
865  }
866 
867  for (i=0; i<n; i++) {
868  nQ = sqrt(q0[i]*q0[i] + qx[i]*qx[i] + qy[i]*qy[i] + qz[i]*qz[i]);
869  q0[i] /= nQ;
870  qx[i] /= nQ;
871  qy[i] /= nQ;
872  qz[i] /= nQ;
873  }
874 }
875 
876 /**
877  * Transforms chiA and chiB from the coprecessing frame to the coorbital frame
878  * using the orbital phase.
879  */
881  gsl_vector **chiA, /**< 3 time-dependent components of chiA in the coorbital frame */
882  gsl_vector **chiB, /**< 3 time-dependent components of chiB in the coorbital frame */
883  gsl_vector *phi_vec /**< The time-dependent orbital phase */
884 ) {
885  int i;
886  int n = phi_vec->size;
887  REAL8 sp, cp, tmp;
888  REAL8 *phi = phi_vec->data;
889 
890  REAL8 *chix = chiA[0]->data;
891  REAL8 *chiy = chiA[1]->data;
892  for (i=0; i<n; i++) {
893  cp = cos(phi[i]);
894  sp = sin(phi[i]);
895  tmp = chix[i];
896  chix[i] = tmp*cp + chiy[i]*sp;
897  chiy[i] = -1*tmp*sp + chiy[i]*cp;
898  }
899 
900  chix = chiB[0]->data;
901  chiy = chiB[1]->data;
902  for (i=0; i<n; i++) {
903  cp = cos(phi[i]);
904  sp = sin(phi[i]);
905  tmp = chix[i];
906  chix[i] = tmp*cp + chiy[i]*sp;
907  chiy[i] = -1*tmp*sp + chiy[i]*cp;
908  }
909 }
910 
911 /*
912  * When integrating the dynamics ODE, we need to evaluate fits.
913  * This helper function computes the fit inputs from the current ODE solution.
914  * The spin components of y are in the coprecessing frame, but the spin
915  * components of x are in the coorbital frame.
916  */
918  REAL8 *x, /**< Result, length 7 */
919  REAL8 q, /**< Mass ratio */
920  REAL8 *y /**< [q0, qx, qy, qz, phi, chiAx, chiAy, chiAz, chiBx, chiBy, chiBz] */
921 ) {
922  REAL8 sp = sin(y[4]);
923  REAL8 cp = cos(y[4]);
924 
925  // q
926  x[0] = q;
927 
928  // chiA, transformed to the coorbital frame
929  x[1] = y[5]*cp + y[6]*sp;
930  x[2] = -1*y[5]*sp + y[6]*cp;
931  x[3] = y[7];
932 
933  // chiB, transformed to the coorbital frame
934  x[4] = y[8]*cp + y[9]*sp;
935  x[5] = -1*y[8]*sp + y[9]*cp;
936  x[6] = y[10];
937 }
938 
939 /*
940  * After evaluating fits, we have many time derivatives, with components in
941  * the coorbital frame. dydt has components in the coprecessing frame, so we
942  * transform them, and we also compute the time derivative of the unit
943  * quaternions using the coprecessing components of Omega
944  */
946  REAL8 *dydt, /**< Result, length 11 */
947  REAL8 *y, /**< ODE solution at the current time step */
948  REAL8 *Omega_coorb_xy, /**< a form of time derivative of the coprecessing frame */
949  REAL8 omega, /**< orbital angular frequency in the coprecessing frame */
950  REAL8 *chiA_dot, /**< chiA time derivative */
951  REAL8 *chiB_dot /**< chiB time derivative */
952 ) {
953  REAL8 omega_quat_x, omega_quat_y;
954 
955  REAL8 cp = cos(y[4]);
956  REAL8 sp = sin(y[4]);
957 
958  // Quaternion derivative:
959  // Omega = 2 * quat^{-1} * dqdt -> dqdt = 0.5 * quat * omegaQuat, where
960  // omegaQuat = [0, Omega_copr_x, Omega_copr_y, 0]
961  omega_quat_x = Omega_coorb_xy[0]*cp - Omega_coorb_xy[1]*sp;
962  omega_quat_y = Omega_coorb_xy[0]*sp + Omega_coorb_xy[1]*cp;
963  dydt[0] = (-0.5)*y[1]*omega_quat_x - 0.5*y[2]*omega_quat_y;
964  dydt[1] = (-0.5)*y[3]*omega_quat_y + 0.5*y[0]*omega_quat_x;
965  dydt[2] = 0.5*y[3]*omega_quat_x + 0.5*y[0]*omega_quat_y;
966  dydt[3] = 0.5*y[1]*omega_quat_y - 0.5*y[2]*omega_quat_x;
967 
968  // Orbital phase derivative
969  dydt[4] = omega;
970 
971  // Spin derivatives
972  dydt[5] = chiA_dot[0]*cp - chiA_dot[1]*sp;
973  dydt[6] = chiA_dot[0]*sp + chiA_dot[1]*cp;
974  dydt[7] = chiA_dot[2];
975  dydt[8] = chiB_dot[0]*cp - chiB_dot[1]*sp;
976  dydt[9] = chiB_dot[0]*sp + chiB_dot[1]*cp;
977  dydt[10] = chiB_dot[2];
978 }
979 
980 /**
981  * Cubic interpolation of 4 data points
982  * This gives a much closer result to scipy.interpolate.InterpolatedUnivariateSpline than using gsl_interp_cspline
983  * (see comment in spline_array_interp)
984  */
986  REAL8 xout, /**< The target x value */
987  REAL8 *x, /**< The x values of the points to interpolate. Length 4, must be increasing. */
988  REAL8 *y /**< The y values of the points to interpolate. Length 4. */
989 ) {
990  gsl_interp_accel *acc = gsl_interp_accel_alloc();
991  gsl_spline *interpolant = gsl_spline_alloc(gsl_interp_polynomial, 4);
992  gsl_spline_init(interpolant, x, y, 4);
993  REAL8 res = gsl_spline_eval(interpolant, xout, acc);
994  gsl_spline_free(interpolant);
995  gsl_interp_accel_free(acc);
996  return res;
997 }
998 
999 /**
1000  * Do cubic spline interpolation using a gsl_interp_cspline.
1001  * Results differ from scipy.interpolate.InterpolatedUnivariateSpline due to different boundary conditions.
1002  * This difference leads to small differences between this implementation and the python
1003  * implementation, especially near the start and end of the waveform.
1004  */
1005 static gsl_vector *spline_array_interp(
1006  gsl_vector *xout, /**< The vector of points onto which we want to interpolate. */
1007  gsl_vector *x, /**< The x values of the data to interpolate. */
1008  gsl_vector *y /**< The y values of the data to interpolate. */
1009 ) {
1010  gsl_interp_accel *acc = gsl_interp_accel_alloc();
1011  gsl_spline *spline = gsl_spline_alloc(gsl_interp_cspline, x->size);
1012  gsl_spline_init(spline, x->data, y->data, x->size);
1013 
1014  gsl_vector *res = gsl_vector_alloc(xout->size);
1015  REAL8 tmp;
1016  for (size_t i=0; i<xout->size; i++) {
1017  tmp = gsl_spline_eval(spline, gsl_vector_get(xout, i), acc);
1018  gsl_vector_set(res, i, tmp);
1019  }
1020  gsl_spline_free(spline);
1021  gsl_interp_accel_free(acc);
1022  return res;
1023 }
1024 
1025 /**
1026  * Computes the orbital angular frequency at a dynamics node.
1027  */
1029  size_t node_index, /**< The index of the dynamics node. */
1030  REAL8 q, /**< The mass ratio. */
1031  REAL8 *y0, /**< The value of the ODE state y = [q0, qx, qy, qz, orbphase, chiAx, chiAy, chiAz, chiBx, chiBy, chiBz] */
1032  PrecessingNRSurData *__sur_data /**< Loaded surrogate data */
1033 ){
1034  REAL8 x[7];
1036  FitData *data = __sur_data->ds_node_data[node_index]->omega_data;
1037  REAL8 omega = PrecessingNRSur_eval_fit(data, x, __sur_data);
1038  return omega;
1039 }
1040 
1041 /**
1042  * Computes a reference time from a reference orbital angular frequency.
1043  */
1045  REAL8 omega_ref, /**< reference orbital angular frequency */
1046  REAL8 q, /**< mass ratio */
1047  REAL8 *chiA0, /**< chiA at reference point */
1048  REAL8 *chiB0, /**< chiB at reference point */
1049  REAL8 *init_quat, /**< coprecessing frame quaternion at reference point */
1050  REAL8 init_orbphase, /**< orbital phase at reference point */
1051  PrecessingNRSurData *__sur_data /**< Loaded surrogate data */
1052 ) {
1053 
1054  REAL8 START_TIME = 100;
1055  if (__sur_data->PrecessingNRSurVersion == 0) {
1056  START_TIME = NRSUR7DQ2_START_TIME;
1057  } else if (__sur_data->PrecessingNRSurVersion == 1) {
1058  START_TIME = NRSUR7DQ4_START_TIME;
1059  }
1060 
1061  if (fabs(omega_ref) < 1.e-10) {
1062  XLAL_PRINT_WARNING("WARNING: Treating omega_ref = 0 as a flag to use t_ref = t_0 = %.2f", START_TIME);
1063  return START_TIME;
1064  }
1065 
1066  // omega_ref <= 0.201 guarantees omega is smaller than the final and merger omegas for all masses and spins.
1067  if (omega_ref > 0.201) XLAL_ERROR_REAL8(XLAL_EINVAL,
1068  "Reference frequency omega_ref=%0.4f > 0.2, too large!\n", omega_ref);
1069 
1070  REAL8 y0[11];
1071  int i;
1072  for (i=0; i<4; i++) y0[i] = init_quat[i];
1073  y0[4] = init_orbphase;
1074  for (i=0; i<3; i++) {
1075  y0[i+5] = chiA0[i];
1076  y0[i+8] = chiB0[i];
1077  }
1078 
1079  REAL8 omega_min = PrecessingNRSur_get_omega(0, q, y0, __sur_data);
1080  if (omega_ref < omega_min) {
1082  "Got omega_ref or omega_low=%0.4f smaller than the minimum omega=%0.4f for this configuration!", omega_ref, omega_min);
1083  }
1084 
1085  // i0=0 is a lower bound; find the first index where omega > omega_ref, and the previous index will have omega <= omega_ref.
1086  size_t imax = 1;
1087  REAL8 omega_max = PrecessingNRSur_get_omega(imax, q, y0, __sur_data);
1088  gsl_vector *t_ds = __sur_data->t_ds;
1089  while ((omega_max <= omega_ref) && (imax < t_ds->size)) {
1090  imax += 1;
1091  omega_min = omega_max;
1092  omega_max = PrecessingNRSur_get_omega(imax, q, y0, __sur_data);
1093  }
1094 
1095  if (omega_max <= omega_ref) {
1096  XLAL_ERROR_REAL8(XLAL_FAILURE, "Tried all nodes and still have omega=%0.4f <= omega_ref=%0.4f\n", omega_min, omega_ref);
1097  }
1098 
1099  // Do a linear interpolation between omega_min and omega_max
1100  REAL8 t_min = gsl_vector_get(t_ds, imax-1);
1101  REAL8 t_max = gsl_vector_get(t_ds, imax);
1102  REAL8 t_ref = (t_min * (omega_max - omega_ref) + t_max * (omega_ref - omega_min)) / (omega_max - omega_min);
1103  return t_ref;
1104 }
1105 
1106 /**
1107  * Compute dydt at a given dynamics node, where y is the numerical solution to the dynamics ODE.
1108  */
1110  REAL8 *dydt, /**< Output: dy/dt evaluated at the ODE time node with index i0. Must have space for 11 entries. */
1111  int i0, /**< Time node index. i0=-1, -2, and -3 are used for time nodes 1/2, 3/2, and 5/2 respectively. */
1112  REAL8 q, /**< Mass ratio */
1113  REAL8 *y, /**< Current ODE state: [q0, qx, qy, qz, orbphase, chiAx, chiAy, chiAz, chiBx, chiBy, chiBz] */
1114  PrecessingNRSurData *__sur_data /**< Loaded surrogate data */
1115 ) {
1116  // Setup fit variables
1117  REAL8 x[7];
1119 
1120  // Get fit data
1121  DynamicsNodeFitData *ds_node;
1122  if (i0 >= 0) {
1123  ds_node = __sur_data->ds_node_data[i0];
1124  } else {
1125  ds_node = __sur_data->ds_half_node_data[-1*i0 - 1];
1126  }
1127 
1128  // Evaluate fits
1129  REAL8 omega, Omega_coorb_xy[2], chiA_dot[3], chiB_dot[3];
1130  omega = PrecessingNRSur_eval_fit(ds_node->omega_data, x, __sur_data);
1131  PrecessingNRSur_eval_vector_fit(Omega_coorb_xy, ds_node->omega_copr_data,
1132  x, __sur_data);
1134  x, __sur_data);
1136  x, __sur_data);
1137  PrecessingNRSur_assemble_dydt(dydt, y, Omega_coorb_xy, omega, chiA_dot, chiB_dot);
1138 }
1139 
1140 /**
1141  * Compute dydt at any time by evaluating dydt at 4 nearby dynamics nodes and
1142  * using cubic spline interpolation to evaluate at the desired time.
1143  */
1145  REAL8 *dydt, /**< Output: dy/dt evaluated at time t. Must have space for 11 entries. */
1146  REAL8 t, /**< Time at which the ODE should be evaluated. */
1147  REAL8 q, /**< Mass ratio */
1148  REAL8 *y, /**< Current ODE state */
1149  PrecessingNRSurData *__sur_data /**< Loaded surrogate data */
1150 ) {
1151  // Make sure we are within the valid range
1152  gsl_vector *t_ds = __sur_data->t_ds;
1153  int i1 = 0;
1154  int imax = t_ds->size - 1;
1155  REAL8 t1 = gsl_vector_get(t_ds, i1);
1156  REAL8 tmax = gsl_vector_get(t_ds, imax);
1157  if (t < t1 || t > tmax) {
1158  XLAL_ERROR_VOID(XLAL_FAILURE, "Tried to get dydt at t=%0.12f, not in valid range [%0.12f, %0.12f]\n", t, t1, tmax);
1159  }
1160 
1161  REAL8 t2 = gsl_vector_get(t_ds, i1+1);
1162  // Put t2 slightly larger than t
1163  while (t2 < t) {
1164  i1 += 1;
1165  t1 = t2;
1166  t2 = gsl_vector_get(t_ds, i1+1);
1167  }
1168 
1169  // Do cubic spline interpolation using 4 data points, cenetered if possible
1170  int i0 = i1-1;
1171  if (i0 < 0) i0 = 0;
1172  if (i0 > imax-3) i0 = imax-3;
1173  REAL8 times[4], derivs[4], dydt0[11], dydt1[11], dydt2[11], dydt3[11];
1174  int j;
1175  for (j=0; j<4; j++) times[j] = gsl_vector_get(t_ds, i0+j);
1176  PrecessingNRSur_get_time_deriv_from_index(dydt0, i0, q, y, __sur_data);
1177  PrecessingNRSur_get_time_deriv_from_index(dydt1, i0+1, q, y, __sur_data);
1178  PrecessingNRSur_get_time_deriv_from_index(dydt2, i0+2, q, y, __sur_data);
1179  PrecessingNRSur_get_time_deriv_from_index(dydt3, i0+3, q, y, __sur_data);
1180 
1181  for (j=0; j<11; j++) {
1182  derivs[0] = dydt0[j];
1183  derivs[1] = dydt1[j];
1184  derivs[2] = dydt2[j];
1185  derivs[3] = dydt3[j];
1186  dydt[j] = cubic_interp(t, times, derivs);
1187  }
1188 
1189 }
1190 
1191 /**
1192  * Initialize the dynamics ODE at a dynamics node.
1193  * Given t_ref, finds the nearest dynamics node and integrates the initial conditions at t_ref
1194  * a tiny bit to obtain the ODE state at the dynamics node.
1195  */
1197  REAL8 *dynamics_data, /**< ODE output */
1198  REAL8 t_ref, /**< reference time. t_ds[i0] will be close to t_ref. */
1199  REAL8 q, /**< mass ratio */
1200  REAL8 *chiA0, /**< chiA at t_ref. */
1201  REAL8 *chiB0, /**< chiB at t_ref. */
1202  REAL8 init_orbphase, /**< orbital phase at t_ref. */
1203  REAL8 *init_quat, /**< quaternion at t_ref. */
1204  REAL8 normA, /**< |chiA| */
1205  REAL8 normB, /**< |chiB| */
1206  PrecessingNRSurData *__sur_data /**< Loaded surrogate data */
1207 ) {
1208  int imin, imax, i0, j;
1209  REAL8 tmin, tmax, t0, dt, y_ref[11], dydt_ref[11], *node_data;
1210  gsl_vector *t_ds;
1211 
1212  // First find i0, the closest dynamics time node to time t_ref, with a binary search
1213  t_ds = __sur_data->t_ds;
1214  imin = 0;
1215  imax = t_ds->size - 1;
1216  tmin = gsl_vector_get(t_ds, imin);
1217  tmax = gsl_vector_get(t_ds, imax);
1218  while (imax - imin > 1) {
1219  i0 = (imax + imin)/2;
1220  t0 = gsl_vector_get(t_ds, i0);
1221  if (t0 > t_ref) {
1222  imax = i0;
1223  tmax = t0;
1224  } else {
1225  imin = i0;
1226  tmin = t0;
1227  }
1228  }
1229 
1230  // Now we have tmin <= t_ref < tmax, and imax = imin+1.
1231  // Step towards the closest of tmin, tmax.
1232  if (fabs(t_ref - tmin) < fabs(tmax - t_ref)) {
1233  i0 = imin;
1234  dt = tmin - t_ref;
1235  } else {
1236  i0 = imax;
1237  dt = tmax - t_ref;
1238  }
1239 
1240  // Setup y at t_ref
1241  for (j=0; j<4; j++) y_ref[j] = init_quat[j];
1242  y_ref[4] = init_orbphase;
1243  for (j=0; j<3; j++) {
1244  y_ref[5+j] = chiA0[j];
1245  y_ref[8+j] = chiB0[j];
1246  }
1247 
1248  // Compute dydt at t_ref
1249  PrecessingNRSur_get_time_deriv(dydt_ref, t_ref, q, y_ref, __sur_data);
1250 
1251  // modify y_ref to be y at t0 = t_ds[i0]
1252  for (j=0; j<11; j++) y_ref[j] += dt * dydt_ref[j];
1253  PrecessingNRSur_normalize_y(normA, normB, y_ref);
1254 
1255  // transfer results to ODE output data structures
1256  node_data = dynamics_data + i0*11; // Point to the output at i0
1257  for (j=0; j<11; j++) node_data[j] = y_ref[j];
1258 
1259  return i0;
1260 }
1261 
1262 /**
1263  * Initializes the AB4 ODE system from the surrogate start time by taking 3 RK4 steps,
1264  * making use of the three additional half-time-step nodes for the RK4 substeps.
1265  * This is the recommended way to initialize the AB4 ODE system - the additional nodes
1266  * are there to increase accuracy during initialization.
1267  * The drawback is that we are forced to accept fRef to be the
1268  * surrogate start frequency, which will depend on masses and spins.
1269  * The ODE system is for the vector of 11 quantities:
1270  * [q0, qx, qy, qz, varphi, chiAx, chiAy, chiAz, chiBx, chiBy, chiBz]
1271  * where (q0, qx, qy, qz) is the coprecessing frame quaternion, varphi is the
1272  * orbital phase, and the spin components are in the coprecessing frame.
1273  */
1275  REAL8 *dynamics_data, /**< A pointer to the start of the ODE output; the first 11 entries
1276  (for node 0) should have already been set.*/
1277  REAL8 *time_steps, /**< Output: first three time steps. Should have size 3. */
1278  REAL8 *dydt0, /**< Output: dydt at node 0. Should have size 11. */
1279  REAL8 *dydt1, /**< Output: dydt at node 1. Should have size 11. */
1280  REAL8 *dydt2, /**< Output: dydt at node 2. Should have size 11. */
1281  REAL8 *dydt3, /**< Output: dydt at node 3. Should have size 11. */
1282  REAL8 normA, /**< |chiA| */
1283  REAL8 normB, /**< |chiB| */
1284  REAL8 q, /**< mass ratio */
1285  PrecessingNRSurData *__sur_data /**< Loaded surrogate data */
1286 ) {
1287  gsl_vector *t_ds = __sur_data->t_ds;
1288  REAL8 t1, t2;
1289  int i, j;
1290  REAL8 k2[11], k3[11], k4[11]; // dydt for all but the first RK4 substep (which will be one of dydt0, dydt1, dydt2)
1291  REAL8 *node_data;
1292  REAL8 y_tmp[11];
1293  REAL8 *dydt[3];
1294  dydt[0] = dydt0;
1295  dydt[1] = dydt1;
1296  dydt[2] = dydt2;
1297 
1298  // Setup time steps
1299  t1 = gsl_vector_get(t_ds, 0);
1300  for (i=0; i<3; i++) {
1301  t2 = gsl_vector_get(t_ds, i+1);
1302  time_steps[i] = t2 - t1;
1303  t1 = t2;
1304  }
1305 
1306  // Three steps of RK4.
1307  for (i=0; i<3; i++ ) {
1308 
1309  // Initial substep
1310  node_data = dynamics_data + 11*i;
1311  PrecessingNRSur_get_time_deriv_from_index(dydt[i], i, q, node_data, __sur_data);
1312 
1313  // Next substep: evaluate at the half node (by using -1, -2, or -3 for i=0, 1, or 2)
1314  for (j=0; j<11; j++) {
1315  y_tmp[j] = node_data[j] + 0.5*time_steps[i]*dydt[i][j];
1316  }
1317  PrecessingNRSur_get_time_deriv_from_index(k2, -1-i, q, y_tmp, __sur_data);
1318 
1319  // Next substep: also at the half node, but update y_tmp
1320  for (j=0; j<11; j++) {
1321  y_tmp[j] = node_data[j] + 0.5*time_steps[i]*k2[j];
1322  }
1323  PrecessingNRSur_get_time_deriv_from_index(k3, -1-i, q, y_tmp, __sur_data);
1324 
1325  // Final substep: evaluate at the next node
1326  for (j=0; j<11; j++) {
1327  y_tmp[j] = node_data[j] + time_steps[i]*k3[j];
1328  }
1329  PrecessingNRSur_get_time_deriv_from_index(k4, i+1, q, y_tmp, __sur_data);
1330 
1331  // Compute the RK4 expression for the next node, normalize, and store the data
1332  for (j=0; j<11; j++) {
1333  y_tmp[j] = node_data[j] + (time_steps[i]/6.0)*(dydt[i][j] + 2*k2[j] + 2*k3[j] + k4[j]);
1334  }
1335  PrecessingNRSur_normalize_y(normA, normB, y_tmp);
1336  node_data = node_data + 11;
1337  for (j=0; j<11; j++) {
1338  node_data[j] = y_tmp[j];
1339  }
1340  }
1341 
1342  // Finally, we need to compute dydt3;
1343  node_data = dynamics_data + 33;
1344  PrecessingNRSur_get_time_deriv_from_index(dydt3, 3, q, node_data, __sur_data);
1345 }
1346 
1347 /**
1348  * Initializes the AB4 ODE system from a single arbitrary dynamics node by taking 3 RK4 steps.
1349  * Compared to PrecessingNRSur_initialize_RK4_with_half_nodes, this may be slightly less accurate
1350  * and slightly more time consuming as we must interpolate many dynamics node fit evaluations,
1351  * but this is more flexible as we can initialize from any node instead of just node 0.
1352  */
1354  REAL8 *dynamics_data, /**< A pointer to the start of the ODE output.
1355  Entries with i0*11 <= i < (i0+1)*11 should have already be computed. */
1356  REAL8 *time_steps, /**< Output: first three time steps. Should have size 3. */
1357  REAL8 *dydt0, /**< Output: dydt at node i0 + 0. Should have size 11. */
1358  REAL8 *dydt1, /**< Output: dydt at node i0 + 1. Should have size 11. */
1359  REAL8 *dydt2, /**< Output: dydt at node i0 + 2. Should have size 11. */
1360  REAL8 *dydt3, /**< Output: dydt at node i0 + 3. Should have size 11. */
1361  REAL8 normA, /**< |chiA| */
1362  REAL8 normB, /**< |chiB| */
1363  REAL8 q, /**< mass ratio */
1364  int i0, /**< the node that is already initialized */
1365  PrecessingNRSurData *__sur_data /**< Loaded surrogate data */
1366 ) {
1367 
1368  gsl_vector *t_ds = __sur_data->t_ds;
1369  REAL8 t1, t2;
1370  int i, j;
1371  REAL8 k2[11], k3[11], k4[11]; // dydt for all but the first RK4 substep (which will be one of dydt0, dydt1, dydt2, dydt3)
1372  REAL8 *node_data;
1373  REAL8 y_tmp[11];
1374  REAL8 *dydt[4];
1375  dydt[0] = dydt0;
1376  dydt[1] = dydt1;
1377  dydt[2] = dydt2;
1378  dydt[3] = dydt3;
1379  int i_start;
1380 
1381  // If possible, do 3 steps backwards with RK4. Otherwise, do 3 steps forwards.
1382  if (i0 > 2) {
1383  i_start = i0 - 3;
1384 
1385  // Setup time steps
1386  t1 = gsl_vector_get(t_ds, i_start);
1387  for (i=0; i<3; i++) {
1388  t2 = gsl_vector_get(t_ds, i_start+i+1);
1389  time_steps[i] = t2 - t1;
1390  t1 = t2;
1391  }
1392 
1393  // Three steps of RK4 BACKWARDS, so include a minus sign beside all dydts.
1394  for (i=0; i<3; i++ ) {
1395 
1396  // Initial substep
1397  node_data = dynamics_data + 11*(i0-i);
1398  PrecessingNRSur_get_time_deriv_from_index(dydt[3-i], i0-i, q, node_data, __sur_data);
1399 
1400  // Next substep: evaluate halfway to the next timestep
1401  t1 = 0.5*(gsl_vector_get(t_ds, i0-i) + gsl_vector_get(t_ds, i0-i-1));
1402  for (j=0; j<11; j++) {
1403  y_tmp[j] = node_data[j] - 0.5*time_steps[2-i]*dydt[3-i][j];
1404  }
1405  PrecessingNRSur_get_time_deriv(k2, t1, q, y_tmp, __sur_data);
1406 
1407  // Next substep: also halfway, but update y_tmp
1408  for (j=0; j<11; j++) {
1409  y_tmp[j] = node_data[j] - 0.5*time_steps[2-i]*k2[j];
1410  }
1411  PrecessingNRSur_get_time_deriv(k3, t1, q, y_tmp, __sur_data);
1412 
1413  // Final substep: evaluate at the next node
1414  for (j=0; j<11; j++) {
1415  y_tmp[j] = node_data[j] - time_steps[2-i]*k3[j];
1416  }
1417  PrecessingNRSur_get_time_deriv_from_index(k4, i0-i-1, q, y_tmp, __sur_data);
1418 
1419  // Compute the RK4 expression for the next node, normalize, and store the data
1420  for (j=0; j<11; j++) {
1421  y_tmp[j] = node_data[j] - (time_steps[2-i]/6.0)*(dydt[3-i][j] + 2*k2[j] + 2*k3[j] + k4[j]);
1422  }
1423  PrecessingNRSur_normalize_y(normA, normB, y_tmp);
1424  node_data = node_data - 11;
1425  for (j=0; j<11; j++) {
1426  node_data[j] = y_tmp[j];
1427  }
1428  }
1429  // Finally, we need to compute dydt0;
1430  node_data = dynamics_data + (i0 - 3)*11;
1431  PrecessingNRSur_get_time_deriv_from_index(dydt0, i0-3, q, node_data, __sur_data);
1432  } else {
1433  i_start = i0;
1434 
1435  // Setup time steps
1436  t1 = gsl_vector_get(t_ds, i0);
1437  for (i=0; i<3; i++) {
1438  t2 = gsl_vector_get(t_ds, i0+i+1);
1439  time_steps[i] = t2 - t1;
1440  t1 = t2;
1441  }
1442 
1443  // Three steps of RK4.
1444  for (i=0; i<3; i++ ) {
1445 
1446  // Initial substep
1447  node_data = dynamics_data + 11*(i+i0);
1448  PrecessingNRSur_get_time_deriv_from_index(dydt[i], i0+i, q, node_data, __sur_data);
1449 
1450  // Next substep: evaluate halfway to the next timestep
1451  t1 = 0.5*(gsl_vector_get(t_ds, i+i0) + gsl_vector_get(t_ds, i+i0+1));
1452  for (j=0; j<11; j++) {
1453  y_tmp[j] = node_data[j] + 0.5*time_steps[i]*dydt[i][j];
1454  }
1455  PrecessingNRSur_get_time_deriv(k2, t1, q, y_tmp, __sur_data);
1456 
1457  // Next substep: also halfway, but update y_tmp
1458  for (j=0; j<11; j++) {
1459  y_tmp[j] = node_data[j] + 0.5*time_steps[i]*k2[j];
1460  }
1461  PrecessingNRSur_get_time_deriv(k3, t1, q, y_tmp, __sur_data);
1462 
1463  // Final substep: evaluate at the next node
1464  for (j=0; j<11; j++) {
1465  y_tmp[j] = node_data[j] + time_steps[i]*k3[j];
1466  }
1467  PrecessingNRSur_get_time_deriv_from_index(k4, i0+i+1, q, y_tmp, __sur_data);
1468 
1469  // Compute the RK4 expression for the next node, normalize, and store the data
1470  for (j=0; j<11; j++) {
1471  y_tmp[j] = node_data[j] + (time_steps[i]/6.0)*(dydt[i][j] + 2*k2[j] + 2*k3[j] + k4[j]);
1472  }
1473  PrecessingNRSur_normalize_y(normA, normB, y_tmp);
1474  node_data = node_data + 11;
1475  for (j=0; j<11; j++) {
1476  node_data[j] = y_tmp[j];
1477  }
1478  }
1479  // Finally, we need to compute dydt3;
1480  node_data = dynamics_data + (i0 + 3)*11;
1481  PrecessingNRSur_get_time_deriv_from_index(dydt3, i0+3, q, node_data, __sur_data);
1482  }
1483 
1484  return i_start;
1485 }
1486 
1487 /**
1488  * Integrates the AB4 ODE system in time forwards, and backwards if needed.
1489  * The system should have already been initialized at 4 adjacent dynamics nodes.
1490  * Output is a flattened array where entries i*11 <= j < (i+1)*11 are
1491  * [q0, qx, qy, qz, varphi, chiAx, chiAy, chiAz, chiBx, chiBy, chiBz] at dynamics node i,
1492  * where (q0, qx, qy, qz) is the coprecessing frame quaternion, varphi is the
1493  * orbital phase, and the spin components are in the coprecessing frame.
1494  */
1496  REAL8 *dynamics_data, /**< ODE output */
1497  REAL8 *time_steps, /**< The first three time steps beginning at i_start. */
1498  REAL8 *dydt0, /**< dydt at node i_start */
1499  REAL8 *dydt1, /**< dydt at node i_start + 1 */
1500  REAL8 *dydt2, /**< dydt at node i_start + 2 */
1501  REAL8 *dydt3, /**< dydt at node i_start + 3 */
1502  REAL8 normA, /**< |chiA| */
1503  REAL8 normB, /**< |chiB| */
1504  REAL8 q, /**< mass ratio */
1505  int i_start, /**< nodes i_start through i_start+3 are already initialized */
1506  PrecessingNRSurData *__sur_data /**< Loaded surrogate data */
1507 ) {
1508  // At each step, we will use the 4 most recent nodes where we have the solution, and integrate up to the first new node.
1509  REAL8 dt1, dt2, dt3, dt4; // Time steps between the 4 known nodes, as well as up to the first new node (dt4)
1510  REAL8 k1[11], k2[11], k3[11], k4[11]; // dydt(y(t); t) evaluated at the 4 known nodes
1511  REAL8 ynext[11]; // Temporary solution
1512  REAL8 *node_data; // Pointer for writing output
1513  int i, j;
1514 
1515  gsl_vector *t_ds = __sur_data->t_ds;
1516 
1517  // Initialize to integrate forwards
1518  dt1 = time_steps[0];
1519  dt2 = time_steps[1];
1520  dt3 = time_steps[2];
1521  for (j=0; j<11; j++) {
1522  k1[j] = dydt0[j];
1523  k2[j] = dydt1[j];
1524  k3[j] = dydt2[j];
1525  k4[j] = dydt3[j];
1526  }
1527 
1528  // Integrate forwards
1529  node_data = dynamics_data + 11*(i_start + 3); // Point to latest known solution
1530  for (i=i_start+4; i< (int)t_ds->size; i++) { // i indexes the output
1531 
1532  // Compute dy = dydt*dt, write it in ynext
1533  dt4 = gsl_vector_get(t_ds, i) - gsl_vector_get(t_ds, i-1);
1534  NRSur_ab4_dy(ynext, k1, k2, k3, k4, dt1, dt2, dt3, dt4, 11);
1535 
1536  // Add the latest known node, to obtain y at the next node
1537  for (j=0; j<11; j++) ynext[j] += node_data[j];
1538 
1539  // Normalize and write output
1540  PrecessingNRSur_normalize_y(normA, normB, ynext);
1541  node_data += 11;
1542  for (j=0; j<11; j++) node_data[j] = ynext[j];
1543 
1544  // Setup for next step
1545  if (i < (int) t_ds->size - 1) {
1546  dt1 = dt2;
1547  dt2 = dt3;
1548  dt3 = dt4;
1549  for (j=0; j<11; j++) {
1550  k1[j] = k2[j];
1551  k2[j] = k3[j];
1552  k3[j] = k4[j];
1553  }
1554  PrecessingNRSur_get_time_deriv_from_index(k4, i, q, node_data, __sur_data);
1555  }
1556  }
1557 
1558  // Initialize to integrate backwards
1559  dt1 = -1 * time_steps[2];
1560  dt2 = -1 * time_steps[1];
1561  dt3 = -1 * time_steps[0];
1562  for (j=0; j<11; j++) {
1563  k1[j] = dydt3[j];
1564  k2[j] = dydt2[j];
1565  k3[j] = dydt1[j];
1566  k4[j] = dydt0[j];
1567  }
1568 
1569  // Integrate backwards
1570  node_data = dynamics_data + 11*(i_start); // Point to earliest known solution
1571  for (i=i_start-1; i>=0; i--) { // i indexes the output
1572 
1573  // Compute dy = dydt*dt, write it in ynext
1574  dt4 = gsl_vector_get(t_ds, i) - gsl_vector_get(t_ds, i+1);
1575  NRSur_ab4_dy(ynext, k1, k2, k3, k4, dt1, dt2, dt3, dt4, 11);
1576 
1577  // Add the earliest known node, to obtain y at the previous
1578  for (j=0; j<11; j++) ynext[j] += node_data[j];
1579 
1580  // Normalize and write output
1581  PrecessingNRSur_normalize_y(normA, normB, ynext);
1582  node_data -= 11;
1583  for (j=0; j<11; j++) node_data[j] = ynext[j];
1584 
1585  // Setup for next step
1586  if (i > 0) {
1587  dt1 = dt2;
1588  dt2 = dt3;
1589  dt3 = dt4;
1590  for (j=0; j<11; j++) {
1591  k1[j] = k2[j];
1592  k2[j] = k3[j];
1593  k3[j] = k4[j];
1594  }
1595  PrecessingNRSur_get_time_deriv_from_index(k4, i, q, node_data, __sur_data);
1596  }
1597  }
1598 }
1599 
1600 /**
1601  * Evaluates a single NRSur coorbital waveoform data piece.
1602  * The dynamics ODE must have already been solved, since this requires the
1603  * spins evaluated at all of the empirical nodes for this waveform data piece.
1604  */
1606  gsl_vector *result, /**< Output: Should have already been assigned space */
1607  REAL8 q, /**< Mass ratio */
1608  gsl_vector **chiA, /**< 3 gsl_vector *s, one for each (coorbital) component */
1609  gsl_vector **chiB, /**< similar to chiA */
1610  WaveformDataPiece *data, /**< The data piece to evaluate */
1611  PrecessingNRSurData *__sur_data /**< Loaded surrogate data */
1612 ) {
1613 
1614  gsl_vector *nodes = gsl_vector_alloc(data->n_nodes);
1615  REAL8 x[7];
1616  int i, j, node_index;
1617 
1618  // Evaluate the fits at the empirical nodes, using the spins at the empirical node times
1619  x[0] = q;
1620  for (i=0; i<data->n_nodes; i++) {
1621  node_index = gsl_vector_long_get(data->empirical_node_indices, i);
1622  for (j=0; j<3; j++) {
1623  x[1+j] = gsl_vector_get(chiA[j], node_index);
1624  x[4+j] = gsl_vector_get(chiB[j], node_index);
1625  }
1626  gsl_vector_set(nodes, i, PrecessingNRSur_eval_fit(data->fit_data[i], x, __sur_data));
1627  }
1628 
1629  // Evaluate the empirical interpolant
1630  gsl_blas_dgemv(CblasTrans, 1.0, data->empirical_interpolant_basis, nodes, 0.0, result);
1631 
1632  gsl_vector_free(nodes);
1633 }
1634 
1635 /************************ Main Waveform Generation Routines ***********/
1636 
1637 /**
1638  * This is the main NRSur dynamics surrogate integration routine.
1639  * Given omega_ref and the system parameters at omega_ref, we find the corresponding t_ref
1640  * and first (if needed) take a tiny time step and evolve the system to the nearest
1641  * dynamics node.
1642  * We then initialize the AB4 ODE system at 4 consecutive dynamics nodes by taking 3 RK4 steps.
1643  * Finally, we integrate forwards (and backwards if needed) to obtain the solution at all
1644  * dynamics nodes.
1645  * Output is a flattened array where entries i*11 <= j < (i+1)*11 are
1646  * [q0, qx, qy, qz, varphi, chiAx, chiAy, chiAz, chiBx, chiBy, chiBz] at dynamics node i,
1647  * where (q0, qx, qy, qz) is the coprecessing frame quaternion, varphi is the
1648  * orbital phase, and the spin components are in the coprecessing frame.
1649  */
1651  REAL8 *dynamics_data, /**< Output: length (n_dynamics_nodes * 11) */
1652  REAL8 q, /**< Mass ratio mA / mB */
1653  REAL8 *chiA0, /**< chiA at the reference point */
1654  REAL8 *chiB0, /**< chiB at the reference point */
1655  REAL8 omega_ref, /**< orbital angular frequency at reference point */
1656  REAL8 init_orbphase, /**< orbital phase at the reference point */
1657  REAL8 *init_quat, /**< coprecessing quaterion at the reference point */
1658  LALDict* LALparams, /**< Dict with extra parameters */
1659  UINT4 PrecessingNRSurVersion /**< 0 for NRSur7dq2, 1 for NRSur7dq4 */
1660 ) {
1661 
1662  // By default we do not allow unlimited_extrapolation
1663  UINT4 unlim_extrap = 0;
1664  if (LALparams != NULL &&
1665  XLALDictContains(LALparams, "unlimited_extrapolation")) {
1666  // Unless the user asks for it
1667  unlim_extrap
1668  = XLALDictLookupUINT4Value(LALparams, "unlimited_extrapolation");
1669  }
1670 
1671  REAL8 normA = sqrt(chiA0[0]*chiA0[0] + chiA0[1]*chiA0[1] + chiA0[2]*chiA0[2]);
1672  REAL8 normB = sqrt(chiB0[0]*chiB0[0] + chiB0[1]*chiB0[1] + chiB0[2]*chiB0[2]);
1673 
1674  // Sanity checks and warnings
1675  REAL8 Q_MAX = 1;
1676  REAL8 CHI_MAX = 0;
1677  REAL8 Q_MAX_WARN = 1;
1678  REAL8 CHI_MAX_WARN = 0;
1679  PrecessingNRSurData *__sur_data;
1680  if (PrecessingNRSurVersion == 0) {
1682  CHI_MAX = NRSUR7DQ2_CHI_MAX;
1683  Q_MAX_WARN = NRSUR7DQ2_Q_MAX_WARN;
1684  CHI_MAX_WARN = NRSUR7DQ2_CHI_MAX_WARN;
1685  __sur_data = &__lalsim_NRSur7dq2_data;
1686  } else if (PrecessingNRSurVersion == 1) {
1688  CHI_MAX = NRSUR7DQ4_CHI_MAX;
1689  Q_MAX_WARN = NRSUR7DQ4_Q_MAX_WARN;
1690  CHI_MAX_WARN = NRSUR7DQ4_CHI_MAX_WARN;
1691  __sur_data = &__lalsim_NRSur7dq4_data;
1692  } else {
1694  "Only 0 or 1 are currently allowed for PrecessingNRSurVersion\n");
1695  }
1696 
1697  if (q < 0.999) {
1699  "Invalid mass ratio q = %0.4f < 1\n", q);
1700  }
1701  if ((q > Q_MAX) && (unlim_extrap == 0)) {
1703  "Too much extrapolation. Mass ratio q = %0.4f > %0.4f, the "
1704  "maximum allowed value.\n", q, Q_MAX);
1705  }
1706  if (q > Q_MAX_WARN) {
1708  "Extrapolating to mass ratio q = %0.4f > %0.4f, the maximum "
1709  "mass ratio used to train the surrogate.\n", q, Q_MAX_WARN);
1710  }
1711  if ((normA > CHI_MAX) && (unlim_extrap == 0)) {
1713  "Too much extrapolation. Spin magnitude |chiA| = %0.4f > %.04f "
1714  "the maximum allowed .\n", normA, CHI_MAX);
1715  }
1716  if ((normB > CHI_MAX) && (unlim_extrap == 0)) {
1718  "Too much extrapolation. Spin magnitude |chiB| = %0.4f > %.04f "
1719  "the maximum allowed value.\n", normB, CHI_MAX);
1720  }
1721  if (normA > CHI_MAX_WARN) {
1723  "Extrapolating to spin magnitude |chiA| = %0.4f > %0.4f, the "
1724  "maximum spin magnitude used to train the surrogate.\n",
1725  normA, CHI_MAX_WARN);
1726  }
1727  if (normB > CHI_MAX_WARN) {
1729  "Extrapolating to spin magnitude |chiB| = %0.4f > %0.4f, the "
1730  "maximum spin magnitude used to train the surrogate.\n",
1731  normB, CHI_MAX_WARN);
1732  }
1733 
1734  REAL8 t_ref = PrecessingNRSur_get_t_ref(omega_ref, q, chiA0, chiB0,
1735  init_quat, init_orbphase, __sur_data);
1736  if (XLAL_IS_REAL8_FAIL_NAN(t_ref)) {
1737  XLAL_ERROR_REAL8(XLAL_FAILURE, "Failed to determine t_ref");
1738  }
1739 
1740  // Initialize the ODE system by stepping to a dynamics node indexed by i0
1741  int i0 = PrecessingNRSur_initialize_at_dynamics_node(dynamics_data, t_ref, q, chiA0, chiB0, init_orbphase, init_quat, normA, normB, __sur_data);
1742 
1743  // To initialize the RK4 ODE integration scheme, we need y to be evaluated at 4 consecutive nodes. Currently we have 1.
1744  // The original method assumed we always use t_ref=-4500M and took 3 steps using RK4, making use of the time nodes 1/2, 3/2, and 5/2.
1745  // If i0=0, we can use that scheme and obtain y at nodes 0, 1, 2, and 3. Otherwise, obtain 4 consecutive nodes starting at i_start,
1746  // where i_start = max(0, i0 - 3).
1747  int i_start;
1748  REAL8 time_steps[3];
1749  REAL8 dydt0[11], dydt1[11], dydt2[11], dydt3[11];
1750 
1751  if (i0 == 0) {
1752  PrecessingNRSur_initialize_RK4_with_half_nodes(dynamics_data, time_steps, dydt0, dydt1, dydt2, dydt3, normA, normB, q, __sur_data);
1753  i_start = 0;
1754  } else {
1755  i_start = PrecessingNRSur_initialize_RK4(dynamics_data, time_steps, dydt0, dydt1, dydt2, dydt3, normA, normB, q, i0, __sur_data);
1756  }
1757 
1758  // Now that we have 4 consecutive evaluated nodes beginning at i_start, we can integrate from i_start backwards to 0
1759  // with AB4 as well as from i_start+3 to the end
1760  PrecessingNRSur_integrate_AB4(dynamics_data, time_steps, dydt0, dydt1, dydt2, dydt3, normA, normB, q, i_start, __sur_data);
1761  return XLAL_SUCCESS;
1762 }
1763 
1764 
1765 /**
1766  * Wrapper for Loading NRSur7dq2 and NRSur7dq4 data.
1767  *
1768  * Calls NRSur7dq2_Init_LALDATA() or NRSur7dq4_Init_LALDATA() if surrogate data
1769  * is not already loaded, and returns loaded surrogate data. If surrogate data
1770  * is already loaded, just returns the loaded data.
1771  */
1773  Approximant approximant /**< approximant (NRSur7dq2 or NRSur7dq4) */
1774 ) {
1775 
1776  PrecessingNRSurData *__sur_data;
1777 
1778  switch (approximant) {
1779  case NRSur7dq2:
1780  #ifdef LAL_PTHREAD_LOCK
1781  (void) pthread_once(&NRSur7dq2_is_initialized,
1783  #else
1785  #endif
1786  __sur_data = &__lalsim_NRSur7dq2_data;
1787  // 0 is for NRSur7dq2
1788  __sur_data->PrecessingNRSurVersion = 0;
1789  break;
1790 
1791  case NRSur7dq4:
1792  #ifdef LAL_PTHREAD_LOCK
1793  (void) pthread_once(&NRSur7dq4_is_initialized,
1795  #else
1797  #endif
1798  __sur_data = &__lalsim_NRSur7dq4_data;
1799  // 1 is for NRSur7dq2
1800  __sur_data->PrecessingNRSurVersion = 1;
1801  break;
1802 
1803  default:
1804  XLAL_ERROR_NULL(XLAL_EINVAL, "Invalid approximant, only NRSur7dq2"
1805  " and NRSur7dq4 are allowed.\n");
1806  }
1807 
1808  return __sur_data;
1809 }
1810 
1811 
1812 
1813 
1814 /* This is the core function of the NRSur7dq2 and NRSur7dq4 models.
1815  * It evaluates the model, and returns waveform modes in the inertial frame, sampled on t_coorb.
1816  * When using a custom ModeArray the user must explicitly supply all modes, -ve and +ve m-modes.
1817  * Note that the co-orbital frame modes are specified in ModeArray, not the inertial frame modes.
1818  */
1820  MultiModalWaveform **h, /**< Output. Dimensionless waveform modes sampled on t_coorb */
1821  REAL8 q, /**< Mass ratio mA / mB */
1822  REAL8 *chiA0, /**< chiA at the reference point */
1823  REAL8 *chiB0, /**< chiB at the reference point */
1824  REAL8 omega_ref, /**< orbital angular frequency at the reference point */
1825  REAL8 init_orbphase, /**< orbital phase at the reference point */
1826  REAL8 *init_quat, /**< coprecessing quaterion at the reference point */
1827  LALValue* ModeArray, /**< Container for the ell and m co-orbital modes to generate */
1828  LALDict* LALparams, /**< Dict with extra parameters */
1829  Approximant approximant /**< approximant (NRSur7dq2 or NRSur7dq4) */
1830 ) {
1831 
1832  // Load surrogate data if needed. If not, just access the loaded data.
1834  if (!__sur_data->setup) {
1835  XLAL_ERROR_NULL(XLAL_EFAILED, "Error loading surrogate data.\n");
1836  }
1837 
1838  // Make sure we didn't request any unavailable modes
1839  int ell, m;
1840  for (ell=NRSUR_LMAX+1; ell <= LAL_SIM_L_MAX_MODE_ARRAY; ell++) {
1841  for (m=-ell; m<=ell; m++) {
1842  if (XLALSimInspiralModeArrayIsModeActive(ModeArray, ell, m) == 1) {
1843  XLAL_ERROR_NULL(XLAL_EINVAL, "A requested mode has ell larger than available");
1844  }
1845  }
1846  }
1847 
1848  // time arrays
1849  gsl_vector *t_ds = __sur_data->t_ds;
1850  gsl_vector *t_coorb = __sur_data->t_coorb;
1851  int n_ds = t_ds->size;
1852  int n_coorb = t_coorb->size;
1853 
1854  // Get dynamics
1855  REAL8 *dynamics_data = XLALCalloc(n_ds * 11, sizeof(REAL8));
1856 
1857  int ret = PrecessingNRSur_IntegrateDynamics(dynamics_data, q, chiA0, chiB0,
1858  omega_ref, init_orbphase, init_quat, LALparams,
1859  __sur_data->PrecessingNRSurVersion);
1860  if(ret != XLAL_SUCCESS) {
1861  XLALFree(dynamics_data);
1862  XLAL_ERROR_NULL(XLAL_FAILURE, "Failed to integrate dynamics");
1863  }
1864 
1865  // Put output into appropriate vectors
1866  int i, j;
1867  gsl_vector *dynamics[11];
1868  for (j=0; j<11; j++) {
1869  dynamics[j] = gsl_vector_alloc(n_ds);
1870  }
1871  for (i=0; i<n_ds; i++) {
1872  for (j=0; j<11; j++) {
1873  gsl_vector_set(dynamics[j], i, dynamics_data[11*i + j]);
1874  }
1875  }
1876 
1877  XLALFree(dynamics_data);
1878 
1879  // Interpolate onto the coorbital time grid
1880  // NOTE: The spins are currently in the coprecessing frame, they are called
1881  // chiA_coorb because they will be overwritten below with the coorbital frame
1882  // spins.
1883  gsl_vector *quat_coorb[4], *chiA_coorb[3], *chiB_coorb[3];
1884  quat_coorb[0] = spline_array_interp(t_coorb, t_ds, dynamics[0]);
1885  quat_coorb[1] = spline_array_interp(t_coorb, t_ds, dynamics[1]);
1886  quat_coorb[2] = spline_array_interp(t_coorb, t_ds, dynamics[2]);
1887  quat_coorb[3] = spline_array_interp(t_coorb, t_ds, dynamics[3]);
1888  gsl_vector *phi_coorb = spline_array_interp(t_coorb, t_ds, dynamics[4]);
1889  chiA_coorb[0] = spline_array_interp(t_coorb, t_ds, dynamics[5]);
1890  chiA_coorb[1] = spline_array_interp(t_coorb, t_ds, dynamics[6]);
1891  chiA_coorb[2] = spline_array_interp(t_coorb, t_ds, dynamics[7]);
1892  chiB_coorb[0] = spline_array_interp(t_coorb, t_ds, dynamics[8]);
1893  chiB_coorb[1] = spline_array_interp(t_coorb, t_ds, dynamics[9]);
1894  chiB_coorb[2] = spline_array_interp(t_coorb, t_ds, dynamics[10]);
1895 
1896  for (j=0; j<11; j++) {
1897  gsl_vector_free(dynamics[j]);
1898  }
1899 
1900  // Normalize spins after interpolation
1901  REAL8 normA = sqrt(chiA0[0]*chiA0[0] + chiA0[1]*chiA0[1] + chiA0[2]*chiA0[2]);
1902  REAL8 normB = sqrt(chiB0[0]*chiB0[0] + chiB0[1]*chiB0[1] + chiB0[2]*chiB0[2]);
1903  PrecessingNRSur_normalize_results(normA, normB, quat_coorb, chiA_coorb, chiB_coorb);
1904 
1905  // Transform spins from coprecessing frame to coorbital frame for use in coorbital waveform surrogate
1906  PrecessingNRSur_rotate_spins(chiA_coorb, chiB_coorb, phi_coorb);
1907 
1908  // Evaluate the coorbital waveform surrogate
1909  MultiModalWaveform *h_coorb = NULL;
1910  MultiModalWaveform_Init(&h_coorb, NRSUR_LMAX, n_coorb);
1911  gsl_vector *data_piece_eval = gsl_vector_alloc(n_coorb);
1912  WaveformDataPiece *data_piece_data;
1913  int i0; // for indexing the (ell, m=0) mode, such that the (ell, m) mode is index (i0 + m).
1914  WaveformFixedEllModeData *ell_data;
1915  for (ell=2; ell<=NRSUR_LMAX; ell++) {
1916  ell_data = __sur_data->coorbital_mode_data[ell - 2];
1917  i0 = ell*(ell+1) - 4;
1918 
1919  // m=0
1920  if (XLALSimInspiralModeArrayIsModeActive(ModeArray, ell, 0) != 1) {
1921  }
1922  else {
1923  data_piece_data = ell_data->m0_real_data;
1924  PrecessingNRSur_eval_data_piece(data_piece_eval, q, chiA_coorb, chiB_coorb, data_piece_data, __sur_data);
1925  gsl_vector_add(h_coorb->modes_real_part[i0], data_piece_eval);
1926 
1927  data_piece_data = ell_data->m0_imag_data;
1928  PrecessingNRSur_eval_data_piece(data_piece_eval, q, chiA_coorb, chiB_coorb, data_piece_data, __sur_data);
1929  gsl_vector_add(h_coorb->modes_imag_part[i0], data_piece_eval);
1930  }
1931 
1932  // Other modes
1933  for (m=1; m<=ell; m++) {
1934  if ((XLALSimInspiralModeArrayIsModeActive(ModeArray, ell, m) != 1) &&
1935  (XLALSimInspiralModeArrayIsModeActive(ModeArray, ell, -m) != 1)) {
1936  continue;
1937  }
1938 
1939  // h^{ell, m} = X_plus + X_minus
1940  // h^{ell, -m} = (X_plus - X_minus)* <- complex conjugate
1941 
1942  // Re[X_plus] gets added to both Re[h^{ell, m}] and Re[h^{ell, -m}]
1943  data_piece_data = ell_data->X_real_plus_data[m-1];
1944  PrecessingNRSur_eval_data_piece(data_piece_eval, q, chiA_coorb, chiB_coorb, data_piece_data, __sur_data);
1945  gsl_vector_add(h_coorb->modes_real_part[i0+m], data_piece_eval);
1946  gsl_vector_add(h_coorb->modes_real_part[i0-m], data_piece_eval);
1947 
1948  // Re[X_minus] gets added to Re[h^{ell, m}] and subtracted from Re[h^{ell, -m}]
1949  data_piece_data = ell_data->X_real_minus_data[m-1];
1950  PrecessingNRSur_eval_data_piece(data_piece_eval, q, chiA_coorb, chiB_coorb, data_piece_data, __sur_data);
1951  gsl_vector_add(h_coorb->modes_real_part[i0+m], data_piece_eval);
1952  gsl_vector_sub(h_coorb->modes_real_part[i0-m], data_piece_eval);
1953 
1954  // Im[X_plus] gets added to Re[h^{ell, m}] and subtracted from Re[h^{ell, -m}]
1955  data_piece_data = ell_data->X_imag_plus_data[m-1];
1956  PrecessingNRSur_eval_data_piece(data_piece_eval, q, chiA_coorb, chiB_coorb, data_piece_data, __sur_data);
1957  gsl_vector_add(h_coorb->modes_imag_part[i0+m], data_piece_eval);
1958  gsl_vector_sub(h_coorb->modes_imag_part[i0-m], data_piece_eval);
1959 
1960  // Im[X_minus] gets added to both Re[h^{ell, m}] and Re[h^{ell, -m}]
1961  data_piece_data = ell_data->X_imag_minus_data[m-1];
1962  PrecessingNRSur_eval_data_piece(data_piece_eval, q, chiA_coorb, chiB_coorb, data_piece_data, __sur_data);
1963  gsl_vector_add(h_coorb->modes_imag_part[i0+m], data_piece_eval);
1964  gsl_vector_add(h_coorb->modes_imag_part[i0-m], data_piece_eval);
1965  }
1966  }
1967 
1968  // Rotate to the inertial frame, write results in h
1969  MultiModalWaveform_Init(h, NRSUR_LMAX, n_coorb);
1970  TransformModesCoorbitalToInertial(*h, h_coorb, quat_coorb, phi_coorb);
1971 
1972  // Cleanup
1973  MultiModalWaveform_Destroy(h_coorb);
1974  for (i=0; i<3; i++) {
1975  gsl_vector_free(chiA_coorb[i]);
1976  gsl_vector_free(chiB_coorb[i]);
1977  gsl_vector_free(quat_coorb[i]);
1978  }
1979  gsl_vector_free(quat_coorb[3]);
1980  gsl_vector_free(phi_coorb);
1981  gsl_vector_free(data_piece_eval);
1982 
1983  return __sur_data;
1984 }
1985 
1986 /**
1987  * Computes the starting frequency of the NRSur7dq2 or NRSur7dq4 waveform
1988  * approximant when evaluated using fRef=0 (which uses the entire surrogate
1989  * waveform starting 4500M for NRSur7dq2 and 4300M for NRSur7dq4, before the
1990  * peak amplitude).
1991  */
1993  REAL8 m1, /**< mass of companion 1 (kg) */
1994  REAL8 m2, /**< mass of companion 2 (kg) */
1995  REAL8 s1x, /**< initial value of S1x */
1996  REAL8 s1y, /**< initial value of S1y */
1997  REAL8 s1z, /**< initial value of S1z */
1998  REAL8 s2x, /**< initial value of S2x */
1999  REAL8 s2y, /**< initial value of S2y */
2000  REAL8 s2z, /**< initial value of S2z */
2001  PrecessingNRSurData *__sur_data /**< Loaded surrogate data */
2002 
2003 ) {
2004  REAL8 Mtot = (m1 + m2) / LAL_MSUN_SI;
2005  REAL8 Mtot_sec = Mtot * LAL_MTSUN_SI; /* Total mass in seconds */
2006  REAL8 q = m1 / m2;
2007  REAL8 y0[11];
2008  y0[0] = 1.0; // Scalar component of reference quaternion
2009  int i;
2010  for (i=1; i<4; i++) y0[i] = 0.0; // Vector components of quaternion
2011  y0[4] = 0.0; // Reference orbital phase - doesn't affect frequency
2012  y0[5] = s1x;
2013  y0[6] = s1y;
2014  y0[7] = s1z;
2015  y0[8] = s2x;
2016  y0[9] = s2y;
2017  y0[10] = s2z;
2018  REAL8 omega_dimless_start = PrecessingNRSur_get_omega(0, q, y0, __sur_data);
2019  // GW freq is roughly twice the orbital freq
2020  REAL8 f_start_Hz = omega_dimless_start / (LAL_PI * Mtot_sec);
2021  return f_start_Hz;
2022 }
2023 
2024 /**
2025  * If m1<m2, swaps the labels for the two BHs such that Bh1 is always heavier.
2026  * Then rotates the in-plane spins by pi. These two together are the same
2027  * as a rigid rotation of the system by pi. This rotation will be undone by
2028  * multiplying the odd-m modes by a minus sign after the modes are generated.
2029  */
2031  REAL8 *m1, /**< Input and Output. mass of companion 1 (kg) */
2032  REAL8 *m2, /**< Input and Output. mass of companion 2 (kg) */
2033  REAL8 *s1x, /**< Input and Output. S1x at reference epoch */
2034  REAL8 *s1y, /**< Input and Output. S1y at reference epoch */
2035  REAL8 *s1z, /**< Input and Output. S1z at reference epoch */
2036  REAL8 *s2x, /**< Input and Output. S2x at reference epoch */
2037  REAL8 *s2y, /**< Input and Output. S2y at reference epoch */
2038  REAL8 *s2z /**< Input and Output. S2z at reference epoch */
2039 ) {
2040 
2041  // BH A is defined to be the one with the larger mass, BH B with the
2042  // smaller mass.
2043  if (*m1 < *m2) {
2044  // Switch the labels for the two BHs
2045  REAL8 tmp = *m1;
2046  *m1 = *m2;
2047  *m2 = tmp;
2048  // For the in-plane spins, also change the signs after swapping. This
2049  // is the same as roting the spins about the z-axis by pi.
2050  tmp = *s1x;
2051  *s1x = -*s2x;
2052  *s2x = -tmp;
2053  tmp = *s1y;
2054  *s1y = -*s2y;
2055  *s2y = -tmp;
2056  // No sign change for z-spins
2057  tmp = *s1z;
2058  *s1z = *s2z;
2059  *s2z = tmp;
2060 
2061  return true;
2062  } else {
2063  return false;
2064  }
2065 }
2066 
2067 
2068 /**
2069  * This function evaluates the NRSur7dq2 or NRSur7dq4 surrogate model and sums
2070  * over all ell <= 4 modes to obtain the + and x polarizations.
2071  * Papers:
2072  * NRSur7dq2: https://arxiv.org/abs/1705.07089
2073  * NRSur7dq4: https://arxiv.org/abs/1905.09300
2074  *
2075  * The system is initialized at a time where the orbital frequency of the
2076  * waveform in the coprecessing frame (Eq.11 of
2077  * https://arxiv.org/abs/1705.07089) is pi * fRef. At this reference point, the
2078  * system is initialized in the coorbital frame, with the z-axis along the
2079  * principal eigenvector of the angular momentum operator acting on the
2080  * waveform, the x-axis along the line of separation from the lighter BH to the
2081  * heavier BH, and the y-axis completing the triad. The polarizations of the
2082  * waveform are returned in the sky of this frame at
2083  * (inclination, pi/2 - phiRef). This agrees with the LAL convention.
2084  *
2085  * Extra options are passed through a LALDict:
2086  * dictParams = lal.CreateDict()
2087  *
2088  * 1. Mode options: When using a custom ModeArray the user must explicitly
2089  * supply all modes, -ve and +ve m-modes. Note that these are modes in the
2090  * co-orbital frame, not the inertial frame.
2091  * Examples:
2092  * # First, create the 'empty' mode array
2093  * ma=lalsim.SimInspiralCreateModeArray()
2094  *
2095  * 1a. If you want all models upto ellMax:
2096  * # add the modes (these are in the coorbital frame)
2097  * for ell in range(2, ellMax+1):
2098  * lalsim.SimInspiralModeArrayActivateAllModesAtL(ma, ell)
2099  * # then insert your ModeArray into the LALDict params with
2100  * lalsim.SimInspiralWaveformParamsInsertModeArray(dictParams, ma)
2101  *
2102  * 1b. If you want only a single mode:
2103  * # add (l,m) and (l,-m) modes (these are in the coorbital frame)
2104  * lalsim.SimInspiralModeArrayActivateMode(ma, l, m)
2105  * lalsim.SimInspiralModeArrayActivateMode(ma, l, -m)
2106  * # then insert your ModeArray into the LALDict params with
2107  * lalsim.SimInspiralWaveformParamsInsertModeArray(params, ma)
2108  *
2109  * 2. Extrapolation: The surrogate models are trained on the following ranges.
2110  * NRSur7dq2: q <= 2.01, |chi_1|, |chi_2| <= 0.81.
2111  * NRSur7dq4: q <= 4.01, |chi_1|, |chi_2| <= 0.81.
2112  * If you want a guarantee of accuracy you should stick to the above ranges.
2113  *
2114  * We allow extrapolation to the following ranges, but with a warning:
2115  * NRSur7dq2: q <= 3.01, |chi_1|, |chi_2| <= 1.
2116  * NRSur7dq4: q <= 6.01, |chi_1|, |chi_2| <= 1.
2117  * We expect the models to be reasonable when extrapolated to these ranges.
2118  *
2119  * Beyond the above ranges, we raise an error. If you want to extrapolate
2120  * beyond these limits you can specify unlimited_extrapolation = 1 in your
2121  * dictParams as follows:
2122  * # USE AT YOUR OWN RISK!!
2123  * lal.DictInsertUINT4Value(dictParams, "unlimited_extrapolation", 1)
2124  */
2126  REAL8TimeSeries **hplus, /**< OUTPUT h_+ vector */
2127  REAL8TimeSeries **hcross, /**< OUTPUT h_x vector */
2128  REAL8 phiRef, /**< azimuthal angle for Ylms */
2129  REAL8 inclination, /**< inclination angle */
2130  REAL8 deltaT, /**< sampling interval (s) */
2131  REAL8 m1, /**< mass of companion 1 (kg) */
2132  REAL8 m2, /**< mass of companion 2 (kg) */
2133  REAL8 distance, /**< distance of source (m) */
2134  REAL8 fMin, /**< start GW frequency (Hz) */
2135  REAL8 fRef, /**< reference GW frequency (Hz) */
2136  REAL8 s1x, /**< initial value of S1x */
2137  REAL8 s1y, /**< initial value of S1y */
2138  REAL8 s1z, /**< initial value of S1z */
2139  REAL8 s2x, /**< initial value of S2x */
2140  REAL8 s2y, /**< initial value of S2y */
2141  REAL8 s2z, /**< initial value of S2z */
2142  LALDict* LALparams, /**< Dict with extra parameters */
2143  Approximant approximant /**< approximant (NRSur7dq2 or NRSur7dq4) */
2144 ) {
2145 
2146  LALValue* ModeArray
2148 
2149  int ell, m;
2150  if ( ModeArray == NULL ) {
2151  // Use all available modes
2152  ModeArray = XLALSimInspiralCreateModeArray();
2153  for (ell=2; ell <= NRSUR_LMAX; ell++) {
2155  }
2156  } // Otherwise, use the specified modes.
2157 
2158  // If m1 < m2, switch the BH labels so that q = m1/m2 >= 1 always.
2159  bool labels_switched = PrecessingNRSur_switch_labels_if_needed(&m1, &m2,
2160  &s1x, &s1y, &s1z, &s2x, &s2y, &s2z);
2161 
2162  // Parameters
2163  REAL8 massA = m1 / LAL_MSUN_SI;
2164  REAL8 massB = m2 / LAL_MSUN_SI;
2165  REAL8 Mtot = massA+massB;
2166  REAL8 Mtot_sec = Mtot * LAL_MTSUN_SI; /* Total mass in seconds */
2167  REAL8 q = massA / massB;
2168  REAL8 chiA0[3], chiB0[3];
2169  chiA0[0] = s1x;
2170  chiA0[1] = s1y;
2171  chiA0[2] = s1z;
2172  chiB0[0] = s2x;
2173  chiB0[1] = s2y;
2174  chiB0[2] = s2z;
2175 
2176  REAL8 omegaMin_dimless;
2177  REAL8 omegaRef_dimless;
2178  int ret = get_dimless_omega(&omegaMin_dimless, &omegaRef_dimless,
2179  fMin, fRef, Mtot_sec);
2180  if(ret != XLAL_SUCCESS) {
2181  if(ModeArray) XLALDestroyValue(ModeArray);
2182  XLAL_ERROR(XLAL_EFUNC, "Failed to process fMin/fRef");
2183  }
2184 
2185  // When generating the modes take the inertial frame to be aligned with
2186  // the coorbital frame at the reference point. This means that the
2187  // quaternion from inertial frame to coprecessing frame is identity, and
2188  // the init_orbphase = 0 between the coprecessing frame and the coorbital
2189  // frame at the reference point. This agrees with the LAL convention, and
2190  // phiRef should only used in the Ylms when computing the polarizations.
2191  REAL8 init_quat[4];
2192  init_quat[0] = 1.0;
2193  init_quat[1] = 0.0;
2194  init_quat[2] = 0.0;
2195  init_quat[3] = 0.0;
2196  REAL8 init_orbphase = 0;
2197 
2198  // Evaluate the model modes
2199  MultiModalWaveform *h_inertial_modes = NULL;
2200  PrecessingNRSurData *__sur_data = PrecessingNRSur_core(&h_inertial_modes,
2201  q, chiA0, chiB0, omegaRef_dimless, init_orbphase, init_quat,
2202  ModeArray, LALparams, approximant);
2203 
2204  if (!h_inertial_modes) {
2205  if(ModeArray) XLALDestroyValue(ModeArray);
2206  return XLAL_FAILURE;
2207  }
2208 
2209  // Setup the time series for h_+ and h_x evaluated on the surrogate time array
2210  size_t length = h_inertial_modes->n_times;
2211  gsl_vector *hplus_model_times = gsl_vector_calloc(length);
2212  gsl_vector *hcross_model_times = gsl_vector_calloc(length);
2213 
2214  // Sum up the modes
2215  COMPLEX16 swsh_coef;// = XLALSpinWeightedSphericalHarmonic(spheretheta, spherephi, -2, swsh_L, swsh_m);
2216  REAL8 c_re, c_im;// = creal(swsh_coef);
2217  REAL8 hmre, hmim;
2218  int i, j;
2219  bool skip_ell_modes;
2220  i=0;
2221  REAL8 prefactor = 1;
2222  for (ell=2; ell <= NRSUR_LMAX; ell++) {
2223  /* If *any* co-orbital frame mode of this ell is non-zero we need to
2224  * sum over all modes of this ell, as the frame rotation will have
2225  * mixed the modes. */
2226  skip_ell_modes = true;
2227  for (m=-ell; m<=ell; m++) {
2228  if (XLALSimInspiralModeArrayIsModeActive(ModeArray, ell, m) == 1) {
2229  skip_ell_modes = false;
2230  }
2231  }
2232  if (skip_ell_modes) {
2233  i += 2*ell + 1;
2234  continue;
2235  }
2236  for (m=-ell; m<=ell; m++) {
2237  /* See Harald Pfeiffer, T18002260-v1 for frame diagram. The
2238  * surrogate frame has its z direction aligned with the orbital
2239  * angular momentum, which agrees with the lowercase source frame z
2240  * in the diagram. The surrogate x direction, however, points along
2241  * the line of ascending nodes (the funny omega with circles on the
2242  * ends). The uppercase Z, which is the direction along which we
2243  * want to evaluate the waveform, is always in the surrogate frame
2244  * (y, z) plane. Z is rotated from z towards the +ive y surrogate
2245  * axis, so we should always evaluate at
2246  * (inclination, pi/2-phiRef). */
2247  swsh_coef = XLALSpinWeightedSphericalHarmonic(inclination,
2248  0.5 * LAL_PI - phiRef, -2, ell, m);
2249  c_re = creal(swsh_coef);
2250  c_im = cimag(swsh_coef);
2251 
2252  // If m1 < m2, the BH labels were switched so that m1 > m2. This
2253  // generates the requested waveform but with a rotation of pi about
2254  // the z-axis. We now undo this by multiplying all odd-m modes by
2255  // -1, which is the same as a pi rotation; the even-m modes don't
2256  // need this fix as they are insensitive to pi rotations about
2257  // z-axis.
2258  if ((m%2 != 0) && (labels_switched)) {
2259  prefactor = -1;
2260  } else {
2261  prefactor = 1;
2262  }
2263  for (j=0; j < (int)length; j++) {
2264  hmre = gsl_vector_get(h_inertial_modes->modes_real_part[i], j);
2265  hmim = gsl_vector_get(h_inertial_modes->modes_imag_part[i], j);
2266  hplus_model_times->data[j] += prefactor
2267  * (c_re * hmre - c_im * hmim);
2268  hcross_model_times->data[j] -= prefactor
2269  * (c_im * hmre + c_re * hmim); // - sign as h = h_+ - i*h_x
2270  }
2271  i += 1;
2272  }
2273  }
2274 
2275  // Scale the amplitude based on the distance
2276  REAL8 a0 = Mtot_sec * LAL_C_SI / distance;
2277  gsl_vector_scale(hplus_model_times, a0);
2278  gsl_vector_scale(hcross_model_times, a0);
2279 
2280  // Setup output times
2281  gsl_vector *model_times = __sur_data->t_coorb;
2282  REAL8 deltaT_dimless = deltaT / Mtot_sec;
2283  REAL8 t0 = gsl_vector_get(model_times, 0);
2284  REAL8 start_freq = PrecessingNRSur_StartFrequency(m1, m2, s1x, s1y, s1z, s2x, s2y, s2z, __sur_data);
2285  if (fMin >= start_freq) {
2286  t0 = PrecessingNRSur_get_t_ref(omegaMin_dimless, q, chiA0, chiB0,
2287  init_quat, init_orbphase, __sur_data);
2288  } else if (fMin > 0) {
2289  // Cleanup and exit
2290  if(ModeArray) XLALDestroyValue(ModeArray);
2291  gsl_vector_free(hplus_model_times);
2292  gsl_vector_free(hcross_model_times);
2293  MultiModalWaveform_Destroy(h_inertial_modes);
2294  XLAL_ERROR_REAL8(XLAL_EDOM, "fMin should be 0 or >= %0.8f for this configuration, got %0.8f", start_freq, fMin);
2295  }
2296  REAL8 tf = gsl_vector_get(model_times, length-1);
2297  int nt = (int) ceil((tf - t0) / deltaT_dimless);
2298  gsl_vector *output_times = gsl_vector_alloc(nt);
2299  for (j=0; j<nt; j++) {
2300  gsl_vector_set(output_times, j, t0 + deltaT_dimless*j);
2301  }
2302 
2303  // Interpolate onto output times
2304  REAL8 t;
2306  XLALGPSAdd( &epoch, Mtot_sec * t0);
2307  *hplus = XLALCreateREAL8TimeSeries("hp: TD waveform", &epoch, 0.0, deltaT, &lalStrainUnit, nt);
2308  *hcross = XLALCreateREAL8TimeSeries("hc: TD waveform", &epoch, 0.0, deltaT, &lalStrainUnit, nt);
2309  gsl_interp_accel *acc = gsl_interp_accel_alloc();
2310  gsl_spline *spl_hplus = gsl_spline_alloc(gsl_interp_cspline, length);
2311  gsl_spline *spl_hcross = gsl_spline_alloc(gsl_interp_cspline, length);
2312  gsl_spline_init(spl_hplus, model_times->data, hplus_model_times->data, length);
2313  gsl_spline_init(spl_hcross, model_times->data, hcross_model_times->data, length);
2314  for (j=0; j<nt; j++) {
2315  t = gsl_vector_get(output_times, j);
2316  (*hplus)->data->data[j] = gsl_spline_eval(spl_hplus, t, acc);
2317  (*hcross)->data->data[j] = gsl_spline_eval(spl_hcross, t, acc);
2318  }
2319 
2320  // Cleanup
2321  gsl_vector_free(hplus_model_times);
2322  gsl_vector_free(hcross_model_times);
2323  gsl_vector_free(output_times);
2324  gsl_interp_accel_free(acc);
2325  gsl_spline_free(spl_hplus);
2326  gsl_spline_free(spl_hcross);
2327  MultiModalWaveform_Destroy(h_inertial_modes);
2328 
2329  if(ModeArray) {
2330  XLALDestroyValue(ModeArray);
2331  }
2332 
2333  return XLAL_SUCCESS;
2334 }
2335 
2336 /**
2337  * This function evaluates the NRSur7dq2 or NRSur7dq4 surrogate model and
2338  * returns the inertial frame modes in the form of a SphHarmTimeSeries.
2339  * Papers:
2340  * NRSur7dq2: https://arxiv.org/abs/1705.07089
2341  * NRSur7dq4: https://arxiv.org/abs/1905.09300
2342  *
2343  * The system is initialized at a time where the orbital frequency of the
2344  * waveform in the coprecessing frame (Eq.11 of
2345  * https://arxiv.org/abs/1705.07089) is pi * fRef. At this reference point, the
2346  * system is initialized in the coorbital frame, with the z-axis along the
2347  * principal eigenvector of the angular momentum operator acting on the
2348  * waveform, the x-axis along the line of separation from the lighter BH to the
2349  * heavier BH, and the y-axis completing the triad. The modes are returned in
2350  * this frame, which agrees with the LAL convention. When combining the modes
2351  * to get the polarizations, the Ylms should be evaluated at
2352  * (inclination, pi/2 - phiRef), following the LAL convention.
2353  *
2354  * Extra options are passed through a LALDict:
2355  * dictParams = lal.CreateDict()
2356  *
2357  * 1. Mode options: When using a custom ModeArray the user must explicitly
2358  * supply all modes, -ve and +ve m-modes. Note that these are modes in the
2359  * co-orbital frame, not the inertial frame.
2360  * Examples:
2361  * # First, create the 'empty' mode array
2362  * ma=lalsim.SimInspiralCreateModeArray()
2363  *
2364  * 1a. If you want all models upto ellMax:
2365  * # add the modes (these are in the coorbital frame)
2366  * for ell in range(2, ellMax+1):
2367  * lalsim.SimInspiralModeArrayActivateAllModesAtL(ma, ell)
2368  * # then insert your ModeArray into the LALDict params with
2369  * lalsim.SimInspiralWaveformParamsInsertModeArray(dictParams, ma)
2370  *
2371  * 1b. If you want only a single mode:
2372  * # add (l,m) and (l,-m) modes (these are in the coorbital frame)
2373  * lalsim.SimInspiralModeArrayActivateMode(ma, l, m)
2374  * lalsim.SimInspiralModeArrayActivateMode(ma, l, -m)
2375  * # then insert your ModeArray into the LALDict params with
2376  * lalsim.SimInspiralWaveformParamsInsertModeArray(params, ma)
2377  *
2378  * 2. Extrapolation: The surrogate models are trained on the following ranges.
2379  * NRSur7dq2: q <= 2.01, |chi_1|, |chi_2| <= 0.81.
2380  * NRSur7dq4: q <= 4.01, |chi_1|, |chi_2| <= 0.81.
2381  * If you want a guarantee of accuracy you should stick to the above ranges.
2382  *
2383  * We allow extrapolation to the following ranges, but with a warning:
2384  * NRSur7dq2: q <= 3.01, |chi_1|, |chi_2| <= 1.
2385  * NRSur7dq4: q <= 6.01, |chi_1|, |chi_2| <= 1.
2386  * We expect the models to be reasonable when extrapolated to these ranges.
2387  *
2388  * Beyond the above ranges, we raise an error. If you want to extrapolate
2389  * beyond these limits you can specify unlimited_extrapolation = 1 in your
2390  * dictParams as follows:
2391  * # USE AT YOUR OWN RISK!!
2392  * lal.DictInsertUINT4Value(dictParams, "unlimited_extrapolation", 1)
2393  */
2395  REAL8 deltaT, /**< sampling interval (s) */
2396  REAL8 m1, /**< mass of companion 1 (kg) */
2397  REAL8 m2, /**< mass of companion 2 (kg) */
2398  REAL8 S1x, /**< x-component of the dimensionless spin of object 1 */
2399  REAL8 S1y, /**< y-component of the dimensionless spin of object 1 */
2400  REAL8 S1z, /**< z-component of the dimensionless spin of object 1 */
2401  REAL8 S2x, /**< x-component of the dimensionless spin of object 2 */
2402  REAL8 S2y, /**< y-component of the dimensionless spin of object 2 */
2403  REAL8 S2z, /**< z-component of the dimensionless spin of object 2 */
2404  REAL8 fMin, /**< start GW frequency (Hz) */
2405  REAL8 fRef, /**< reference GW frequency (Hz) */
2406  REAL8 distance, /**< distance of source (m) */
2407  LALDict* LALparams, /**< Dict with extra parameters */
2408  Approximant approximant /**< approximant (NRSur7dq2 or NRSur7dq4) */
2409 ) {
2410  SphHarmTimeSeries *hlms = NULL;
2411 
2412  LALValue* ModeArray
2414 
2415  int ell, m;
2416  if ( ModeArray == NULL ) {
2417  // Use all available modes
2418  ModeArray = XLALSimInspiralCreateModeArray();
2419  for (ell=2; ell <= NRSUR_LMAX; ell++) {
2421  }
2422  } // Otherwise, use the specified modes.
2423 
2424  // If m1 < m2, switch the BH labels so that q = m1/m2 >= 1 always.
2425  bool labels_switched = PrecessingNRSur_switch_labels_if_needed(&m1, &m2,
2426  &S1x, &S1y, &S1z, &S2x, &S2y, &S2z);
2427 
2428  // Parameters
2429  REAL8 massA = m1 / LAL_MSUN_SI;
2430  REAL8 massB = m2 / LAL_MSUN_SI;
2431  REAL8 Mtot = massA+massB;
2432  REAL8 Mtot_sec = Mtot * LAL_MTSUN_SI; /* Total mass in seconds */
2433  REAL8 q = massA / massB;
2434  REAL8 chiA0[3], chiB0[3];
2435  chiA0[0] = S1x;
2436  chiA0[1] = S1y;
2437  chiA0[2] = S1z;
2438  chiB0[0] = S2x;
2439  chiB0[1] = S2y;
2440  chiB0[2] = S2z;
2441 
2442  REAL8 omegaMin_dimless;
2443  REAL8 omegaRef_dimless;
2444  int ret = get_dimless_omega(&omegaMin_dimless, &omegaRef_dimless,
2445  fMin, fRef, Mtot_sec);
2446  if(ret != XLAL_SUCCESS) {
2447  if(ModeArray) XLALDestroyValue(ModeArray);
2448  XLAL_ERROR_NULL(XLAL_EFUNC, "Failed to process fMin/fRef");
2449  }
2450 
2451  // When generating the modes take the inertial frame to be aligned with
2452  // the coorbital frame at the reference point. This means that the
2453  // quaternion from inertial frame to coprecessing frame is identity, and
2454  // the init_orbphase = 0 between the coprecessing frame and the coorbital
2455  // frame at the reference point. This agrees with the LAL convention, and
2456  // phiRef should only used in the Ylms when computing the polarizations.
2457  REAL8 init_quat[4];
2458  init_quat[0] = 1.0;
2459  init_quat[1] = 0.0;
2460  init_quat[2] = 0.0;
2461  init_quat[3] = 0.0;
2462  REAL8 init_orbphase = 0;
2463 
2464  // Evaluate the model modes
2465  MultiModalWaveform *h_inertial = NULL;
2466  PrecessingNRSurData *__sur_data = PrecessingNRSur_core(&h_inertial, q,
2467  chiA0, chiB0, omegaRef_dimless, init_orbphase, init_quat,
2468  ModeArray, LALparams, approximant);
2469  if (!h_inertial) {
2470  if(ModeArray) XLALDestroyValue(ModeArray);
2471  XLAL_PRINT_INFO("PrecessingNRSur_core failed!");
2472  return hlms;
2473  }
2474 
2475  // Scale the amplitude based on the distance
2476  REAL8 a0 = Mtot_sec * LAL_C_SI / distance;
2477 
2478  // Setup dimensionless output times
2479  gsl_vector *model_times = __sur_data->t_coorb;
2480  REAL8 deltaT_dimless = deltaT / Mtot_sec;
2481  REAL8 t0 = gsl_vector_get(model_times, 0);
2482  REAL8 start_freq = PrecessingNRSur_StartFrequency(m1, m2, S1x, S1y, S1z, S2x, S2y, S2z, __sur_data);
2483  if (fMin >= start_freq) {
2484  t0 = PrecessingNRSur_get_t_ref(omegaMin_dimless, q, chiA0, chiB0,
2485  init_quat, init_orbphase, __sur_data);
2486  } else if (fMin > 0) {
2487  // Cleanup and exit
2488  if(ModeArray) XLALDestroyValue(ModeArray);
2489  MultiModalWaveform_Destroy(h_inertial);
2490  XLAL_ERROR_NULL(XLAL_EDOM, "fMin should be 0 or >= %0.8f for this configuration, got %0.8f", start_freq, fMin);
2491  }
2492  size_t length = model_times->size;
2493  REAL8 tf = gsl_vector_get(model_times, length-1);
2494  int nt = (int) ceil((tf - t0) / deltaT_dimless);
2495  gsl_vector *output_times = gsl_vector_alloc(nt);
2496  int j;
2497  for (j=0; j<nt; j++) {
2498  gsl_vector_set(output_times, j, t0 + deltaT_dimless*j);
2499  }
2500 
2501  // Setup interpolation onto dimensionless output times
2502  REAL8 t;
2504  XLALGPSAdd( &epoch, Mtot_sec * t0);
2505  gsl_interp_accel *acc = gsl_interp_accel_alloc();
2506 
2507  // Create LAL modes
2508  COMPLEX16TimeSeries *tmp_mode;
2509  int i=0;
2510  bool skip_ell_modes;
2511  char mode_name[32];
2512  REAL8 prefactor = 1;
2513  for (ell=2; ell <= NRSUR_LMAX; ell++) {
2514  /* If *any* co-orbital frame mode of this ell is non-zero we need to sum over all modes of this ell,
2515  * as the frame rotation will have mixed the modes. */
2516  skip_ell_modes = true;
2517  for (m=-ell; m<=ell; m++) {
2518  if (XLALSimInspiralModeArrayIsModeActive(ModeArray, ell, m) == 1) skip_ell_modes = false;
2519  }
2520  if (skip_ell_modes) {
2521  i += 2*ell + 1;
2522  continue;
2523  }
2524  for (m=-ell; m<=ell; m++) {
2525 
2526  // If m1 < m2, the BH labels were switched so that m1 > m2. This
2527  // generates the requested waveform but with a rotation of pi about
2528  // the z-axis. We now undo this by multiplying all odd-m modes by
2529  // -1, which is the same as a pi rotation; the even-m modes don't
2530  // need this fix as they are insensitive to pi rotations about
2531  // z-axis.
2532  if ((m%2 != 0) && (labels_switched)) {
2533  prefactor = -1;
2534  } else {
2535  prefactor = 1;
2536  }
2537 
2538  gsl_vector_scale(h_inertial->modes_real_part[i], prefactor*a0);
2539  gsl_vector_scale(h_inertial->modes_imag_part[i], prefactor*a0);
2540  snprintf(mode_name, sizeof(mode_name), "(%d, %d) mode", ell, m);
2541  tmp_mode = XLALCreateCOMPLEX16TimeSeries(mode_name, &epoch, 0.0,
2542  deltaT, &lalStrainUnit, nt);
2543  gsl_spline *spl_re = gsl_spline_alloc(gsl_interp_cspline, length);
2544  gsl_spline *spl_im = gsl_spline_alloc(gsl_interp_cspline, length);
2545  gsl_spline_init(spl_re, model_times->data,
2546  h_inertial->modes_real_part[i]->data, length);
2547  gsl_spline_init(spl_im, model_times->data,
2548  h_inertial->modes_imag_part[i]->data, length);
2549  for (j=0; j<nt; j++) {
2550  t = gsl_vector_get(output_times, j);
2551  tmp_mode->data->data[j] = gsl_spline_eval(spl_re, t, acc);
2552  tmp_mode->data->data[j] += I * gsl_spline_eval(spl_im, t, acc);
2553  }
2554  gsl_spline_free(spl_re);
2555  gsl_spline_free(spl_im);
2556  hlms = XLALSphHarmTimeSeriesAddMode(hlms, tmp_mode, ell, m);
2558  i += 1;
2559  }
2560  }
2561 
2562  // Cleanup
2563  MultiModalWaveform_Destroy(h_inertial);
2564  gsl_vector_free(output_times);
2565  gsl_interp_accel_free(acc);
2566  if(ModeArray) {
2567  XLALDestroyValue(ModeArray);
2568  }
2569 
2570  return hlms;
2571 }
2572 
2573 /**
2574  * This function evaluates the NRSur7dq2 or NRSur7dq4 surrogate model and
2575  * returns the precessing frame dynamics.
2576  * Papers:
2577  * NRSur7dq2: https://arxiv.org/abs/1705.07089
2578  * NRSur7dq4: https://arxiv.org/abs/1905.09300
2579  *
2580  * Example usage:
2581  * t_dynamics, quat0, quat1, quat2, quat3, orbphase, chiAx, chiAy, chiAz,
2582  * chiBx, chiBy, chiBz = lalsim.PrecessingNRSurDynamics(
2583  * q, chiA0x, chiA0y, chiA0z, chiB0x, chiB0y, chiB0z,
2584  * omegaRef_dimless, init_quat0, init_quat1, init_quat2, init_quat3,
2585  * init_orbphase, LALparams, approxTag)
2586  *
2587  * INPUTS:
2588  * q
2589  * mass ratio, mA/mB >= 1.
2590  * chiA0x
2591  * chiA0y
2592  * chiA0z
2593  * spin of the heavier BH at the reference epoch, in the coorbital
2594  * frame, as defined in the above papers. This agrees with the LAL
2595  * convention, see LALSimInspiral.h for definition of the LAL frame.
2596  * chiB0x
2597  * chiB0y
2598  * chiB0z
2599  * spin of the lighter BH at the reference epoch, in the coorbital
2600  * frame.
2601  * omegaRef_dimless
2602  * reference dimensionless orbital frequency in the coprecessing
2603  * frame, this is used to set the reference epoch. The coprecessing
2604  * frame is defined in the above papers.
2605  * init_quat0
2606  * init_quat1
2607  * init_quat2
2608  * init_quat3
2609  * coprecessing frame quaternion at the reference epoch. To follow
2610  * the LAL convention this should be [1, 0, 0, 0], but the surrogate
2611  * allows arbitrary unit quaternions.
2612  * init_orbphase
2613  * orbital phase in the coprecessing frame at the reference epoch. To
2614  * follow the LAL convention this should be 0, but the surrogate
2615  * allows arbitrary initial orbital phase.
2616  * LALparams
2617  * LAL dictionary containing additional parameters, if required.
2618  * Initialized with: LALparams = lal.CreateDict()
2619  *
2620  * Extrapolation options:
2621  * The surrogate models are trained on the following ranges:
2622  * NRSur7dq2: q <= 2.01, |chi_1|, |chi_2| <= 0.81.
2623  * NRSur7dq4: q <= 4.01, |chi_1|, |chi_2| <= 0.81.
2624  * If you want a guarantee of accuracy you should stick to the above
2625  * ranges.
2626  *
2627  * We allow extrapolation to the following ranges, but with a warning:
2628  * NRSur7dq2: q <= 3.01, |chi_1|, |chi_2| <= 1.
2629  * NRSur7dq4: q <= 6.01, |chi_1|, |chi_2| <= 1.
2630  * We expect the models to be reasonable when extrapolated to these
2631  * ranges.
2632  *
2633  * Beyond the above ranges, we raise an error. If you want to
2634  * extrapolate beyond these limits you can specify
2635  * unlimited_extrapolation = 1 in your LALparams as follows:
2636  * USE AT YOUR OWN RISK!!
2637  * lal.DictInsertUINT4Value(LALparams,"unlimited_extrapolation",1)
2638  * approxTag
2639  * LAL object specifying the surrogate model to use. Example:
2640  * approxTag = lalsim.SimInspiralGetApproximantFromString('NRSur7dq4')
2641  * OUTPUTS:
2642  * t_dynamics
2643  * time values at which the dynamics are returned. These are the
2644  * dynamics time nodes as described in the above papers and are
2645  * nonuniform and sparse.
2646  * quat0
2647  * quat1
2648  * quat2
2649  * quat3
2650  * time series of coprecessing frame quaternions. These are used to
2651  * transform between the inertial frame and the coprecessing frames.
2652  * orbphase
2653  * orbital phase time series in the coprecessing frame. This is used
2654  * to transform between the coprecessing and the coorbital frames.
2655  * chiAx
2656  * chiAy
2657  * chiAz
2658  * time series of spin of the heavier BH in the coprecessing frame.
2659  * For convenience for the main use case with surrogate remnant fits,
2660  * these are in coprecessing frame rather then the reference frame in
2661  * which the input spins are specified.
2662  * chiBx
2663  * chiBy
2664  * chiBz
2665  * time series of spin of the heavier BH in the coprecessing frame.
2666  */
2668  gsl_vector **t_dynamics, /**< Output: Time array at which the dynamics are returned. */
2669  gsl_vector **quat0, /**< Output: Time series of 0th index of coprecessing frame quaternion. */
2670  gsl_vector **quat1, /**< Output: Time series of 1st index of coprecessing frame quaternion. */
2671  gsl_vector **quat2, /**< Output: Time series of 2nd index of coprecessing frame quaternion. */
2672  gsl_vector **quat3, /**< Output: Time series of 3rd index of coprecessing frame quaternion. */
2673  gsl_vector **orbphase, /**< Output: Time series of orbital phase in the coprecessing frame. */
2674  gsl_vector **chiAx, /**< Output: Time series of x-comp of dimensionless spin of BhA in the coprecessing frame. */
2675  gsl_vector **chiAy, /**< Output: Time series of y-comp of dimensionless spin of BhA in the coprecessing frame. */
2676  gsl_vector **chiAz, /**< Output: Time series of z-comp of dimensionless spin of BhA in the coprecessing frame. */
2677  gsl_vector **chiBx, /**< Output: Time series of x-comp of dimensionless spin of BhB in the coprecessing frame. */
2678  gsl_vector **chiBy, /**< Output: Time series of y-comp of dimensionless spin of BhB in the coprecessing frame. */
2679  gsl_vector **chiBz, /**< Output: Time series of z-comp of dimensionless spin of BhB in the coprecessing frame. */
2680  REAL8 q, /**< mass ratio m1/m2 >= 1. */
2681  REAL8 chiA0x, /**< x-comp of dimensionless spin of BhA in the coorbital frame at the reference epoch. */
2682  REAL8 chiA0y, /**< y-comp of dimensionless spin of BhA in the coorbital frame at the reference epoch. */
2683  REAL8 chiA0z, /**< z-comp of dimensionless spin of BhA in the coorbital frame at the reference epoch. */
2684  REAL8 chiB0x, /**< x-comp of dimensionless spin of BhB in the coorbital frame at the reference epoch. */
2685  REAL8 chiB0y, /**< y-comp of dimensionless spin of BhB in the coorbital frame at the reference epoch. */
2686  REAL8 chiB0z, /**< z-comp of dimensionless spin of BhB in the coorbital frame at the reference epoch. */
2687  REAL8 omegaRef_dimless, /**< Dimensionless orbital frequency (rad/M) in the coprecessing frame at the reference epoch.*/
2688  REAL8 init_quat0, /**< 0th comp of the coprecessing frame quaternion at the reference epoch.*/
2689  REAL8 init_quat1, /**< 1st comp of the coprecessing frame quaternion at the reference epoch.*/
2690  REAL8 init_quat2, /**< 2nd comp of the coprecessing frame quaternion at the reference epoch.*/
2691  REAL8 init_quat3, /**< 3rd comp of the coprecessing frame quaternion at the reference epoch.*/
2692  REAL8 init_orbphase, /**< orbital phase in the coprecessing frame at the reference epoch. */
2693  LALDict* LALparams, /**< Dict with extra parameters. */
2694  Approximant approximant /**< approximant (NRSur7dq2 or NRSur7dq4). */
2695 ) {
2696 
2697  // Load surrogate data if needed. If not, just access the loaded data.
2699  if (!__sur_data->setup) {
2700  XLAL_ERROR(XLAL_EFAILED, "Error loading surrogate data.\n");
2701  }
2702 
2703  // Input spins at reference epoch
2704  // The input values are in the coorbital frame at omegaRef_dimless, but
2705  // the surrogate wants them in the coprecessing frame, so we transform
2706  // from the coorbital frame to the coprecessing frame here.
2707  REAL8 chiA0[3], chiB0[3];
2708  chiA0[0] = chiA0x * cos(init_orbphase) - chiA0y * sin(init_orbphase);
2709  chiA0[1] = chiA0y * cos(init_orbphase) + chiA0x * sin(init_orbphase);
2710  chiA0[2] = chiA0z;
2711  chiB0[0] = chiB0x * cos(init_orbphase) - chiB0y * sin(init_orbphase);
2712  chiB0[1] = chiB0y * cos(init_orbphase) + chiB0x * sin(init_orbphase);
2713  chiB0[2] = chiB0z;
2714 
2715  // Coprecessing frame quaternion at reference epoch
2716  REAL8 init_quat[4];
2717  init_quat[0] = init_quat0;
2718  init_quat[1] = init_quat1;
2719  init_quat[2] = init_quat2;
2720  init_quat[3] = init_quat3;
2721 
2722  // Get surrogate dynamics. The spins are returned in the coprecessing
2723  // frame.
2724  int n_ds = (__sur_data->t_ds)->size;
2725  REAL8 *dynamics_data = XLALCalloc(n_ds * 11, sizeof(REAL8));
2726  int ret = PrecessingNRSur_IntegrateDynamics(dynamics_data, q,
2727  chiA0, chiB0, omegaRef_dimless, init_orbphase, init_quat,
2728  LALparams, __sur_data->PrecessingNRSurVersion);
2729  if(ret != XLAL_SUCCESS) {
2730  XLALFree(dynamics_data);
2731  XLAL_ERROR(XLAL_FAILURE, "Failed to integrate dynamics");
2732  }
2733 
2734  // Put output into appropriate vectors
2735  int i, j;
2736  gsl_vector *dynamics[11];
2737  for (j=0; j<11; j++) {
2738  dynamics[j] = gsl_vector_alloc(n_ds);
2739  }
2740  for (i=0; i<n_ds; i++) {
2741  for (j=0; j<11; j++) {
2742  gsl_vector_set(dynamics[j], i, dynamics_data[11*i + j]);
2743  }
2744  }
2745 
2746  // We want to make a copy as *t_dynamics gets destroyed after it is
2747  // returned by SWIG. So, if we just pass __sur_data->t_ds, it fails
2748  // upon a second call.
2749  *t_dynamics = gsl_vector_alloc(n_ds);
2750  gsl_vector_memcpy(*t_dynamics, __sur_data->t_ds);
2751 
2752  *quat0 = dynamics[0];
2753  *quat1 = dynamics[1];
2754  *quat2 = dynamics[2];
2755  *quat3 = dynamics[3];
2756  *orbphase = dynamics[4];
2757  *chiAx = dynamics[5]; // These are in the coprecessing frame
2758  *chiAy = dynamics[6];
2759  *chiAz = dynamics[7];
2760  *chiBx = dynamics[8];
2761  *chiBy = dynamics[9];
2762  *chiBz = dynamics[10];
2763 
2764  XLALFree(dynamics_data);
2765 
2766  return XLAL_SUCCESS;
2767 }
2768 
2769 /** @} */
2770 /** @} */
2771 
struct tagLALH5File LALH5File
struct tagLALH5Dataset LALH5Dataset
int XLALDictContains(const LALDict *dict, const char *key)
UINT4 XLALDictLookupUINT4Value(LALDict *dict, const char *key)
static const double Q_MAX
static PrecessingNRSurData __lalsim_NRSur7dq2_data
Global surrogate data.
static PrecessingNRSurData __lalsim_NRSur7dq4_data
static const double NRSUR7DQ4_Q_MAX
static const double NRSUR7DQ4_Q_FIT_OFFSET
static const char NRSUR7DQ4_DATAFILE[]
static const double NRSUR7DQ4_START_TIME
static const double NRSUR7DQ4_Q_FIT_SLOPE
static const double NRSUR7DQ2_Q_FIT_SLOPE
static const double NRSUR7DQ2_Q_MAX
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 const double NRSUR7DQ4_CHI_MAX
static const double NRSUR7DQ2_START_TIME
static const char NRSUR7DQ2_DATAFILE[]
Auxiliary functions for SEOBNRv1/v2 reduced order modeling codes in LALSimIMRSEOBNRv2ROMDoubleSpin....
LALValue * XLALSimInspiralWaveformParamsLookupModeArray(LALDict *params)
const char * name
static void NRSur_ab4_dy(REAL8 *dy, REAL8 *k1, REAL8 *k2, REAL8 *k3, REAL8 *k4, REAL8 t12, REAL8 t23, REAL8 t34, REAL8 t45, int dim)
static void MultiModalWaveform_Destroy(MultiModalWaveform *wave)
Destroy a MultiModalWaveform structure; free all associated memory.
static void TransformModesCoorbitalToInertial(MultiModalWaveform *h, MultiModalWaveform *h_coorb, gsl_vector **quat, gsl_vector *orbphase)
Transforms modes from the coorbital frame to the inertial frame.
static void MultiModalWaveform_Init(MultiModalWaveform **wave, int LMax, int n_times)
Initialize a MultiModalWaveforme.
static UNUSED int get_dimless_omega(REAL8 *omegaMin_dimless, REAL8 *omegaRef_dimless, const REAL8 fMin, const REAL8 fRef, const REAL8 Mtot_sec)
Wrapper to get dimensionless omegaMin/omegaRef from fMin/fRef.
void XLALDestroyValue(LALValue *value)
static REAL8 sum(REAL8 *arr1, REAL8 *arr2, int k)
double dt
Definition: bh_ringdown.c:113
double i
Definition: bh_ringdown.c:118
sigmaKerr data[0]
#define XLAL_FILE_RESOLVE_PATH(fname)
UINT4Vector * XLALH5DatasetQueryDims(LALH5Dataset UNUSED *dset)
LALH5File * XLALH5GroupOpen(LALH5File UNUSED *file, const char UNUSED *name)
LALH5File * XLALH5FileOpen(const char UNUSED *path, const char UNUSED *mode)
LALH5Dataset * XLALH5DatasetRead(LALH5File UNUSED *file, const char UNUSED *name)
#define LAL_C_SI
#define LAL_MSUN_SI
#define LAL_PI
#define LAL_MTSUN_SI
#define LIGOTIMEGPSZERO
double complex COMPLEX16
double REAL8
uint32_t UINT4
void * XLALCalloc(size_t m, size_t n)
void * XLALMalloc(size_t n)
void XLALFree(void *p)
#define LAL_SIM_L_MAX_MODE_ARRAY
Maximum L spherical harmonic mode that is supported in Mode Array.
Approximant
Enum that specifies the PN approximant to be used in computing the waveform.
@ NRSur7dq4
q=4 extension of NRSur7dq2, arxiv: 1905.09300
@ NRSur7dq2
Time domain, fully precessing NR surrogate model with up to ell=4 modes, arxiv: 1705....
LALValue * XLALSimInspiralCreateModeArray(void)
Create a LALValue pointer to store the mode array.
int XLALSimInspiralModeArrayIsModeActive(LALValue *modes, unsigned l, int m)
LALValue * XLALSimInspiralModeArrayActivateAllModesAtL(LALValue *modes, unsigned l)
static PrecessingNRSurData * PrecessingNRSur_LoadData(Approximant approximant)
Wrapper for Loading NRSur7dq2 and NRSur7dq4 data.
static REAL8 PrecessingNRSur_get_omega(size_t node_index, REAL8 q, REAL8 *y0, PrecessingNRSurData *__sur_data)
Computes the orbital angular frequency at a dynamics node.
static bool NRSur7dq4_IsSetup(void)
Helper function which returns whether or not the global NRSur7dq4 surrogate data has been initialized...
static PrecessingNRSurData * PrecessingNRSur_core(MultiModalWaveform **h, REAL8 q, REAL8 *chiA0, REAL8 *chiB0, REAL8 omega_ref, REAL8 init_orbphase, REAL8 *init_quat, LALValue *ModeArray, LALDict *LALparams, Approximant approximant)
static void PrecessingNRSur_normalize_y(REAL8 chiANorm, REAL8 chiBNorm, REAL8 *y)
static int PrecessingNRSur_initialize_RK4(REAL8 *dynamics_data, REAL8 *time_steps, REAL8 *dydt0, REAL8 *dydt1, REAL8 *dydt2, REAL8 *dydt3, REAL8 normA, REAL8 normB, REAL8 q, int i0, PrecessingNRSurData *__sur_data)
Initializes the AB4 ODE system from a single arbitrary dynamics node by taking 3 RK4 steps.
static int PrecessingNRSur_Init(PrecessingNRSurData *data, LALH5File *file, UINT4 PrecessingNRSurVersion)
Initialize a PrecessingNRSurData structure from an open hdf5 file.
static void NRSur7dq2_Init_LALDATA(void)
This needs to be called once, before __lalsim_NRSur7dq2_data is used.
static void PrecessingNRSur_rotate_spins(gsl_vector **chiA, gsl_vector **chiB, gsl_vector *phi_vec)
Transforms chiA and chiB from the coprecessing frame to the coorbital frame using the orbital phase.
SphHarmTimeSeries * XLALSimInspiralPrecessingNRSurModes(REAL8 deltaT, REAL8 m1, REAL8 m2, REAL8 S1x, REAL8 S1y, REAL8 S1z, REAL8 S2x, REAL8 S2y, REAL8 S2z, REAL8 fMin, REAL8 fRef, REAL8 distance, LALDict *LALparams, Approximant approximant)
This function evaluates the NRSur7dq2 or NRSur7dq4 surrogate model and returns the inertial frame mod...
static void PrecessingNRSur_LoadWaveformDataPiece(LALH5File *sub, WaveformDataPiece **data, bool invert_sign)
Loads a single NRSur coorbital waveform data piece from file into a WaveformDataPiece.
static void PrecessingNRSur_get_time_deriv_from_index(REAL8 *dydt, int i0, REAL8 q, REAL8 *y, PrecessingNRSurData *__sur_data)
Compute dydt at a given dynamics node, where y is the numerical solution to the dynamics ODE.
int XLALSimInspiralPrecessingNRSurPolarizations(REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, REAL8 phiRef, REAL8 inclination, REAL8 deltaT, REAL8 m1, REAL8 m2, REAL8 distance, REAL8 fMin, REAL8 fRef, REAL8 s1x, REAL8 s1y, REAL8 s1z, REAL8 s2x, REAL8 s2y, REAL8 s2z, LALDict *LALparams, Approximant approximant)
This function evaluates the NRSur7dq2 or NRSur7dq4 surrogate model and sums over all ell <= 4 modes t...
static void PrecessingNRSur_get_time_deriv(REAL8 *dydt, REAL8 t, REAL8 q, REAL8 *y, PrecessingNRSurData *__sur_data)
Compute dydt at any time by evaluating dydt at 4 nearby dynamics nodes and using cubic spline interpo...
static bool NRSur7dq2_IsSetup(void)
Helper function which returns whether or not the global NRSur7dq2 surrogate data has been initialized...
static void PrecessingNRSur_assemble_dydt(REAL8 *dydt, REAL8 *y, REAL8 *Omega_coorb_xy, REAL8 omega, REAL8 *chiA_dot, REAL8 *chiB_dot)
static void PrecessingNRSur_LoadCoorbitalEllModes(WaveformFixedEllModeData **coorbital_mode_data, LALH5File *file, int i)
Load the WaveformFixedEllModeData from file for a single value of ell.
static void PrecessingNRSur_integrate_AB4(REAL8 *dynamics_data, REAL8 *time_steps, REAL8 *dydt0, REAL8 *dydt1, REAL8 *dydt2, REAL8 *dydt3, REAL8 normA, REAL8 normB, REAL8 q, int i_start, PrecessingNRSurData *__sur_data)
Integrates the AB4 ODE system in time forwards, and backwards if needed.
static void NRSur7dq4_Init_LALDATA(void)
This needs to be called once, before __lalsim_NRSur7dq4_data is used.
static int NRSur7dq4_effective_spins(REAL8 *chiHat, REAL8 *chi_a, const REAL8 q, const REAL8 chi1z, const REAL8 chi2z)
int XLALPrecessingNRSurDynamics(gsl_vector **t_dynamics, gsl_vector **quat0, gsl_vector **quat1, gsl_vector **quat2, gsl_vector **quat3, gsl_vector **orbphase, gsl_vector **chiAx, gsl_vector **chiAy, gsl_vector **chiAz, gsl_vector **chiBx, gsl_vector **chiBy, gsl_vector **chiBz, REAL8 q, REAL8 chiA0x, REAL8 chiA0y, REAL8 chiA0z, REAL8 chiB0x, REAL8 chiB0y, REAL8 chiB0z, REAL8 omegaRef_dimless, REAL8 init_quat0, REAL8 init_quat1, REAL8 init_quat2, REAL8 init_quat3, REAL8 init_orbphase, LALDict *LALparams, Approximant approximant)
This function evaluates the NRSur7dq2 or NRSur7dq4 surrogate model and returns the precessing frame d...
static int PrecessingNRSur_initialize_at_dynamics_node(REAL8 *dynamics_data, REAL8 t_ref, REAL8 q, REAL8 *chiA0, REAL8 *chiB0, REAL8 init_orbphase, REAL8 *init_quat, REAL8 normA, REAL8 normB, PrecessingNRSurData *__sur_data)
Initialize the dynamics ODE at a dynamics node.
static void PrecessingNRSur_LoadFitData(FitData **fit_data, LALH5File *sub, const char *name)
Loads a single fit for NRSur7dq2 or NRSur7dq4.
static void PrecessingNRSur_normalize_results(REAL8 normA, REAL8 normB, gsl_vector **quat, gsl_vector **chiA, gsl_vector **chiB)
static int PrecessingNRSur_IntegrateDynamics(REAL8 *dynamics_data, REAL8 q, REAL8 *chiA0, REAL8 *chiB0, REAL8 omega_ref, REAL8 init_orbphase, REAL8 *init_quat, LALDict *LALparams, UINT4 PrecessingNRSurVersion)
This is the main NRSur dynamics surrogate integration routine.
static void NRSur7dq4_eval_vector_fit(REAL8 *res, VectorFitData *data, REAL8 *x)
REAL8 NRSur7dq4_eval_fit(FitData *data, REAL8 *x)
static void NRSur7dq2_eval_vector_fit(REAL8 *res, VectorFitData *data, REAL8 *x)
static bool PrecessingNRSur_switch_labels_if_needed(REAL8 *m1, REAL8 *m2, REAL8 *s1x, REAL8 *s1y, REAL8 *s1z, REAL8 *s2x, REAL8 *s2y, REAL8 *s2z)
If m1<m2, swaps the labels for the two BHs such that Bh1 is always heavier.
static void PrecessingNRSur_LoadDynamicsNode(DynamicsNodeFitData **ds_node_data, LALH5File *sub, int i, UINT4 PrecessingNRSurVersion)
Loads the data for a single dynamics node into a DynamicsNodeFitData struct.
REAL8 ipow(REAL8 base, int exponent)
Helper function for integer powers.
static void PrecessingNRSur_eval_vector_fit(REAL8 *res, VectorFitData *data, REAL8 *x, PrecessingNRSurData *__sur_data)
static gsl_vector * spline_array_interp(gsl_vector *xout, gsl_vector *x, gsl_vector *y)
Do cubic spline interpolation using a gsl_interp_cspline.
static REAL8 PrecessingNRSur_get_t_ref(REAL8 omega_ref, REAL8 q, REAL8 *chiA0, REAL8 *chiB0, REAL8 *init_quat, REAL8 init_orbphase, PrecessingNRSurData *__sur_data)
Computes a reference time from a reference orbital angular frequency.
REAL8 PrecessingNRSur_eval_fit(FitData *data, REAL8 *x, PrecessingNRSurData *__sur_data)
static void PrecessingNRSur_initialize_RK4_with_half_nodes(REAL8 *dynamics_data, REAL8 *time_steps, REAL8 *dydt0, REAL8 *dydt1, REAL8 *dydt2, REAL8 *dydt3, REAL8 normA, REAL8 normB, REAL8 q, PrecessingNRSurData *__sur_data)
Initializes the AB4 ODE system from the surrogate start time by taking 3 RK4 steps,...
static void PrecessingNRSur_eval_data_piece(gsl_vector *result, REAL8 q, gsl_vector **chiA, gsl_vector **chiB, WaveformDataPiece *data, PrecessingNRSurData *__sur_data)
Evaluates a single NRSur coorbital waveoform data piece.
static void PrecessingNRSur_ds_fit_x(REAL8 *x, REAL8 q, REAL8 *y)
REAL8 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...
REAL8 NRSur7dq2_eval_fit(FitData *data, REAL8 *x)
static REAL8 cubic_interp(REAL8 xout, REAL8 *x, REAL8 *y)
Cubic interpolation of 4 data points This gives a much closer result to scipy.interpolate....
static void NRSur7dq4_LoadVectorFitData(VectorFitData **vector_fit_data, LALH5File *sub, const char *name, const size_t size)
Loads a vector fit for NRSur7dq4.
SphHarmTimeSeries * XLALSphHarmTimeSeriesAddMode(SphHarmTimeSeries *appended, const COMPLEX16TimeSeries *inmode, UINT4 l, INT4 m)
Prepend a node to a linked list of SphHarmTimeSeries, or create a new head.
static const INT4 m
static const INT4 q
COMPLEX16 XLALSpinWeightedSphericalHarmonic(REAL8 theta, REAL8 phi, int s, int l, int m)
REAL8TimeSeries * XLALCreateREAL8TimeSeries(const CHAR *name, const LIGOTimeGPS *epoch, REAL8 f0, REAL8 deltaT, const LALUnit *sampleUnits, size_t length)
void XLALDestroyCOMPLEX16TimeSeries(COMPLEX16TimeSeries *series)
COMPLEX16TimeSeries * XLALCreateCOMPLEX16TimeSeries(const CHAR *name, const LIGOTimeGPS *epoch, REAL8 f0, REAL8 deltaT, const LALUnit *sampleUnits, size_t length)
const LALUnit lalStrainUnit
#define XLAL_ERROR_VOID(...)
#define XLAL_ERROR_REAL8(...)
#define XLAL_PRINT_INFO(...)
#define XLAL_ERROR_NULL(...)
#define XLAL_ERROR(...)
#define XLAL_PRINT_WARNING(...)
#define XLAL_CHECK_ABORT(assertion)
int XLALPrintError(const char *fmt,...) _LAL_GCC_PRINTF_FORMAT_(1
#define XLAL_IS_REAL8_FAIL_NAN(val)
XLAL_SUCCESS
XLAL_EFUNC
XLAL_EDOM
XLAL_EIO
XLAL_EINVAL
XLAL_EFAILED
XLAL_FAILURE
LIGOTimeGPS * XLALGPSAdd(LIGOTimeGPS *epoch, REAL8 dt)
list x
list y
path
int pthread_once(pthread_once_t *once_control, void(*init_routine)(void))
COMPLEX16Sequence * data
COMPLEX16 * data
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.
gsl_vector ** modes_real_part
The real part of each mode - n_modes vectors of length n_times.
gsl_vector ** modes_imag_part
The imaginary part of each mode - n_modes vectors of length n_times.
int n_times
Number of time samples in each mode.
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.
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
Structure to carry a collection of spherical harmonic modes in COMPLEX16 time series.
UINT4 * data
Data used in a single vector fit NOTE: basisFunctionOrders, coefs, componentIndices,...
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),...
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.
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.
LIGOTimeGPS epoch
Definition: unicorn.c:20
double deltaT
Definition: unicorn.c:24