Loading [MathJax]/extensions/TeX/AMSsymbols.js
LAL 7.7.0.1-5e288d3
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
TriggerInterpolate.h
Go to the documentation of this file.
1
2/*
3 * Interpolation of time and SNR of triggers
4 *
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with with program; see the file COPYING. If not, write to the
17 * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
18 * MA 02110-1301 USA
19 *
20 * Copyright (C) 2012-2020 Leo Singer
21 */
22
23#ifndef _TRIGGERINTERPOLATE_H
24#define _TRIGGERINTERPOLATE_H
25
26#include <lal/LALAtomicDatatypes.h>
27
28#if defined(__cplusplus)
29extern "C" {
30#elif 0
31} /* so that editors will match preceding brace */
32#endif
33
34
35/**
36 * \defgroup TriggerInterpolation Trigger Interpolation
37 * \ingroup lal_tools
38 *
39 * This module implements several algorithms for performing sub-sample
40 * interpolation of peaks in the output of a matched filter. The interface for
41 * each algorithm is the same. To illustrate, here is a description of the
42 * interface for the 'QuadraticFit' algorithm.
43 *
44 * In this example, the interpolation will take into account a given sample,
45 * the 5 samples before, and the 5 samples after, for a total of 11 samples.
46 *
47 * To apply the interpolant to a time series, call
48 * \c XLALTriggerInterpolateQuadraticFit. It is requried that |y[0]| >= |y[i]|
49 * for all i from -\c window to +\c window. Suppose that you have a complex
50 * array \c y:
51 *
52 * \code{.c}
53 * COMPLEX16 y[] = {...};
54 * \endcode
55 *
56 * and suppose that the maximum of |y[i]| occurs at i = 8.
57 * To interpolate the time and complex value at which the maximum occurs, do
58 * the following:
59 *
60 * \code{.c}
61 * double tmax;
62 * COMPLEX16 ymax;
63 * int result = XLALTriggerInterpolateQuadraticFit(interp, &tmax, &ymax, &y[8], 5);
64 * \endcode
65 *
66 * Upon success, the return value is zero, the interpolated index of the
67 * maximum sample is 8 + \c *tmax, and the interpolated value is in \c ymax.
68 * Upon failure, the return value is nonzero and neither \c *tmax nor \c *ymax
69 * are modified.
70 *
71 * The recommended algorithms are \c XLALTriggerInterpolateCubicSplineAmpPhase.
72 * and \c XLALTriggerInterpolateQuadraticFit. The
73 * \c XLALTriggerInterpolateCubicSpline and \c XLALTriggerInterpolateLanczos
74 * methods should not be used because they tend to exhibit oscillations between
75 * samples and may perform no better than the nearest-neighbor method; these
76 * two functions are preserved only for historical interest. The
77 * \c XLALTriggerInterpolateCubicSplineAmpPhase method is preferred over the
78 * \c XLALTriggerInterpolateQuadraticFit method because the former provides a
79 * continuous and smooth model for the signal at all times, whereas the latter
80 * does not.
81 *
82 * \{
83 */
84
85
86/**
87 * Catmull-Rom cubic spline interpolation on amplitude phase
88 *
89 * Find the Catmull-Rom cubic spline interpolants that pass through the
90 * amplitude and phase of the data.
91 *
92 * \warning It is an error if the window size is less than 2.
93 */
95 double *tmax,
96 COMPLEX16 *ymax,
97 const COMPLEX16 *y,
98 unsigned int window);
99
100
101/**
102 * Catmull-Rom cubic spline interpolation
103 *
104 * Find the Catmull-Rom cubic spline interpolants that pass through the real and
105 * imaginary parts of the data. Find the interpolated time by solving for the
106 * maximum of the sum of the squares of the two spline interpolants.
107 *
108 * \warning It is an error if the window size is less than 2.
109 */
111 double *tmax,
112 COMPLEX16 *ymax,
113 const COMPLEX16 *y,
114 unsigned int window);
115
116
117/**
118 * Lanczos interpolation
119 *
120 * Convolve the data with a Lanczos reconstruction kernel and find the maximum
121 * of the absolute value using a one-dimensional numerical method.
122 *
123 * \warning It is an error if the window size is less than 1.
124 */
126 double *tmax,
127 COMPLEX16 *ymax,
128 const COMPLEX16 *y,
129 unsigned int window);
130
131
132/**
133 * Nearest-neighbor interpolation
134 *
135 * Perform trivial, nearest-neighbor interpolation.
136 * Always sets <tt>*tmax = 0, *ymax = data[0]</tt>.
137 */
139 double *tmax,
140 COMPLEX16 *ymax,
141 const COMPLEX16 *y,
142 unsigned int window);
143
144
145/**
146 * Quadratic fit
147 *
148 * Fit a quadratic function to the absolute value of the data. Report the sample
149 * index of the vertex if the vertex is a maximum. Note that this is not
150 * actually an interpolant because it is not guaranteed to agree with the input
151 * data points, and that it always sets *ymax = data[0].
152 *
153 * \warning It is an error if the window size is less than 1.
154 */
156 double *tmax,
157 COMPLEX16 *ymax,
158 const COMPLEX16 *y,
159 unsigned int window);
160
161
162/** \} */
163
164
165#if 0
166{ /* so that editors will match succeeding brace */
167#elif defined(__cplusplus)
168} /* extern "C" */
169#endif
170
171#endif /* _TRIGGERINTERPOLATE_h */
double complex COMPLEX16
Double-precision floating-point complex number (16 bytes total)
int XLALTriggerInterpolateQuadraticFit(double *tmax, COMPLEX16 *ymax, const COMPLEX16 *y, unsigned int window)
Quadratic fit.
int XLALTriggerInterpolateCubicSplineAmpPhase(double *tmax, COMPLEX16 *ymax, const COMPLEX16 *y, unsigned int window)
Catmull-Rom cubic spline interpolation on amplitude phase.
int XLALTriggerInterpolateCubicSpline(double *tmax, COMPLEX16 *ymax, const COMPLEX16 *y, unsigned int window)
Catmull-Rom cubic spline interpolation.
int XLALTriggerInterpolateLanczos(double *tmax, COMPLEX16 *ymax, const COMPLEX16 *y, unsigned int window)
Lanczos interpolation.
int XLALTriggerInterpolateNearestNeighbor(double *tmax, COMPLEX16 *ymax, const COMPLEX16 *y, unsigned int window)
Nearest-neighbor interpolation.