Coverage for /home/andrew/Documents/Research/gwent/gwent/binary.py : 87%

Hot-keys on this page
r m x p toggle line displays
j k next/prev highlighted chunk
0 (zero) top of page
1 (one) first highlighted chunk
1import numpy as np
2import os
3import astropy.constants as const
4import astropy.units as u
5import scipy.interpolate as interp
6from astropy.cosmology import z_at_value
7from astropy.cosmology import WMAP9 as cosmo
9import gwent
10from .waveform import Get_Waveform
11from . import utils
13current_path = os.path.abspath(gwent.__path__[0])
14load_directory = os.path.join(current_path, "LoadFiles/")
17class BinaryBlackHole:
18 """Base Class for frequency domain strains from Binary Black Holes.
20 Parameters
21 ----------
22 M : float
23 Total mass of the black hole binary (m1+m2)
24 q : float
25 Mass ratio of the black hole binary (m1/m2, m1<m2)
26 z : float
27 Redshift of the black hole binary
29 load_location : string, optional
30 the directory of the loaded file, (ie. '/path/to/file')
32 Notes
33 -----
34 IMRPhenomD waveforms calibrated for q = m1/m2 < 18
36 """
38 def __init__(self, *args, **kwargs):
39 if len(args) == 3:
40 [M, q, z] = args
41 elif len(args) == 5:
42 [M, q, z, _, _] = args
43 else:
44 raise ValueError(
45 "args must be a list of 3 ([M,q,z]) or 5 ([M,q,z,chi1,chi2])"
46 )
47 self.M = M
48 self.q = q
49 self.z = z
51 for keys, value in kwargs.items():
52 if keys == "load_location":
53 self.load_location = value
55 if hasattr(self, "load_location"):
56 self.Load_Data()
58 @property
59 def M(self):
60 self._M = utils.make_quant(self._M, "M_sun")
61 return self._M
63 @M.setter
64 def M(self, value):
65 self.var_dict = ["M", value]
66 self._M = self._return_value
68 @property
69 def q(self):
70 return self._q
72 @q.setter
73 def q(self, value):
74 self.var_dict = ["q", value]
75 self._q = self._return_value
77 @property
78 def z(self):
79 return self._z
81 @z.setter
82 def z(self, value):
83 self.var_dict = ["z", value]
84 self._z = self._return_value
86 @property
87 def h_f(self):
88 if not hasattr(self, "_h_f"):
89 raise NotImplementedError(
90 "The strain must be defined inside BBHFrequencyDomain or BBHTimeDomain classes."
91 )
92 return self._h_f
94 @h_f.setter
95 def h_f(self, value):
96 self._h_f = value
98 @property
99 def f(self):
100 if not hasattr(self, "_f"):
101 raise NotImplementedError(
102 "The frequency must be defined inside BBHFrequencyDomain or BBHTimeDomain classes."
103 )
104 return self._f
106 @f.setter
107 def f(self, value):
108 self._f = value
110 @property
111 def var_dict(self):
112 return self._var_dict
114 @var_dict.setter
115 def var_dict(self, value):
116 utils.Get_Var_Dict(self, value)
118 def Load_Data(self):
119 if hasattr(self, "load_location"):
120 if os.path.exists(self.load_location):
121 self._load_data = np.loadtxt(self.load_location)
122 else:
123 raise IOError(
124 "File %s does not exist, please assign load_location a correct filepath."
125 % self.load_location
126 )
127 else:
128 raise ValueError(
129 'load_location is not assigned, please set with name_of_BBH.load_location="path/to/file".'
130 )
133class BBHFrequencyDomain(BinaryBlackHole):
134 """Subclass of BinaryBlackHole for BBH GWs generated in the frequency domain.
136 Parameters
137 ----------
138 chi1 : float
139 The dimensionless spin parameter abs(a/m) for black hole m1.
140 chi2 : float
141 The dimensionless spin parameter abs(a/m) for black hole m2
143 f_low : float, optional
144 The lowest frequency in natural units (Mf, G=c=1) at which the BBH waveform is calculated
145 nfreqs : int, optional
146 The number of frequencies at which the BBH waveform is calculated
147 instrument: object, optional
148 If assigned, the optimal frequency (ie. most sensitive frequency) of the detector is used as
149 the binary's GW frequency
151 Notes
152 -----
153 IMRPhenomD waveforms calibrated for aligned spins chi_1, chi_2 = abs(a/m) <= 0.85 or if q=1 abs(a/m)<0.98
155 """
157 def __init__(self, *args, **kwargs):
158 super().__init__(*args, **kwargs)
159 [_, _, _, chi1, chi2] = args
160 self.chi1 = chi1
161 self.chi2 = chi2
163 for keys, value in kwargs.items():
164 if keys == "f_low":
165 self.f_low = value
166 elif keys == "f_high":
167 self.f_high = value
168 elif keys == "nfreqs":
169 self.nfreqs = value
170 elif keys == "instrument":
171 self.instrument = value
172 self.Check_Freq_Evol()
173 elif keys == "f_gw":
174 self.f_gw = value
175 if not hasattr(self, "nfreqs"):
176 self.nfreqs = int(1e3)
177 if not hasattr(self, "f_low"):
178 self.f_low = 1e-5
180 self.Get_Fitcoeffs()
182 @property
183 def chi1(self):
184 return self._chi1
186 @chi1.setter
187 def chi1(self, value):
188 self.var_dict = ["chi1", value]
189 self._chi1 = self._return_value
191 @property
192 def chi2(self):
193 return self._chi2
195 @chi2.setter
196 def chi2(self, value):
197 self.var_dict = ["chi2", value]
198 self._chi2 = self._return_value
200 @property
201 def instrument(self):
202 return self._instrument
204 @instrument.setter
205 def instrument(self, value):
206 self._instrument = value
208 @property
209 def h_gw(self):
210 if not hasattr(self, "_h_gw"):
211 if hasattr(self, "_instrument"):
212 self._h_gw = Get_Mono_Strain(self, frame="source").to("")
213 else:
214 raise ValueError(
215 "No instrument assigned, please fix it. "
216 'Try: "source.instrument = instrument".'
217 )
218 return self._h_gw
220 @h_gw.setter
221 def h_gw(self, value):
222 self._h_gw = value
224 @h_gw.deleter
225 def h_gw(self):
226 del self._h_gw
228 @property
229 def f_gw(self):
230 if not hasattr(self, "_f_gw"):
231 if hasattr(self, "_instrument"):
232 self._f_gw = self.instrument.f_opt
233 else:
234 raise ValueError(
235 "No GW frequency or instrument assigned, please fix it. "
236 "Try: assigning source.f_gw or source.instrument."
237 )
238 return self._f_gw
240 @f_gw.setter
241 def f_gw(self, value):
242 self._f_gw = value
244 @f_gw.deleter
245 def f_gw(self):
246 del self._f_gw
248 @property
249 def h_f(self):
250 if not hasattr(self, "_h_f"):
251 if not (hasattr(self, "_phenomD_f") and hasattr(self, "_phenomD_h")):
252 self.Get_PhenomD_Strain()
253 [_, self._h_f] = Strain_Conv(self, self._phenomD_f, self._phenomD_h)
254 return self._h_f
256 @h_f.deleter
257 def h_f(self):
258 del self._h_f
260 @property
261 def f(self):
262 if not hasattr(self, "_f"):
263 if not (hasattr(self, "_phenomD_f") and hasattr(self, "_phenomD_h")):
264 self.Get_PhenomD_Strain()
265 [self._f, _] = Strain_Conv(self, self._phenomD_f, self._phenomD_h)
266 return self._f
268 @f.deleter
269 def f(self):
270 del self._f
272 def Get_Fitcoeffs(self):
273 """Loads Quasi-Normal Mode fitting files for speed later."""
274 fit_coeffs_filedirectory = os.path.join(
275 load_directory, "PhenomDFiles/fitcoeffsWEB.dat"
276 )
277 self._fitcoeffs = np.loadtxt(fit_coeffs_filedirectory)
279 def Get_PhenomD_Strain(self):
280 """Gets the BBH's frequency and waveform from IMRPhenomD."""
281 if not hasattr(self, "_fitcoeffs"):
282 self.Get_Fitcoeffs()
283 [self._phenomD_f, self._phenomD_h] = Get_Waveform(self)
285 def Get_F_Dot(self, freq, frame="observer"):
286 """Calculates the change in frequency of a binary black hole at a given frequency.
288 Parameters
289 ----------
290 freq : float
291 the binary GW frequency in the corresponding frame.
292 frame : str, {'observer','source'}
293 Determines whether the given frequency is in the source or observer frame.
295 """
296 freq = utils.make_quant(freq, "1/s")
297 m_conv = const.G / const.c ** 3 # Converts M = [M] to M = [sec]
298 eta = self.q / (1 + self.q) ** 2
300 M_time = self.M.to("kg") * m_conv
301 M_chirp = eta ** (3 / 5) * M_time
303 if frame == "observer":
304 f_obs_source = freq * (1 + self.z)
305 elif frame == "source":
306 f_obs_source = freq
307 else:
308 raise ValueError("The reference frame can only be observer or source.")
310 return (
311 96
312 / 5
313 * np.pi ** (8 / 3)
314 * (M_chirp ** (5 / 3))
315 * (f_obs_source ** (11 / 3))
316 )
318 def Get_Time_From_Merger(self, freq, frame="observer"):
319 """Calculates the time from merger of a binary black hole given a frequency.
321 Parameters
322 ----------
323 freq : float
324 the binary GW frequency in the corresponding frame.
325 frame : str, {'observer','source'}
326 Determines whether the given frequency is in the source or observer frame.
328 """
329 freq = utils.make_quant(freq, "1/s")
330 m_conv = const.G / const.c ** 3 # Converts M = [M] to M = [sec]
331 eta = self.q / (1 + self.q) ** 2
333 M_time = self.M.to("kg") * m_conv
334 M_chirp = eta ** (3 / 5) * M_time
336 if frame == "observer":
337 f_obs_source = freq * (1 + self.z)
338 elif frame == "source":
339 f_obs_source = freq
340 else:
341 raise ValueError("The reference frame can only be observer or source.")
343 return 5 * (M_chirp) ** (-5 / 3) * (8 * np.pi * f_obs_source) ** (-8 / 3)
345 def Get_Source_Freq(self, tau, frame="observer"):
346 """Calculates the binary black hole's gravitational wave frequency given a time from merger
348 Parameters
349 ----------
350 tau : int, float, Quantity
351 the time from merger in the respective frame.
352 If not an astropy quantity, assumed to be a time in seconds
353 frame : str, {'observer','source'}
354 Determines whether the given frequency is in the source or observer frame.
356 """
357 tau = utils.make_quant(tau, "s")
358 if frame == "observer":
359 tau_source = tau / (1 + self.z)
360 elif frame == "source":
361 tau_source = tau
362 else:
363 raise ValueError("The reference frame can only be observer or source.")
365 m_conv = const.G / const.c ** 3 # Converts M = [M] to M = [sec]
366 eta = self.q / (1 + self.q) ** 2
368 M_time = self.M.to("kg") * m_conv
369 M_chirp = eta ** (3 / 5) * M_time
371 return 1.0 / 8.0 / np.pi / M_chirp * (5 * M_chirp / tau_source) ** (3.0 / 8.0)
373 def Check_Freq_Evol(
374 self, T_evol=None, T_evol_frame="observer", f_gw_frame="source"
375 ):
376 """Checks the frequency evolution of the black hole binary.
378 Parameters
379 ----------
380 T_evol : int,float, Quantity
381 The length of time the binary may evolve
382 T_evol_frame : str, {'observer','source'}
383 Determines whether the given T_evol is in the source or observer frame.
384 f_gw_frame : str, {'observer','source'}
385 Determines whether the frequency is in the source or observer frame.
386 May not be used if source has an instrument assigned.
388 Notes
389 -----
390 If the frequency of the binary does evolve over more than one bin,
391 (ie f(T_evol)-f(t_init) = delf_obs < 1/T_evol), it is monochromatic, so we set the frequency
392 to the optimal frequency of the detector
394 Otherwise it is chirping and evolves over the observation and we
395 set the starting frequency we observe it at to f(Tobs), which is the
396 frequency at an observation time before merger
398 To get the change in frequency, we use eqn 41 from Hazboun,Romano, and Smith (2019) <https://arxiv.org/abs/1907.04341>
399 which uses binomial expansion of f_T_evol_inst - f_init_inst and thus will never be imaginary
401 """
403 m_conv = const.G / const.c ** 3 # Converts M = [M] to M = [sec]
404 eta = self.q / (1 + self.q) ** 2
406 M_time = self.M.to("kg") * m_conv
407 M_chirp_source = eta ** (3 / 5) * M_time
409 if T_evol is not None:
410 T_evol = utils.make_quant(T_evol, "s")
411 if hasattr(self, "_instrument"):
412 if hasattr(self, "_f_gw"):
413 t_init_source = self.Get_Time_From_Merger(
414 self.f_gw, frame=f_gw_frame
415 )
416 else:
417 # Assumes f_init is the optimal frequency in the instrument frame to get t_init_source
418 t_init_source = self.Get_Time_From_Merger(
419 self.instrument.f_opt, frame="observer"
420 )
421 # f(T_evol), the frequency of the source at T_evol before merger
422 f_T_evol_source = self.Get_Source_Freq(T_evol, frame=T_evol_frame)
423 elif hasattr(self, "_f_gw") and not hasattr(self, "_instrument"):
424 t_init_source = self.Get_Time_From_Merger(self.f_gw, frame=f_gw_frame)
425 # f(T_evol), the frequency of the source at T_evol before merger
426 f_T_evol_source = self.Get_Source_Freq(T_evol, frame=T_evol_frame)
427 else:
428 raise ValueError(
429 "Must either assign T_evol a value, or assign the source an instrument."
430 )
431 else:
432 if hasattr(self, "instrument"):
433 T_evol_frame == "observer"
434 T_evol = utils.make_quant(np.max(np.unique(self.instrument.T_obs)), "s")
435 # Assumes f_init is the optimal frequency in the instrument frame to get t_init_source
436 t_init_source = self.Get_Time_From_Merger(
437 self.instrument.f_opt, frame="observer"
438 )
439 # f(T_evol), the frequency of the source at T_evol before merger
440 f_T_evol_source = self.Get_Source_Freq(T_evol, frame="observer")
441 else:
442 raise ValueError(
443 "Must either assign T_evol a value, or assign the source an instrument."
444 )
446 # Assumes t_init is in source frame, can either be randomly drawn
447 # t_init_source = np.random.uniform(0,100)*u.yr
449 if T_evol_frame == "observer":
450 T_obs = T_evol
451 T_evol_source = T_evol / (1 + self.z)
452 elif T_evol_frame == "source":
453 T_obs = T_evol * (1 + self.z)
454 T_evol_source = T_evol
455 else:
456 raise ValueError("The reference frame can only be observer or source.")
458 # f(T_evol) in the instrument frame (ie. the observed frequency)
459 self.f_T_obs = f_T_evol_source / (1 + self.z)
461 # t_init_source = make_quant(t_init_source,'s')
462 # f_init_source = self.Get_Source_Freq(t_init_source)
463 # self.f_init = f_init_source/(1+self.z)
464 # f_after_T_obs_source = self.Get_Source_Freq((t_init_source-T_obs_source))
465 # self.f_T_obs = f_after_T_obs_source/(1+self.z)
466 # delf_obs_source_exact = f_after_T_obs_source-f_init_source
468 delf_source = (
469 1.0
470 / 8.0
471 / np.pi
472 / M_chirp_source
473 * (5 * M_chirp_source / t_init_source) ** (3.0 / 8.0)
474 * (3 * T_evol_source / 8 / t_init_source)
475 )
476 delf_obs = delf_source / (1 + self.z)
477 # print('delf_obs: ',delf_obs)
478 # print('1/T_obs: ',1/T_obs)
479 if delf_obs < (1 / T_obs):
480 self.ismono = True
481 else:
482 self.ismono = False
485class BBHTimeDomain(BinaryBlackHole):
486 """Subclass of BinaryBlackHole for input in the time domain"""
488 def __init__(self, *args, **kwargs):
489 super().__init__(*args, **kwargs)
490 self.Get_hf_from_hcross_hplus()
492 @property
493 def t(self):
494 if not hasattr(self, "_t"):
495 self._t = self._load_data[:, 0]
496 self._t = utils.make_quant(self._t, "s")
497 return self._t
499 @property
500 def h_plus_t(self):
501 if not hasattr(self, "_h_plus_t"):
502 self._h_plus_t = self._load_data[:, 1]
503 return self._h_plus_t
505 @property
506 def h_cross_t(self):
507 if not hasattr(self, "_h_cross_t"):
508 self._h_cross_t = self._load_data[:, 1]
509 return self._h_cross_t
511 @property
512 def h_f(self):
513 if not hasattr(self, "_h_f"):
514 [natural_f, natural_h] = self.Get_hf_from_hcross_hplus()
515 [_, self._h_f] = Strain_Conv(self, natural_f, natural_h)
516 return self._h_f
518 @h_f.deleter
519 def h_f(self):
520 del self._h_f
522 @property
523 def f(self):
524 if not hasattr(self, "_f"):
525 [natural_f, natural_h] = self.Get_hf_from_hcross_hplus()
526 [self._f, _] = Strain_Conv(self, natural_f, natural_h)
527 return self._f
529 @f.deleter
530 def f(self):
531 del self._f
533 def Get_hf_from_hcross_hplus(self, interp_res="coarse", windowing="left"):
534 """Converts dimensionless, time domain strain to frequency space using a windowed fft
536 Parameters
537 ----------
538 interp_res : {'coarse','fine'}, optional
539 'coarse' uses maximum difference between subsequent time steps for interpolation
540 'fine' uses minimum difference between subsequent time steps for interpolation
541 windowing : {'left','right','all'}, optional
542 'left' windows the left side of the time data
543 'right' windows the right side of the time data
544 'all' windows the both the left and right side of the time data
546 Returns
547 -------
548 natural_f : array
549 The frequency of the input source in natural units (G=c=1)
550 natural_h : array
551 The strain of the input source in natural units (G=c=1)
553 """
555 # Interpolate time to evenly sampled data, can be fine or coarse
556 diff_t = np.diff(self.t.value)
557 if interp_res == "fine":
558 dt = min(diff_t)
559 elif interp_res == "coarse":
560 dt = max(diff_t)
562 interp_t = np.arange(self.t[0].value, self.t[-1].value, dt)
563 # interpolate strain to evenly sampled data for FFT
564 h_cross_t = interp.interp1d(self.t, self.h_cross_t, kind="cubic")
565 h_plus_t = interp.interp1d(self.t, self.h_plus_t, kind="cubic")
566 interp_h_cross_t = h_cross_t(interp_t)
567 interp_h_plus_t = h_plus_t(interp_t)
569 # Filter/Window
570 hann_window = np.hanning(len(interp_t)) # Two sided
571 if windowing == "left":
572 #########################
573 """Applies window to first (left) half"""
574 first_half = hann_window[
575 : int(len(interp_t) / 2)
576 ] # Only need tapering on first half of waveform
577 second_half = np.ones(
578 len(interp_t) - len(first_half)
579 ) # no windowing on second half of waveform
580 #########################
581 window = np.append(
582 first_half, second_half
583 ) # Only apply window to first half of waveform
584 elif windowing == "right":
585 #########################
586 """Applies window to second (right) half"""
587 second_half = hann_window[
588 int(len(interp_t) / 2) :
589 ] # Only need tapering on second half of waveform
590 first_half = np.ones(
591 len(interp_t) - len(second_half)
592 ) # no windowing on first half of waveform
593 #########################
594 window = np.append(first_half, second_half)
595 elif windowing == "all":
596 window = hann_window
597 # Window!
598 win_h_cross_t = np.multiply(interp_h_cross_t, window)
599 win_h_plus_t = np.multiply(interp_h_plus_t, window)
601 # FFT the two polarizations
602 h_cross_f = np.fft.fft(win_h_cross_t)
603 h_plus_f = np.fft.fft(win_h_plus_t)
604 freqs = np.fft.fftfreq(len(interp_t), d=dt)
606 # cut = np.abs(freqs).argmax() #Cut off the negative frequencies
607 f_cut_low = 3e-3 # Low Cutoff frequency
608 f_cut_high = 1.5e-1 # High Cutoff frequency
609 cut_low = np.abs(
610 freqs - f_cut_low
611 ).argmin() # Cut off frequencies lower than a frequency
612 cut_high = np.abs(
613 freqs - f_cut_high
614 ).argmin() # Cut off frequencies higher than a frequency
615 # cut=int(len(freqs)*0.9) #Cut off percentage of frequencies
616 h_cross_f = h_cross_f[cut_low:cut_high]
617 h_plus_f = h_plus_f[cut_low:cut_high]
618 natural_f = freqs[cut_low:cut_high]
620 # Combine them for raw spectral power
621 natural_h_f = np.sqrt((np.abs(h_cross_f)) ** 2 + (np.abs(h_plus_f)) ** 2)
622 return [natural_f, natural_h_f]
625def Strain_Conv(source, natural_f, natural_h):
626 """Converts frequency and strain in natural units (G=c=1) to Hertz and dimensionless, respectively.
628 Parameters
629 ----------
630 source
631 Instance of gravitational wave source class
632 natural_f : array [Mf]
633 the frequency of the source in natural units (G=c=1)
634 natural_h : array [Mf]
635 the strain of the source in natural units (G=c=1)
637 """
638 DL = cosmo.luminosity_distance(source.z)
639 DL = DL.to("m")
641 m_conv = const.G / const.c ** 3 # Converts M = [M] to M = [sec]
642 M_redshifted_time = source.M.to("kg") * (1 + source.z) * m_conv
644 # frequency and strain of source in detector frame
645 freq_conv = 1 / M_redshifted_time
646 # Normalized factor to match Stationary phase approx at low frequencies
647 strain_conv = np.sqrt(5 / 24 / np.pi) * (const.c / DL) * M_redshifted_time ** 2
649 f = natural_f * freq_conv
650 h_f = natural_h * strain_conv
651 return [f, h_f]
654def Get_Char_Strain(source):
655 """Converts source strain to characteristic strain
657 Parameters
658 ----------
659 source : object
660 Instance of gravitational wave source class
662 """
663 h_char = np.sqrt(4 * source.f ** 2 * source.h_f ** 2)
664 return h_char
667def Get_Mono_Char_Strain(
668 source, T_obs, inc=None, T_evol_frame="observer", f_gw_frame="observer"
669):
670 """Calculates the characteristic strain from a binary black hole observed for some observation time.
672 Parameters
673 ----------
674 source : object
675 Instance of gravitational wave source class
676 T_obs : int,float, Quantity
677 The length of time for which the binary is observed
678 inc : float, optional
679 The inclination of the source in radians. If inc is None, the strain is \
680 sky and inclination averaged strain from Robson et al. 2019 (eqn 27) <https://arxiv.org/pdf/1803.01944.pdf> \
681 T_evol_frame : str, {'observer','source'}
682 Determines whether the given T_evol is in the source or observer frame.
683 f_gw_frame : str, {'observer','source'}
684 Determines whether the frequency is in the source or observer frame.
685 May not be used if source has an instrument assigned.
687 Returns
688 -------
689 float
690 The characteristic strain of a monochromatic source in the source frame.
691 """
692 strain_amp = Get_Mono_Strain(source, inc=inc, frame=f_gw_frame)
694 f_gw = utils.make_quant(source.f_gw, "Hz")
695 if f_gw_frame == "observer":
696 f_obs_source = f_gw * (1 + source.z)
697 elif f_gw_frame == "source":
698 f_obs_source = f_gw
699 else:
700 raise ValueError("The reference frame can only be observer or source.")
702 f_dot_source = source.Get_F_Dot(f_obs_source, frame="source")
703 # print(f_dot_source)
704 """
705 T_obs = utils.make_quant(T_obs, "s")
706 if T_evol_frame == "observer":
707 T_obs_source = T_obs / (1 + source.z)
708 elif T_evol_frame == "source":
709 T_obs_source = T_obs
710 else:
711 raise ValueError("The reference frame can only be observer or source.")
713 return np.sqrt(f_obs_source*T_obs_source)*strain_amp
714 """
715 return f_obs_source * np.sqrt(2 / f_dot_source) * strain_amp
718def Get_Mono_Strain(source, inc=None, frame="source"):
719 """Calculates the fourier domain strain from a monochromatic binary black hole.
721 Parameters
722 ----------
723 source : object
724 Instance of gravitational wave source class
725 inc : float, optional
726 The inclination of the source in radians. If inc is None, the strain is \
727 sky and inclination averaged strain from Robson et al. 2019 (eqn 27) <https://arxiv.org/pdf/1803.01944.pdf> \
728 frame : str, {'source','observer'}
729 Determines whether the frequency is in the source or observer frame. \
730 May not be used if source has an instrument assigned.
732 Returns
733 -------
734 float
735 The strain of a monochromatic source in the source frame.
737 """
738 f_gw = utils.make_quant(source.f_gw, "Hz")
739 if frame == "observer":
740 f_obs_source = f_gw * (1 + source.z)
741 elif frame == "source":
742 f_obs_source = f_gw
743 else:
744 raise ValueError("The reference frame can only be observer or source.")
746 DL = cosmo.luminosity_distance(source.z)
747 DL = DL.to("m")
749 # Converts M = [M] to M = [sec]
750 m_conv = const.G / const.c ** 3
752 eta = source.q / (1 + source.q) ** 2
753 M_redshifted_time = source.M.to("kg") * (1 + source.z) * m_conv
754 M_chirp = eta ** (3 / 5) * M_redshifted_time
756 if inc is not None:
757 if inc > np.pi or inc < -np.pi:
758 raise ValueError("Inclination must be between -pi and pi.")
759 a = (1 + np.cos(inc) ** 2) / 2
760 b = np.cos(inc)
761 const_val = 4.0 * np.sqrt(a ** 2 + b ** 2)
762 else:
763 const_val = 8 / np.sqrt(5)
765 return (
766 const_val
767 * (const.c / DL)
768 * (np.pi * f_obs_source) ** (2.0 / 3.0)
769 * M_chirp ** (5.0 / 3.0)
770 )