LAL  7.5.0.1-89842e6
CubicSplineTriggerInterpolantTest.c
Go to the documentation of this file.
1 /*
2  * This program is free software; you can redistribute it and/or modify
3  * it under the terms of the GNU General Public License as published by
4  * the Free Software Foundation; either version 2 of the License, or
5  * (at your option) any later version.
6  *
7  * This program is distributed in the hope that it will be useful,
8  * but WITHOUT ANY WARRANTY; without even the implied warranty of
9  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
10  * GNU General Public License for more details.
11  *
12  * You should have received a copy of the GNU General Public License
13  * along with with program; see the file COPYING. If not, write to the
14  * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
15  * MA 02110-1301 USA
16  *
17  * Copyright (C) 2012 Leo Singer
18  */
19 
20 
21 #include <complex.h>
22 #include <math.h>
23 #include <stdlib.h>
24 
25 #include <lal/TriggerInterpolation.h>
26 
27 int main(__attribute__ ((unused)) int argc, __attribute__ ((unused)) char **argv)
28 {
29  int result;
30  double tmax;
31 
32  CubicSplineTriggerInterpolant *interp = XLALCreateCubicSplineTriggerInterpolant(2);
33  if (!interp)
34  exit(EXIT_FAILURE);
35 
36  {
37  /* Concave-up, increasing series: maximum should occur at upper sample boundary */
38  const COMPLEX16 y[] = {0, 1, 8, 27, 256};
39  COMPLEX16 ymax;
40 
41  result = XLALCOMPLEX16ApplyCubicSplineTriggerInterpolant(interp, &tmax, &ymax, &y[2]);
42  if (result)
43  exit(EXIT_FAILURE);
44 
45  if (tmax != 1)
46  exit(EXIT_FAILURE);
47  if (creal(ymax) != 27)
48  exit(EXIT_FAILURE);
49  if (cimag(ymax) != 0)
50  exit(EXIT_FAILURE);
51  }
52 
53  {
54  /* Concave-up, increasing series: maximum should occur at upper sample boundary */
55  const COMPLEX8 y[] = {0, 1, 8, 27, 256};
56  COMPLEX8 ymax;
57 
58  result = XLALCOMPLEX8ApplyCubicSplineTriggerInterpolant(interp, &tmax, &ymax, &y[2]);
59  if (result)
60  exit(EXIT_FAILURE);
61 
62  if (tmax != 1)
63  exit(EXIT_FAILURE);
64  if (crealf(ymax) != 27)
65  exit(EXIT_FAILURE);
66  if (cimagf(ymax) != 0)
67  exit(EXIT_FAILURE);
68  }
69 
70  {
71  /* Concave-up, increasing series: maximum should occur at upper sample boundary */
72  const REAL8 y[] = {0, 1, 8, 27, 256};
73  REAL8 ymax;
74 
75  result = XLALREAL8ApplyCubicSplineTriggerInterpolant(interp, &tmax, &ymax, &y[2]);
76  if (result)
77  exit(EXIT_FAILURE);
78 
79  if (tmax != 1)
80  exit(EXIT_FAILURE);
81  if (ymax != 27)
82  exit(EXIT_FAILURE);
83  }
84 
85  {
86  /* Concave-up, decreasing series: maximum should occur at lower sample boundary */
87  const REAL4 y[] = {256, 27, 8, 1, 0};
88  REAL4 ymax;
89 
90  result = XLALREAL4ApplyCubicSplineTriggerInterpolant(interp, &tmax, &ymax, &y[2]);
91  if (result)
92  exit(EXIT_FAILURE);
93 
94  if (tmax != -1)
95  exit(EXIT_FAILURE);
96  if (ymax != 27)
97  exit(EXIT_FAILURE);
98  }
99 
100  {
101  /* Test some random sample data. */
102  const COMPLEX16 y[] = {0, 1.23412285-1.56122275*I, -2.06157986+0.78298714*I, -0.52811449+0.47859189*I, 0.42450302+2.06877714*I};
103  COMPLEX16 ymax;
104 
105  result = XLALCOMPLEX16ApplyCubicSplineTriggerInterpolant(interp, &tmax, &ymax, &y[2]);
106  if (result)
107  exit(EXIT_FAILURE);
108 
109  if (fabs(0.108106552865 - tmax) > 1e-8)
110  exit(EXIT_FAILURE);
111  if (fabs(-2.10041938439 - creal(ymax)) > 1e-8)
112  exit(EXIT_FAILURE);
113  if (fabs(0.85409051094 - cimag(ymax)) > 1e-8)
114  exit(EXIT_FAILURE);
115  }
116 
117  {
118  /* Test some random sample data. */
119  const COMPLEX8 y[] = {0, 1.23412285-1.56122275*I, -2.06157986+0.78298714*I, -0.52811449+0.47859189*I, 0.42450302+2.06877714*I};
120  COMPLEX8 ymax;
121 
122  result = XLALCOMPLEX8ApplyCubicSplineTriggerInterpolant(interp, &tmax, &ymax, &y[2]);
123  if (result)
124  exit(EXIT_FAILURE);
125 
126  if (fabs(0.108106552865 - tmax) > 1e-6)
127  exit(EXIT_FAILURE);
128  if (fabs(-2.10041938439 - crealf(ymax)) > 1e-6)
129  exit(EXIT_FAILURE);
130  if (fabs(0.85409051094 - cimagf(ymax)) > 1e-6)
131  exit(EXIT_FAILURE);
132  }
133 
135  exit(EXIT_SUCCESS);
136 }
int main(__attribute__((unused)) int argc, __attribute__((unused)) char **argv)
#define __attribute__(x)
Definition: getdate.c:146
int XLALCOMPLEX8ApplyCubicSplineTriggerInterpolant(CubicSplineTriggerInterpolant *interp, double *tmax, COMPLEX8 *ymax, const COMPLEX8 *y)
CubicSplineTriggerInterpolant * XLALCreateCubicSplineTriggerInterpolant(unsigned int window)
Constructor.
int XLALCOMPLEX16ApplyCubicSplineTriggerInterpolant(CubicSplineTriggerInterpolant *interp, double *t, COMPLEX16 *y, const COMPLEX16 *data)
Perform interpolation using the allocated workspace.
void XLALDestroyCubicSplineTriggerInterpolant(CubicSplineTriggerInterpolant *interp)
Destructor.
int XLALREAL4ApplyCubicSplineTriggerInterpolant(CubicSplineTriggerInterpolant *interp, double *tmax, REAL4 *ymax, const REAL4 *y)
int XLALREAL8ApplyCubicSplineTriggerInterpolant(CubicSplineTriggerInterpolant *interp, double *tmax, REAL8 *ymax, const REAL8 *y)
double complex COMPLEX16
Double-precision floating-point complex number (16 bytes total)
double REAL8
Double precision real floating-point number (8 bytes).
float complex COMPLEX8
Single-precision floating-point complex number (8 bytes total)
float REAL4
Single precision real floating-point number (4 bytes).