Skip to content

sgnligo.sources.sim_inspiral_source

A source element to generate gravitational wave waveforms from injection files.

This module provides the SimInspiralSource class which reads injection parameters from XML (LIGOLW) or HDF5 files and generates time-domain waveforms projected onto each detector with proper time delays, antenna patterns, and phase corrections.

All injection parameters are stored internally as LAL dictionaries (lal.Dict) in SI units, which allows pass-through of all LAL-supported parameters including tidal deformability, higher-order modes, and approximant-specific settings.

CachedWaveform dataclass

Cached waveform data for an injection.

Each detector has its own LAL REAL8TimeSeries which includes: - epoch: GPS time of first sample (accounts for light travel time) - deltaT: Sample spacing (1/sample_rate) - data: The strain values

Using LAL TimeSeries preserves sub-sample timing information needed for proper interpolation when injecting into buffers.

Source code in sgnligo/sources/sim_inspiral_source.py
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
@dataclass
class CachedWaveform:
    """Cached waveform data for an injection.

    Each detector has its own LAL REAL8TimeSeries which includes:
    - epoch: GPS time of first sample (accounts for light travel time)
    - deltaT: Sample spacing (1/sample_rate)
    - data: The strain values

    Using LAL TimeSeries preserves sub-sample timing information needed
    for proper interpolation when injecting into buffers.
    """

    injection_id: int
    strain: Dict[str, lal.REAL8TimeSeries]  # Per-IFO strain TimeSeries {ifo: series}

    def get_end_gps(self, ifo: str) -> float:
        """Get the GPS end time for a specific IFO."""
        ts = self.strain[ifo]
        epoch = float(ts.epoch)
        duration = ts.data.length * ts.deltaT
        return epoch + duration

    def get_max_end_gps(self) -> float:
        """Get the latest end GPS time across all IFOs."""
        return max(self.get_end_gps(ifo) for ifo in self.strain)

get_end_gps(ifo)

Get the GPS end time for a specific IFO.

Source code in sgnligo/sources/sim_inspiral_source.py
230
231
232
233
234
235
def get_end_gps(self, ifo: str) -> float:
    """Get the GPS end time for a specific IFO."""
    ts = self.strain[ifo]
    epoch = float(ts.epoch)
    duration = ts.data.length * ts.deltaT
    return epoch + duration

get_max_end_gps()

Get the latest end GPS time across all IFOs.

Source code in sgnligo/sources/sim_inspiral_source.py
237
238
239
def get_max_end_gps(self) -> float:
    """Get the latest end GPS time across all IFOs."""
    return max(self.get_end_gps(ifo) for ifo in self.strain)

SimInspiralSource dataclass

Bases: TSSource


              flowchart TD
              sgnligo.sources.sim_inspiral_source.SimInspiralSource[SimInspiralSource]

              

              click sgnligo.sources.sim_inspiral_source.SimInspiralSource href "" "sgnligo.sources.sim_inspiral_source.SimInspiralSource"
            

Source element that generates GW waveforms from an injection file or test mode.

Reads injection parameters from XML (LIGOLW SimInspiralTable) or HDF5 files and generates time-domain waveforms projected onto each detector with proper time delays, antenna patterns, and phase corrections using LALSimulation.

Alternatively, test mode can be used to generate periodic test injections (BNS, NSBH, or BBH) every 30 seconds, positioned directly overhead of State College, PA.

All injection parameters are stored as InjectionParams objects, which use LAL dictionaries internally for pass-through of any LAL-supported parameters including tidal deformability, higher-order modes, and approximant-specific settings.

The source outputs zeros + summed injection signals on separate pads for each interferometer. Combine with a noise source using an Adder transform to create realistic data with injections.

Parameters:

Name Type Description Default
injection_file Optional[str]

Path to injection file (XML or HDF5). Mutually exclusive with test_mode.

None
test_mode Optional[str]

Generate test injections instead of loading from file. Valid values: "bns" (100 Mpc), "nsbh" (200 Mpc), "bbh" (500 Mpc). Mutually exclusive with injection_file.

None
ifos Optional[List[str]]

List of interferometer prefixes, e.g., ["H1", "L1", "V1"]

None
sample_rate int

Output sample rate in Hz (default: 16384)

16384
f_min float

Minimum frequency for waveform generation in Hz (default: 10)

10.0
approximant_override Optional[str]

Override approximant for all injections (optional)

None
Example

Using injection file

source = SimInspiralSource( ... name="Injections", ... injection_file="injections.xml", ... ifos=["H1", "L1"], ... t0=1234567890, ... end=1234567900, ... )

Using test mode

source = SimInspiralSource( ... name="TestInjections", ... test_mode="bbh", ... ifos=["H1", "L1"], ... t0=1234567890, ... end=1234567900, ... )

Source code in sgnligo/sources/sim_inspiral_source.py
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
@dataclass
class SimInspiralSource(TSSource):
    """Source element that generates GW waveforms from an injection file or test mode.

    Reads injection parameters from XML (LIGOLW SimInspiralTable) or HDF5 files
    and generates time-domain waveforms projected onto each detector with proper
    time delays, antenna patterns, and phase corrections using LALSimulation.

    Alternatively, test mode can be used to generate periodic test injections
    (BNS, NSBH, or BBH) every 30 seconds, positioned directly overhead of
    State College, PA.

    All injection parameters are stored as InjectionParams objects, which use
    LAL dictionaries internally for pass-through of any LAL-supported parameters
    including tidal deformability, higher-order modes, and approximant-specific
    settings.

    The source outputs zeros + summed injection signals on separate pads for
    each interferometer. Combine with a noise source using an Adder transform
    to create realistic data with injections.

    Args:
        injection_file: Path to injection file (XML or HDF5). Mutually exclusive
            with test_mode.
        test_mode: Generate test injections instead of loading from file.
            Valid values: "bns" (100 Mpc), "nsbh" (200 Mpc), "bbh" (500 Mpc).
            Mutually exclusive with injection_file.
        ifos: List of interferometer prefixes, e.g., ["H1", "L1", "V1"]
        sample_rate: Output sample rate in Hz (default: 16384)
        f_min: Minimum frequency for waveform generation in Hz (default: 10)
        approximant_override: Override approximant for all injections (optional)

    Example:
        >>> # Using injection file
        >>> source = SimInspiralSource(
        ...     name="Injections",
        ...     injection_file="injections.xml",
        ...     ifos=["H1", "L1"],
        ...     t0=1234567890,
        ...     end=1234567900,
        ... )

        >>> # Using test mode
        >>> source = SimInspiralSource(
        ...     name="TestInjections",
        ...     test_mode="bbh",
        ...     ifos=["H1", "L1"],
        ...     t0=1234567890,
        ...     end=1234567900,
        ... )
    """

    injection_file: Optional[str] = None
    test_mode: Optional[str] = None
    ifos: Optional[List[str]] = None
    sample_rate: int = 16384
    f_min: float = 10.0
    approximant_override: Optional[str] = None

    # Internal state (not user-configurable)
    _waveform_cache: WaveformCache = field(init=False, repr=False)
    _injections: List[InjectionParams] = field(
        init=False, repr=False, default_factory=list
    )
    _channel_dict: Dict[str, str] = field(init=False, repr=False, default_factory=dict)

    def __post_init__(self):
        """Initialize the source."""
        # Validate mutual exclusivity of injection_file and test_mode
        if self.injection_file and self.test_mode:
            raise ValueError("Cannot specify both injection_file and test_mode")
        if not self.injection_file and not self.test_mode:
            raise ValueError("Must specify either injection_file or test_mode")

        # Validate test_mode value
        if self.test_mode:
            test_mode_lower = self.test_mode.lower()
            if test_mode_lower not in TEST_MODE_PARAMS:
                raise ValueError(
                    f"Invalid test_mode '{self.test_mode}'. "
                    f"Must be one of: {list(TEST_MODE_PARAMS.keys())}"
                )
            # Normalize to lowercase
            self.test_mode = test_mode_lower

        if self.ifos is None:
            self.ifos = ["H1", "L1"]

        # Create channel names for source pads
        self._channel_dict = {ifo: f"{ifo}:INJ-STRAIN" for ifo in self.ifos}
        self.source_pad_names = list(self._channel_dict.values())

        # Load injections from file or start with empty list for test mode
        if self.injection_file:
            self._injections = load_injections(self.injection_file)
        else:
            # Test mode: injections will be generated dynamically
            self._injections = []

        # Initialize waveform cache
        self._waveform_cache = WaveformCache(
            injections=self._injections,
            ifos=self.ifos,
            sample_rate=self.sample_rate,
            f_min=self.f_min,
            approximant_override=self.approximant_override,
            test_mode=self.test_mode,
        )

        # Call parent's post_init
        super().__post_init__()

        # Set buffer params for each pad
        for _ifo, channel in self._channel_dict.items():
            pad = self.srcs[channel]
            self.set_pad_buffer_params(
                pad=pad,
                sample_shape=(),
                rate=self.sample_rate,
            )

    def _ifo_from_pad(self, pad: SourcePad) -> str:
        """Get IFO prefix from pad."""
        for ifo, channel in self._channel_dict.items():
            if self.srcs[channel] == pad:
                return ifo
        raise ValueError(f"Unknown pad: {pad}")

    def new(self, pad: SourcePad) -> TSFrame:
        """Generate a new frame with injection signals.

        Args:
            pad: Source pad requesting new data

        Returns:
            TSFrame containing injection signals (zeros + waveforms)
        """
        # Get frame prepared by base class
        frame = self.prepare_frame(pad)
        buffer = frame.buffers[0]

        # Determine IFO from pad
        ifo = self._ifo_from_pad(pad)

        # Get buffer time window
        # Convert from Offset to GPS seconds
        # SGNTS Offset is normalized to Offset.MAX_RATE internally, regardless of
        # the source's sample rate. This is a core SGNTS design choice.
        stride_max = Offset.sample_stride(Offset.MAX_RATE)
        buf_start = buffer.offset / stride_max
        buf_end = buffer.end_offset / stride_max

        # Get number of samples from the buffer's expected shape
        num_samples = buffer.shape[0]

        # Create a LAL REAL8TimeSeries as target for injections
        # This preserves exact GPS epoch for proper sub-sample interpolation
        target = lal.CreateREAL8TimeSeries(
            f"{ifo}:STRAIN",
            lal.LIGOTimeGPS(buf_start),
            0.0,  # f0
            1.0 / self.sample_rate,  # deltaT
            lal.StrainUnit,
            num_samples,
        )
        target.data.data[:] = 0.0

        # Find and add all overlapping injections
        # SimAddInjectionREAL8TimeSeries handles sub-sample interpolation
        overlapping = self._waveform_cache.get_overlapping_injections(
            buf_start, buf_end
        )
        for inj_id in overlapping:
            self._waveform_cache.add_injection_to_target(inj_id, ifo, target)

        # Copy result to output buffer
        buffer.set_data(target.data.data)
        return frame

    def internal(self) -> None:
        """Internal processing - cleanup expired waveforms."""
        super().internal()

        # Cleanup expired waveforms from cache
        # Convert from Offset to GPS seconds
        # SGNTS Offset is normalized to Offset.MAX_RATE internally
        stride_max = Offset.sample_stride(Offset.MAX_RATE)
        current_gps = self.current_end / stride_max
        self._waveform_cache.cleanup_expired(current_gps)

__post_init__()

Initialize the source.

Source code in sgnligo/sources/sim_inspiral_source.py
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
def __post_init__(self):
    """Initialize the source."""
    # Validate mutual exclusivity of injection_file and test_mode
    if self.injection_file and self.test_mode:
        raise ValueError("Cannot specify both injection_file and test_mode")
    if not self.injection_file and not self.test_mode:
        raise ValueError("Must specify either injection_file or test_mode")

    # Validate test_mode value
    if self.test_mode:
        test_mode_lower = self.test_mode.lower()
        if test_mode_lower not in TEST_MODE_PARAMS:
            raise ValueError(
                f"Invalid test_mode '{self.test_mode}'. "
                f"Must be one of: {list(TEST_MODE_PARAMS.keys())}"
            )
        # Normalize to lowercase
        self.test_mode = test_mode_lower

    if self.ifos is None:
        self.ifos = ["H1", "L1"]

    # Create channel names for source pads
    self._channel_dict = {ifo: f"{ifo}:INJ-STRAIN" for ifo in self.ifos}
    self.source_pad_names = list(self._channel_dict.values())

    # Load injections from file or start with empty list for test mode
    if self.injection_file:
        self._injections = load_injections(self.injection_file)
    else:
        # Test mode: injections will be generated dynamically
        self._injections = []

    # Initialize waveform cache
    self._waveform_cache = WaveformCache(
        injections=self._injections,
        ifos=self.ifos,
        sample_rate=self.sample_rate,
        f_min=self.f_min,
        approximant_override=self.approximant_override,
        test_mode=self.test_mode,
    )

    # Call parent's post_init
    super().__post_init__()

    # Set buffer params for each pad
    for _ifo, channel in self._channel_dict.items():
        pad = self.srcs[channel]
        self.set_pad_buffer_params(
            pad=pad,
            sample_shape=(),
            rate=self.sample_rate,
        )

internal()

Internal processing - cleanup expired waveforms.

Source code in sgnligo/sources/sim_inspiral_source.py
624
625
626
627
628
629
630
631
632
633
def internal(self) -> None:
    """Internal processing - cleanup expired waveforms."""
    super().internal()

    # Cleanup expired waveforms from cache
    # Convert from Offset to GPS seconds
    # SGNTS Offset is normalized to Offset.MAX_RATE internally
    stride_max = Offset.sample_stride(Offset.MAX_RATE)
    current_gps = self.current_end / stride_max
    self._waveform_cache.cleanup_expired(current_gps)

new(pad)

Generate a new frame with injection signals.

Parameters:

Name Type Description Default
pad SourcePad

Source pad requesting new data

required

Returns:

Type Description
TSFrame

TSFrame containing injection signals (zeros + waveforms)

Source code in sgnligo/sources/sim_inspiral_source.py
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
def new(self, pad: SourcePad) -> TSFrame:
    """Generate a new frame with injection signals.

    Args:
        pad: Source pad requesting new data

    Returns:
        TSFrame containing injection signals (zeros + waveforms)
    """
    # Get frame prepared by base class
    frame = self.prepare_frame(pad)
    buffer = frame.buffers[0]

    # Determine IFO from pad
    ifo = self._ifo_from_pad(pad)

    # Get buffer time window
    # Convert from Offset to GPS seconds
    # SGNTS Offset is normalized to Offset.MAX_RATE internally, regardless of
    # the source's sample rate. This is a core SGNTS design choice.
    stride_max = Offset.sample_stride(Offset.MAX_RATE)
    buf_start = buffer.offset / stride_max
    buf_end = buffer.end_offset / stride_max

    # Get number of samples from the buffer's expected shape
    num_samples = buffer.shape[0]

    # Create a LAL REAL8TimeSeries as target for injections
    # This preserves exact GPS epoch for proper sub-sample interpolation
    target = lal.CreateREAL8TimeSeries(
        f"{ifo}:STRAIN",
        lal.LIGOTimeGPS(buf_start),
        0.0,  # f0
        1.0 / self.sample_rate,  # deltaT
        lal.StrainUnit,
        num_samples,
    )
    target.data.data[:] = 0.0

    # Find and add all overlapping injections
    # SimAddInjectionREAL8TimeSeries handles sub-sample interpolation
    overlapping = self._waveform_cache.get_overlapping_injections(
        buf_start, buf_end
    )
    for inj_id in overlapping:
        self._waveform_cache.add_injection_to_target(inj_id, ifo, target)

    # Copy result to output buffer
    buffer.set_data(target.data.data)
    return frame

WaveformCache

Manages waveform generation and caching.

Waveforms are generated on-demand and cached until they're fully consumed (i.e., the pipeline has moved past the waveform's end time).

Source code in sgnligo/sources/sim_inspiral_source.py
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
class WaveformCache:
    """Manages waveform generation and caching.

    Waveforms are generated on-demand and cached until they're fully
    consumed (i.e., the pipeline has moved past the waveform's end time).
    """

    def __init__(
        self,
        injections: List[InjectionParams],
        ifos: List[str],
        sample_rate: int,
        f_min: float,
        approximant_override: Optional[str] = None,
        test_mode: Optional[str] = None,
    ):
        """Initialize the waveform cache.

        Args:
            injections: List of InjectionParams with injection parameters
            ifos: List of interferometer prefixes
            sample_rate: Output sample rate in Hz
            f_min: Minimum frequency for waveform generation
            approximant_override: Override approximant for all injections
            test_mode: If set, generate test injections ("bns", "nsbh", "bbh")
        """
        self.injections = injections
        self.ifos = ifos
        self.sample_rate = sample_rate
        self.f_min = f_min
        self.approximant_override = approximant_override
        self.cache: Dict[int, CachedWaveform] = {}

        # Test mode configuration
        self._test_mode = test_mode
        self._generated_test_times: Set[float] = set()

        # Pre-compute injection time windows for fast overlap queries
        self._injection_windows: List[tuple[float, float]] = []
        for inj in injections:
            lal_d = inj.to_lal_dict()
            pre_dur, post_dur = estimate_waveform_duration(lal_d, f_min)
            t_co_gps = inj.t_co_gps
            start = t_co_gps - pre_dur
            end = t_co_gps + post_dur
            self._injection_windows.append((start, end))

    def _ensure_test_injections_for_range(
        self, start_gps: float, end_gps: float
    ) -> None:
        """Generate test injections whose waveforms could overlap the range.

        For test mode, this method dynamically creates injections at every
        TEST_INJECTION_INTERVAL seconds. We need to generate injections with
        coalescence times that could produce waveforms overlapping [start_gps, end_gps].

        Since waveforms extend before (inspiral) and after (ringdown) the coalescence
        time, we need to look for t_co values where:
        - t_co - pre_dur < end_gps  (waveform starts before buffer ends)
        - t_co + post_dur > start_gps  (waveform ends after buffer starts)

        We use a conservative lookahead to ensure we don't miss any injections.

        Args:
            start_gps: Start of the time range (GPS seconds)
            end_gps: End of the time range (GPS seconds)
        """
        if not self._test_mode:
            return

        # Use a conservative estimate for maximum waveform duration before coalescence
        # BNS from 10-20 Hz can have inspirals of 100+ seconds, BBH is shorter
        # We'll use 300 seconds as a safe upper bound for any test mode
        max_pre_duration_estimate = 300.0

        # Expand the range to catch all injections whose waveforms could overlap
        # Look back by max_pre_duration (for late coalescence times with long inspirals)
        # Look ahead by one interval (for early coalescence with long post-merger)
        search_start = start_gps - max_pre_duration_estimate
        search_end = end_gps + TEST_INJECTION_INTERVAL

        # Find the first 30-second boundary at or after search_start
        first_t_co = (
            math.ceil(search_start / TEST_INJECTION_INTERVAL) * TEST_INJECTION_INTERVAL
        )

        t_co = first_t_co
        while t_co <= search_end:
            if t_co not in self._generated_test_times:
                self._generated_test_times.add(t_co)
                # Generate the injection using shared factory
                inj = create_test_injection(self._test_mode, t_co)
                # Add to injections list
                self.injections.append(inj)
                # Compute and add time window
                lal_d = inj.to_lal_dict()
                pre_dur, post_dur = estimate_waveform_duration(lal_d, self.f_min)
                self._injection_windows.append((t_co - pre_dur, t_co + post_dur))
            t_co += TEST_INJECTION_INTERVAL

    def get_overlapping_injections(self, buf_start: float, buf_end: float) -> List[int]:
        """Find injections that overlap the buffer time window.

        For test mode, this also ensures test injections exist for the query range.

        Args:
            buf_start: Buffer start GPS time
            buf_end: Buffer end GPS time

        Returns:
            List of injection indices that overlap the buffer
        """
        # For test mode, ensure we have injections for this time range
        self._ensure_test_injections_for_range(buf_start, buf_end)

        overlapping = []
        for i, (inj_start, inj_end) in enumerate(self._injection_windows):
            if inj_start < buf_end and inj_end > buf_start:
                overlapping.append(i)
        return overlapping

    def _generate_and_cache(self, inj_id: int) -> None:
        """Generate waveform for an injection and cache it.

        Args:
            inj_id: Index of injection in self.injections
        """
        inj = self.injections[inj_id]
        params = inj.to_lal_dict()

        # Generate h+, hx waveforms
        hp, hc = generate_waveform_td(
            params, self.sample_rate, self.f_min, self.approximant_override
        )

        # Get coalescence time and shift waveform epochs to absolute GPS time
        # before projection. The generated waveforms have epochs relative to
        # coalescence (typically negative).
        t_co_gps = inj.t_co_gps
        hp.epoch += t_co_gps
        hc.epoch += t_co_gps

        # Project onto each detector
        # Each detector has different timing due to light travel time delays
        # applied by SimDetectorStrainREAL8TimeSeries
        strain_dict: Dict[str, lal.REAL8TimeSeries] = {}

        for ifo in self.ifos:
            strain_dict[ifo] = project_to_detector(hp, hc, params, ifo)

        # Cache the waveform
        self.cache[inj_id] = CachedWaveform(
            injection_id=inj_id,
            strain=strain_dict,
        )

    def add_injection_to_target(
        self,
        inj_id: int,
        ifo: str,
        target: lal.REAL8TimeSeries,
    ) -> None:
        """Add an injection waveform to a target buffer using proper interpolation.

        Uses SimAddInjectionREAL8TimeSeries which performs sub-sample
        re-interpolation in the frequency domain to properly align the
        injection epoch to integer sample boundaries in the target.

        Args:
            inj_id: Injection index
            ifo: Interferometer prefix
            target: Target LAL REAL8TimeSeries to add injection into (modified in place)
        """
        if inj_id not in self.cache:
            self._generate_and_cache(inj_id)

        cached = self.cache[inj_id]
        source_strain = cached.strain[ifo]

        # SimAddInjectionREAL8TimeSeries handles:
        # - Finding overlapping region based on epochs
        # - Sub-sample interpolation for proper alignment
        # - Adding only the overlapping portion to target
        # The source may be modified (padded for alignment) but remains reusable
        lalsim.SimAddInjectionREAL8TimeSeries(target, source_strain, None)

    def cleanup_expired(self, current_gps: float) -> None:
        """Remove waveforms that are fully consumed.

        A waveform is expired when all IFOs have moved past it.

        Args:
            current_gps: Current GPS time of the pipeline
        """
        expired = []
        for k, v in self.cache.items():
            # Use the latest end time across all IFOs
            if v.get_max_end_gps() < current_gps:
                expired.append(k)
        for k in expired:
            del self.cache[k]

__init__(injections, ifos, sample_rate, f_min, approximant_override=None, test_mode=None)

Initialize the waveform cache.

Parameters:

Name Type Description Default
injections List[InjectionParams]

List of InjectionParams with injection parameters

required
ifos List[str]

List of interferometer prefixes

required
sample_rate int

Output sample rate in Hz

required
f_min float

Minimum frequency for waveform generation

required
approximant_override Optional[str]

Override approximant for all injections

None
test_mode Optional[str]

If set, generate test injections ("bns", "nsbh", "bbh")

None
Source code in sgnligo/sources/sim_inspiral_source.py
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
def __init__(
    self,
    injections: List[InjectionParams],
    ifos: List[str],
    sample_rate: int,
    f_min: float,
    approximant_override: Optional[str] = None,
    test_mode: Optional[str] = None,
):
    """Initialize the waveform cache.

    Args:
        injections: List of InjectionParams with injection parameters
        ifos: List of interferometer prefixes
        sample_rate: Output sample rate in Hz
        f_min: Minimum frequency for waveform generation
        approximant_override: Override approximant for all injections
        test_mode: If set, generate test injections ("bns", "nsbh", "bbh")
    """
    self.injections = injections
    self.ifos = ifos
    self.sample_rate = sample_rate
    self.f_min = f_min
    self.approximant_override = approximant_override
    self.cache: Dict[int, CachedWaveform] = {}

    # Test mode configuration
    self._test_mode = test_mode
    self._generated_test_times: Set[float] = set()

    # Pre-compute injection time windows for fast overlap queries
    self._injection_windows: List[tuple[float, float]] = []
    for inj in injections:
        lal_d = inj.to_lal_dict()
        pre_dur, post_dur = estimate_waveform_duration(lal_d, f_min)
        t_co_gps = inj.t_co_gps
        start = t_co_gps - pre_dur
        end = t_co_gps + post_dur
        self._injection_windows.append((start, end))

add_injection_to_target(inj_id, ifo, target)

Add an injection waveform to a target buffer using proper interpolation.

Uses SimAddInjectionREAL8TimeSeries which performs sub-sample re-interpolation in the frequency domain to properly align the injection epoch to integer sample boundaries in the target.

Parameters:

Name Type Description Default
inj_id int

Injection index

required
ifo str

Interferometer prefix

required
target REAL8TimeSeries

Target LAL REAL8TimeSeries to add injection into (modified in place)

required
Source code in sgnligo/sources/sim_inspiral_source.py
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
def add_injection_to_target(
    self,
    inj_id: int,
    ifo: str,
    target: lal.REAL8TimeSeries,
) -> None:
    """Add an injection waveform to a target buffer using proper interpolation.

    Uses SimAddInjectionREAL8TimeSeries which performs sub-sample
    re-interpolation in the frequency domain to properly align the
    injection epoch to integer sample boundaries in the target.

    Args:
        inj_id: Injection index
        ifo: Interferometer prefix
        target: Target LAL REAL8TimeSeries to add injection into (modified in place)
    """
    if inj_id not in self.cache:
        self._generate_and_cache(inj_id)

    cached = self.cache[inj_id]
    source_strain = cached.strain[ifo]

    # SimAddInjectionREAL8TimeSeries handles:
    # - Finding overlapping region based on epochs
    # - Sub-sample interpolation for proper alignment
    # - Adding only the overlapping portion to target
    # The source may be modified (padded for alignment) but remains reusable
    lalsim.SimAddInjectionREAL8TimeSeries(target, source_strain, None)

cleanup_expired(current_gps)

Remove waveforms that are fully consumed.

A waveform is expired when all IFOs have moved past it.

Parameters:

Name Type Description Default
current_gps float

Current GPS time of the pipeline

required
Source code in sgnligo/sources/sim_inspiral_source.py
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
def cleanup_expired(self, current_gps: float) -> None:
    """Remove waveforms that are fully consumed.

    A waveform is expired when all IFOs have moved past it.

    Args:
        current_gps: Current GPS time of the pipeline
    """
    expired = []
    for k, v in self.cache.items():
        # Use the latest end time across all IFOs
        if v.get_max_end_gps() < current_gps:
            expired.append(k)
    for k in expired:
        del self.cache[k]

get_overlapping_injections(buf_start, buf_end)

Find injections that overlap the buffer time window.

For test mode, this also ensures test injections exist for the query range.

Parameters:

Name Type Description Default
buf_start float

Buffer start GPS time

required
buf_end float

Buffer end GPS time

required

Returns:

Type Description
List[int]

List of injection indices that overlap the buffer

Source code in sgnligo/sources/sim_inspiral_source.py
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
def get_overlapping_injections(self, buf_start: float, buf_end: float) -> List[int]:
    """Find injections that overlap the buffer time window.

    For test mode, this also ensures test injections exist for the query range.

    Args:
        buf_start: Buffer start GPS time
        buf_end: Buffer end GPS time

    Returns:
        List of injection indices that overlap the buffer
    """
    # For test mode, ensure we have injections for this time range
    self._ensure_test_injections_for_range(buf_start, buf_end)

    overlapping = []
    for i, (inj_start, inj_end) in enumerate(self._injection_windows):
        if inj_start < buf_end and inj_end > buf_start:
            overlapping.append(i)
    return overlapping

estimate_waveform_duration(params, f_min)

Estimate the duration of a waveform before and after merger.

Parameters:

Name Type Description Default
params Dict

LAL dict with injection parameters (masses in SI units)

required
f_min float

Minimum frequency for waveform generation (Hz)

required

Returns:

Type Description
tuple[float, float]

Tuple of (pre_merger_duration, post_merger_duration) in seconds

Source code in sgnligo/sources/sim_inspiral_source.py
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
def estimate_waveform_duration(params: lal.Dict, f_min: float) -> tuple[float, float]:
    """Estimate the duration of a waveform before and after merger.

    Args:
        params: LAL dict with injection parameters (masses in SI units)
        f_min: Minimum frequency for waveform generation (Hz)

    Returns:
        Tuple of (pre_merger_duration, post_merger_duration) in seconds
    """
    # Masses are already in SI units from the dict
    m1_si = get_lal_dict_value(params, "mass1")
    m2_si = get_lal_dict_value(params, "mass2")
    spin1z = get_lal_dict_value(params, "spin1z")
    spin2z = get_lal_dict_value(params, "spin2z")

    mtotal_si = m1_si + m2_si

    # Chirp time (inspiral duration)
    try:
        tchirp = lalsim.SimInspiralChirpTimeBound(f_min, m1_si, m2_si, spin1z, spin2z)
    except Exception:
        # Fallback estimate using leading-order chirp time
        mchirp = (m1_si * m2_si) ** (3.0 / 5.0) / mtotal_si ** (1.0 / 5.0)
        tchirp = (
            5.0
            / 256.0
            * (lal.C_SI**3 / (lal.G_SI * mchirp)) ** (5.0 / 3.0)
            * (np.pi * f_min) ** (-8.0 / 3.0)
        )

    # Merge time bound
    try:
        tmerge = lalsim.SimInspiralMergeTimeBound(m1_si, m2_si)
    except Exception:
        tmerge = 0.1  # Conservative estimate

    # Ringdown time bound
    try:
        # Final spin estimate (simple approximation)
        final_spin = min(abs(spin1z) + abs(spin2z), 0.998)
        tring = lalsim.SimInspiralRingdownTimeBound(mtotal_si, final_spin)
    except Exception:
        tring = 0.5  # Conservative estimate

    # Add safety margins
    pre_merger = tchirp + 1.0  # Extra second before
    post_merger = tmerge + tring + 0.5  # Extra half second after

    return pre_merger, post_merger

generate_waveform_td(params, sample_rate, f_min, approximant_override=None)

Generate time-domain h+, hx waveforms from LAL parameter dict.

Uses LALSimulation's unified generator interface which automatically handles both time-domain and frequency-domain approximants. All parameters from the input dict are passed through to the generator.

Parameters:

Name Type Description Default
params Dict

LAL dict with injection parameters (masses/distance in SI units)

required
sample_rate int

Sample rate in Hz

required
f_min float

Minimum frequency in Hz

required
approximant_override Optional[str]

Override the approximant from the dict

None

Returns:

Type Description
tuple[REAL8TimeSeries, REAL8TimeSeries]

Tuple of (hp, hc) as LAL REAL8TimeSeries objects

Source code in sgnligo/sources/sim_inspiral_source.py
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
def generate_waveform_td(
    params: lal.Dict,
    sample_rate: int,
    f_min: float,
    approximant_override: Optional[str] = None,
) -> tuple[lal.REAL8TimeSeries, lal.REAL8TimeSeries]:
    """Generate time-domain h+, hx waveforms from LAL parameter dict.

    Uses LALSimulation's unified generator interface which automatically
    handles both time-domain and frequency-domain approximants. All
    parameters from the input dict are passed through to the generator.

    Args:
        params: LAL dict with injection parameters (masses/distance in SI units)
        sample_rate: Sample rate in Hz
        f_min: Minimum frequency in Hz
        approximant_override: Override the approximant from the dict

    Returns:
        Tuple of (hp, hc) as LAL REAL8TimeSeries objects
    """
    # Get or override approximant (validated at load time)
    approximant_str = approximant_override or get_lal_dict_string(params, "approximant")
    approximant = lalsim.GetApproximantFromString(approximant_str)

    # Add runtime parameters to dict (these are not stored in injection files)
    lalsim.SimInspiralWaveformParamsInsertDeltaT(params, 1.0 / sample_rate)
    lalsim.SimInspiralWaveformParamsInsertF22Start(params, f_min)

    # Use f_ref from dict if present, otherwise use f_min
    f_ref = get_lal_dict_value(params, "f_ref", f_min)
    lalsim.SimInspiralWaveformParamsInsertF22Ref(params, f_ref)

    # Create generator and add conditioning for FD->TD conversion
    gen = lalsim.SimInspiralChooseGenerator(approximant, None)
    lalsim.SimInspiralGeneratorAddStandardConditioning(gen)

    # Generate TD waveform - all parameters from dict are used
    hp, hc = lalsim.SimInspiralGenerateTDWaveform(params, gen)

    return hp, hc

project_to_detector(hp, hc, params, ifo)

Project h+, hx waveforms onto a detector.

Uses LALSimulation's accurate projection which handles: - Antenna response (F+, Fx) - Light travel time delays - Phase corrections

Parameters:

Name Type Description Default
hp REAL8TimeSeries

Plus polarization time series

required
hc REAL8TimeSeries

Cross polarization time series

required
params Dict

LAL dict with sky position and polarization

required
ifo str

Interferometer prefix (H1, L1, V1, etc.)

required

Returns:

Type Description
REAL8TimeSeries

Detector strain as REAL8TimeSeries

Source code in sgnligo/sources/sim_inspiral_source.py
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
def project_to_detector(
    hp: lal.REAL8TimeSeries,
    hc: lal.REAL8TimeSeries,
    params: lal.Dict,
    ifo: str,
) -> lal.REAL8TimeSeries:
    """Project h+, hx waveforms onto a detector.

    Uses LALSimulation's accurate projection which handles:
    - Antenna response (F+, Fx)
    - Light travel time delays
    - Phase corrections

    Args:
        hp: Plus polarization time series
        hc: Cross polarization time series
        params: LAL dict with sky position and polarization
        ifo: Interferometer prefix (H1, L1, V1, etc.)

    Returns:
        Detector strain as REAL8TimeSeries
    """
    # Get detector
    detector = lal.cached_detector_by_prefix[ifo]

    # Extract sky position and polarization (validated at load time)
    ra = get_lal_dict_value(params, "ra")
    dec = get_lal_dict_value(params, "dec")
    psi = get_lal_dict_value(params, "psi")

    # Use LALSimulation's accurate detector strain function
    strain = lalsim.SimDetectorStrainREAL8TimeSeries(
        hp,
        hc,
        ra,
        dec,
        psi,
        detector,
    )

    return strain