Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALInference 4.1.9.1-5e288d3
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
LALInferenceGenerateROQ.h
Go to the documentation of this file.
1/*
2 * LALInferenceROQ.h: Reduced order quadrature basis and interpolant generation
3 *
4 * Copyright (C) 2014 Matthew Pitkin, Rory Smith
5 *
6 *
7 * This program is free software; you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation; either version 2 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with with program; see the file COPYING. If not, write to the
19 * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
20 * MA 02110-1301 USA
21 */
22#ifndef LALINFERENCEROQ_H
23#define LALINFERENCEROQ_H
24
25#include <lal/LALStdlib.h>
26#include <lal/LALConstants.h>
27#include <lal/LALDatatypes.h>
28#include <lal/LALInference.h>
29#include <lal/XLALError.h>
30
31#include <gsl/gsl_vector.h>
32#include <gsl/gsl_matrix.h>
33#include <gsl/gsl_blas.h>
34#include <gsl/gsl_cblas.h>
35#include <gsl/gsl_linalg.h>
36#include <gsl/gsl_complex.h>
37
38#ifdef __cplusplus
39extern "C" {
40#endif
41
42
43/** A structure to hold a real (double precision) interpolant matrix and interpolation node indices */
44typedef struct tagLALInferenceREALROQInterpolant{
45 REAL8Array *B; /**< The interpolant matrix */
46 UINT4 *nodes; /**< The nodes (indices) for the interpolation */
48
49/** A structure to hold a complex (double precision) interpolant matrix and interpolation node indices */
50typedef struct tagLALInferenceCOMPLEXROQInterpolant{
51 COMPLEX16Array *B; /**< The interpolant matrix */
52 UINT4 *nodes; /**< The nodes (indices) for the interpolation */
54
55/* function to create or enrich a real orthonormal basis set from a training set of models */
57 const REAL8Vector *delta,
58 REAL8 tolerance,
59 REAL8Array **TS,
60 UINT4Vector **greedypoints);
61
63 const REAL8Vector *delta,
64 REAL8 tolerance,
65 COMPLEX16Array **TS,
66 UINT4Vector **greedypoints);
67
68/* functions to test the basis */
70 const REAL8Vector *delta,
71 const REAL8Array *RB,
72 REAL8Array **testmodels);
73
75 const REAL8Vector *delta,
76 const COMPLEX16Array *RB,
77 COMPLEX16Array **testmodels);
78
80 REAL8 tolerance,
81 const REAL8Array *RB,
82 REAL8Array **testmodels);
83
85 REAL8 tolerance,
86 const COMPLEX16Array *RB,
87 COMPLEX16Array **testmodels);
88
89/* functions to enrich the training model set and basis set */
91 const REAL8 tolerance,
92 REAL8Array **RB,
93 UINT4Vector **greedypoints,
94 const REAL8Array *testmodels,
95 REAL8Array **testmodelsnew);
96
98 const REAL8 tolerance,
99 COMPLEX16Array **RB,
100 UINT4Vector **greedypoints,
101 const COMPLEX16Array *testmodels,
102 COMPLEX16Array **testmodelsnew);
103
104/* function to create the empirical interpolant */
106
108
109/* create ROQ weights for interpolant to calculate the linear <d|h> terms */
111
112/* create ROQ weights for interpolant to calculate the quadratic model terms real(<h|h>) */
114
115/* create ROQ weights for interpolant to calculate the data dot model terms */
117
118/* calculate ROQ version of the dot product (where the model vector is the model just computed at the interpolant points) */
121
122/* memory destruction */
125
126#ifdef __cplusplus
127}
128#endif
129
130#endif
INT4 LALInferenceTestCOMPLEX16OrthonormalBasis(const REAL8Vector *delta, REAL8 tolerance, const COMPLEX16Array *RB, COMPLEX16Array **testmodels)
Test the reduced basis against another set of waveforms.
COMPLEX16 LALInferenceROQCOMPLEX16DotProduct(COMPLEX16Vector *weights, COMPLEX16Vector *model)
Calculate the dot product of two complex vectors using the ROQ iterpolant.
REAL8 LALInferenceEnrichCOMPLEX16Basis(const REAL8Vector *delta, const REAL8 tolerance, COMPLEX16Array **RB, UINT4Vector **greedypoints, const COMPLEX16Array *testmodels, COMPLEX16Array **testmodelsnew)
Expand the complex training waveforms with ones from a set of new training waveforms.
REAL8 LALInferenceGenerateCOMPLEX16OrthonormalBasis(COMPLEX16Array **RB, const REAL8Vector *delta, REAL8 tolerance, COMPLEX16Array **TS, UINT4Vector **greedypoints)
Create a orthonormal basis set from a training set of complex waveforms.
INT4 LALInferenceTestREAL8OrthonormalBasis(const REAL8Vector *delta, REAL8 tolerance, const REAL8Array *RB, REAL8Array **testmodels)
Test the reduced basis against another set of waveforms.
void LALInferenceRemoveCOMPLEXROQInterpolant(LALInferenceCOMPLEXROQInterpolant *a)
Free memory for a LALInferenceCOMPLEXROQInterpolant.
REAL8 LALInferenceGenerateREAL8OrthonormalBasis(REAL8Array **RB, const REAL8Vector *delta, REAL8 tolerance, REAL8Array **TS, UINT4Vector **greedypoints)
Create a orthonormal basis set from a training set of real waveforms.
COMPLEX16Vector * LALInferenceGenerateCOMPLEX16LinearWeights(COMPLEX16Array *B, COMPLEX16Vector *data, REAL8Vector *vars)
Create the weights for the ROQ interpolant for the complex data and model dot product <d|h>
REAL8Vector * LALInferenceGenerateREAL8LinearWeights(REAL8Array *B, REAL8Vector *data, REAL8Vector *vars)
Create the weights for the ROQ interpolant for the linear data and model dot product <d|h>
REAL8 LALInferenceROQREAL8DotProduct(REAL8Vector *weights, REAL8Vector *model)
Calculate the dot product of two vectors using the ROQ iterpolant.
REAL8Vector * LALInferenceGenerateQuadraticWeights(REAL8Array *B, REAL8Vector *vars)
Create the weights for the ROQ interpolant for the model quadratic model term real(<h|h>).
REAL8 LALInferenceEnrichREAL8Basis(const REAL8Vector *delta, const REAL8 tolerance, REAL8Array **RB, UINT4Vector **greedypoints, const REAL8Array *testmodels, REAL8Array **testmodelsnew)
Expand the real training waveforms with ones from a set of new training waveforms.
LALInferenceCOMPLEXROQInterpolant * LALInferenceGenerateCOMPLEXROQInterpolant(COMPLEX16Array *RB)
Create a complex empirical interpolant from a set of orthonormal basis functions.
void LALInferenceValidateCOMPLEX16OrthonormalBasis(REAL8Vector **projerr, const REAL8Vector *delta, const COMPLEX16Array *RB, COMPLEX16Array **testmodels)
Validate the complex reduced basis against another set of waveforms.
void LALInferenceValidateREAL8OrthonormalBasis(REAL8Vector **projerr, const REAL8Vector *delta, const REAL8Array *RB, REAL8Array **testmodels)
Validate the real reduced basis against another set of waveforms.
void LALInferenceRemoveREALROQInterpolant(LALInferenceREALROQInterpolant *a)
Free memory for a LALInferenceREALROQInterpolant.
LALInferenceREALROQInterpolant * LALInferenceGenerateREALROQInterpolant(REAL8Array *RB)
Create a real empirical interpolant from a set of orthonormal basis functions.
double complex COMPLEX16
double REAL8
uint32_t UINT4
int32_t INT4
static const INT4 a
A structure to hold a complex (double precision) interpolant matrix and interpolation node indices.
COMPLEX16Array * B
The interpolant matrix.
UINT4 * nodes
The nodes (indices) for the interpolation.
A structure to hold a real (double precision) interpolant matrix and interpolation node indices.
UINT4 * nodes
The nodes (indices) for the interpolation.
REAL8Array * B
The interpolant matrix.