Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALSimulation 6.2.0.1-5e288d3
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
LALSimFindAttachTime.h
Go to the documentation of this file.
1/*
2 * Copyright (C) 2015 S. Babak
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#ifndef _LALSIMFINDATTACHTIME_H
21#define _LALSIMFINDATTACHTIME_H
22
23#include <math.h>
24
25#include <lal/LALGSL.h>
26
27#include <lal/LALDatatypes.h>
28#include <lal/TimeSeries.h>
29#include <lal/LALStdlib.h>
30#include <lal/LALConstants.h>
31#include <lal/TimeSeries.h>
32#include <lal/Units.h>
33
34
35#include <lal/LALSimIMR.h>
36#include <lal/Date.h>
37#include <lal/Units.h>
38#include <lal/SeqFactories.h>
39
40//#include "LALSimIMRSpinEOB.h"
41
42
43#include <gsl/gsl_linalg.h>
44#include <gsl/gsl_interp.h>
45#include <gsl/gsl_spline.h>
46#include <gsl/gsl_vector.h>
47#include <gsl/gsl_matrix.h>
48
49
50
51
52#if defined(__cplusplus)
53extern "C" {
54#elif 0
55} /* so that editors will match preceding brace */
56#endif
57
59 REAL8Array *dynamicsHi,
60 unsigned int numdynvars,
61 unsigned int retLenHi,
62 SpinEOBParams seobParams,
63 SpinEOBHCoeffs seobCoeffs,
64 REAL8 m1,
65 REAL8 m2,
66 REAL8Vector *radiusVec,
67 int *found,
68 REAL8* tMaxOmega,
69 INT4 use_optimized
70);
71
73 REAL8Vector *timeHi,
74 COMPLEX16Vector *hP22,
75 int *found);
76
77
79 REAL8Vector *timeHi,
80 COMPLEX16Vector *hP22,
81 REAL8Vector *radiusVec,
82 int *found,
83 REAL8* tMaxAmp
84);
85
87 REAL8Vector * signal1,
88 REAL8Vector * signal2,
89 REAL8* ratio,
90 const REAL8 tAtt,
91 const INT4 l,
92 const INT4 m,
93 const REAL8 dt,
94 const REAL8 mass1,
95 const REAL8 mass2,
96 const REAL8 spin1x,
97 const REAL8 spin1y,
98 const REAL8 spin1z,
99 const REAL8 spin2x,
100 const REAL8 spin2y,
101 const REAL8 spin2z,
102 REAL8Vector * timeVec,
103 REAL8Vector * matchrange,
105 const REAL8 JLN,
106 REAL8 * timediff
107);
108
110 REAL8Vector * signal1,
111 REAL8Vector * signal2,
114 REAL8* ratio22,
115 REAL8* ratio2m2,
116 REAL8* tAtt,
117 const REAL8 thr,
118 const REAL8 dt,
119 const REAL8 mass1,
120 const REAL8 mass2,
121 const REAL8 spin1x,
122 const REAL8 spin1y,
123 const REAL8 spin1z,
124 const REAL8 spin2x,
125 const REAL8 spin2y,
126 const REAL8 spin2z,
127 REAL8Vector * timeVec,
128 REAL8Vector * matchrange,
130 const REAL8 JLN,
131 const REAL8 combsize,
132 const REAL8 tMaxOmega,
133 const REAL8 tMaxAmp
134);
135
136
137
138#endif
double XLALSimLocateAmplTime(REAL8Vector *timeHi, COMPLEX16Vector *hP22, REAL8Vector *radiusVec, int *found, REAL8 *tMaxAmp)
double XLALSimLocateMaxAmplTime(REAL8Vector *timeHi, COMPLEX16Vector *hP22, int *found)
int XLALSimAdjustRDattachmentTime(REAL8Vector *signal1, REAL8Vector *signal2, COMPLEX16TimeSeries *h22, COMPLEX16TimeSeries *h2m2, REAL8 *ratio22, REAL8 *ratio2m2, REAL8 *tAtt, const REAL8 thr, const REAL8 dt, const REAL8 mass1, const REAL8 mass2, const REAL8 spin1x, const REAL8 spin1y, const REAL8 spin1z, const REAL8 spin2x, const REAL8 spin2y, const REAL8 spin2z, REAL8Vector *timeVec, REAL8Vector *matchrange, Approximant approximant, const REAL8 JLN, const REAL8 combsize, const REAL8 tMaxOmega, const REAL8 tMaxAmp)
double XLALSimLocateOmegaTime(REAL8Array *dynamicsHi, unsigned int numdynvars, unsigned int retLenHi, SpinEOBParams seobParams, SpinEOBHCoeffs seobCoeffs, REAL8 m1, REAL8 m2, REAL8Vector *radiusVec, int *found, REAL8 *tMaxOmega, INT4 use_optimized)
INT4 XLALSimCheckRDattachment(REAL8Vector *signal1, REAL8Vector *signal2, REAL8 *ratio, const REAL8 tAtt, const INT4 l, const INT4 m, const REAL8 dt, const REAL8 mass1, const REAL8 mass2, const REAL8 spin1x, const REAL8 spin1y, const REAL8 spin1z, const REAL8 spin2x, const REAL8 spin2y, const REAL8 spin2z, REAL8Vector *timeVec, REAL8Vector *matchrange, Approximant approximant, const REAL8 JLN, REAL8 *timediff)
int l
Definition: bh_qnmode.c:135
double dt
Definition: bh_ringdown.c:113
double REAL8
int32_t INT4
Approximant
Enum that specifies the PN approximant to be used in computing the waveform.
static const INT4 m
Parameters for the spinning EOB model, used in calculating the Hamiltonian.
Parameters for the spinning EOB model.