Hide keyboard shortcuts

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, sys 

3import astropy.constants as const 

4import astropy.units as u 

5import scipy.interpolate as interp 

6import scipy.stats as stats 

7 

8from astropy.cosmology import z_at_value 

9from astropy.cosmology import WMAP9 as cosmo 

10 

11import gwent 

12from . import utils 

13 

14import hasasia.sensitivity as hassens 

15import hasasia.sim as hassim 

16 

17current_path = os.path.abspath(gwent.__path__[0]) 

18sys.path.insert(0, current_path + "/vendor/pygwinc_clone") 

19import gwinc 

20 

21load_directory = os.path.join(current_path, "LoadFiles/") 

22 

23 

24class PTA: 

25 r""" 

26 Class to make a PTA instrument using the methods of Hazboun, Romano, Smith (2019) <https://arxiv.org/abs/1907.04341> 

27 

28 Parameters 

29 ---------- 

30 

31 name : string 

32 name of the instrument 

33 

34 n_p : int 

35 the number of pulsars in the PTA 

36 

37 T_obs : float,array, list, optional 

38 the observation time of the PTA in [years] 

39 If T_obs is a float, the observation time is used as the observation time for all pulsars 

40 If T_obs is a list of length of n_p, the observation times are used as the corresponding pulsar observation times. 

41 If T_obs is a list of length 2, it is assumed the values are the minimum and maximum observation time values 

42 (ie. [min,max]) from which individual pulsar observation times are uniformly sampled. 

43 sigma : float, array, list, optional 

44 the rms error on the pulsar TOAs in [sec] 

45 If sigma is a float, the given rms error is used for all pulsars. 

46 If sigma is a list of length of n_p, the amplitudes are used as the corresponding pulsar rms error. 

47 If sigma is a list of length 2, it is assumed the values are the minimum and maximum rms errors values 

48 (ie. [min,max]) from which individual pulsar rms errors are uniformly sampled. 

49 cadence : float, array, list, optional 

50 How often the pulsars are observed in [num/year] 

51 If cadence is a float, the given cadence is used for all pulsars. 

52 If cadence is a list of length of n_p, the amplitudes are used as the corresponding pulsar cadence. 

53 If cadence is a list of length 2, it is assumed the values are the minimum and maximum cadence values 

54 (ie. [min,max]) from which individual pulsar cadences are uniformly sampled. 

55 sb_amp : float, optional 

56 Amplitude of the Stochastic background added as red noise 

57 sb_alpha : float, optional 

58 the Stochastic background power law, if empty and sb_amp is set, it is assumed to be -2/3 

59 rn_amp : array, list, optional 

60 Individual pulsar red noise amplitude. 

61 If rn_amp is a list of length of n_p, the amplitudes are used as the corresponding pulsar RN injection. 

62 If rn_amp is a list of length 2, it is assumed the values are the minimum and maximum RN amplitude values 

63 (ie. [min,max]) from which individual pulsar RN amplitudes are uniformly sampled. 

64 rn_alpha : array, list, optional 

65 Individual pulsar red noise alpha (power law spectral index). 

66 If rn_alpha is a list of length of n_p, the alpha indices are used as the corresponding pulsar RN injection. 

67 If rn_alpha is a list of length 2, it is assumed the values are the minimum and maximum RN alpha values 

68 (ie. [min,max]) from which individual pulsar RN alphas are uniformly sampled. 

69 phi : list, array, optional 

70 Individual pulsar longitude in ecliptic coordinates. 

71 If not defined, NANOGrav 11yr pulsar locations are used. 

72 If n_p > 34 (the number of pulsars in the 11yr dataset), 

73 it draws more pulsars from distributions based on the NANOGrav 11yr pulsars. 

74 theta : array, list, optional 

75 Individual pulsar colatitude in ecliptic coordinates. 

76 If not defined, NANOGrav 11yr pulsar locations are used. 

77 If n_p > 34 (the number of pulsars in the 11yr dataset), 

78 it draws more pulsars from distributions based on the NANOGrav 11yr pulsars. 

79 use_11yr : bool, optional 

80 Uses the NANOGrav 11yr noise as the individual pulsar noises, 

81 if n_p > 34 (the number of pulsars in the 11yr dataset), 

82 it draws more pulsars from distributions based on the NANOGrav 11yr pulsar noise 

83 use_rn : bool, optional 

84 If no rn_amp assigned, uses the NANOGrav 11yr noise as the individual pulsar RN noises, 

85 if n_p > 34 (the number of pulsars in the 11yr dataset), 

86 it draws more pulsars from distributions based on the NANOGrav 11yr pulsar noise 

87 load_location : string, optional 

88 If you want to load a PTA curve from a file, it's the file path 

89 I_type : string, {'E','A','h'} 

90 Sets the type of input data. 

91 'E' is the effective strain spectral density $S_{n}(f)$ ('ENSD'), 

92 'A' is the amplitude spectral density, $\sqrt{S_{n}(f)}$ ('ASD'), 

93 'h' is the characteristic strain $h_{n}(f)$ ('h') 

94 f_low : float, optional 

95 Assigned lowest frequency of PTA (default assigns 1/(5*T_obs)) 

96 f_high : float, optional 

97 Assigned highest frequency of PTA (default is Nyquist freq cadence/2) 

98 nfreqs : int, optional 

99 Number of frequencies in logspace the sensitivity is calculated 

100 nbins : int, optional 

101 Used to add values to every bin for sampled parameters. Default is 8 for smooth, non-zero distributions. 

102 Changing this could change distribution, so be wary, not sure how much it affects anything. 

103 """ 

104 

105 def __init__(self, name, *args, **kwargs): 

106 self.name = name 

107 if len(args) != 0: 

108 if len(args) == 1: 

109 self.n_p = args[0] 

110 else: 

111 raise ValueError("Too many args, not enough kwargs!") 

112 

113 for keys, value in kwargs.items(): 

114 if keys == "T_obs": 

115 self.T_obs = value 

116 elif keys == "sigma": 

117 self.sigma = value 

118 elif keys == "cadence": 

119 self.cadence = value 

120 elif keys == "rn_amp": 

121 self.rn_amp = value 

122 elif keys == "rn_alpha": 

123 self.rn_alpha = value 

124 elif keys == "phi": 

125 self.phi = value 

126 elif keys == "theta": 

127 self.theta = value 

128 elif keys == "sb_amp": 

129 self.sb_amp = value 

130 elif keys == "sb_alpha": 

131 self.sb_alpha = value 

132 elif keys == "use_11yr": 

133 self.use_11yr = value 

134 elif keys == "use_rn": 

135 self.use_rn = value 

136 elif keys == "load_location": 

137 self.load_location = value 

138 elif keys == "I_type": 

139 self.I_type = value 

140 elif keys == "f_low": 

141 self.f_low = utils.make_quant(value, "Hz") 

142 elif keys == "f_high": 

143 self.f_high = utils.make_quant(value, "Hz") 

144 elif keys == "nfreqs": 

145 self.nfreqs = value 

146 elif keys == "nbins": 

147 self.nbins = value 

148 else: 

149 raise ValueError("%s is not an accepted input option." % keys) 

150 

151 if not hasattr(self, "nfreqs"): 

152 self.nfreqs = int(1e3) 

153 if not hasattr(self, "nbins"): 

154 self.nbins = 8 

155 if hasattr(self, "load_location"): 

156 Load_Data(self) 

157 

158 if not hasattr(self, "use_11yr"): 

159 self.use_11yr = False 

160 if not hasattr(self, "use_rn"): 

161 self.use_rn = False 

162 

163 if hasattr(self, "f_low") and hasattr(self, "f_high"): 

164 self.fT = ( 

165 np.logspace( 

166 np.log10(self.f_low.value), np.log10(self.f_high.value), self.nfreqs 

167 ) 

168 * u.Hz 

169 ) 

170 

171 @property 

172 def n_p(self): 

173 return self._n_p 

174 

175 @n_p.setter 

176 def n_p(self, value): 

177 self.var_dict = ["n_p", value] 

178 self._n_p = self._return_value 

179 

180 @property 

181 def T_obs(self): 

182 return self._T_obs 

183 

184 @T_obs.setter 

185 def T_obs(self, value): 

186 self.var_dict = ["T_obs", value] 

187 if not isinstance(self._return_value, u.Quantity): 

188 self._return_value = utils.make_quant(self._return_value, "yr") 

189 self._T_obs = self._return_value 

190 

191 @property 

192 def cadence(self): 

193 return self._cadence 

194 

195 @cadence.setter 

196 def cadence(self, value): 

197 self.var_dict = ["cadence", value] 

198 if not isinstance(self._return_value, u.Quantity): 

199 self._return_value = utils.make_quant(self._return_value, "1/yr") 

200 self._cadence = self._return_value 

201 

202 @property 

203 def sigma(self): 

204 return self._sigma 

205 

206 @sigma.setter 

207 def sigma(self, value): 

208 self.var_dict = ["sigma", value] 

209 if not isinstance(self._return_value, u.Quantity): 

210 self._return_value = utils.make_quant(self._return_value, "s") 

211 self._sigma = self._return_value 

212 

213 @property 

214 def phi(self): 

215 return self._phi 

216 

217 @phi.setter 

218 def phi(self, value): 

219 self.var_dict = ["phi", value] 

220 self._phi = self._return_value 

221 

222 @property 

223 def theta(self): 

224 return self._theta 

225 

226 @theta.setter 

227 def theta(self, value): 

228 self.var_dict = ["theta", value] 

229 self._theta = self._return_value 

230 

231 @property 

232 def rn_amp(self): 

233 return self._rn_amp 

234 

235 @rn_amp.setter 

236 def rn_amp(self, value): 

237 self.var_dict = ["rn_amp", value] 

238 self._rn_amp = self._return_value 

239 

240 @property 

241 def rn_alpha(self): 

242 return self._rn_alpha 

243 

244 @rn_alpha.setter 

245 def rn_alpha(self, value): 

246 self.var_dict = ["rn_alpha", value] 

247 self._rn_alpha = self._return_value 

248 

249 @property 

250 def var_dict(self): 

251 return self._var_dict 

252 

253 @var_dict.setter 

254 def var_dict(self, value): 

255 utils.Get_Var_Dict(self, value) 

256 

257 @property 

258 def fT(self): 

259 if not hasattr(self, "_fT"): 

260 # frequency sampled from 1/observation time to nyquist frequency (c/2) 

261 # 5 is the default value for now (from Hazboun et al. 2019) 

262 if not hasattr(self, "_T_obs") or not hasattr(self, "_cadence"): 

263 self.Init_PTA() 

264 T_obs_sec = np.max(self.T_obs).to("s").value 

265 cadence_sec = np.max(self.cadence).to("1/s").value 

266 self._fT = ( 

267 np.logspace( 

268 np.log10(1 / (5 * T_obs_sec)), 

269 np.log10(cadence_sec / 2), 

270 self.nfreqs, 

271 ) 

272 * u.Hz 

273 ) 

274 return self._fT 

275 

276 @fT.setter 

277 def fT(self, value): 

278 self._fT = value 

279 

280 @fT.deleter 

281 def fT(self): 

282 del self._fT 

283 

284 @property 

285 def h_n_f(self): 

286 """Effective Strain Noise Amplitude""" 

287 if not hasattr(self, "_h_n_f"): 

288 if hasattr(self, "_I_data"): 

289 if self._I_Type == "h": 

290 self._h_n_f = self._I_data[:, 1] 

291 elif self._I_Type == "ENSD": 

292 self._h_n_f = np.sqrt(self.S_n_f * self.fT) 

293 elif self._I_Type == "ASD": 

294 S_n_f_sqrt = self._I_data[:, 1] 

295 self._h_n_f = S_n_f_sqrt * np.sqrt(self.fT.value) 

296 else: 

297 if not hasattr(self, "_sensitivitycurve"): 

298 self.Init_PTA() 

299 self._h_n_f = self._sensitivitycurve.h_c 

300 return self._h_n_f 

301 

302 @h_n_f.setter 

303 def h_n_f(self, value): 

304 self._h_n_f = value 

305 

306 @h_n_f.deleter 

307 def h_n_f(self): 

308 del self._h_n_f 

309 

310 @property 

311 def S_n_f(self): 

312 """Effective noise power amplitude""" 

313 if not hasattr(self, "_S_n_f"): 

314 if hasattr(self, "_I_data"): 

315 if self._I_Type == "ASD": 

316 S_n_f_sqrt = self._I_data[:, 1] 

317 self._S_n_f = S_n_f_sqrt ** 2 / u.Hz 

318 elif self._I_Type == "ENSD": 

319 self._S_n_f = self._I_data[:, 1] / u.Hz 

320 elif self._I_Type == "h": 

321 self._S_n_f = self.h_n_f ** 2 / self.fT 

322 else: 

323 if not hasattr(self, "_sensitivitycurve"): 

324 self.Init_PTA() 

325 self._S_n_f = self._sensitivitycurve.S_eff 

326 self._S_n_f = utils.make_quant(self._S_n_f, "1/Hz") 

327 return self._S_n_f 

328 

329 @S_n_f.setter 

330 def S_n_f(self, value): 

331 self._S_n_f = value 

332 

333 @S_n_f.deleter 

334 def S_n_f(self): 

335 del self._S_n_f 

336 

337 @property 

338 def f_opt(self): 

339 """The optimal frequency of the instrument ie. the frequecy at the lowest strain noise PSD""" 

340 self._f_opt = self.fT[np.argmin(self.S_n_f)] 

341 return self._f_opt 

342 

343 def Load_NANOGrav_11yr_Params(self): 

344 """Loads in NANOGrav 11yr data 

345 

346 Notes 

347 ----- 

348 The file is in the form of observation times (T_obs) in the first column, 

349 sky locations (phi,theta) in the second and third columns, 

350 Individual Pulsar cadences and WN RMS (sigmas) in the fourth and fifth, 

351 RN Amplitudes, and RN Alphas in the last two columns. 

352 """ 

353 NANOGrav_11yr_params_filedirectory = os.path.join( 

354 load_directory, "InstrumentFiles/NANOGrav/NANOGrav_11yr_params.txt" 

355 ) 

356 self._NANOGrav_11yr_params = np.loadtxt(NANOGrav_11yr_params_filedirectory) 

357 

358 def Get_Param_Distributions(self, var_name, NG_11yr_idx): 

359 """Gets the noise parameter values (sigma, Rn_amplitudes,RN alphas) and sky locations (phis, thetas) 

360 and generates populated arrays from which distributions can be made. If no user values for a param are given, 

361 it uses the NANOGrav 11yr parameters. 

362 

363 var_name : string 

364 The name of the noise parameter to assign sampled parameters 

365 NG_11yr_idx : int 

366 Index of corresponding value in NANOGrav 11yr params 

367 """ 

368 

369 if not hasattr(self, var_name): 

370 samp_var = self._NANOGrav_11yr_params[NG_11yr_idx] 

371 if var_name == "phi": 

372 # Add non-zero probability of picking 0 and 2pi 

373 return np.append(samp_var, np.linspace(0.0, 2 * np.pi, self.nbins)) 

374 elif var_name == "theta": 

375 # Add non-zero probability of picking 0 and pi 

376 return np.append(samp_var, np.linspace(0.0, np.pi, self.nbins)) 

377 elif var_name == "rn_amp": 

378 return np.append( 

379 samp_var, 

380 np.logspace( 

381 min(np.log10(samp_var)), 

382 max(np.log10(samp_var)), 

383 self.nbins, 

384 ), 

385 ) 

386 else: 

387 return np.append( 

388 samp_var, np.linspace(min(samp_var), max(samp_var), self.nbins) 

389 ) 

390 else: 

391 var = getattr(self, var_name) 

392 if isinstance(var, u.Quantity): 

393 var = var.value 

394 if isinstance(var, (list, np.ndarray)): 

395 if self.var_dict[var_name]["sampled"] == False: 

396 if len(var) == self.n_p: 

397 return var 

398 elif len(var) == 1: 

399 return np.ones(self.n_p) * var 

400 else: 

401 if self.var_dict["n_p"]["sampled"] == True: 

402 unique_vals = np.unique(var) 

403 if len(unique_vals) == 1: 

404 return unique_vals[0] 

405 else: 

406 if var_name == "rn_amp": 

407 return np.append( 

408 var, 

409 np.logspace( 

410 min(np.log10(unique_vals)), 

411 max(np.log10(unique_vals)), 

412 self.nbins, 

413 ), 

414 ) 

415 else: 

416 return np.append( 

417 var, 

418 np.linspace( 

419 min(unique_vals), 

420 max(unique_vals), 

421 self.nbins, 

422 ), 

423 ) 

424 else: 

425 raise ValueError( 

426 "{} must be a single value, or the same length as n_p: {}".format( 

427 var_name, self.n_p 

428 ) 

429 ) 

430 else: 

431 if len(var) == 2: 

432 # Uniformly sample in logspace 

433 if var_name == "rn_amp": 

434 samp_var = np.random.uniform( 

435 np.log10(var[0]), np.log10(var[1]), size=self.n_p 

436 ) 

437 else: 

438 # Uniformly Sample in linspace 

439 samp_var = np.random.uniform(var[0], var[1], size=self.n_p) 

440 elif len(var) == self.n_p: 

441 samp_var = var 

442 else: 

443 raise ValueError( 

444 "To sample {}, it must be either [min,max], or an array of individual pulsar {} of length n_p: {}".format( 

445 var_name, var_name, self.n_p 

446 ) 

447 ) 

448 

449 if var_name == "rn_amp": 

450 add_var = np.logspace(min(samp_var), max(samp_var), self.nbins) 

451 samp_var = 10 ** samp_var 

452 return np.append(samp_var, add_var) 

453 else: 

454 return np.append( 

455 samp_var, 

456 np.linspace(min(samp_var), max(samp_var), self.nbins), 

457 ) 

458 else: 

459 if var_name in self.var_dict.keys(): 

460 if self.var_dict[var_name]["sampled"] == False: 

461 return np.ones(self.n_p) * var 

462 else: 

463 self.var_dict[var_name]["sampled"] == True 

464 samp_var = self._NANOGrav_11yr_params[NG_11yr_idx] 

465 if var_name == "phi": 

466 # Add non-zero probability of picking 0 and 2pi 

467 return np.append( 

468 samp_var, np.linspace(0.0, 2 * np.pi, self.nbins) 

469 ) 

470 elif var_name == "theta": 

471 # Add non-zero probability of picking 0 and pi 

472 return np.append(samp_var, np.linspace(0.0, np.pi, self.nbins)) 

473 elif var_name == "rn_amp": 

474 return np.append( 

475 samp_var, 

476 np.logspace( 

477 min(np.log10(samp_var)), 

478 max(np.log10(samp_var)), 

479 self.nbins, 

480 ), 

481 ) 

482 else: 

483 return np.append( 

484 samp_var, 

485 np.linspace(min(samp_var), max(samp_var), self.nbins), 

486 ) 

487 

488 def Get_Sample_Draws(self, var_name, num_draws): 

489 """For observation times, all noise parameters (sigma, Rn_amplitudes,RN alphas), cadence, and sky locations (phis, thetas), 

490 uses the individual parameter value ranges and generates distributions from which to draw new parameters. 

491 

492 Notes 

493 ----- 

494 To draw from the generated distributions, one does draws = self._distribution.rvs(size=sample_size) 

495 """ 

496 var_list = ["T_obs", "phi", "theta", "cadence", "sigma", "rn_amp", "rn_alpha"] 

497 

498 for i in range(len(var_list)): 

499 if var_name == var_list[i]: 

500 NG_11yr_idx = i 

501 

502 samp_var = self.Get_Param_Distributions(var_name, NG_11yr_idx) 

503 if isinstance(samp_var, (list, np.ndarray)): 

504 if len(samp_var) > 1: 

505 if var_name in ["rn_amp"]: 

506 var_hist = np.histogram( 

507 samp_var, 

508 bins=np.logspace( 

509 min(np.log10(samp_var)), max(np.log10(samp_var)), self.nbins 

510 ), 

511 density=True, 

512 ) 

513 else: 

514 var_hist = np.histogram(samp_var, bins=self.nbins, density=True) 

515 var_dist = stats.rv_histogram(var_hist) 

516 return var_dist.rvs(size=num_draws) 

517 else: 

518 return np.ones(num_draws) * samp_var[0] 

519 else: 

520 return np.ones(num_draws) * samp_var 

521 

522 def Init_PTA(self): 

523 """Initializes a PTA in hasasia 

524 

525 Notes 

526 ----- 

527 Assigns pulsar parameters based on what the initial values were given per parameter. 

528 If necessary parameters are left unassigned, it uses 11yr values for n_p <= 34, and samples from the 11yr parameters if n_p > 34 

529 If a range of values were given for a parameter, the per pulsar parameters are drawn from a uniform distribution 

530 assigns the new pulsar parameters to the corresponding PTA class parameter. 

531 See Hazboun, Romano, Smith (2019) <https://arxiv.org/abs/1907.04341> for details 

532 

533 """ 

534 if not hasattr(self, "_NANOGrav_11yr_params"): 

535 self.Load_NANOGrav_11yr_Params() 

536 

537 [ 

538 NG_T_obs, 

539 NG_phis, 

540 NG_thetas, 

541 NG_cadences, 

542 NG_sigmas, 

543 NG_rn_amps, 

544 NG_rn_alphas, 

545 ] = self._NANOGrav_11yr_params 

546 var_list = ["T_obs", "phi", "theta", "cadence", "sigma", "rn_amp", "rn_alpha"] 

547 

548 for i, var in enumerate(var_list): 

549 if var in self.var_dict.keys(): 

550 if self.var_dict[var]["sampled"] == True: 

551 setattr(self, var, self.Get_Sample_Draws(var, self.n_p)) 

552 else: 

553 if self.var_dict["n_p"]["sampled"] == True: 

554 prev_var = getattr(self, var) 

555 if isinstance(prev_var, u.Quantity): 

556 prev_var = prev_var.value 

557 if isinstance(prev_var, (list, np.ndarray)): 

558 if len(prev_var) > self.n_p: 

559 setattr(self, var, prev_var[: self.n_p]) 

560 elif len(prev_var) < self.n_p: 

561 n_added_p = self.n_p - len(prev_var) 

562 var_draw = self.Get_Sample_Draws(var, n_added_p) 

563 setattr( 

564 self, 

565 var, 

566 np.append(prev_var, var_draw), 

567 ) 

568 else: 

569 pass 

570 else: 

571 # Constant values for all pulsars 

572 setattr(self, var, self.Get_Param_Distributions(var, i)) 

573 else: 

574 # Constant values for all pulsars 

575 setattr(self, var, self.Get_Param_Distributions(var, i)) 

576 else: 

577 # Assign/sample values for values needed to make a sensitivity curve 

578 if var in ["T_obs", "phi", "theta", "cadence", "sigma"]: 

579 # 34 pulsars in the 11yr dataset (ie. len(phis)) 

580 if self.use_11yr: 

581 if self.n_p <= len(self._NANOGrav_11yr_params[i]): 

582 setattr( 

583 self, var, self._NANOGrav_11yr_params[i][: self.n_p] 

584 ) 

585 else: 

586 n_added_p = self.n_p - len(self._NANOGrav_11yr_params[i]) 

587 var_draw = self.Get_Sample_Draws(var, n_added_p) 

588 setattr( 

589 self, 

590 var, 

591 np.append(self._NANOGrav_11yr_params[i], var_draw), 

592 ) 

593 else: 

594 setattr(self, var, self.Get_Sample_Draws(var, self.n_p)) 

595 else: 

596 if self.use_rn: 

597 if self.use_11yr: 

598 if self.n_p <= len(self._NANOGrav_11yr_params[i]): 

599 setattr( 

600 self, var, self._NANOGrav_11yr_params[i][: self.n_p] 

601 ) 

602 else: 

603 n_added_p = self.n_p - len( 

604 self._NANOGrav_11yr_params[i] 

605 ) 

606 var_draw = self.Get_Sample_Draws(var, n_added_p) 

607 setattr( 

608 self, 

609 var, 

610 np.append(self._NANOGrav_11yr_params[i], var_draw), 

611 ) 

612 else: 

613 setattr(self, var, self.Get_Sample_Draws(var, self.n_p)) 

614 

615 if hasattr(self, "rn_amp"): 

616 if hasattr(self, "sb_amp"): 

617 psrs = hassim.sim_pta( 

618 timespan=self.T_obs.value, 

619 cad=self.cadence.value, 

620 sigma=self.sigma.value, 

621 phi=self.phi, 

622 theta=self.theta, 

623 Npsrs=self.n_p, 

624 A_rn=self.rn_amp, 

625 alpha=self.rn_alpha, 

626 A_gwb=self.sb_amp, 

627 freqs=self.fT.value, 

628 ) 

629 else: 

630 psrs = hassim.sim_pta( 

631 timespan=self.T_obs.value, 

632 cad=self.cadence.value, 

633 sigma=self.sigma.value, 

634 phi=self.phi, 

635 theta=self.theta, 

636 Npsrs=self.n_p, 

637 A_rn=self.rn_amp, 

638 alpha=self.rn_alpha, 

639 freqs=self.fT.value, 

640 ) 

641 elif hasattr(self, "sb_amp"): 

642 if not hasattr(self, "sb_alpha"): 

643 self.sb_alpha = -2 / 3.0 

644 # Make a set of psrs with the same parameters with a sb as red noise 

645 psrs = hassim.sim_pta( 

646 timespan=self.T_obs.value, 

647 cad=self.cadence.value, 

648 sigma=self.sigma.value, 

649 phi=self.phi, 

650 theta=self.theta, 

651 Npsrs=self.n_p, 

652 A_rn=self.sb_amp, 

653 alpha=self.sb_alpha, 

654 freqs=self.fT.value, 

655 ) 

656 else: 

657 psrs = hassim.sim_pta( 

658 timespan=self.T_obs.value, 

659 cad=self.cadence.value, 

660 sigma=self.sigma.value, 

661 phi=self.phi, 

662 theta=self.theta, 

663 Npsrs=self.n_p, 

664 freqs=self.fT.value, 

665 ) 

666 # Turn of sampling for an initialized PTA (except if sampling n_p) 

667 for key, value in self.var_dict.items(): 

668 if key == "n_p": 

669 pass 

670 else: 

671 value["sampled"] = False 

672 # Get Spectra of pulsars 

673 spectra = [] 

674 for p in psrs: 

675 sp = hassens.Spectrum(p, freqs=self.fT.value) 

676 spectra.append(sp) 

677 

678 self._sensitivitycurve = hassens.DeterSensitivityCurve(spectra) 

679 

680 

681class Interferometer: 

682 r""" 

683 Class to make an interferometer 

684 

685 Parameters 

686 ---------- 

687 

688 name : string 

689 name of the instrument 

690 

691 T_obs : float 

692 the observation time of the PTA in [years] 

693 

694 load_location : string, optional 

695 If you want to load an instrument curve from a file, it's the file path 

696 I_type : string, {'E','A','h'} 

697 Sets the type of input data. 

698 'E' is the effective strain spectral density $S_{n}(f)$ ('ENSD'), 

699 'A' is the amplitude spectral density, $\sqrt{S_{n}(f)}$ ('ASD'), 

700 'h' is the characteristic strain $h_{n}(f)$ ('h') 

701 f_low : float, optional 

702 Assigned lowest frequency of instrument (default is assigned in particular child classes) 

703 f_high : float, optional 

704 Assigned highest frequency of instrument (default is assigned in particular child classes) 

705 nfreqs : int, optional 

706 Number of frequencies in logspace the sensitivity is calculated (default is 1e3) 

707 

708 """ 

709 

710 def __init__(self, name, T_obs, **kwargs): 

711 self.name = name 

712 self.T_obs = T_obs 

713 for keys, value in kwargs.items(): 

714 if keys == "load_location": 

715 self.load_location = value 

716 elif keys == "I_type": 

717 self.I_type = value 

718 elif keys == "f_low": 

719 self.f_low = utils.make_quant(value, "Hz") 

720 elif keys == "f_high": 

721 self.f_high = utils.make_quant(value, "Hz") 

722 elif keys == "nfreqs": 

723 self.nfreqs = value 

724 

725 if hasattr(self, "load_location"): 

726 Load_Data(self) 

727 

728 @property 

729 def T_obs(self): 

730 return self._T_obs 

731 

732 @T_obs.setter 

733 def T_obs(self, value): 

734 self.var_dict = ["T_obs", value] 

735 if not isinstance(self._return_value, u.Quantity): 

736 self._return_value = utils.make_quant(self._return_value, "yr") 

737 self._T_obs = self._return_value 

738 

739 @property 

740 def var_dict(self): 

741 return self._var_dict 

742 

743 @var_dict.setter 

744 def var_dict(self, value): 

745 utils.Get_Var_Dict(self, value) 

746 

747 @property 

748 def fT(self): 

749 if not hasattr(self, "_fT"): 

750 if hasattr(self, "_I_data"): 

751 self._fT = self._I_data[:, 0] * u.Hz 

752 if isinstance(self, SpaceBased): 

753 self.Set_T_Function_Type() 

754 if isinstance(self, GroundBased): 

755 self._fT = ( 

756 np.logspace( 

757 np.log10(self.f_low.value), 

758 np.log10(self.f_high.value), 

759 self.nfreqs, 

760 ) 

761 * u.Hz 

762 ) 

763 return self._fT 

764 

765 @fT.setter 

766 def fT(self, value): 

767 self._fT = value 

768 

769 @fT.deleter 

770 def fT(self): 

771 del self._fT 

772 

773 @property 

774 def f_opt(self): 

775 """The optimal frequency of the instrument ie. the frequecy at the lowest strain noise PSD""" 

776 self._f_opt = self.fT[np.argmin(self.S_n_f)] 

777 return self._f_opt 

778 

779 @property 

780 def P_n_f(self): 

781 """Strain power sensitivity. """ 

782 raise NotImplementedError( 

783 "Power Spectral Density method must be defined inside SpaceBased or GroundBased classes." 

784 ) 

785 

786 @property 

787 def S_n_f(self): 

788 """Effective Noise Power Specral Density""" 

789 if not hasattr(self, "_S_n_f"): 

790 if hasattr(self, "_I_data"): 

791 if self._I_Type == "ASD": 

792 S_n_f_sqrt = self._I_data[:, 1] 

793 self._S_n_f = S_n_f_sqrt ** 2 / u.Hz 

794 elif self._I_Type == "ENSD": 

795 self._S_n_f = self._I_data[:, 1] / u.Hz 

796 elif self._I_Type == "h": 

797 self._S_n_f = self.h_n_f ** 2 / self.fT 

798 else: 

799 raise NotImplementedError( 

800 "Effective Noise Power Spectral Density method must be defined inside SpaceBased or GroundBased classes." 

801 ) 

802 return self._S_n_f 

803 

804 @S_n_f.deleter 

805 def S_n_f(self): 

806 del self._S_n_f 

807 

808 @property 

809 def h_n_f(self): 

810 """Characteristic Strain/effective strain noise amplitude""" 

811 if not hasattr(self, "_h_n_f"): 

812 if hasattr(self, "_I_data") and self._I_Type == "h": 

813 self._h_n_f = self._I_data[:, 1] 

814 else: 

815 self._h_n_f = np.sqrt(self.fT * self.S_n_f) 

816 return self._h_n_f 

817 

818 @h_n_f.deleter 

819 def h_n_f(self): 

820 del self._h_n_f 

821 

822 

823class GroundBased(Interferometer): 

824 """ 

825 Class to make a Ground Based Instrument using the Interferometer base class 

826 

827 Parameters 

828 ---------- 

829 noise_dict : dictionary, optional 

830 A nested noise dictionary that has the main variable parameter name(s) in the top level, 

831 the next level of the dictionary contains the subparameter variable name(s) and the desired value 

832 to which the subparameter will be changed. The subparameter value can also be an array/list of the 

833 [value,min,max] if one wishes to vary the instrument over then min/max range. 

834 

835 """ 

836 

837 def __init__(self, name, T_obs, **kwargs): 

838 super().__init__(name, T_obs, **kwargs) 

839 

840 for keys, value in kwargs.items(): 

841 if keys == "noise_dict": 

842 if isinstance(value, dict): 

843 self.noise_dict = value 

844 else: 

845 raise ValueError(keys + " must be a dictionary of noise sources.") 

846 

847 if not hasattr(self, "nfreqs"): 

848 self.nfreqs = int(1e3) 

849 if not hasattr(self, "f_low"): 

850 self.f_low = 1.0 * u.Hz 

851 if not hasattr(self, "f_high"): 

852 self.f_high = 1e4 * u.Hz 

853 

854 if not hasattr(self, "load_location"): 

855 if not hasattr(self, "noise_dict"): 

856 self.Init_GroundBased() 

857 else: 

858 self.Set_Noise_Dict(self.noise_dict) 

859 

860 @property 

861 def P_n_f(self): 

862 """Power Spectral Density. """ 

863 err_mssg = "Currently we only calculate the Effective Noise Power Spectral Density for Ground Based detectors.\n" 

864 err_mssg += "i.e. We do not separate the transfer function from the Power Spectral Density" 

865 raise NotImplementedError(err_mssg) 

866 

867 @property 

868 def S_n_f(self): 

869 """Effective Noise Power Spectral Density""" 

870 if not hasattr(self, "_S_n_f"): 

871 if hasattr(self, "_I_data"): 

872 if self._I_Type == "ASD": 

873 S_n_f_sqrt = self._I_data[:, 1] 

874 self._S_n_f = S_n_f_sqrt ** 2 / u.Hz 

875 elif self._I_Type == "ENSD": 

876 self._S_n_f = self._I_data[:, 1] / u.Hz 

877 elif self._I_Type == "h": 

878 self._S_n_f = self.h_n_f ** 2 / self.fT 

879 else: 

880 if not any( 

881 hasattr(self, attr) 

882 for attr in ["_noise_budget", "_ifo", "_base_inst"] 

883 ): 

884 self.Init_GroundBased() 

885 self._S_n_f = ( 

886 self._noise_budget(self.fT.value, ifo=self._ifo).calc() / u.Hz 

887 ) 

888 return self._S_n_f 

889 

890 @S_n_f.deleter 

891 def S_n_f(self): 

892 del self._S_n_f 

893 

894 def Init_GroundBased(self): 

895 """Initialized the Ground Based detector in gwinc""" 

896 base_inst = [ 

897 name for name in self.name.split() if name in gwinc.available_ifos() 

898 ] 

899 if len(base_inst) == 1: 

900 self._base_inst = base_inst[0] 

901 else: 

902 print( 

903 "You must select a base instrument model from ", 

904 [model for model in gwinc.available_ifos()], 

905 ) 

906 print( 

907 "Setting base instrument to aLIGO. To change base instrument, include different model in class name and reinitialize." 

908 ) 

909 self._base_inst = "aLIGO" 

910 

911 self._noise_budget, self._init_ifo, _, _ = gwinc.load_ifo(self._base_inst) 

912 self._ifo = gwinc.precompIFO(self.fT.value, self._init_ifo) 

913 

914 def Set_Noise_Dict(self, noise_dict): 

915 """Sets new values in the nested dictionary of variable noise values 

916 

917 Parameters 

918 ---------- 

919 

920 noise_dict : dictionary 

921 A nested noise dictionary that has the main variable parameter name(s) in the top level, 

922 the next level of the dictionary contains the subparameter variable name(s) and the desired value 

923 to which the subparameter will be changed. The subparameter value can also be an array/list of the 

924 [value,min,max] if one wishes to vary the instrument over then min/max range. 

925 

926 Examples 

927 -------- 

928 obj.Set_Noise_Dict({'Infrastructure':{'Length':[3000,1000,5000],'Temp':500},'Laser':{'Wavelength':1e-5,'Power':130}}) 

929 

930 """ 

931 if not hasattr(self, "_ifo"): 

932 self.Init_GroundBased() 

933 if isinstance(noise_dict, dict): 

934 for base_noise, inner_noise_dict in noise_dict.items(): 

935 if base_noise in self._ifo.keys(): 

936 for sub_noise, sub_noise_val in inner_noise_dict.items(): 

937 if sub_noise in self._ifo[base_noise].keys(): 

938 if isinstance(sub_noise_val, dict): 

939 for ( 

940 sub_sub_noise, 

941 sub_sub_noise_val, 

942 ) in sub_noise_val.items(): 

943 self.var_dict = [ 

944 base_noise 

945 + " " 

946 + sub_noise 

947 + " " 

948 + sub_sub_noise, 

949 sub_sub_noise_val, 

950 ] 

951 setattr( 

952 getattr(self._ifo, base_noise)[sub_noise], 

953 sub_sub_noise, 

954 self._return_value, 

955 ) 

956 self._ifo = gwinc.precompIFO( 

957 self.fT.value, self._ifo 

958 ) 

959 else: 

960 self.var_dict = [ 

961 base_noise + " " + sub_noise, 

962 sub_noise_val, 

963 ] 

964 setattr( 

965 getattr(self._ifo, base_noise), 

966 sub_noise, 

967 self._return_value, 

968 ) 

969 self._ifo = gwinc.precompIFO(self.fT.value, self._ifo) 

970 else: 

971 raise ValueError( 

972 sub_noise 

973 + " is not a subparameter variable noise source.\ 

974 Try calling Get_Noise_Dict on your GroundBased object to find acceptable variables." 

975 ) 

976 else: 

977 err_mssg = ( 

978 base_noise 

979 + " is not a valid parameter variable noise source.\n" 

980 ) 

981 err_mssg += "Try calling Get_Noise_Dict on your GroundBased object to find acceptable variables." 

982 raise ValueError(err_mssg) 

983 else: 

984 raise ValueError("Input must be a dictionary of noise sources.") 

985 

986 def Get_Noise_Dict(self): 

987 """Gets and prints the available variable noises in the detector design""" 

988 i = 0 

989 for key_1, item_1 in self._ifo.items(): 

990 print(key_1, "Parameters:") 

991 for key_2, item_2 in item_1.items(): 

992 if isinstance(item_2, np.ndarray): 

993 i += 1 

994 print(" ", key_2, ": array of shape", item_2.shape) 

995 elif isinstance(item_2, list): 

996 i += 1 

997 print(" ", key_2, ": array of shape", len(item_2)) 

998 elif isinstance(item_2, (int, float)): 

999 i += 1 

1000 print(" ", key_2, ":", item_2) 

1001 elif isinstance(item_2, gwinc.struct.Struct): 

1002 print(" ", key_2, "Subparameters:") 

1003 for key_3, item_3 in item_2.items(): 

1004 if isinstance(item_3, np.ndarray): 

1005 i += 1 

1006 print( 

1007 " ", " ", key_3, ": array of shape", item_3.shape 

1008 ) 

1009 elif isinstance(item_3, list): 

1010 i += 1 

1011 print( 

1012 " ", " ", key_3, ": array of shape", len(item_3) 

1013 ) 

1014 elif isinstance(item_3, (int, float)): 

1015 i += 1 

1016 print(" ", " ", key_3, ":", item_3) 

1017 else: 

1018 i += 1 

1019 print(" ", key_2, ":", item_2) 

1020 

1021 print(" ") 

1022 print("Number of Variables: ", i) 

1023 

1024 

1025class SpaceBased(Interferometer): 

1026 """ 

1027 Class to make a Space Based Instrument using the Interferometer base class 

1028 

1029 Parameters 

1030 ---------- 

1031 L : float 

1032 the armlength the of detector in [meters] 

1033 A_acc : float 

1034 the Amplitude of the Acceleration Noise in [meters/second^2] 

1035 f_acc_break_low : float 

1036 the lower break frequency of the acceleration noise in [Hz] 

1037 f_acc_break_high : float 

1038 the higher break frequency of the acceleration noise in [Hz] 

1039 A_IFO : float 

1040 the amplitude of the interferometer 

1041 

1042 T_type : string, {'N','A'} 

1043 Picks the transfer function generation method 

1044 'N' uses the numerically approximated method in Robson, Cornish, and Liu, 2019 

1045 'A' uses the analytic fit in Larson, Hiscock, and Hellings, 2000 

1046 Background : Boolean 

1047 Add in a Galactic Binary Confusion Noise 

1048 

1049 """ 

1050 

1051 def __init__(self, name, T_obs, *args, **kwargs): 

1052 super().__init__(name, T_obs, **kwargs) 

1053 

1054 for keys, value in kwargs.items(): 

1055 if keys == "T_type": 

1056 self.T_type = value 

1057 elif keys == "Background": 

1058 self.Background = value 

1059 

1060 if not hasattr(self, "nfreqs"): 

1061 self.nfreqs = int(1e3) 

1062 if not hasattr(self, "f_low"): 

1063 self.f_low = 1e-5 * u.Hz 

1064 if not hasattr(self, "f_high"): 

1065 self.f_high = 1.0 * u.Hz 

1066 if not hasattr(self, "Background"): 

1067 self.Background = False 

1068 

1069 if len(args) != 0: 

1070 [L, A_acc, f_acc_break_low, f_acc_break_high, A_IFO, f_IFO_break] = args 

1071 self.L = L 

1072 self.A_acc = A_acc 

1073 self.f_acc_break_low = f_acc_break_low 

1074 self.f_acc_break_high = f_acc_break_high 

1075 self.A_IFO = A_IFO 

1076 self.f_IFO_break = f_IFO_break 

1077 

1078 if not hasattr(self, "load_location"): 

1079 if not hasattr(self, "T_type"): 

1080 self.T_type = "N" 

1081 self.Set_T_Function_Type() 

1082 

1083 @property 

1084 def L(self): 

1085 return self._L 

1086 

1087 @L.setter 

1088 def L(self, value): 

1089 self.var_dict = ["L", value] 

1090 if not isinstance(self._return_value, u.Quantity): 

1091 self._return_value = utils.make_quant(self._return_value, "m") 

1092 self._L = self._return_value 

1093 

1094 @property 

1095 def A_acc(self): 

1096 return self._A_acc 

1097 

1098 @A_acc.setter 

1099 def A_acc(self, value): 

1100 self.var_dict = ["A_acc", value] 

1101 if not isinstance(self._return_value, u.Quantity): 

1102 self._return_value = utils.make_quant(self._return_value, "m/s2") 

1103 self._A_acc = self._return_value 

1104 

1105 @property 

1106 def f_acc_break_low(self): 

1107 return self._f_acc_break_low 

1108 

1109 @f_acc_break_low.setter 

1110 def f_acc_break_low(self, value): 

1111 self.var_dict = ["f_acc_break_low", value] 

1112 if not isinstance(self._return_value, u.Quantity): 

1113 self._return_value = utils.make_quant(self._return_value, "Hz") 

1114 self._f_acc_break_low = self._return_value 

1115 

1116 @property 

1117 def f_acc_break_high(self): 

1118 return self._f_acc_break_high 

1119 

1120 @f_acc_break_high.setter 

1121 def f_acc_break_high(self, value): 

1122 self.var_dict = ["f_acc_break_high", value] 

1123 if not isinstance(self._return_value, u.Quantity): 

1124 self._return_value = utils.make_quant(self._return_value, "Hz") 

1125 self._f_acc_break_high = self._return_value 

1126 

1127 @property 

1128 def A_IFO(self): 

1129 return self._A_IFO 

1130 

1131 @A_IFO.setter 

1132 def A_IFO(self, value): 

1133 self.var_dict = ["A_IFO", value] 

1134 if not isinstance(self._return_value, u.Quantity): 

1135 self._return_value = utils.make_quant(self._return_value, "m") 

1136 self._A_IFO = self._return_value 

1137 

1138 @property 

1139 def f_IFO_break(self): 

1140 return self._f_IFO_break 

1141 

1142 @f_IFO_break.setter 

1143 def f_IFO_break(self, value): 

1144 self.var_dict = ["f_IFO_break", value] 

1145 if not isinstance(self._return_value, u.Quantity): 

1146 self._return_value = utils.make_quant(self._return_value, "Hz") 

1147 self._f_IFO_break = self._return_value 

1148 

1149 @property 

1150 def P_n_f(self): 

1151 """Power Spectral Density""" 

1152 if not hasattr(self, "_P_n_f"): 

1153 if not hasattr(self, "_T_Function_Type"): 

1154 self.Set_T_Function_Type() 

1155 

1156 P_acc = ( 

1157 self.A_acc ** 2 

1158 * (1 + (self.f_acc_break_low / self.fT) ** 2) 

1159 * (1 + (self.fT / (self.f_acc_break_high)) ** 4) 

1160 / (2 * np.pi * self.fT) ** 4 

1161 ) # Acceleration Noise 

1162 P_IMS = self.A_IFO ** 2 * ( 

1163 1 + (self.f_IFO_break / self.fT) ** 4 

1164 ) # Displacement noise of the interferometric TM--to-TM 

1165 

1166 f_trans = const.c / 2 / np.pi / self.L # Transfer frequency 

1167 self._P_n_f = ( 

1168 (P_IMS + 2 * (1 + np.cos(self.fT.value / f_trans.value) ** 2) * P_acc) 

1169 / self.L ** 2 

1170 / u.Hz 

1171 ) 

1172 return self._P_n_f 

1173 

1174 @P_n_f.deleter 

1175 def P_n_f(self): 

1176 del self._P_n_f 

1177 

1178 @property 

1179 def S_n_f(self): 

1180 """Effective Noise Power Specral Density""" 

1181 if not hasattr(self, "_S_n_f"): 

1182 if hasattr(self, "_I_data"): 

1183 if self._I_Type == "ASD": 

1184 S_n_f_sqrt = self._I_data[:, 1] 

1185 self._S_n_f = S_n_f_sqrt ** 2 / u.Hz 

1186 elif self._I_Type == "ENSD": 

1187 self._S_n_f = self._I_data[:, 1] / u.Hz 

1188 elif self._I_Type == "h": 

1189 self._S_n_f = self.h_n_f ** 2 / self.fT 

1190 else: 

1191 if self.Background: 

1192 self._S_n_f = ( 

1193 self.P_n_f + self.Add_Background() 

1194 ) / self.transferfunction ** 2 

1195 else: 

1196 self._S_n_f = self.P_n_f / self.transferfunction ** 2 

1197 

1198 return self._S_n_f 

1199 

1200 @S_n_f.deleter 

1201 def S_n_f(self): 

1202 del self._S_n_f 

1203 

1204 def Load_Transfer_Function(self): 

1205 # Numerical transfer function 

1206 Numerical_Transfer_Function_filedirectory = os.path.join( 

1207 load_directory, "NumericalTransferFunction/transfer.dat" 

1208 ) 

1209 Numerical_Transfer_Function_data = np.loadtxt( 

1210 Numerical_Transfer_Function_filedirectory 

1211 ) 

1212 self._transferfunctiondata = Numerical_Transfer_Function_data 

1213 

1214 def Get_Numeric_Transfer_Function(self): 

1215 if not hasattr(self, "_transferfunctiondata"): 

1216 self.Load_Transfer_Function() 

1217 

1218 fc = const.c / (2 * self.L) # light round trip freq 

1219 LISA_Transfer_Function_f = fc * self._transferfunctiondata[:, 0] 

1220 

1221 idx_f_5 = np.abs(LISA_Transfer_Function_f - self.f_low).argmin() 

1222 idx_f_1 = np.abs(LISA_Transfer_Function_f - self.f_high).argmin() 

1223 

1224 # 3/10 is normalization 2/5sin(openingangle) 

1225 # Some papers use 3/20, not summing over 2 independent low-freq data channels 

1226 self.transferfunction = ( 

1227 np.sqrt(3 / 10) * self._transferfunctiondata[idx_f_5:idx_f_1, 1] 

1228 ) 

1229 self.fT = LISA_Transfer_Function_f[idx_f_5:idx_f_1] 

1230 

1231 def Get_Analytic_Transfer_Function(self, openingangle=None): 

1232 # Response function approximation from Calculation described by Cornish, Robson, Liu 2019 

1233 self.fT = ( 

1234 np.logspace( 

1235 np.log10(self.f_low.value), np.log10(self.f_high.value), self.nfreqs 

1236 ) 

1237 * u.Hz 

1238 ) 

1239 f_L = const.c / 2 / np.pi / self.L # Transfer frequency 

1240 # 3/10 is normalization 2/5sin(openingangle) 

1241 if isinstance(openingangle, (int, float, u.Quantity)): 

1242 R_f = (2 * np.sin(openingangle) ** 2) / 5 / (1 + 0.6 * (self.fT / f_L) ** 2) 

1243 else: 

1244 R_f = 3 / 10 / (1 + 0.6 * (self.fT / f_L) ** 2) 

1245 self.transferfunction = np.sqrt(R_f) 

1246 

1247 def Set_T_Function_Type(self): 

1248 if self.T_type == "n" or self.T_type == "N": 

1249 self._T_type = "numeric" 

1250 elif self.T_type == "a" or self.T_type == "A": 

1251 self._T_type = "analytic" 

1252 else: 

1253 print("\nYou can get the transfer function via 2 methods:") 

1254 print( 

1255 ' *To use the numerically approximated method in Robson, Cornish, and Liu, 2019, input "N".' 

1256 ) 

1257 print( 

1258 ' *To use the analytic fit in Larson, Hiscock, and Hellings, 2000, input "A".' 

1259 ) 

1260 calc_type = input("Please select the calculation type: ") 

1261 self.Set_T_Function_Type(calc_type) 

1262 

1263 if self._T_type == "numeric": 

1264 self.Get_Numeric_Transfer_Function() 

1265 if self._T_type == "analytic": 

1266 self.Get_Analytic_Transfer_Function() 

1267 

1268 def Add_Background(self): 

1269 """ 

1270 Galactic confusions noise parameters as a function of T_obs 

1271 """ 

1272 A = 1.4e-44 

1273 agam = 1100 

1274 bgam = 3 / 10 

1275 afk = 0.0016 

1276 bfk = -2 / 9 

1277 f_k = afk * self.T_obs ** bfk 

1278 gamma = agam * self.T_obs ** bgam 

1279 

1280 f = self.fT.value 

1281 S_c_f = ( 

1282 A 

1283 * (f ** (-7 / 3)) 

1284 * (1 + np.tanh(gamma.value * (f_k.value - f))) 

1285 * (1 / u.Hz) 

1286 ) # White Dwarf Background Noise 

1287 return S_c_f 

1288 

1289 

1290def Load_Data(detector): 

1291 """ 

1292 Function to load in a file to initialize any detector. 

1293 

1294 Parameters 

1295 ---------- 

1296 detector : object 

1297 Instance of a detector class 

1298 

1299 """ 

1300 if not hasattr(detector, "I_type"): 

1301 print("Is the data:") 

1302 print(' *Effective Noise Spectral Density - "E"') 

1303 print(' *Amplitude Spectral Density- "A"') 

1304 print(' *Effective Strain - "h"') 

1305 detector.I_type = input("Please enter one of the answers in quotations: ") 

1306 Load_Data(detector) 

1307 

1308 if detector.I_type == "E" or detector.I_type == "e": 

1309 detector._I_Type = "ENSD" 

1310 elif detector.I_type == "A" or detector.I_type == "a": 

1311 detector._I_Type = "ASD" 

1312 elif detector.I_type == "h" or detector.I_type == "H": 

1313 detector._I_Type = "h" 

1314 else: 

1315 print("Is the data:") 

1316 print(' *Effective Noise Spectral Density - "E"') 

1317 print(' *Amplitude Spectral Density- "A"') 

1318 print(' *Effective Strain - "h"') 

1319 detector.I_type = input("Please enter one of the answers in quotations: ") 

1320 Load_Data(detector) 

1321 

1322 detector._I_data = np.loadtxt(detector.load_location) 

1323 detector.fT = detector._I_data[:, 0] * u.Hz