LALSimIMRPhenomX_internals.c 87.7 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
/*
 * Copyright (C) 2018 Geraint Pratten
 *
 * 	This code builds on:
 * 		LALSimIMRPhenomD_internals.c
 *
 *  This program is free software; you can redistribute it and/or modify
 *  it under the terms of the GNU General Public License as published by
 *  the Free Software Foundation; either version 2 of the License, or
 *  (at your option) any later version.
 *
 *  This program is distributed in the hope that it will be useful,
 *  but WITHOUT ANY WARRANTY; without even the implied warranty of
 *  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 *  GNU General Public License for more details.
 *
 *  You should have received a copy of the GNU General Public License
 *  along with with program; see the file COPYING. If not, write to the
19
20
 *  Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
 *  MA  02110-1301  USA
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
 */


/**
 * \author Geraint Pratten
 *
 * Internal functions for IMRPhenomX, arXiv:2001.11412
 */

/* IMRPhenomX Header Files */
#include "LALSimIMRPhenomXUtilities.h"
#include "LALSimIMRPhenomX_internals.h"

#include "LALSimIMRPhenomX_inspiral.c"
#include "LALSimIMRPhenomX_intermediate.c"
#include "LALSimIMRPhenomX_ringdown.c"
#include "LALSimIMRPhenomX_qnm.c"
//#include "LALSimIMRPhenomX_tidal.c"
//#include "LALSimIMRPhenomX_precessing.c"

/* LAL Header Files */
#include <lal/LALSimIMR.h>
#include <lal/LALConstants.h>
#include <lal/Date.h>
#include <lal/FrequencySeries.h>
#include <lal/Units.h>

/* GSL Header Files */
#include <gsl/gsl_linalg.h>

/* This struct is used to pre-cache useful powers of frequency, avoiding numerous expensive operations */
int IMRPhenomX_Initialize_Powers(IMRPhenomX_UsefulPowers *p, REAL8 number)
{
	XLAL_CHECK(0 != p, XLAL_EFAULT, "p is NULL");
	XLAL_CHECK(number >= 0, XLAL_EDOM, "number must be non-negative");

	double sixth      = pow(number, 1.0 / 6.0);
	double m_sixth    = 1.0 / sixth;

	p->one_sixth      = sixth;
	p->m_one_sixth    = m_sixth;

	p->m_one          = 1.0 / number;
	p->itself         = number;

	p->one_third      = sixth * sixth;
	p->m_one_third    = 1.0 / (p->one_third);

	p->two_thirds     = p->one_third * p->one_third;
	p->m_two_thirds   = 1.0 / p->two_thirds;

	p->four_thirds    = p->two_thirds * p->two_thirds;
	p->m_four_thirds  = 1.0 / p->four_thirds;

	p->five_thirds    = p->four_thirds * p->one_third;
	p->m_five_thirds  = 1.0 / p->five_thirds;

	p->seven_thirds   = p->four_thirds * number;
	p->m_seven_thirds = 1.0 / p->seven_thirds;

	p->eight_thirds   = p->seven_thirds * p->one_third;
	p->m_eight_thirds = 1.0 / p->eight_thirds;

	p->two            = number   * number;
	p->three          = p->two   * number;
	p->four           = p->three * number;
	p->five           = p->four  * number;

	p->m_two          = 1.0 / p->two;
	p->m_three        = 1.0 / p->three;
	p->m_four         = 1.0 / p->four;

	p->seven_sixths   = p->one_sixth   * p->itself;
	p->m_seven_sixths = p->m_one_sixth * p->m_one;

	p->log            = log(number);
	p->sqrt           = p->one_sixth*p->one_sixth*p->one_sixth;

	return XLAL_SUCCESS;
}

/* A stripped down version of IMRPhenomX_Initialize_Powers for main production loop */
int IMRPhenomX_Initialize_Powers_Light(IMRPhenomX_UsefulPowers *p, REAL8 number)
{
	XLAL_CHECK(0 != p, XLAL_EFAULT, "p is NULL");
	XLAL_CHECK(number >= 0, XLAL_EDOM, "number must be non-negative");

	double sixth      = pow(number, 1.0 / 6.0);
	double m_sixth    = 1.0 / sixth;

	p->one_sixth      = sixth;
	p->m_one_sixth    = m_sixth;

	p->m_one          = 1.0 / number;
	p->itself         = number;

	p->one_third      = sixth * sixth;
	p->two_thirds     = p->one_third * p->one_third;

	p->m_two          = p->m_one   * p->m_one;
	p->m_three        = p->m_two   * p->m_one;

	p->seven_sixths   = p->one_sixth   * p->itself;
	p->m_seven_sixths = p->m_one_sixth * p->m_one;

	p->log            = log(number);

	return XLAL_SUCCESS;
}


int IMRPhenomXSetWaveformVariables(
	IMRPhenomXWaveformStruct *wf,
  const REAL8 m1_SI,
  const REAL8 m2_SI,
  const REAL8 chi1L_In,
  const REAL8 chi2L_In,
  const REAL8 deltaF,
  const REAL8 fRef,
  const REAL8 phi0,
  const REAL8 f_min,
  const REAL8 f_max,
  const REAL8 distance,
  const REAL8 inclination,
  LALDict *LALParams,
	UNUSED const UINT4 debug
)
{
	/* Copy model version to struct */
	wf->IMRPhenomXInspiralPhaseVersion      = XLALSimInspiralWaveformParamsLookupPhenomXInspiralPhaseVersion(LALParams);
151
	wf->LALparams = LALParams;
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
	wf->IMRPhenomXIntermediatePhaseVersion  = XLALSimInspiralWaveformParamsLookupPhenomXIntermediatePhaseVersion(LALParams);
	wf->IMRPhenomXRingdownPhaseVersion      = XLALSimInspiralWaveformParamsLookupPhenomXRingdownPhaseVersion(LALParams);

	wf->IMRPhenomXInspiralAmpVersion        = XLALSimInspiralWaveformParamsLookupPhenomXInspiralAmpVersion(LALParams);
	wf->IMRPhenomXIntermediateAmpVersion    = XLALSimInspiralWaveformParamsLookupPhenomXIntermediateAmpVersion(LALParams);
	wf->IMRPhenomXRingdownAmpVersion        = XLALSimInspiralWaveformParamsLookupPhenomXRingdownAmpVersion(LALParams);

	wf->debug = PHENOMXDEBUG;

	if(PHENOMXDEBUG)
	{
		printf("\n");
		printf("Inspiral Amp Version		: %d\n",wf->IMRPhenomXInspiralAmpVersion);
		printf("Intermediate Amp Version	: %d\n",wf->IMRPhenomXIntermediateAmpVersion);
		printf("Ringdown Amp Version		: %d\n",wf->IMRPhenomXRingdownAmpVersion);
		printf("\n");
		printf("Inspiral Phase Version		: %d\n",wf->IMRPhenomXInspiralPhaseVersion);
		printf("Intermediate Phase Version	: %d\n",wf->IMRPhenomXIntermediatePhaseVersion);
		printf("Ringdown Phase Version		: %d\n",wf->IMRPhenomXRingdownPhaseVersion);
		printf("\n");
	}

	/*
	First perform a sanity check to make sure that a legitimate waveform version has been called. If not, fail immediately.
	Every time a new version for one of the regions is added, user must add case below.
	*/
	int InspPhaseFlag = wf->IMRPhenomXInspiralPhaseVersion;

	/* Phase : Inspiral */
	switch( InspPhaseFlag )
	{
		// Canonical TaylorF2 up to 3.5PN with 4 pseudo-PN coefficients
		case 104:
		{
			break;
		}
		// Canonical TaylorF2 up to 3.5PN with 5 pseudo-PN coefficients
		case 105:
		{
			break;
		}
		// Extended TaylorF2 up to 4.5PN with 4 pseudo-PN coefficients
		case 114:
		{
			break;
		}
		// Extended TaylorF2 up to 4.5PN with 5 pseudo-PN coefficients
		case 115:
		{
			break;
		}
		// We must pass a recognised inspiral phase version. Otherwise fail.
		default:
		{
			XLAL_ERROR(XLAL_EINVAL, "Error in IMRPhenomX_SetWaveformVariables: IMRPhenomXInspiralPhaseVersion not recognized. Recommended flag is 104.\n");
		}
	}

	int IntPhaseFlag = wf->IMRPhenomXIntermediatePhaseVersion;
	/* Phase : Intermediate */
	switch( IntPhaseFlag )
	{
		// 4 intermediate coefficients
		case 104:
		{
			break;
		}
		// 5 intermediate coefficients
		case 105:
		{
			break;
		}
		default:
		{
			XLAL_ERROR(XLAL_EINVAL,"Error in IMRPhenomX_SetWaveformVariables: IMRPhenomXIntermediatePhaseVersion not recognized. Recommended flag is 105.\n");
		}
	}

	int RDPhaseFlag = wf->IMRPhenomXRingdownPhaseVersion;
	/* Phase : Ringdown */
	switch( RDPhaseFlag )
	{
		// 5 ringdown coefficients
		case 105:
		{
			break;
		}
		default:
		{
			XLAL_ERROR(XLAL_EINVAL,"Error in IMRPhenomX_SetWaveformVariables: IMRPhenomXRingdownPhaseVersion not recognized. Recommended flag is 105.\n");
		}
	}

	int InsAmpFlag = wf->IMRPhenomXInspiralAmpVersion;
	/* Amplitude : Inspiral */
	switch( InsAmpFlag )
	{
		// Canonical TaylorF2 with 4 pseudo-PN coefficients
		case 103:
		{
			break;
		}
		default:
		{
			XLAL_ERROR(XLAL_EINVAL,"Error in IMRPhenomX_SetWaveformVariables: IMRPhenomXInspiralAmpVersion not recognized. Recommended flag is 103.\n");
		}
	}

	int IntAmpFlag = wf->IMRPhenomXIntermediateAmpVersion;
	/* Amplitude : Intermediate */
	switch( IntAmpFlag )
	{
		// Intermediate region constructed from 4th order polynomial
		case 1043:
		{
			break;
		}
		// Intermediate region constructed from a 4th order polynomial
		case 104:
		{
			break;
		}
		// Intermediate region constructed from a 5th order polynomial
		case 105:
		{
			break;
		}
		default:
		{
			XLAL_ERROR(XLAL_EINVAL,"Error in IMRPhenomX_SetWaveformVariables: IMRPhenomXIntermediateAmpVersion not recognized. Recommended flag is 104.\n");
		}
	}

	int RDAmpFlag = wf->IMRPhenomXRingdownAmpVersion;
	/* Amplitude : Ringdown */
	switch( RDAmpFlag )
	{
		case 103:
		{
			break;
		}
		default:
		{
			XLAL_ERROR(XLAL_EINVAL,"Error in IMRPhenomX_SetWaveformVariables: IMRPhenomXRingdownAmpVersion not recognized. Recommended flag is 103.\n");
		}
	}

	/* Rescale to mass in solar masses */
	REAL8 m1_In      = m1_SI / LAL_MSUN_SI; // Mass 1 in solar masses
	REAL8 m2_In      = m2_SI / LAL_MSUN_SI; // Mass 2 in solar masses

	REAL8 m1, m2, chi1L, chi2L;

	/* Check if m1 >= m2, if not then swap masses/spins */
	if(m1_In >= m2_In)
	{
		chi1L = chi1L_In;
		chi2L = chi2L_In;
		m1    = m1_In;
		m2    = m2_In;
	}
	else
	{
		XLAL_PRINT_WARNING("Warning: m1 < m2, swapping the masses and spins.\n");
		chi1L = chi2L_In;
		chi2L = chi1L_In;
		m1    = m2_In;
		m2    = m1_In;
	}

	/* Check that physical spins have been called. Perform internal nudge to check for numerical round-off issues. */
	if(chi1L > 1.0)
	{
		IMRPhenomX_InternalNudge(chi1L,1.0,1e-6);
	}
	if(chi2L > 1.0)
	{
		IMRPhenomX_InternalNudge(chi2L,1.0,1e-6);
	}
	if(chi1L < -1.0)
	{
		IMRPhenomX_InternalNudge(chi1L,-1.0,1e-6);
	}
	if(chi2L < -1.0)
	{
		IMRPhenomX_InternalNudge(chi2L,-1.0,1e-6);
	}
	/* If spins are still unphysical after checking for small round-off errors, fail. */
	if( chi1L > 1.0 || chi1L < -1.0 || chi2L > 1.0 || chi2L < -1.0)
	{
		//XLALPrintError("Unphyiscal spins: must be in the range [-1,1].\n");
		XLAL_ERROR(XLAL_EDOM, "Unphysical spins requested: must obey the Kerr bound [-1,1].\n");
	}

	// Symmetric mass ratio
	REAL8 delta = fabs((m1 - m2) / (m1+m2));
	REAL8 eta   = fabs(0.25 * (1.0 - delta*delta) ); // use fabs to prevent negative sign due to roundoff
349
	REAL8 q   = ((m1 > m2) ? (m1 / m2) : (m2 / m1));
350

351
	/* If eta > 0.25, e.g. roundoff, then set to 0.25 */
352
353
	if(eta > 0.25)
	{
354
		eta = 0.25;
355
356
357
358
359
	}
	if(eta > 0.25 || eta < 0.0)
	{
		XLAL_ERROR(XLAL_EDOM,"Unphysical eta: must be between 0 and 0.25\n");
	}
360
	if(eta == 0.25) q = 1.;
361
362
363
364
365
366
367
368
369
370
	//if(eta < MAX_ALLOWED_ETA)
	//{
	//	XLAL_PRINT_WARNING("Warning: The model is not calibrated against numerical-relativity mass-ratios greater than 18.\n");
	//}

	/* Check that masses are correctly ordered. Swap if necessary.
	return_code = IMRPhenomX_CheckMassOrdering(m1,m2,chi1L,chi2L);
	XLAL_CHECK(XLAL_SUCCESS == return_code, XLAL_EFUNC, "Failure: IMRPhenomX_CheckMassOrdering");
	*/

371

372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397

	/* Masses definitions. Note that m1 and m2 are the component masses in solar masses. */
	wf->m1_SI     = m1 * LAL_MSUN_SI;																						// Mass 1 (larger) in SI units
	wf->m2_SI     = m2 * LAL_MSUN_SI;																						// Mass 2 (smaller) in SI units
	wf->q         = q;																								    			// Mass ratio >= 1
	wf->eta       = eta;                                                        // Symmetric mass ratio
	wf->Mtot_SI   = wf->m1_SI + wf->m2_SI;                                      // Total mass in SI units
	wf->Mtot      = m1 + m2;																										// Total mass in solar masss
	wf->m1        = m1 / (wf->Mtot);                         						        // Mass 1 (larger), dimensionless (i.e. m1 \in [0,1])
	wf->m2        = m2 / (wf->Mtot);						                                // Mass 2 (smaller), dimensionless
	wf->M_sec     = wf->Mtot * LAL_MTSUN_SI;                                    // Conversion factor Hz to geometric frequency, total mass in seconds
	wf->delta     = delta;                                                      // PN symmetry coefficient

	wf->eta2 = eta * eta;
	wf->eta3 = eta * wf->eta2;

	if(wf->debug)
	{
		printf("m1 (Msun) = %.6f\n",m1);
		printf("m2 (Msun) = %.6f\n",m2);
		printf("m1_SI     = %.6f\n",wf->m1_SI);
		printf("m2_SI     = %.6f\n",wf->m2_SI);
		printf("m1        = %.6f\n",wf->m1);
		printf("m2        = %.6f\n",wf->m2);
	}

398
	if(wf->q > 1000.0)
399
	{
400
		XLAL_PRINT_WARNING("Warning: The model is not supported for mass ratios > 1000.\n");
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
	}

	/* Spins */
	wf->chi1L     = chi1L;
	wf->chi2L     = chi2L;

	wf->chi1L2L   = chi1L * chi2L;

	/* Useful powers of spin */
	wf->chi1L2    = chi1L*chi1L;
	wf->chi1L3    = chi1L*chi1L*chi1L;

	wf->chi2L2    = chi2L*chi2L;
	wf->chi2L3    = chi2L*chi2L*chi2L;

	/* Spin parameterisations */
	wf->chiEff    = XLALSimIMRPhenomXchiEff(eta,chi1L,chi2L);
	wf->chiPNHat  = XLALSimIMRPhenomXchiPNHat(eta,chi1L,chi2L);
	wf->STotR     = XLALSimIMRPhenomXSTotR(eta,chi1L,chi2L);
	wf->dchi      = XLALSimIMRPhenomXdchi(chi1L,chi2L);

	wf->SigmaL    = (wf->chi2L * wf->m2) - (wf->chi1L * wf->m1); 										// SigmaL = (M/m2)*(S2.L) - (M/m2)*(S1.L)
	wf->SL        = wf->chi1L * (wf->m1 * wf->m1) + wf->chi2L * (wf->m2 * wf->m2);  // SL = S1.L + S2.L

	if(wf->debug)
	{
427
428
429
430
431
432
433
		printf("eta     : %.16e\n",wf->eta);
		printf("q       : %.16e\n",wf->q);
		printf("chi1L   : %.16e\n",wf->chi1L);
		printf("chi2L   : %.16e\n",wf->chi2L);
		printf("chi_eff : %.16e\n",wf->chiEff);
		printf("chi_hat : %.16e\n",wf->chiPNHat);
		printf("STotR   : %.16e\n",wf->STotR);
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
	}

	/* If no reference frequency is passed, set it to the starting GW frequency */
	wf->fRef      = (fRef == 0.0) ? f_min : fRef;
	wf->phiRef_In = phi0;
	wf->phi0      = phi0;							// Orbital phase at reference frequency (as passed from lalsimulation)
	wf->beta      = LAL_PI*0.5 - phi0; 				// Azimuthal angle of binary at reference frequency
	wf->phifRef   = 0.0;							// This is calculated later

	/* Geometric reference frequency */
	wf->MfRef     = XLALSimIMRPhenomXUtilsHztoMf(wf->fRef,wf->Mtot);
	wf->piM       = LAL_PI * wf->M_sec;
	wf->v_ref     = cbrt(wf->piM * wf->fRef);

	wf->deltaF    = deltaF;

	/* Define the default end of the waveform as: 0.3 Mf. This value is chosen such that the 44 mode also shows the ringdown part. */
	wf->fCutDef   = 0.3;
	// If chieff is very high the ringdown of the 44 is almost cut it out when using 0.3, so we increase a little bit the cut of freq up 0.33.
	if(wf->chiEff > 0.99) wf->fCutDef = 0.33;

	/* Minimum and maximum frequency */
	wf->fMin      = f_min;
	wf->fMax      = f_max;

	/* Convert fCut to physical cut-off frequency */
	wf->fCut      = wf->fCutDef / wf->M_sec;

	/* Sanity check that minimum start frequency is less than cut-off frequency */
	if (wf->fCut <= wf->fMin)
	{
		XLALPrintError("(fCut = %g Hz) <= f_min = %g\n", wf->fCut, wf->fMin);
	}

	if(wf->debug)
	{
		printf("fRef : %.6f\n",wf->fRef);
		printf("phi0 : %.6f\n",wf->phi0);
		printf("fCut : %.6f\n",wf->fCut);
		printf("fMin : %.6f\n",wf->fMin);
		printf("fMax : %.6f\n",wf->fMax);
	}

	/* By default f_max_prime is f_max. If fCut < fMax, then use fCut, i.e. waveform up to fMax will be zeros */
	wf->f_max_prime   = wf->fMax;
	wf->f_max_prime   = wf->fMax ? wf->fMax : wf->fCut;
	wf->f_max_prime   = (wf->f_max_prime > wf->fCut) ? wf->fCut : wf->f_max_prime;

	if (wf->f_max_prime <= wf->fMin)
	{
		XLALPrintError("f_max <= f_min\n");
	}

	if(wf->debug)
	{
		printf("fMin        = %.6f\n",wf->fMin);
		printf("fMax        = %.6f\n",wf->fMax);
		printf("f_max_prime = %.6f\n",wf->f_max_prime);
	}

	/* Final Mass and Spin */
	wf->Mfinal    = XLALSimIMRPhenomXFinalMass2017(wf->eta,wf->chi1L,wf->chi2L);
	wf->afinal    = XLALSimIMRPhenomXFinalSpin2017(wf->eta,wf->chi1L,wf->chi2L);

	/* Ringdown and damping frequency of final BH */
#if QNMfits == 1
	wf->fRING     = evaluate_QNMfit_fring22(wf->afinal) / (wf->Mfinal);
	wf->fDAMP     = evaluate_QNMfit_fdamp22(wf->afinal) / (wf->Mfinal);
#else
	wf->fRING     = interpolateQNMData_fring_22(wf->afinal) / (wf->Mfinal);
	wf->fDAMP     = interpolateQNMData_fdamp_22(wf->afinal) / (wf->Mfinal);
#endif


	if(wf->debug)
	{
		printf("Mf  = %.6f\n",wf->Mfinal);
		printf("af  = %.6f\n",wf->afinal);
		printf("frd = %.6f\n",wf->fRING);
		printf("fda = %.6f\n",wf->fDAMP);
	}

	if(wf->Mfinal > 1.0)
	{
		XLAL_ERROR(XLAL_EDOM, "IMRPhenomX_FinalMass2018: Final mass > 1.0 not physical.");
	}
	if(fabs(wf->afinal) > 1.0)
	{
		XLAL_ERROR(XLAL_EDOM, "IMRPhenomX_FinalSpin2018: Final spin > 1.0 is not physical.");
	}

	/* Fit to the hybrid minimum energy circular orbit (MECO), Cabero et al, Phys.Rev. D95 (2017) */
	wf->fMECO       = XLALSimIMRPhenomXfMECO(wf->eta,wf->chi1L,wf->chi2L);

	/* Innermost stable circular orbit (ISCO), e.g. Ori et al, Phys.Rev. D62 (2000) 124022 */
	wf->fISCO       = XLALSimIMRPhenomXfISCO(wf->afinal);

	if(wf->debug)
	{
		printf("fMECO = %.6f\n",wf->fMECO);
		printf("fISCO = %.6f\n",wf->fISCO);
	}

	if(wf->fMECO > wf->fISCO)
	{
		/* If MECO > fISCO, throw an error - this may be the case for very corner cases in the parameter space (e.g. q ~ 1000, chi ~ 0.999) */
		XLAL_ERROR(XLAL_EDOM, "Error: f_MECO cannot be greater than f_ISCO.");
	}

	/* Distance and inclination */
	wf->distance    = distance;
	wf->inclination = inclination;

	/* Amplitude normalization */
	wf->amp0        = wf->Mtot * LAL_MRSUN_SI * wf->Mtot * LAL_MTSUN_SI / wf->distance;
	wf->ampNorm     = sqrt(2.0/3.0) * sqrt(wf->eta) * powers_of_lalpi.m_one_sixth;

	if(wf->debug)
	{
		printf("\n\neta     = %.6f\n",wf->eta);
		printf("\nampNorm = %e\n",wf->ampNorm);
		printf("amp0 : %e",wf->amp0);
	}

	/* Phase normalization */
	wf->dphase0     = 5.0 / (128.0 * pow(LAL_PI,5.0/3.0));

	if(wf->debug)
	{
		printf("\n\n **** Sanity checks complete. Waveform struct has been initialized. **** \n\n");
	}

	return XLAL_SUCCESS;
}

int IMRPhenomXGetAmplitudeCoefficients(
	IMRPhenomXWaveformStruct *pWF,
	IMRPhenomXAmpCoefficients *pAmp
)
{
	INT4 debug = PHENOMXDEBUG;

	REAL8 V1, V2, V3, V4;
	REAL8 F1, F2, F3, F4;
	REAL8 d1, d4;

	/* Declare local spin variables for convenience */
	REAL8 chi1L         = pWF->chi1L;
	REAL8 chi2L         = pWF->chi2L;
	REAL8 chi1L2        = pWF->chi1L2;
	REAL8 chi1L3        = pWF->chi1L3;
	REAL8 chi2L2        = pWF->chi2L2;

	// Local PN symmetry parameter
	REAL8 delta         = pWF->delta;

	// Local powers of the symmetric mass ratio
	REAL8 eta           = pWF->eta;
	REAL8 eta2          = pWF->eta2;
	REAL8 eta3          = pWF->eta3;

	if(debug)
	{
		printf("chi1L  = %.6f\n",chi1L);
		printf("chi1L2 = %.6f\n",chi1L2);
		printf("chi1L3 = %.6f\n",chi1L3);
		printf("chi2L2 = %.6f\n",chi2L2);
		printf("delta  = %.6f\n",delta);
		printf("eta2   = %.6f\n",eta2);
		printf("eta3   = %.6f\n",eta3);
	}

606
	// Phenomenological ringdown coefficients: note that \gamma2 = \lambda in arXiv:2001.11412 and \gamma3 = \sigma in arXiv:2001.11412
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
	pAmp->gamma2        = IMRPhenomX_Ringdown_Amp_22_gamma2(pWF->eta,pWF->STotR,pWF->dchi,pWF->delta,pWF->IMRPhenomXRingdownAmpVersion);
	pAmp->gamma3        = IMRPhenomX_Ringdown_Amp_22_gamma3(pWF->eta,pWF->STotR,pWF->dchi,pWF->delta,pWF->IMRPhenomXRingdownAmpVersion);

	/* Get peak ringdown frequency, Eq. 5.14 in arXiv:2001.11412: Abs[fring + fdamp * gamma3 * (Sqrt[1 - gamma2^2] - 1)/gamma2 ] */
	pAmp->fAmpRDMin     = IMRPhenomX_Ringdown_Amp_22_PeakFrequency(pAmp->gamma2,pAmp->gamma3,pWF->fRING,pWF->fDAMP,pWF->IMRPhenomXRingdownAmpVersion);

	// Value of v1RD at fAmpRDMin is used to calcualte gamma1 \propto the amplitude of the deformed Lorentzian
	pAmp->v1RD          = IMRPhenomX_Ringdown_Amp_22_v1(pWF->eta,pWF->STotR,pWF->dchi,pWF->delta,pWF->IMRPhenomXRingdownAmpVersion);
	F1                  = pAmp->fAmpRDMin;

	// Solving linear equation: ansatzRD(F1) == v1
	pAmp->gamma1 = ( pAmp->v1RD / (pWF->fDAMP * pAmp->gamma3) ) * (F1*F1 - 2.0*F1*pWF->fRING + pWF->fRING*pWF->fRING + pWF->fDAMP*pWF->fDAMP*pAmp->gamma3*pAmp->gamma3)
													* exp( ( (F1 - pWF->fRING) * pAmp->gamma2 ) / (pWF->fDAMP * pAmp->gamma3) );

	/* Pre-cache these constants for the ringdown amplitude */
	pAmp->gammaR   = pAmp->gamma2 / (pWF->fDAMP * pAmp->gamma3);
	pAmp->gammaD2  = (pAmp->gamma3 * pWF->fDAMP) * (pAmp->gamma3 * pWF->fDAMP);
	pAmp->gammaD13 = pWF->fDAMP * pAmp->gamma1 * pAmp->gamma3;

	if(debug)
	{
		printf("V1     = %.6f\n",pAmp->v1RD);
		printf("gamma1 = %.6f\n",pAmp->gamma1);
		printf("gamma2 = %.6f\n",pAmp->gamma2);
		printf("gamma3 = %.6f\n",pAmp->gamma3);
		printf("fmax   = %.6f\n",pAmp->fAmpRDMin);
	}

	/* Inspiral-Intermediate transition frequencies, see Eq. 5.7 in arXiv:2001.11412 */
	pAmp->fAmpInsMin    = 0.0026;
	pAmp->fAmpInsMax    = pWF->fMECO + 0.25*(pWF->fISCO - pWF->fMECO);

	/* Inspiral-Intermediate matching frequency, Eq. 5.16 in arXiv:2001.11412 */
	pAmp->fAmpMatchIN   = pAmp->fAmpInsMax;

	/* Intermediate-Ringdown matching frequency */
	//pAmp->fAmpMatchIM   = pAmp->fAmpRDMin;

	/* End of intermerdiate region, Eq. 5.16 in arXiv:2001.11412 */
	pAmp->fAmpIntMax    = pAmp->fAmpRDMin;

	/* TaylorF2 PN Amplitude Coefficients from usual PN re-expansion of Hlm's */
	pAmp->pnInitial     = 1.0;
	pAmp->pnOneThird    = 0.0;
	pAmp->pnTwoThirds   = ( (-969 + 1804*eta)/672. ) * powers_of_lalpi.two_thirds;
	pAmp->pnThreeThirds = ( (81*(chi1L + chi2L) + 81*chi1L*delta - 81*chi2L*delta - 44*(chi1L + chi2L)*eta)/48. ) * powers_of_lalpi.itself;
	pAmp->pnFourThirds  = ( (-27312085 - 10287648*chi1L2*(1 + delta)
													+ 24*(428652*chi2L2*(-1 + delta) + (-1975055 + 10584*(81*chi1L2 - 94*chi1L*chi2L + 81*chi2L2))*eta
													+ 1473794*eta2))/8.128512e6 ) * powers_of_lalpi.one_third * powers_of_lalpi.itself;
	pAmp->pnFiveThirds  = ( (-6048*chi1L3*(-1 - delta + (3 + delta)*eta) + chi2L*(-((287213 + 6048*chi2L2)*(-1 + delta))
													+ 4*(-93414 + 1512*chi2L2*(-3 + delta) + 2083*delta)*eta - 35632*eta2)
													+ chi1L*(287213*(1 + delta) - 4*eta*(93414 + 2083*delta + 8908*eta))
													+ 42840*(-1 + 4*eta)*LAL_PI)/32256. ) * powers_of_lalpi.five_thirds;
	pAmp->pnSixThirds   = ( (-1242641879927 + 12.0*(28.0*(-3248849057.0
													+ 11088*(163199*chi1L2 - 266498*chi1L*chi2L + 163199*chi2L2))*eta2
													+ 27026893936*eta3 - 116424*(147117*(-(chi2L2*(-1.0 + delta)) + chi1L2*(1.0 + delta))
													+ 60928*(chi1L + chi2L + chi1L*delta - chi2L*delta)*LAL_PI)
													+ eta*(545384828789.0 - 77616*(638642*chi1L*chi2L + chi1L2*(-158633 + 282718*delta)
													- chi2L2*(158633.0 + 282718.0*delta) - 107520.0*(chi1L + chi2L)*LAL_PI
													+ 275520*powers_of_lalpi.two))))/6.0085960704e10 ) * powers_of_lalpi.two;


	/* These are the same order as the canonical pseudo-PN terms but we can add higher order PN information as and when available here */
	pAmp->pnSevenThirds = 0.0;
	pAmp->pnEightThirds = 0.0;
	pAmp->pnNineThirds  = 0.0;

	/* Generate Pseudo-PN Coefficients for Amplitude. */
	switch( pWF->IMRPhenomXInspiralAmpVersion)
	{
		case 103:
		{
			pAmp->CollocationValuesAmpIns[0] = 0.0; // This is not currently used but may be invoked in future recalibration. Leave here, its harmless.
			pAmp->CollocationValuesAmpIns[1] = IMRPhenomX_Inspiral_Amp_22_v2(pWF->eta,pWF->chiPNHat,pWF->dchi,pWF->delta,pWF->IMRPhenomXInspiralAmpVersion);
			pAmp->CollocationValuesAmpIns[2] = IMRPhenomX_Inspiral_Amp_22_v3(pWF->eta,pWF->chiPNHat,pWF->dchi,pWF->delta,pWF->IMRPhenomXInspiralAmpVersion);
			pAmp->CollocationValuesAmpIns[3] = IMRPhenomX_Inspiral_Amp_22_v4(pWF->eta,pWF->chiPNHat,pWF->dchi,pWF->delta,pWF->IMRPhenomXInspiralAmpVersion);

			pAmp->CollocationPointsAmpIns[0] = 0.25 * pAmp->fAmpMatchIN;
			pAmp->CollocationPointsAmpIns[1] = 0.50 * pAmp->fAmpMatchIN;
			pAmp->CollocationPointsAmpIns[2] = 0.75 * pAmp->fAmpMatchIN;
			pAmp->CollocationPointsAmpIns[3] = 1.00 * pAmp->fAmpMatchIN;
			break;
		}
		default:
		{
			XLAL_ERROR(XLAL_EINVAL, "Error in IMRPhenomXGetAmplitudeCoefficients: IMRPhenomXInspiralAmpVersion is not valid.\n");
		}
	}

	// Set collocation points and values
	V1 = pAmp->CollocationValuesAmpIns[1];
	V2 = pAmp->CollocationValuesAmpIns[2];
	V3 = pAmp->CollocationValuesAmpIns[3];
	F1 = pAmp->CollocationPointsAmpIns[1];
	F2 = pAmp->CollocationPointsAmpIns[2];
	F3 = pAmp->CollocationPointsAmpIns[3];

	if(debug)
	{
		printf("\nAmplitude pseudo PN coeffcients\n");
		printf("fAmpMatchIN = %.6f\n",pAmp->fAmpMatchIN);
		printf("V1   = %.6f\n",V1);
		printf("V2   = %.6f\n",V2);
		printf("V3   = %.6f\n",V3);
		printf("F1   = %e\n",F1);
		printf("F2   = %e\n",F2);
		printf("F3   = %e\n",F3);
	}

	/* Using the collocation points and their values, we solve a linear system of equations to recover the pseudo-PN coefficients */
	pAmp->rho1 = IMRPhenomX_Inspiral_Amp_22_rho1(V1,V2,V3,F1,F2,F3,pWF->IMRPhenomXInspiralAmpVersion);
	pAmp->rho2 = IMRPhenomX_Inspiral_Amp_22_rho2(V1,V2,V3,F1,F2,F3,pWF->IMRPhenomXInspiralAmpVersion);
	pAmp->rho3 = IMRPhenomX_Inspiral_Amp_22_rho3(V1,V2,V3,F1,F2,F3,pWF->IMRPhenomXInspiralAmpVersion);

	/* Set TaylorF2 pseudo-PN coefficients */
	pAmp->pnSevenThirds = pAmp->rho1;
	pAmp->pnEightThirds = pAmp->rho2;
	pAmp->pnNineThirds  = pAmp->rho3;

	if(debug)
	{
		printf("\nTaylorF2 PN Amplitude Coefficients\n");
		printf("pnTwoThirds   = %.6f\n",pAmp->pnTwoThirds);
		printf("pnThreeThirds = %.6f\n",pAmp->pnThreeThirds);
		printf("pnFourThirds  = %.6f\n",pAmp->pnFourThirds);
		printf("pnFiveThirds  = %.6f\n",pAmp->pnFiveThirds);
		printf("pnSixThirds   = %.6f\n",pAmp->pnSixThirds);
		printf("\n");
		printf("powers_of_lalpi.itself = %.6f\n",powers_of_lalpi.itself);
		printf("powers_of_lalpi.four_thirds = %.6f\n",powers_of_lalpi.four_thirds);
		printf("powers_of_lalpi.five_thirds = %.6f\n",powers_of_lalpi.five_thirds);
		printf("\n");
		printf("Pseudo-PN Amplitude Coefficients (Agrees with MMA).\n");
		printf("alpha1 = %.6f\n",pAmp->pnSevenThirds);
		printf("alpha2 = %.6f\n",pAmp->pnEightThirds);
		printf("alpha3  = %.6f\n",pAmp->pnNineThirds);
	}

	/* Now we can reconstruct the intermediate region, which depends on the derivative of the inspiral and ringdown fits at the boundaries */
	F1     = pAmp->fAmpMatchIN;
	F4     = pAmp->fAmpRDMin;

	// Initialise useful powers for the boundary points
	IMRPhenomX_UsefulPowers powers_of_F1;
	IMRPhenomX_Initialize_Powers(&powers_of_F1,F1);

	IMRPhenomX_UsefulPowers powers_of_F4;
	IMRPhenomX_Initialize_Powers(&powers_of_F4,F4);

	/* Calculate derivative of amplitude on boundary points */
	d1     = IMRPhenomX_Inspiral_Amp_22_DAnsatz(F1,pWF,pAmp);
	d4     = IMRPhenomX_Ringdown_Amp_22_DAnsatz(F4,pWF,pAmp);

	double inspF1 = IMRPhenomX_Inspiral_Amp_22_Ansatz(F1,&powers_of_F1,pWF,pAmp);
	double rdF4   = IMRPhenomX_Ringdown_Amp_22_Ansatz(F4,pWF,pAmp);

763
	/*
764
		Use d1 and d4 calculated above to get the derivative of the amplitude on the boundaries:
765

766
		D[ f^{7/6} Ah^{-1}[f] , f] = ( (7/6) * f^{1/6} / Ah[f] ) - ( f^{7/6} Ah'[f] ) / (Ah[f]^2)
767

768
		where d1 = Ah'[F1], inspF1 = Ah[F1] etc.
769

770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
	*/
	d1     = ((7.0/6.0) * powers_of_F1.one_sixth / inspF1) - ( powers_of_F1.seven_sixths * d1 / (inspF1*inspF1) );
	d4     = ((7.0/6.0) * powers_of_F4.one_sixth / rdF4)   - ( powers_of_F4.seven_sixths * d4 / (rdF4*rdF4) );

	if(debug)
	{
		printf("d1 = %.6f\n",d1);
		printf("d4 = %.6f\n",d4);
	}

	/* Pass inverse of points to reconstruction function... */
	switch ( pWF->IMRPhenomXIntermediateAmpVersion )
	{
		case 1043:
		{
			// Use a 5th order polynomial in intermediate - great agreement to calibration points but poor extrapolation
			F2     = F1 + (1.0/3.0) * (F4-F1);
			F3     = F1 + (2.0/3.0) * (F4-F1);

			V1     = powers_of_F1.m_seven_sixths * inspF1;
			V2     = IMRPhenomX_Intermediate_Amp_22_v2(pWF->eta,pWF->STotR,pWF->dchi,pWF->delta,pWF->IMRPhenomXIntermediateAmpVersion);
			V3     = IMRPhenomX_Intermediate_Amp_22_v3(pWF->eta,pWF->STotR,pWF->dchi,pWF->delta,pWF->IMRPhenomXIntermediateAmpVersion);
			V4     = powers_of_F4.m_seven_sixths * rdF4;

			V1 		 = 1.0 / V1;
			V2 		 = 1.0 / V2;
			V3 		 = 1.0 / V3;
			V4 		 = 1.0 / V4;
			break;
		}
		case 104:
		{
			// Use a 4th order polynomial in intermediate - good extrapolation, recommended default fit
			F2     = F1 + (1.0/2.0) * (F4-F1);
			F3     = 0.0;

			V1     = powers_of_F1.m_seven_sixths * inspF1;
			V2     = IMRPhenomX_Intermediate_Amp_22_vA(pWF->eta,pWF->STotR,pWF->dchi,pWF->delta,pWF->IMRPhenomXIntermediateAmpVersion);
			V4     = powers_of_F4.m_seven_sixths * rdF4;

			V1 		 = 1.0 / V1;
			V2 		 = 1.0 / V2;
			V3 		 = 0.0;
			V4 		 = 1.0 / V4;

			break;
		}
		case 105:
		{
			// Use a 5th order polynomial in intermediate - great agreement to calibration points but poor extrapolation
			F2     = F1 + (1.0/3.0) * (F4-F1);
			F3     = F1 + (2.0/3.0) * (F4-F1);

			V1     = powers_of_F1.m_seven_sixths * inspF1;
			V2     = IMRPhenomX_Intermediate_Amp_22_v2(pWF->eta,pWF->STotR,pWF->dchi,pWF->delta,pWF->IMRPhenomXIntermediateAmpVersion);
			V3     = IMRPhenomX_Intermediate_Amp_22_v3(pWF->eta,pWF->STotR,pWF->dchi,pWF->delta,pWF->IMRPhenomXIntermediateAmpVersion);
			V4     = powers_of_F4.m_seven_sixths * rdF4;

			V1 		 = 1.0 / V1;
			V2 		 = 1.0 / V2;
			V3 		 = 1.0 / V3;
			V4 		 = 1.0 / V4;
			break;
		}
		default:
		{
			XLAL_ERROR(XLAL_EINVAL, "Error: IMRPhenomXIntermediateAmpVersion is not valid.\n");
		}
	}

	if(debug)
	{
		printf("\nIntermediate Region: \n");
		printf("F1 = %.6f\n",F1);
		printf("F2 = %.6f\n",F2);
		printf("F3 = %.6f\n",F3);
		printf("F4 = %.6f\n",F4);
		printf("\n");
		printf("Insp@F1 = %.6f\n",IMRPhenomX_Inspiral_Amp_22_Ansatz(F1,&powers_of_F1,pWF,pAmp));
		printf("d1 = %.6f\n",d1);
		printf("d4 = %.6f\n",d4);
		printf("V1 = %.6f\n",V1);
		printf("V2 = %.6f\n",V2);
		printf("V3 = %.6f\n",V3);
		printf("V4 = %.6f\n",V4);
	}

	/*
	Reconstruct the phenomenological coefficients for the intermediate ansatz:
	*/
	pAmp->delta0 = IMRPhenomX_Intermediate_Amp_22_delta0(d1,d4,V1,V2,V3,V4,F1,F2,F3,F4,pWF->IMRPhenomXIntermediateAmpVersion);
	pAmp->delta1 = IMRPhenomX_Intermediate_Amp_22_delta1(d1,d4,V1,V2,V3,V4,F1,F2,F3,F4,pWF->IMRPhenomXIntermediateAmpVersion);
	pAmp->delta2 = IMRPhenomX_Intermediate_Amp_22_delta2(d1,d4,V1,V2,V3,V4,F1,F2,F3,F4,pWF->IMRPhenomXIntermediateAmpVersion);
	pAmp->delta3 = IMRPhenomX_Intermediate_Amp_22_delta3(d1,d4,V1,V2,V3,V4,F1,F2,F3,F4,pWF->IMRPhenomXIntermediateAmpVersion);
	pAmp->delta4 = IMRPhenomX_Intermediate_Amp_22_delta4(d1,d4,V1,V2,V3,V4,F1,F2,F3,F4,pWF->IMRPhenomXIntermediateAmpVersion);
	pAmp->delta5 = IMRPhenomX_Intermediate_Amp_22_delta5(d1,d4,V1,V2,V3,V4,F1,F2,F3,F4,pWF->IMRPhenomXIntermediateAmpVersion);

	if(debug)
	{
  	    printf("delta0 = %.6f\n",pAmp->delta0);
		printf("delta1 = %.6f\n",pAmp->delta1);
		printf("delta2 = %.6f\n",pAmp->delta2);
		printf("delta3 = %.6f\n",pAmp->delta3);
		printf("delta4 = %.6f\n",pAmp->delta4);
		printf("delta5 = %.6f\n",pAmp->delta5);
	}

	return XLAL_SUCCESS;
}

/*
 * Function to populate the IMRPhenomXPhaseCoefficients struct:
*/
int IMRPhenomXGetPhaseCoefficients(
	IMRPhenomXWaveformStruct *pWF,
	IMRPhenomXPhaseCoefficients *pPhase
)
{
	//IMRPhenomXPhaseCoefficients *pPhase = (IMRPhenomXPhaseCoefficients *) XLALMalloc(sizeof(IMRPhenomXPhaseCoefficients));

	const INT4 debug = PHENOMXDEBUG; // pWF->debug;
	/* GSL objects for solving system of equations via LU decomposition */
	gsl_vector *b, *x;
	gsl_matrix *A;
	gsl_permutation *p;
	int s;

	REAL8 deltax;
	REAL8 xmin;
	REAL8 fi;
	//REAL8 gpoints4[4]     = {1.0, 3.0/4, 1.0/4, 0.0};
	//REAL8 gpoints5[5]     = {1.0, 1.0/2 + 1.0/(2.0*sqrt(2.0)), 1.0/2, 1.0/2 - 1.0/(2*sqrt(2.0)), 0.};

	REAL8 gpoints4[4]     = {0.0, 1.0/4.0, 3.0/4.0, 1.0};
	REAL8 gpoints5[5]     = {0.0, 1.0/2 - 1.0/(2*sqrt(2.0)), 1.0/2.0, 1.0/2 + 1.0/(2.0*sqrt(2.0)), 1.0};

	//pPhase->fPhaseMatchIN = pWF->fMECO;
	//pPhase->fPhaseMatchIM = 0.6 * (0.5 * pWF->fRING + pWF->fISCO);

	// Matching regions
910

911
912
	/* This is Eq. 5.11 in the paper */
	double fIMmatch = 0.6 * (0.5 * pWF->fRING + pWF->fISCO);
913

914
915
	/* This is the MECO frequency */
	double fINmatch = pWF->fMECO;
916

917
918
919
920
	/* This is Eq. 5.10 in the paper */
	double deltaf   = (fIMmatch - fINmatch) * 0.03;

	// Transition frequency is just below the MECO frequency and just above the RD fitting region
921

922
923
924
925
926
927
	/* These are defined in Eq. 7.7 and the text just below, f_H = fPhaseMatchIM and f_L = fPhaseMatchIN */
	pPhase->fPhaseMatchIN  = fINmatch - 1.0*deltaf;
	pPhase->fPhaseMatchIM  = fIMmatch + 0.5*deltaf;

	/* Defined in Eq. 7.4, this is f_L */
	pPhase->fPhaseInsMin  = 0.0026;
928

929
930
931
932
933
	/* Defined in Eq. 7.4, this is f_H */
	pPhase->fPhaseInsMax  = 1.020 * pWF->fMECO;

	/* Defined in Eq. 7.12, this is f_L */
	pPhase->fPhaseRDMin   = fIMmatch;
934

935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
	/* Defined in Eq. 7.12, this is f_L */
	pPhase->fPhaseRDMax   = pWF->fRING + 1.25*pWF->fDAMP;

	pPhase->phiNorm    		= -(3.0 * powers_of_lalpi.m_five_thirds) / (128.0);

	/* For convenience, define some variables here */
	REAL8 chi1L           = pWF->chi1L;
	REAL8 chi2L           = pWF->chi2L;

	REAL8 chi1L2L         = chi1L * chi2L;

	REAL8 chi1L2          = pWF->chi1L * pWF->chi1L;
	REAL8 chi1L3          = pWF->chi1L * chi1L2;

	REAL8 chi2L2          = pWF->chi2L * pWF->chi2L;
	REAL8 chi2L3          = pWF->chi2L * chi2L2;

	REAL8 eta             = pWF->eta;
	REAL8 eta2            = eta*eta;
	REAL8 eta3            = eta*eta2;

	REAL8 delta           = pWF->delta;

958
959
	LALDict *LALparams    = pWF->LALparams;

960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
	/* Pre-initialize all phenomenological coefficients */
	pPhase->a0 = 0.0;
	pPhase->a1 = 0.0;
	pPhase->a2 = 0.0;
	pPhase->a3 = 0.0;
	pPhase->a4 = 0.0;

	pPhase->b0 = 0.0;
	pPhase->b1 = 0.0;
	pPhase->b2 = 0.0;
	pPhase->b3 = 0.0;
	pPhase->b4 = 0.0;

	pPhase->c0 = 0.0;
	pPhase->c1 = 0.0;
	pPhase->c2 = 0.0;
	pPhase->c3 = 0.0;
	pPhase->c4 = 0.0;
	pPhase->cL = 0.0;

	pPhase->sigma0 = 0.0;
	pPhase->sigma1 = 0.0;
	pPhase->sigma2 = 0.0;
	pPhase->sigma3 = 0.0;
	pPhase->sigma4 = 0.0;
	pPhase->sigma5 = 0.0;
986

987
988
	/*
		The general strategy is to initialize a linear system of equations:
989

990
991
		A.x = b

992
		- A is a matrix with the coefficients of the ansatz evaluated at the collocation nodes.
993
		- b is a vector of the value of the collocation points
994
		- x is the solution vector (i.e. the coefficients) that we must solve for.
995

996
		We choose to do this using a standard LU decomposition.
997
998
999
1000
	*/

	/* Generate list of collocation points */
	/*
For faster browsing, not all history is shown. View entire blame