Coverage for /home/andrew/Documents/Research/gwent/gwent/detector.py : 86%

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
6import scipy.stats as stats
8from astropy.cosmology import z_at_value
9from astropy.cosmology import WMAP9 as cosmo
11import gwent
12from . import utils
14import gwinc
15import hasasia.sensitivity as hassens
16import hasasia.sim as hassim
18current_path = os.path.abspath(gwent.__path__[0])
19load_directory = os.path.join(current_path, "LoadFiles/")
22class PTA:
23 r"""
24 Class to make a PTA instrument using the methods of Hazboun, Romano, Smith (2019) <https://arxiv.org/abs/1907.04341>
26 Parameters
27 ----------
29 name : string
30 name of the instrument
32 n_p : int
33 the number of pulsars in the PTA
35 T_obs : float,array, list, optional
36 the observation time of the PTA in [years]
37 If T_obs is a float, the observation time is used as the observation time for all pulsars
38 If T_obs is a list of length of n_p, the observation times are used as the corresponding pulsar observation times.
39 If T_obs is a list of length 2, it is assumed the values are the minimum and maximum observation time values
40 (ie. [min,max]) from which individual pulsar observation times are uniformly sampled.
41 sigma : float, array, list, optional
42 the rms error on the pulsar TOAs in [sec]
43 If sigma is a float, the given rms error is used for all pulsars.
44 If sigma is a list of length of n_p, the amplitudes are used as the corresponding pulsar rms error.
45 If sigma is a list of length 2, it is assumed the values are the minimum and maximum rms errors values
46 (ie. [min,max]) from which individual pulsar rms errors are uniformly sampled.
47 cadence : float, array, list, optional
48 How often the pulsars are observed in [num/year]
49 If cadence is a float, the given cadence is used for all pulsars.
50 If cadence is a list of length of n_p, the amplitudes are used as the corresponding pulsar cadence.
51 If cadence is a list of length 2, it is assumed the values are the minimum and maximum cadence values
52 (ie. [min,max]) from which individual pulsar cadences are uniformly sampled.
53 sb_amp : float, optional
54 Amplitude of the Stochastic background added as red noise
55 sb_alpha : float, optional
56 the Stochastic background power law, if empty and sb_amp is set, it is assumed to be -2/3
57 rn_amp : array, list, optional
58 Individual pulsar red noise amplitude.
59 If rn_amp is a list of length of n_p, the amplitudes are used as the corresponding pulsar RN injection.
60 If rn_amp is a list of length 2, it is assumed the values are the minimum and maximum RN amplitude values
61 (ie. [min,max]) from which individual pulsar RN amplitudes are uniformly sampled.
62 rn_alpha : array, list, optional
63 Individual pulsar red noise alpha (power law spectral index).
64 If rn_alpha is a list of length of n_p, the alpha indices are used as the corresponding pulsar RN injection.
65 If rn_alpha is a list of length 2, it is assumed the values are the minimum and maximum RN alpha values
66 (ie. [min,max]) from which individual pulsar RN alphas are uniformly sampled.
67 phi : list, array, optional
68 Individual pulsar longitude in ecliptic coordinates.
69 If not defined, NANOGrav 11yr pulsar locations are used.
70 If n_p > 34 (the number of pulsars in the 11yr dataset),
71 it draws more pulsars from distributions based on the NANOGrav 11yr pulsars.
72 theta : array, list, optional
73 Individual pulsar colatitude in ecliptic coordinates.
74 If not defined, NANOGrav 11yr pulsar locations are used.
75 If n_p > 34 (the number of pulsars in the 11yr dataset),
76 it draws more pulsars from distributions based on the NANOGrav 11yr pulsars.
77 use_11yr : bool, optional
78 Uses the NANOGrav 11yr noise as the individual pulsar noises,
79 if n_p > 34 (the number of pulsars in the 11yr dataset),
80 it draws more pulsars from distributions based on the NANOGrav 11yr pulsar noise
81 use_rn : bool, optional
82 If no rn_amp assigned, uses the NANOGrav 11yr noise as the individual pulsar RN noises,
83 if n_p > 34 (the number of pulsars in the 11yr dataset),
84 it draws more pulsars from distributions based on the NANOGrav 11yr pulsar noise
85 load_location : string, optional
86 If you want to load a PTA curve from a file, it's the file path
87 I_type : string, {'E','A','h'}
88 Sets the type of input data.
89 'E' is the effective strain spectral density $S_{n}(f)$ ('ENSD'),
90 'A' is the amplitude spectral density, $\sqrt{S_{n}(f)}$ ('ASD'),
91 'h' is the characteristic strain $h_{n}(f)$ ('h')
92 f_low : float, optional
93 Assigned lowest frequency of PTA (default assigns 1/(5*T_obs))
94 f_high : float, optional
95 Assigned highest frequency of PTA (default is Nyquist freq cadence/2)
96 nfreqs : int, optional
97 Number of frequencies in logspace the sensitivity is calculated
98 nbins : int, optional
99 Used to add values to every bin for sampled parameters. Default is 8 for smooth, non-zero distributions.
100 Changing this could change distribution, so be wary, not sure how much it affects anything.
101 """
103 def __init__(self, name, *args, **kwargs):
104 self.name = name
105 if len(args) != 0:
106 if len(args) == 1:
107 self.n_p = args[0]
108 else:
109 raise ValueError("Too many args, not enough kwargs!")
111 for keys, value in kwargs.items():
112 if keys == "T_obs":
113 self.T_obs = value
114 elif keys == "sigma":
115 self.sigma = value
116 elif keys == "cadence":
117 self.cadence = value
118 elif keys == "rn_amp":
119 self.rn_amp = value
120 elif keys == "rn_alpha":
121 self.rn_alpha = value
122 elif keys == "phi":
123 self.phi = value
124 elif keys == "theta":
125 self.theta = value
126 elif keys == "sb_amp":
127 self.sb_amp = value
128 elif keys == "sb_alpha":
129 self.sb_alpha = value
130 elif keys == "use_11yr":
131 self.use_11yr = value
132 elif keys == "use_rn":
133 self.use_rn = value
134 elif keys == "load_location":
135 self.load_location = value
136 elif keys == "I_type":
137 self.I_type = value
138 elif keys == "f_low":
139 self.f_low = utils.make_quant(value, "Hz")
140 elif keys == "f_high":
141 self.f_high = utils.make_quant(value, "Hz")
142 elif keys == "nfreqs":
143 self.nfreqs = value
144 elif keys == "nbins":
145 self.nbins = value
146 else:
147 raise ValueError("%s is not an accepted input option." % keys)
149 if not hasattr(self, "nfreqs"):
150 self.nfreqs = int(1e3)
151 if not hasattr(self, "nbins"):
152 self.nbins = 8
153 if hasattr(self, "load_location"):
154 Load_Data(self)
156 if not hasattr(self, "use_11yr"):
157 self.use_11yr = False
158 if not hasattr(self, "use_rn"):
159 self.use_rn = False
161 if hasattr(self, "f_low") and hasattr(self, "f_high"):
162 self.fT = (
163 np.logspace(
164 np.log10(self.f_low.value), np.log10(self.f_high.value), self.nfreqs
165 )
166 * u.Hz
167 )
169 @property
170 def n_p(self):
171 return self._n_p
173 @n_p.setter
174 def n_p(self, value):
175 self.var_dict = ["n_p", value]
176 self._n_p = self._return_value
178 @property
179 def T_obs(self):
180 return self._T_obs
182 @T_obs.setter
183 def T_obs(self, value):
184 self.var_dict = ["T_obs", value]
185 if not isinstance(self._return_value, u.Quantity):
186 self._return_value = utils.make_quant(self._return_value, "yr")
187 self._T_obs = self._return_value
189 @property
190 def cadence(self):
191 return self._cadence
193 @cadence.setter
194 def cadence(self, value):
195 self.var_dict = ["cadence", value]
196 if not isinstance(self._return_value, u.Quantity):
197 self._return_value = utils.make_quant(self._return_value, "1/yr")
198 self._cadence = self._return_value
200 @property
201 def sigma(self):
202 return self._sigma
204 @sigma.setter
205 def sigma(self, value):
206 self.var_dict = ["sigma", value]
207 if not isinstance(self._return_value, u.Quantity):
208 self._return_value = utils.make_quant(self._return_value, "s")
209 self._sigma = self._return_value
211 @property
212 def phi(self):
213 return self._phi
215 @phi.setter
216 def phi(self, value):
217 self.var_dict = ["phi", value]
218 self._phi = self._return_value
220 @property
221 def theta(self):
222 return self._theta
224 @theta.setter
225 def theta(self, value):
226 self.var_dict = ["theta", value]
227 self._theta = self._return_value
229 @property
230 def rn_amp(self):
231 return self._rn_amp
233 @rn_amp.setter
234 def rn_amp(self, value):
235 self.var_dict = ["rn_amp", value]
236 self._rn_amp = self._return_value
238 @property
239 def rn_alpha(self):
240 return self._rn_alpha
242 @rn_alpha.setter
243 def rn_alpha(self, value):
244 self.var_dict = ["rn_alpha", value]
245 self._rn_alpha = self._return_value
247 @property
248 def var_dict(self):
249 return self._var_dict
251 @var_dict.setter
252 def var_dict(self, value):
253 utils.Get_Var_Dict(self, value)
255 @property
256 def fT(self):
257 if not hasattr(self, "_fT"):
258 # frequency sampled from 1/observation time to nyquist frequency (c/2)
259 # 5 is the default value for now (from Hazboun et al. 2019)
260 if not hasattr(self, "_T_obs") or not hasattr(self, "_cadence"):
261 self.Init_PTA()
262 T_obs_sec = np.max(self.T_obs).to("s").value
263 cadence_sec = np.max(self.cadence).to("1/s").value
264 self._fT = (
265 np.logspace(
266 np.log10(1 / (5 * T_obs_sec)),
267 np.log10(cadence_sec / 2),
268 self.nfreqs,
269 )
270 * u.Hz
271 )
272 return self._fT
274 @fT.setter
275 def fT(self, value):
276 self._fT = value
278 @fT.deleter
279 def fT(self):
280 del self._fT
282 @property
283 def h_n_f(self):
284 """Effective Strain Noise Amplitude"""
285 if not hasattr(self, "_h_n_f"):
286 if hasattr(self, "_I_data"):
287 if self._I_Type == "h":
288 self._h_n_f = self._I_data[:, 1]
289 elif self._I_Type == "ENSD":
290 self._h_n_f = np.sqrt(self.S_n_f * self.fT)
291 elif self._I_Type == "ASD":
292 S_n_f_sqrt = self._I_data[:, 1]
293 self._h_n_f = S_n_f_sqrt * np.sqrt(self.fT.value)
294 else:
295 if not hasattr(self, "_sensitivitycurve"):
296 self.Init_PTA()
297 self._h_n_f = self._sensitivitycurve.h_c
298 return self._h_n_f
300 @h_n_f.setter
301 def h_n_f(self, value):
302 self._h_n_f = value
304 @h_n_f.deleter
305 def h_n_f(self):
306 del self._h_n_f
308 @property
309 def S_n_f(self):
310 """Effective noise power amplitude"""
311 if not hasattr(self, "_S_n_f"):
312 if hasattr(self, "_I_data"):
313 if self._I_Type == "ASD":
314 S_n_f_sqrt = self._I_data[:, 1]
315 self._S_n_f = S_n_f_sqrt ** 2 / u.Hz
316 elif self._I_Type == "ENSD":
317 self._S_n_f = self._I_data[:, 1] / u.Hz
318 elif self._I_Type == "h":
319 self._S_n_f = self.h_n_f ** 2 / self.fT
320 else:
321 if not hasattr(self, "_sensitivitycurve"):
322 self.Init_PTA()
323 self._S_n_f = self._sensitivitycurve.S_eff
324 self._S_n_f = utils.make_quant(self._S_n_f, "1/Hz")
325 return self._S_n_f
327 @S_n_f.setter
328 def S_n_f(self, value):
329 self._S_n_f = value
331 @S_n_f.deleter
332 def S_n_f(self):
333 del self._S_n_f
335 @property
336 def f_opt(self):
337 """The optimal frequency of the instrument ie. the frequecy at the lowest strain noise PSD"""
338 self._f_opt = self.fT[np.argmin(self.S_n_f)]
339 return self._f_opt
341 def Load_NANOGrav_11yr_Params(self):
342 """Loads in NANOGrav 11yr data
344 Notes
345 -----
346 The file is in the form of observation times (T_obs) in the first column,
347 sky locations (phi,theta) in the second and third columns,
348 Individual Pulsar cadences and WN RMS (sigmas) in the fourth and fifth,
349 RN Amplitudes, and RN Alphas in the last two columns.
350 """
351 NANOGrav_11yr_params_filedirectory = os.path.join(
352 load_directory, "InstrumentFiles/NANOGrav/NANOGrav_11yr_params.txt"
353 )
354 self._NANOGrav_11yr_params = np.loadtxt(NANOGrav_11yr_params_filedirectory)
356 def Get_Param_Distributions(self, var_name, NG_11yr_idx):
357 """Gets the noise parameter values (sigma, Rn_amplitudes,RN alphas) and sky locations (phis, thetas)
358 and generates populated arrays from which distributions can be made. If no user values for a param are given,
359 it uses the NANOGrav 11yr parameters.
361 var_name : string
362 The name of the noise parameter to assign sampled parameters
363 NG_11yr_idx : int
364 Index of corresponding value in NANOGrav 11yr params
365 """
367 if not hasattr(self, var_name):
368 samp_var = self._NANOGrav_11yr_params[NG_11yr_idx]
369 if var_name == "phi":
370 # Add non-zero probability of picking 0 and 2pi
371 return np.append(samp_var, np.linspace(0.0, 2 * np.pi, self.nbins))
372 elif var_name == "theta":
373 # Add non-zero probability of picking 0 and pi
374 return np.append(samp_var, np.linspace(0.0, np.pi, self.nbins))
375 elif var_name == "rn_amp":
376 return np.append(
377 samp_var,
378 np.logspace(
379 min(np.log10(samp_var)),
380 max(np.log10(samp_var)),
381 self.nbins,
382 ),
383 )
384 else:
385 return np.append(
386 samp_var, np.linspace(min(samp_var), max(samp_var), self.nbins)
387 )
388 else:
389 var = getattr(self, var_name)
390 if isinstance(var, u.Quantity):
391 var = var.value
392 if isinstance(var, (list, np.ndarray)):
393 if self.var_dict[var_name]["sampled"] == False:
394 if len(var) == self.n_p:
395 return var
396 elif len(var) == 1:
397 return np.ones(self.n_p) * var
398 else:
399 if self.var_dict["n_p"]["sampled"] == True:
400 unique_vals = np.unique(var)
401 if len(unique_vals) == 1:
402 return unique_vals[0]
403 else:
404 if var_name == "rn_amp":
405 return np.append(
406 var,
407 np.logspace(
408 min(np.log10(unique_vals)),
409 max(np.log10(unique_vals)),
410 self.nbins,
411 ),
412 )
413 else:
414 return np.append(
415 var,
416 np.linspace(
417 min(unique_vals),
418 max(unique_vals),
419 self.nbins,
420 ),
421 )
422 else:
423 raise ValueError(
424 "{} must be a single value, or the same length as n_p: {}".format(
425 var_name, self.n_p
426 )
427 )
428 else:
429 if len(var) == 2:
430 # Uniformly sample in logspace
431 if var_name == "rn_amp":
432 samp_var = np.random.uniform(
433 np.log10(var[0]), np.log10(var[1]), size=self.n_p
434 )
435 else:
436 # Uniformly Sample in linspace
437 samp_var = np.random.uniform(var[0], var[1], size=self.n_p)
438 elif len(var) == self.n_p:
439 samp_var = var
440 else:
441 raise ValueError(
442 "To sample {}, it must be either [min,max], or an array of individual pulsar {} of length n_p: {}".format(
443 var_name, var_name, self.n_p
444 )
445 )
447 if var_name == "rn_amp":
448 add_var = np.logspace(min(samp_var), max(samp_var), self.nbins)
449 samp_var = 10 ** samp_var
450 return np.append(samp_var, add_var)
451 else:
452 return np.append(
453 samp_var,
454 np.linspace(min(samp_var), max(samp_var), self.nbins),
455 )
456 else:
457 if var_name in self.var_dict.keys():
458 if self.var_dict[var_name]["sampled"] == False:
459 return np.ones(self.n_p) * var
460 else:
461 self.var_dict[var_name]["sampled"] == True
462 samp_var = self._NANOGrav_11yr_params[NG_11yr_idx]
463 if var_name == "phi":
464 # Add non-zero probability of picking 0 and 2pi
465 return np.append(
466 samp_var, np.linspace(0.0, 2 * np.pi, self.nbins)
467 )
468 elif var_name == "theta":
469 # Add non-zero probability of picking 0 and pi
470 return np.append(samp_var, np.linspace(0.0, np.pi, self.nbins))
471 elif var_name == "rn_amp":
472 return np.append(
473 samp_var,
474 np.logspace(
475 min(np.log10(samp_var)),
476 max(np.log10(samp_var)),
477 self.nbins,
478 ),
479 )
480 else:
481 return np.append(
482 samp_var,
483 np.linspace(min(samp_var), max(samp_var), self.nbins),
484 )
486 def Get_Sample_Draws(self, var_name, num_draws):
487 """For observation times, all noise parameters (sigma, Rn_amplitudes,RN alphas), cadence, and sky locations (phis, thetas),
488 uses the individual parameter value ranges and generates distributions from which to draw new parameters.
490 Notes
491 -----
492 To draw from the generated distributions, one does draws = self._distribution.rvs(size=sample_size)
493 """
494 var_list = ["T_obs", "phi", "theta", "cadence", "sigma", "rn_amp", "rn_alpha"]
496 for i in range(len(var_list)):
497 if var_name == var_list[i]:
498 NG_11yr_idx = i
500 samp_var = self.Get_Param_Distributions(var_name, NG_11yr_idx)
501 if isinstance(samp_var, (list, np.ndarray)):
502 if len(samp_var) > 1:
503 if var_name in ["rn_amp"]:
504 var_hist = np.histogram(
505 samp_var,
506 bins=np.logspace(
507 min(np.log10(samp_var)), max(np.log10(samp_var)), self.nbins
508 ),
509 density=True,
510 )
511 else:
512 var_hist = np.histogram(samp_var, bins=self.nbins, density=True)
513 var_dist = stats.rv_histogram(var_hist)
514 return var_dist.rvs(size=num_draws)
515 else:
516 return np.ones(num_draws) * samp_var[0]
517 else:
518 return np.ones(num_draws) * samp_var
520 def Init_PTA(self):
521 """Initializes a PTA in hasasia
523 Notes
524 -----
525 Assigns pulsar parameters based on what the initial values were given per parameter.
526 If necessary parameters are left unassigned, it uses 11yr values for n_p <= 34, and samples from the 11yr parameters if n_p > 34
527 If a range of values were given for a parameter, the per pulsar parameters are drawn from a uniform distribution
528 assigns the new pulsar parameters to the corresponding PTA class parameter.
529 See Hazboun, Romano, Smith (2019) <https://arxiv.org/abs/1907.04341> for details
531 """
532 if not hasattr(self, "_NANOGrav_11yr_params"):
533 self.Load_NANOGrav_11yr_Params()
535 [
536 NG_T_obs,
537 NG_phis,
538 NG_thetas,
539 NG_cadences,
540 NG_sigmas,
541 NG_rn_amps,
542 NG_rn_alphas,
543 ] = self._NANOGrav_11yr_params
544 var_list = ["T_obs", "phi", "theta", "cadence", "sigma", "rn_amp", "rn_alpha"]
546 for i, var in enumerate(var_list):
547 if var in self.var_dict.keys():
548 if self.var_dict[var]["sampled"] == True:
549 setattr(self, var, self.Get_Sample_Draws(var, self.n_p))
550 else:
551 if self.var_dict["n_p"]["sampled"] == True:
552 prev_var = getattr(self, var)
553 if isinstance(prev_var, u.Quantity):
554 prev_var = prev_var.value
555 if isinstance(prev_var, (list, np.ndarray)):
556 if len(prev_var) > self.n_p:
557 setattr(self, var, prev_var[: self.n_p])
558 elif len(prev_var) < self.n_p:
559 n_added_p = self.n_p - len(prev_var)
560 var_draw = self.Get_Sample_Draws(var, n_added_p)
561 setattr(
562 self,
563 var,
564 np.append(prev_var, var_draw),
565 )
566 else:
567 pass
568 else:
569 # Constant values for all pulsars
570 setattr(self, var, self.Get_Param_Distributions(var, i))
571 else:
572 # Constant values for all pulsars
573 setattr(self, var, self.Get_Param_Distributions(var, i))
574 else:
575 # Assign/sample values for values needed to make a sensitivity curve
576 if var in ["T_obs", "phi", "theta", "cadence", "sigma"]:
577 # 34 pulsars in the 11yr dataset (ie. len(phis))
578 if self.use_11yr:
579 if self.n_p <= len(self._NANOGrav_11yr_params[i]):
580 setattr(
581 self, var, self._NANOGrav_11yr_params[i][: self.n_p]
582 )
583 else:
584 n_added_p = self.n_p - len(self._NANOGrav_11yr_params[i])
585 var_draw = self.Get_Sample_Draws(var, n_added_p)
586 setattr(
587 self,
588 var,
589 np.append(self._NANOGrav_11yr_params[i], var_draw),
590 )
591 else:
592 setattr(self, var, self.Get_Sample_Draws(var, self.n_p))
593 else:
594 if self.use_rn:
595 if self.use_11yr:
596 if self.n_p <= len(self._NANOGrav_11yr_params[i]):
597 setattr(
598 self, var, self._NANOGrav_11yr_params[i][: self.n_p]
599 )
600 else:
601 n_added_p = self.n_p - len(
602 self._NANOGrav_11yr_params[i]
603 )
604 var_draw = self.Get_Sample_Draws(var, n_added_p)
605 setattr(
606 self,
607 var,
608 np.append(self._NANOGrav_11yr_params[i], var_draw),
609 )
610 else:
611 setattr(self, var, self.Get_Sample_Draws(var, self.n_p))
613 if hasattr(self, "rn_amp"):
614 if hasattr(self, "sb_amp"):
615 psrs = hassim.sim_pta(
616 timespan=self.T_obs.value,
617 cad=self.cadence.value,
618 sigma=self.sigma.value,
619 phi=self.phi,
620 theta=self.theta,
621 Npsrs=self.n_p,
622 A_rn=self.rn_amp,
623 alpha=self.rn_alpha,
624 A_gwb=self.sb_amp,
625 freqs=self.fT.value,
626 )
627 else:
628 psrs = hassim.sim_pta(
629 timespan=self.T_obs.value,
630 cad=self.cadence.value,
631 sigma=self.sigma.value,
632 phi=self.phi,
633 theta=self.theta,
634 Npsrs=self.n_p,
635 A_rn=self.rn_amp,
636 alpha=self.rn_alpha,
637 freqs=self.fT.value,
638 )
639 elif hasattr(self, "sb_amp"):
640 if not hasattr(self, "sb_alpha"):
641 self.sb_alpha = -2 / 3.0
642 # Make a set of psrs with the same parameters with a sb as red noise
643 psrs = hassim.sim_pta(
644 timespan=self.T_obs.value,
645 cad=self.cadence.value,
646 sigma=self.sigma.value,
647 phi=self.phi,
648 theta=self.theta,
649 Npsrs=self.n_p,
650 A_rn=self.sb_amp,
651 alpha=self.sb_alpha,
652 freqs=self.fT.value,
653 )
654 else:
655 psrs = hassim.sim_pta(
656 timespan=self.T_obs.value,
657 cad=self.cadence.value,
658 sigma=self.sigma.value,
659 phi=self.phi,
660 theta=self.theta,
661 Npsrs=self.n_p,
662 freqs=self.fT.value,
663 )
664 # Turn of sampling for an initialized PTA (except if sampling n_p)
665 for key, value in self.var_dict.items():
666 if key == "n_p":
667 pass
668 else:
669 value["sampled"] = False
670 # Get Spectra of pulsars
671 spectra = []
672 for p in psrs:
673 sp = hassens.Spectrum(p, freqs=self.fT.value)
674 spectra.append(sp)
676 self._sensitivitycurve = hassens.DeterSensitivityCurve(spectra)
679class Interferometer:
680 r"""
681 Class to make an interferometer
683 Parameters
684 ----------
686 name : string
687 name of the instrument
689 T_obs : float
690 the observation time of the PTA in [years]
692 load_location : string, optional
693 If you want to load an instrument curve from a file, it's the file path
694 I_type : string, {'E','A','h'}
695 Sets the type of input data.
696 'E' is the effective strain spectral density $S_{n}(f)$ ('ENSD'),
697 'A' is the amplitude spectral density, $\sqrt{S_{n}(f)}$ ('ASD'),
698 'h' is the characteristic strain $h_{n}(f)$ ('h')
699 f_low : float, optional
700 Assigned lowest frequency of instrument (default is assigned in particular child classes)
701 f_high : float, optional
702 Assigned highest frequency of instrument (default is assigned in particular child classes)
703 nfreqs : int, optional
704 Number of frequencies in logspace the sensitivity is calculated (default is 1e3)
706 """
708 def __init__(self, name, T_obs, **kwargs):
709 self.name = name
710 self.T_obs = T_obs
711 for keys, value in kwargs.items():
712 if keys == "load_location":
713 self.load_location = value
714 elif keys == "I_type":
715 self.I_type = value
716 elif keys == "f_low":
717 self.f_low = utils.make_quant(value, "Hz")
718 elif keys == "f_high":
719 self.f_high = utils.make_quant(value, "Hz")
720 elif keys == "nfreqs":
721 self.nfreqs = value
723 if hasattr(self, "load_location"):
724 Load_Data(self)
726 @property
727 def T_obs(self):
728 return self._T_obs
730 @T_obs.setter
731 def T_obs(self, value):
732 self.var_dict = ["T_obs", value]
733 if not isinstance(self._return_value, u.Quantity):
734 self._return_value = utils.make_quant(self._return_value, "yr")
735 self._T_obs = self._return_value
737 @property
738 def var_dict(self):
739 return self._var_dict
741 @var_dict.setter
742 def var_dict(self, value):
743 utils.Get_Var_Dict(self, value)
745 @property
746 def fT(self):
747 if not hasattr(self, "_fT"):
748 if hasattr(self, "_I_data"):
749 self._fT = self._I_data[:, 0] * u.Hz
750 if isinstance(self, SpaceBased):
751 self.Set_T_Function_Type()
752 if isinstance(self, GroundBased):
753 self._fT = (
754 np.logspace(
755 np.log10(self.f_low.value),
756 np.log10(self.f_high.value),
757 self.nfreqs,
758 )
759 * u.Hz
760 )
761 return self._fT
763 @fT.setter
764 def fT(self, value):
765 self._fT = value
767 @fT.deleter
768 def fT(self):
769 del self._fT
771 @property
772 def f_opt(self):
773 """The optimal frequency of the instrument ie. the frequecy at the lowest strain noise PSD"""
774 self._f_opt = self.fT[np.argmin(self.S_n_f)]
775 return self._f_opt
777 @property
778 def P_n_f(self):
779 """Strain power sensitivity. """
780 raise NotImplementedError(
781 "Power Spectral Density method must be defined inside SpaceBased or GroundBased classes."
782 )
784 @property
785 def S_n_f(self):
786 """Effective Noise Power Specral Density"""
787 if not hasattr(self, "_S_n_f"):
788 if hasattr(self, "_I_data"):
789 if self._I_Type == "ASD":
790 S_n_f_sqrt = self._I_data[:, 1]
791 self._S_n_f = S_n_f_sqrt ** 2 / u.Hz
792 elif self._I_Type == "ENSD":
793 self._S_n_f = self._I_data[:, 1] / u.Hz
794 elif self._I_Type == "h":
795 self._S_n_f = self.h_n_f ** 2 / self.fT
796 else:
797 raise NotImplementedError(
798 "Effective Noise Power Spectral Density method must be defined inside SpaceBased or GroundBased classes."
799 )
800 return self._S_n_f
802 @S_n_f.deleter
803 def S_n_f(self):
804 del self._S_n_f
806 @property
807 def h_n_f(self):
808 """Characteristic Strain/effective strain noise amplitude"""
809 if not hasattr(self, "_h_n_f"):
810 if hasattr(self, "_I_data") and self._I_Type == "h":
811 self._h_n_f = self._I_data[:, 1]
812 else:
813 self._h_n_f = np.sqrt(self.fT * self.S_n_f)
814 return self._h_n_f
816 @h_n_f.deleter
817 def h_n_f(self):
818 del self._h_n_f
821class GroundBased(Interferometer):
822 """
823 Class to make a Ground Based Instrument using the Interferometer base class
825 Parameters
826 ----------
827 noise_dict : dictionary, optional
828 A nested noise dictionary that has the main variable parameter name(s) in the top level,
829 the next level of the dictionary contains the subparameter variable name(s) and the desired value
830 to which the subparameter will be changed. The subparameter value can also be an array/list of the
831 [value,min,max] if one wishes to vary the instrument over then min/max range.
833 """
835 def __init__(self, name, T_obs, **kwargs):
836 super().__init__(name, T_obs, **kwargs)
838 for keys, value in kwargs.items():
839 if keys == "noise_dict":
840 if isinstance(value, dict):
841 self.noise_dict = value
842 else:
843 raise ValueError(keys + " must be a dictionary of noise sources.")
845 if not hasattr(self, "nfreqs"):
846 self.nfreqs = int(1e3)
847 if not hasattr(self, "f_low"):
848 self.f_low = 1.0 * u.Hz
849 if not hasattr(self, "f_high"):
850 self.f_high = 1e4 * u.Hz
852 if not hasattr(self, "load_location"):
853 if not hasattr(self, "noise_dict"):
854 self.Init_GroundBased()
855 else:
856 self.Set_Noise_Dict(self.noise_dict)
858 @property
859 def P_n_f(self):
860 """Power Spectral Density. """
861 err_mssg = "Currently we only calculate the Effective Noise Power Spectral Density for Ground Based detectors.\n"
862 err_mssg += "i.e. We do not separate the transfer function from the Power Spectral Density"
863 raise NotImplementedError(err_mssg)
865 @property
866 def S_n_f(self):
867 """Effective Noise Power Spectral Density"""
868 if not hasattr(self, "_S_n_f"):
869 if hasattr(self, "_I_data"):
870 if self._I_Type == "ASD":
871 S_n_f_sqrt = self._I_data[:, 1]
872 self._S_n_f = S_n_f_sqrt ** 2 / u.Hz
873 elif self._I_Type == "ENSD":
874 self._S_n_f = self._I_data[:, 1] / u.Hz
875 elif self._I_Type == "h":
876 self._S_n_f = self.h_n_f ** 2 / self.fT
877 else:
878 if not any(
879 hasattr(self, attr)
880 for attr in ["_noise_budget", "_ifo", "_base_inst"]
881 ):
882 self.Init_GroundBased()
883 self._S_n_f = (
884 self._noise_budget(self.fT.value, ifo=self._ifo).calc() / u.Hz
885 )
886 return self._S_n_f
888 @S_n_f.deleter
889 def S_n_f(self):
890 del self._S_n_f
892 def Init_GroundBased(self):
893 """Initialized the Ground Based detector in gwinc"""
894 base_inst = [
895 name for name in self.name.split() if name in gwinc.available_ifos()
896 ]
897 if len(base_inst) == 1:
898 self._base_inst = base_inst[0]
899 else:
900 print(
901 "You must select a base instrument model from ",
902 [model for model in gwinc.available_ifos()],
903 )
904 print(
905 "Setting base instrument to aLIGO. To change base instrument, include different model in class name and reinitialize."
906 )
907 self._base_inst = "aLIGO"
909 self._noise_budget, self._init_ifo, _, _ = gwinc.load_ifo(self._base_inst)
910 self._ifo = gwinc.precompIFO(self.fT.value, self._init_ifo)
912 def Set_Noise_Dict(self, noise_dict):
913 """Sets new values in the nested dictionary of variable noise values
915 Parameters
916 ----------
918 noise_dict : dictionary
919 A nested noise dictionary that has the main variable parameter name(s) in the top level,
920 the next level of the dictionary contains the subparameter variable name(s) and the desired value
921 to which the subparameter will be changed. The subparameter value can also be an array/list of the
922 [value,min,max] if one wishes to vary the instrument over then min/max range.
924 Examples
925 --------
926 obj.Set_Noise_Dict({'Infrastructure':{'Length':[3000,1000,5000],'Temp':500},'Laser':{'Wavelength':1e-5,'Power':130}})
928 """
929 if not hasattr(self, "_ifo"):
930 self.Init_GroundBased()
931 if isinstance(noise_dict, dict):
932 for base_noise, inner_noise_dict in noise_dict.items():
933 if base_noise in self._ifo.keys():
934 for sub_noise, sub_noise_val in inner_noise_dict.items():
935 if sub_noise in self._ifo[base_noise].keys():
936 if isinstance(sub_noise_val, dict):
937 for (
938 sub_sub_noise,
939 sub_sub_noise_val,
940 ) in sub_noise_val.items():
941 self.var_dict = [
942 base_noise
943 + " "
944 + sub_noise
945 + " "
946 + sub_sub_noise,
947 sub_sub_noise_val,
948 ]
949 setattr(
950 getattr(self._ifo, base_noise)[sub_noise],
951 sub_sub_noise,
952 self._return_value,
953 )
954 self._ifo = gwinc.precompIFO(
955 self.fT.value, self._ifo
956 )
957 else:
958 self.var_dict = [
959 base_noise + " " + sub_noise,
960 sub_noise_val,
961 ]
962 setattr(
963 getattr(self._ifo, base_noise),
964 sub_noise,
965 self._return_value,
966 )
967 self._ifo = gwinc.precompIFO(self.fT.value, self._ifo)
968 else:
969 raise ValueError(
970 sub_noise
971 + " is not a subparameter variable noise source.\
972 Try calling Get_Noise_Dict on your GroundBased object to find acceptable variables."
973 )
974 else:
975 err_mssg = (
976 base_noise
977 + " is not a valid parameter variable noise source.\n"
978 )
979 err_mssg += "Try calling Get_Noise_Dict on your GroundBased object to find acceptable variables."
980 raise ValueError(err_mssg)
981 else:
982 raise ValueError("Input must be a dictionary of noise sources.")
984 def Get_Noise_Dict(self):
985 """Gets and prints the available variable noises in the detector design"""
986 i = 0
987 for key_1, item_1 in self._ifo.items():
988 print(key_1, "Parameters:")
989 for key_2, item_2 in item_1.items():
990 if isinstance(item_2, np.ndarray):
991 i += 1
992 print(" ", key_2, ": array of shape", item_2.shape)
993 elif isinstance(item_2, list):
994 i += 1
995 print(" ", key_2, ": array of shape", len(item_2))
996 elif isinstance(item_2, (int, float)):
997 i += 1
998 print(" ", key_2, ":", item_2)
999 elif isinstance(item_2, gwinc.struct.Struct):
1000 print(" ", key_2, "Subparameters:")
1001 for key_3, item_3 in item_2.items():
1002 if isinstance(item_3, np.ndarray):
1003 i += 1
1004 print(
1005 " ", " ", key_3, ": array of shape", item_3.shape
1006 )
1007 elif isinstance(item_3, list):
1008 i += 1
1009 print(
1010 " ", " ", key_3, ": array of shape", len(item_3)
1011 )
1012 elif isinstance(item_3, (int, float)):
1013 i += 1
1014 print(" ", " ", key_3, ":", item_3)
1015 else:
1016 i += 1
1017 print(" ", key_2, ":", item_2)
1019 print(" ")
1020 print("Number of Variables: ", i)
1023class SpaceBased(Interferometer):
1024 """
1025 Class to make a Space Based Instrument using the Interferometer base class
1027 Parameters
1028 ----------
1029 L : float
1030 the armlength the of detector in [meters]
1031 A_acc : float
1032 the Amplitude of the Acceleration Noise in [meters/second^2]
1033 f_acc_break_low : float
1034 the lower break frequency of the acceleration noise in [Hz]
1035 f_acc_break_high : float
1036 the higher break frequency of the acceleration noise in [Hz]
1037 A_IFO : float
1038 the amplitude of the interferometer
1040 T_type : string, {'N','A'}
1041 Picks the transfer function generation method
1042 'N' uses the numerically approximated method in Robson, Cornish, and Liu, 2019
1043 'A' uses the analytic fit in Larson, Hiscock, and Hellings, 2000
1044 Background : Boolean
1045 Add in a Galactic Binary Confusion Noise
1047 """
1049 def __init__(self, name, T_obs, *args, **kwargs):
1050 super().__init__(name, T_obs, **kwargs)
1052 for keys, value in kwargs.items():
1053 if keys == "T_type":
1054 self.T_type = value
1055 elif keys == "Background":
1056 self.Background = value
1058 if not hasattr(self, "nfreqs"):
1059 self.nfreqs = int(1e3)
1060 if not hasattr(self, "f_low"):
1061 self.f_low = 1e-5 * u.Hz
1062 if not hasattr(self, "f_high"):
1063 self.f_high = 1.0 * u.Hz
1064 if not hasattr(self, "Background"):
1065 self.Background = False
1067 if len(args) != 0:
1068 [L, A_acc, f_acc_break_low, f_acc_break_high, A_IFO, f_IFO_break] = args
1069 self.L = L
1070 self.A_acc = A_acc
1071 self.f_acc_break_low = f_acc_break_low
1072 self.f_acc_break_high = f_acc_break_high
1073 self.A_IFO = A_IFO
1074 self.f_IFO_break = f_IFO_break
1076 if not hasattr(self, "load_location"):
1077 if not hasattr(self, "T_type"):
1078 self.T_type = "N"
1079 self.Set_T_Function_Type()
1081 @property
1082 def L(self):
1083 return self._L
1085 @L.setter
1086 def L(self, value):
1087 self.var_dict = ["L", value]
1088 if not isinstance(self._return_value, u.Quantity):
1089 self._return_value = utils.make_quant(self._return_value, "m")
1090 self._L = self._return_value
1092 @property
1093 def A_acc(self):
1094 return self._A_acc
1096 @A_acc.setter
1097 def A_acc(self, value):
1098 self.var_dict = ["A_acc", value]
1099 if not isinstance(self._return_value, u.Quantity):
1100 self._return_value = utils.make_quant(self._return_value, "m/s2")
1101 self._A_acc = self._return_value
1103 @property
1104 def f_acc_break_low(self):
1105 return self._f_acc_break_low
1107 @f_acc_break_low.setter
1108 def f_acc_break_low(self, value):
1109 self.var_dict = ["f_acc_break_low", value]
1110 if not isinstance(self._return_value, u.Quantity):
1111 self._return_value = utils.make_quant(self._return_value, "Hz")
1112 self._f_acc_break_low = self._return_value
1114 @property
1115 def f_acc_break_high(self):
1116 return self._f_acc_break_high
1118 @f_acc_break_high.setter
1119 def f_acc_break_high(self, value):
1120 self.var_dict = ["f_acc_break_high", value]
1121 if not isinstance(self._return_value, u.Quantity):
1122 self._return_value = utils.make_quant(self._return_value, "Hz")
1123 self._f_acc_break_high = self._return_value
1125 @property
1126 def A_IFO(self):
1127 return self._A_IFO
1129 @A_IFO.setter
1130 def A_IFO(self, value):
1131 self.var_dict = ["A_IFO", value]
1132 if not isinstance(self._return_value, u.Quantity):
1133 self._return_value = utils.make_quant(self._return_value, "m")
1134 self._A_IFO = self._return_value
1136 @property
1137 def f_IFO_break(self):
1138 return self._f_IFO_break
1140 @f_IFO_break.setter
1141 def f_IFO_break(self, value):
1142 self.var_dict = ["f_IFO_break", value]
1143 if not isinstance(self._return_value, u.Quantity):
1144 self._return_value = utils.make_quant(self._return_value, "Hz")
1145 self._f_IFO_break = self._return_value
1147 @property
1148 def P_n_f(self):
1149 """Power Spectral Density"""
1150 if not hasattr(self, "_P_n_f"):
1151 if not hasattr(self, "_T_Function_Type"):
1152 self.Set_T_Function_Type()
1154 P_acc = (
1155 self.A_acc ** 2
1156 * (1 + (self.f_acc_break_low / self.fT) ** 2)
1157 * (1 + (self.fT / (self.f_acc_break_high)) ** 4)
1158 / (2 * np.pi * self.fT) ** 4
1159 ) # Acceleration Noise
1160 P_IMS = self.A_IFO ** 2 * (
1161 1 + (self.f_IFO_break / self.fT) ** 4
1162 ) # Displacement noise of the interferometric TM--to-TM
1164 f_trans = const.c / 2 / np.pi / self.L # Transfer frequency
1165 self._P_n_f = (
1166 (P_IMS + 2 * (1 + np.cos(self.fT.value / f_trans.value) ** 2) * P_acc)
1167 / self.L ** 2
1168 / u.Hz
1169 )
1170 return self._P_n_f
1172 @P_n_f.deleter
1173 def P_n_f(self):
1174 del self._P_n_f
1176 @property
1177 def S_n_f(self):
1178 """Effective Noise Power Specral Density"""
1179 if not hasattr(self, "_S_n_f"):
1180 if hasattr(self, "_I_data"):
1181 if self._I_Type == "ASD":
1182 S_n_f_sqrt = self._I_data[:, 1]
1183 self._S_n_f = S_n_f_sqrt ** 2 / u.Hz
1184 elif self._I_Type == "ENSD":
1185 self._S_n_f = self._I_data[:, 1] / u.Hz
1186 elif self._I_Type == "h":
1187 self._S_n_f = self.h_n_f ** 2 / self.fT
1188 else:
1189 if self.Background:
1190 self._S_n_f = (
1191 self.P_n_f + self.Add_Background()
1192 ) / self.transferfunction ** 2
1193 else:
1194 self._S_n_f = self.P_n_f / self.transferfunction ** 2
1196 return self._S_n_f
1198 @S_n_f.deleter
1199 def S_n_f(self):
1200 del self._S_n_f
1202 def Load_Transfer_Function(self):
1203 # Numerical transfer function
1204 Numerical_Transfer_Function_filedirectory = os.path.join(
1205 load_directory, "NumericalTransferFunction/transfer.dat"
1206 )
1207 Numerical_Transfer_Function_data = np.loadtxt(
1208 Numerical_Transfer_Function_filedirectory
1209 )
1210 self._transferfunctiondata = Numerical_Transfer_Function_data
1212 def Get_Numeric_Transfer_Function(self):
1213 if not hasattr(self, "_transferfunctiondata"):
1214 self.Load_Transfer_Function()
1216 fc = const.c / (2 * self.L) # light round trip freq
1217 LISA_Transfer_Function_f = fc * self._transferfunctiondata[:, 0]
1219 idx_f_5 = np.abs(LISA_Transfer_Function_f - self.f_low).argmin()
1220 idx_f_1 = np.abs(LISA_Transfer_Function_f - self.f_high).argmin()
1222 # 3/10 is normalization 2/5sin(openingangle)
1223 # Some papers use 3/20, not summing over 2 independent low-freq data channels
1224 self.transferfunction = (
1225 np.sqrt(3 / 10) * self._transferfunctiondata[idx_f_5:idx_f_1, 1]
1226 )
1227 self.fT = LISA_Transfer_Function_f[idx_f_5:idx_f_1]
1229 def Get_Analytic_Transfer_Function(self, openingangle=None):
1230 # Response function approximation from Calculation described by Cornish, Robson, Liu 2019
1231 self.fT = (
1232 np.logspace(
1233 np.log10(self.f_low.value), np.log10(self.f_high.value), self.nfreqs
1234 )
1235 * u.Hz
1236 )
1237 f_L = const.c / 2 / np.pi / self.L # Transfer frequency
1238 # 3/10 is normalization 2/5sin(openingangle)
1239 if isinstance(openingangle, (int, float, u.Quantity)):
1240 R_f = (2 * np.sin(openingangle) ** 2) / 5 / (1 + 0.6 * (self.fT / f_L) ** 2)
1241 else:
1242 R_f = 3 / 10 / (1 + 0.6 * (self.fT / f_L) ** 2)
1243 self.transferfunction = np.sqrt(R_f)
1245 def Set_T_Function_Type(self):
1246 if self.T_type == "n" or self.T_type == "N":
1247 self._T_type = "numeric"
1248 elif self.T_type == "a" or self.T_type == "A":
1249 self._T_type = "analytic"
1250 else:
1251 print("\nYou can get the transfer function via 2 methods:")
1252 print(
1253 ' *To use the numerically approximated method in Robson, Cornish, and Liu, 2019, input "N".'
1254 )
1255 print(
1256 ' *To use the analytic fit in Larson, Hiscock, and Hellings, 2000, input "A".'
1257 )
1258 calc_type = input("Please select the calculation type: ")
1259 self.Set_T_Function_Type(calc_type)
1261 if self._T_type == "numeric":
1262 self.Get_Numeric_Transfer_Function()
1263 if self._T_type == "analytic":
1264 self.Get_Analytic_Transfer_Function()
1266 def Add_Background(self):
1267 """
1268 Galactic confusions noise parameters as a function of T_obs
1269 """
1270 A = 1.4e-44
1271 agam = 1100
1272 bgam = 3 / 10
1273 afk = 0.0016
1274 bfk = -2 / 9
1275 f_k = afk * self.T_obs ** bfk
1276 gamma = agam * self.T_obs ** bgam
1278 f = self.fT.value
1279 S_c_f = (
1280 A
1281 * (f ** (-7 / 3))
1282 * (1 + np.tanh(gamma.value * (f_k.value - f)))
1283 * (1 / u.Hz)
1284 ) # White Dwarf Background Noise
1285 return S_c_f
1288def Load_Data(detector):
1289 """
1290 Function to load in a file to initialize any detector.
1292 Parameters
1293 ----------
1294 detector : object
1295 Instance of a detector class
1297 """
1298 if not hasattr(detector, "I_type"):
1299 print("Is the data:")
1300 print(' *Effective Noise Spectral Density - "E"')
1301 print(' *Amplitude Spectral Density- "A"')
1302 print(' *Effective Strain - "h"')
1303 detector.I_type = input("Please enter one of the answers in quotations: ")
1304 Load_Data(detector)
1306 if detector.I_type == "E" or detector.I_type == "e":
1307 detector._I_Type = "ENSD"
1308 elif detector.I_type == "A" or detector.I_type == "a":
1309 detector._I_Type = "ASD"
1310 elif detector.I_type == "h" or detector.I_type == "H":
1311 detector._I_Type = "h"
1312 else:
1313 print("Is the data:")
1314 print(' *Effective Noise Spectral Density - "E"')
1315 print(' *Amplitude Spectral Density- "A"')
1316 print(' *Effective Strain - "h"')
1317 detector.I_type = input("Please enter one of the answers in quotations: ")
1318 Load_Data(detector)
1320 detector._I_data = np.loadtxt(detector.load_location)
1321 detector.fT = detector._I_data[:, 0] * u.Hz