LALSimulation  5.4.0.1-fe68b98
LALSimIMRPhenomXHM_internals.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2019 Marta Colleoni, Cecilio Garcia Quiros
3  *
4  * This program is free software; you can redistribute it and/or modify
5  * it under the terms of the GNU General Public License as published by
6  * the Free Software Foundation; either version 2 of the License, or
7  * (at your option) any later version.
8  *
9  * This program is distributed in the hope that it will be useful,
10  * but WITHOUT ANY WARRANTY; without even the implied warranty of
11  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  * GNU General Public License for more details.
13  *
14  * You should have received a copy of the GNU General Public License
15  * along with with program; see the file COPYING. If not, write to the
16  * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
17  * MA 02110-1301 USA
18  */
19 //
20 // LALSimIMRPhenomXHM_internals.h
21 //
22 //
23 // Created by Marta on 06/02/2019.
24 //
25 
26 
27 #ifndef _LALSIM_IMR_PHENOMXHM_INTERNALS_H
28 #define _LALSIM_IMR_PHENOMXHM_INTERNALS_H
29 
30 
31 #ifdef __cplusplus
32 extern "C" {
33  #endif
34  #ifdef __GNUC__
35  #define UNUSED __attribute__ ((unused))
36  #else
37  #define UNUSED
38  #endif
39 
40 
41  #include <stdlib.h>
42  #include <stdio.h>
43  #include <math.h>
44  #include <complex.h>
45 
46  #include <lal/LALStdlib.h>
47  #include <lal/LALConstants.h>
48  #include <lal/Date.h>
49  #include <lal/FrequencySeries.h>
50  #include <lal/Units.h>
51 
54 
55  #include <lal/LALSimIMR.h>
56 
57  #define NMAX_INSPIRAL_COEFFICIENTS 13
58 
59  #define FALSE_ZERO 1.0e-15
60 
61  // You should not declare static functions here, since this file is included in other files apart form the source one.
62 
63  /*********** Useful Powers of pi **************/
65 
66  /**************** QNMs and mixing coefficients ************** */
67  void IMRPhenomXHM_Initialize_QNMs(QNMFits *qnmsFits);
69 
70  /****************** Initialization of higher-mode waveform struct ******************/
71  void IMRPhenomXHM_SetHMWaveformVariables(int ell, int emm, IMRPhenomXHMWaveformStruct *wf, IMRPhenomXWaveformStruct *wf22, QNMFits *qnms,LALDict *LALParams);
72 
73 
74  /*************** AMPLITUDE ****************/
75 
76  /***************** set pointers to par-space fits *******************/
78 
79  /***************** Cutting frequencies for amplitude ******************************/
82 
83  double RescaleFactor(IMRPhenomX_UsefulPowers *powers_of_Mf, IMRPhenomXHMAmpCoefficients *pAmp, UINT2 rescalefactor);
84 
85  /***************** Amplitude coefficients ****************/
89 
90  /**************** Amplitude reconstruction ******************/
91  // Functions that return suitable ansatz at a given frequency
95 
96 
97 
98  /***************** PHASE *******************/
99 
100  /***************** set pointers to par-space fits *******************/
103 
104  /***************** intermediate region reconstruction utilities *******************/
105  void EquidistantNodes(double nodes[], double fmin, double fmax, int npts);
107 
108  /***************** spheroidal->spherical harmonic conversion *******************/
111 
112  /***************** spheroidal-harmonic ringdown reconstruction ********************/
115 
116  /***************** phase coefficients *****************/
119 
120  /**************** Phase reconstruction ******************/
121  // functions that return suitable ansatz at a given frequency
125 
126  // functions that return the phase derivative at a given dimensionless frequency
129 
130 
132  UINT4 ell,
133  UINT4 emm,
135  );
136 
137  // Debugging function
138  int ParametersToFile(
139  IMRPhenomXWaveformStruct *pWF, /**< Wf structure for the 22 mode*/
140  IMRPhenomXHMWaveformStruct *pWFHM, /**< Wf structure for the lm mode*/
141  IMRPhenomXHMAmpCoefficients *pAmp, /**< Coefficients struct of the lm Amplitude */
142  UNUSED IMRPhenomXHMPhaseCoefficients *pPhase /**< Coefficients struct of the lm Phase */
143  );
144 
145 #ifdef __cplusplus
146 }
147 #endif
148 
149 
150 #endif /* LALSimIMRPhenomXHM_internals_h */
int IMRPhenomXHM_PN21AmpSign(double ff, IMRPhenomXWaveformStruct *wf22)
REAL8 IMRPhenomXHM_Amplitude_ModeMixing(IMRPhenomX_UsefulPowers *powers_of_Mf, IMRPhenomXHMAmpCoefficients *pAmp, IMRPhenomXHMPhaseCoefficients *pPhase, IMRPhenomXHMWaveformStruct *pWF, IMRPhenomXAmpCoefficients *pAmp22, IMRPhenomXPhaseCoefficients *pPhase22, IMRPhenomXWaveformStruct *pWF22)
void IMRPhenomXHM_SetHMWaveformVariables(int ell, int emm, IMRPhenomXHMWaveformStruct *wf, IMRPhenomXWaveformStruct *wf22, QNMFits *qnms, LALDict *LALParams)
void IMRPhenomXHM_Initialize_QNMs(QNMFits *qnmsFits)
double IMRPhenomXHM_Amplitude_fcutInsp(IMRPhenomXHMWaveformStruct *pWFHM, IMRPhenomXWaveformStruct *pWF22)
int ParametersToFile(IMRPhenomXWaveformStruct *pWF, IMRPhenomXHMWaveformStruct *pWFHM, IMRPhenomXHMAmpCoefficients *pAmp, UNUSED IMRPhenomXHMPhaseCoefficients *pPhase)
void IMRPhenomXHM_GetPNAmplitudeCoefficients(IMRPhenomXHMAmpCoefficients *pAmp, IMRPhenomXHMWaveformStruct *pWFHM, IMRPhenomXWaveformStruct *pWF22)
double IMRPhenomXHM_dPhase_noModeMixing(double f, IMRPhenomX_UsefulPowers *powers_of_f, IMRPhenomXHMPhaseCoefficients *pPhase, IMRPhenomXHMWaveformStruct *pWF, IMRPhenomXWaveformStruct *pWF22)
void IMRPhenomXHM_Ringdown_CollocPtsFreqs(IMRPhenomXHMPhaseCoefficients *pPhase, IMRPhenomXHMWaveformStruct *pWFHM, IMRPhenomXWaveformStruct *pWF22)
void Get21PNAmplitudeCoefficients(IMRPhenomXHMAmpCoefficients *pAmp, IMRPhenomXWaveformStruct *pWF22)
void IMRPhenomXHM_Initialize_MixingCoeffs(IMRPhenomXHMWaveformStruct *wf, IMRPhenomXWaveformStruct *wf22)
double IMRPhenomXHM_dPhase_ModeMixing(double f, IMRPhenomX_UsefulPowers *powers_of_f, IMRPhenomXHMAmpCoefficients *pAmp, IMRPhenomXHMPhaseCoefficients *pPhase, IMRPhenomXHMWaveformStruct *pWF, IMRPhenomXAmpCoefficients *pAmp22, IMRPhenomXPhaseCoefficients *pPhase22, IMRPhenomXWaveformStruct *pWF22)
REAL8 IMRPhenomXHM_Phase_ModeMixingRecycle(IMRPhenomX_UsefulPowers *powers_of_Mf, COMPLEX16 wf22, IMRPhenomXHMAmpCoefficients *pAmp, IMRPhenomXHMPhaseCoefficients *pPhase, IMRPhenomXHMWaveformStruct *pWF)
void IMRPhenomXHM_GetAmplitudeCoefficients(IMRPhenomXHMAmpCoefficients *pAmp, IMRPhenomXHMPhaseCoefficients *pPhase, IMRPhenomXAmpCoefficients *pAmp22, IMRPhenomXPhaseCoefficients *pPhase22, IMRPhenomXHMWaveformStruct *pWFHM, IMRPhenomXWaveformStruct *pWF22)
COMPLEX16 SpheroidalToSpherical(IMRPhenomX_UsefulPowers *powers_of_Mf, IMRPhenomXAmpCoefficients *pAmp22, IMRPhenomXPhaseCoefficients *pPhase22, IMRPhenomXHMAmpCoefficients *pAmplm, IMRPhenomXHMPhaseCoefficients *pPhaselm, IMRPhenomXHMWaveformStruct *pWFlm, IMRPhenomXWaveformStruct *pWF22)
REAL8 IMRPhenomXHM_Amplitude_ModeMixingRecycle(IMRPhenomX_UsefulPowers *powers_of_Mf, COMPLEX16 wf22, IMRPhenomXHMAmpCoefficients *pAmp, IMRPhenomXHMPhaseCoefficients *pPhase, IMRPhenomXHMWaveformStruct *pWF)
REAL8 IMRPhenomXHM_Amplitude_noModeMixing(IMRPhenomX_UsefulPowers *powers_of_Mf, IMRPhenomXHMAmpCoefficients *pAmp, IMRPhenomXHMWaveformStruct *pWF)
void IMRPhenomXHM_Intermediate_CollocPtsFreqs(IMRPhenomXHMPhaseCoefficients *pPhase, IMRPhenomXHMWaveformStruct *pWFHM, IMRPhenomXWaveformStruct *pWF22)
double GetfcutInsp(IMRPhenomXWaveformStruct *pWF22, IMRPhenomXHMWaveformStruct *pWFHM)
double IMRPhenomXHM_Amplitude_fcutRD(IMRPhenomXHMWaveformStruct *pWFHM, IMRPhenomXWaveformStruct *pWF22)
void GetSpheroidalCoefficients(IMRPhenomXHMPhaseCoefficients *pPhase, IMRPhenomXPhaseCoefficients *pPhase22, IMRPhenomXHMWaveformStruct *pWFHM, IMRPhenomXWaveformStruct *pWF22)
void IMRPhenomXHM_FillAmpFitsArray(IMRPhenomXHMAmpCoefficients *pAmp)
REAL8 IMRPhenomXHM_Phase_noModeMixing(IMRPhenomX_UsefulPowers *powers_of_Mf, IMRPhenomXHMPhaseCoefficients *pPhase, IMRPhenomXHMWaveformStruct *pWF, IMRPhenomXWaveformStruct *pWF22)
void IMRPhenomXHM_FillPhaseFitsArray(IMRPhenomXHMPhaseCoefficients *pPhase)
void EquidistantNodes(double nodes[], double fmin, double fmax, int npts)
REAL8 IMRPhenomXHM_Phase_ModeMixing(IMRPhenomX_UsefulPowers *powers_of_Mf, IMRPhenomXHMAmpCoefficients *pAmp, IMRPhenomXHMPhaseCoefficients *pPhase, IMRPhenomXHMWaveformStruct *pWF, IMRPhenomXAmpCoefficients *pAmp22, IMRPhenomXPhaseCoefficients *pPhase22, IMRPhenomXWaveformStruct *pWF22)
REAL8 IMRPhenomXHM_GenerateRingdownFrequency(UINT4 ell, UINT4 emm, IMRPhenomXWaveformStruct *wf22)
void IMRPhenomXHM_GetPhaseCoefficients(IMRPhenomXHMAmpCoefficients *pAmp, IMRPhenomXHMPhaseCoefficients *pPhase, IMRPhenomXAmpCoefficients *pAmp22, IMRPhenomXPhaseCoefficients *pPhase22, IMRPhenomXHMWaveformStruct *pWFHM, IMRPhenomXWaveformStruct *pWF22, LALDict *lalParams)
COMPLEX16 SpheroidalToSphericalRecycle(IMRPhenomX_UsefulPowers *powers_of_Mf, COMPLEX16 wf22, IMRPhenomXHMAmpCoefficients *pAmplm, IMRPhenomXHMPhaseCoefficients *pPhaselm, IMRPhenomXHMWaveformStruct *pWFlm)
IMRPhenomX_UsefulPowers powers_of_lalpiHM
double RescaleFactor(IMRPhenomX_UsefulPowers *powers_of_Mf, IMRPhenomXHMAmpCoefficients *pAmp, UINT2 rescalefactor)
double complex COMPLEX16
double REAL8
uint16_t UINT2
uint32_t UINT4