LALPulsar  6.1.0.1-89842e6
SimulateHeterodynedCW

Detailed Description

The module provides the HeterodynedCWSimulator() class for simulating a signal from a continuous wave source after application of a heterodyned as described in Equations 7 and 8 of [19] .

An example usage to generate the complex heterodyned signal time series is:

from lalpulsar.simulateHeterodynedCW import HeterodynedCWSimulator
from lalpulsar.PulsarParametersWrapper import PulsarParametersPy
import lal
import numpy as np
# set the pulsar parameters
par = PulsarParametersPy()
par['RAJ'] = lal.TranslateHMStoRAD('01:23:34.5')
par['DECJ'] = lal.TranslateDMStoRAD('-45:01:23.4')
par['F'] = [123.456789, -9.87654321e-12] # frequency and first derivative
pepoch = lal.TranslateStringMJDTTtoGPS('58000') # frequency epoch
par['PEPOCH'] = pepoch.gpsSeconds + 1e-9*pepoch.gpsNanoSeconds
par['H0'] = 5.6e-26 # GW amplitude
par['COSIOTA'] = -0.2 # cosine of inclination angle
par['PSI'] = 0.4 # polarization angle (rads)
par['PHI0'] = 2.3 # initial phase (rads)
# set the GPS times of the data
times = np.arange(1000000000.0, 1000086400., 3600)
# set the detector
det = 'H1' # the LIGO Hanford Observatory
# create the HeterodynedCWSimulator object
het = HeterodynedCWSimulator(par, det, times=times)
# get the model complex strain time series
model = het.model(usephase=True)

An example of getting the time series for a signal that has phase parameters that are not identical to the heterodyned parameters would be:

from lalpulsar.simulateHeterodynedCW import HeterodynedCWSimulator
from lalpulsar.PulsarParametersWrapper import PulsarParametersPy
import lal
import numpy as np
# set the "heterodyne" pulsar parameters
par = PulsarParametersPy()
par['RAJ'] = lal.TranslateHMStoRAD('01:23:34.6')
par['DECJ'] = lal.TranslateDMStoRAD('-45:01:23.5')
par['F'] = [123.4567, -9.876e-12] # frequency and first derivative
pepoch = lal.TranslateStringMJDTTtoGPS('58000') # frequency epoch
par['PEPOCH'] = pepoch.gpsSeconds + 1e-9*pepoch.gpsNanoSeconds
# set the times
times = np.arange(1000000000., 1000000600., 60)
# set the detector
det = 'H1' # the LIGO Hanford Observatory
# create the HeterodynedCWSimulator object
het = HeterodynedCWSimulator(par, det, times=times)
# set the updated parameters
parupdate = PulsarParametersPy()
parupdate['RAJ'] = lal.TranslateHMStoRAD('01:23:34.5')
parupdate['DECJ'] = lal.TranslateDMStoRAD('-45:01:23.4')
parupdate['F'] = [123.456789, -9.87654321e-12] # frequency and first derivative
pepoch = lal.TranslateStringMJDTTtoGPS('58000') # frequency epoch
parupdate['PEPOCH'] = pepoch.gpsSeconds + 1e-9*pepoch.gpsNanoSeconds
parupdate['H0'] = 5.6e-26 # GW amplitude
parupdate['COSIOTA'] = -0.2 # cosine of inclination angle
parupdate['PSI'] = 0.4 # polarization angle (rads)
parupdate['PHI0'] = 2.3 # initial phase (rads)
# get the model complex strain time series
model = het.model(parupdate, usephase=True, updateSSB=True)

Data Structures

class  lalpulsar.simulateHeterodynedCW.HeterodynedCWSimulator
 

Variables

string lalpulsar.simulateHeterodynedCW.DOWNLOAD_URL = "https://git.ligo.org/lscsoft/lalsuite/raw/master/lalpulsar/lib/{}"
 

Variable Documentation

◆ DOWNLOAD_URL

string lalpulsar.simulateHeterodynedCW.DOWNLOAD_URL = "https://git.ligo.org/lscsoft/lalsuite/raw/master/lalpulsar/lib/{}"

Definition at line 129 of file simulateHeterodynedCW.py.