-
Notifications
You must be signed in to change notification settings - Fork 0
Expand file tree
/
Copy pathdirac_to_dirac_approx_function_i.h
More file actions
192 lines (177 loc) · 7.58 KB
/
dirac_to_dirac_approx_function_i.h
File metadata and controls
192 lines (177 loc) · 7.58 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
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
126
127
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
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
#ifndef DIRAC_TO_DIRAC_APPROX_FUNCTION_I_H
#define DIRAC_TO_DIRAC_APPROX_FUNCTION_I_H
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_matrix_float.h>
#include <gsl/gsl_matrix_long_double.h>
#include <type_traits>
#include "approximate_options.h"
#include "gsl_minimizer.h"
#include "gsl_vector_matrix_types.h"
/**
* @brief interface for the gausian mixture to dirac approximation with a custom
* weight function
*
* @tparam T type of the vector (float, double)
*/
template <typename T>
class dirac_to_dirac_approx_function_i {
public:
using GSLVectorType = typename GSLTemplateTypeAlias<T>::VectorType;
using GSLVectorViewType = typename GSLTemplateTypeAlias<T>::VectorViewType;
using GSLMatrixType = typename GSLTemplateTypeAlias<T>::MatrixType;
using wXf = void (*)(const double* x, double* res, size_t L, size_t N);
using wXd = wXf; // (*)(const double* x, double* res, size_t L, size_t N);
virtual ~dirac_to_dirac_approx_function_i() = default;
/**
* @brief reduce the data points using raw pointers
*
* @param y input data points
* @param M number of input data points
* @param L number of data points for reduction
* @param N dimension of the data
* @param x first guess for the reduction and return value
* @param wXcallback callback for the weight function
* @param wXDcallback callback for the gradient of the weight function
* @param result minimizer result
* @param options options for minimizer
* @return true, on success, false otherwise
*/
virtual bool approximate(const T* y, size_t M, size_t L, size_t N, T* x,
wXf wXcallback, wXd wXDcallback,
GslminimizerResult* result,
const ApproximateOptions& options) = 0;
/**
* @brief calculate modified van mises distance based on x and y
*
* @param distance pointer to distance value to be calculated
* @param y input data points
* @param M number of elements in y
* @param L number of elements in x
* @param N dimension of the data
* @param bMax bMax
* @param x input data points
* @param wXcallback callback for the weight function
* @param wXDcallback callback for the gradient of the weight function
* @return true, on success, false otherwise
*/
virtual void modified_van_mises_distance_sq(T* distance, const T* y, size_t M,
size_t L, size_t N, size_t bMax,
T* x, wXf wXcallback,
wXd wXDcallback) = 0;
/**
* @brief calculate modified van mises distance based on x and y
*
* @param gradient pointer to gradient to be calculated
* @param y input data points
* @param M number of elements in y
* @param L number of elements in x
* @param N dimension of the data
* @param bMax bMax
* @param x input data points
* @param wXcallback callback for the weight function
* @param wXDcallback callback for the gradient of the weight function
* @return true, on success, false otherwise
*/
virtual void modified_van_mises_distance_sq_derivative(
T* gradient, const T* y, size_t M, size_t L, size_t N, size_t bMax, T* x,
wXf wXcallback, wXd wXDcallback) = 0;
/**
* @brief reduce the data points using gsl vectors
*
* @param y input data points
* @param L number of data points for reduction
* @param N dimension of the data
* @param x first guess for the reduction and return value
* @param wXcallback callback for the weight function
* @param wXDcallback callback for the gradient of the weight function
* @param result minimizer result
* @param options options for minimizer
* @return true, on success, false otherwise
*/
virtual bool approximate(const GSLVectorType* y, size_t L, size_t N,
GSLVectorType* x, wXf wXcallback, wXd wXDcallback,
GslminimizerResult* result,
const ApproximateOptions& options) = 0;
/**
* @brief calculate modified van mises distance based on x and y
*
* @param distance pointer to distance value to be calculated
* @param y input data points
* @param L number of elements in x
* @param N dimension of the data
* @param bMax bMax
* @param x input data points
* @param wXcallback callback for the weight function
* @param wXDcallback callback for the gradient of the weight function
* @return true, on success, false otherwise
*/
virtual void modified_van_mises_distance_sq(T* distance,
const GSLVectorType* y, size_t L,
size_t N, size_t bMax,
GSLVectorType* x, wXf wXcallback,
wXd wXDcallback) = 0;
/**
* @brief calculate modified van mises distance based on x and y
*
* @param gradient pointer to gradient to be calculated
* @param y input data points
* @param L number of elements in x
* @param N dimension of the data
* @param bMax bMax
* @param x input data points
* @param wXcallback callback for the weight function
* @param wXDcallback callback for the gradient of the weight function
* @return true, on success, false otherwise
*/
virtual void modified_van_mises_distance_sq_derivative(
GSLMatrixType* gradient, const GSLVectorType* y, size_t L, size_t N,
size_t bMax, GSLVectorType* x, wXf wXcallback, wXd wXDcallback) = 0;
/**
* @brief reduce the data points using gsl matricies where possible
*
* @param y input data points
* @param L number of data points for reduction
* @param x first guess for the reduction and return value
* @param wXcallback callback for the weight function
* @param wXDcallback callback for the gradient of the weight function
* @param result minimizer result
* @param options options for minimizer
* @return true, on success, false otherwise
*/
virtual bool approximate(GSLMatrixType* y, size_t L, GSLMatrixType* x,
wXf wXcallback, wXd wXDcallback,
GslminimizerResult* result,
const ApproximateOptions& options) = 0;
/**
* @brief calculate modified van mises distance based on x and y
*
* @param distance pointer to distance value to be calculated
* @param y input data points
* @param L number of elements in x
* @param bMax bMax
* @param x input data points
* @param wXcallback callback for the weight function
* @param wXDcallback callback for the gradient of the weight function
* @return true, on success, false otherwise
*/
virtual void modified_van_mises_distance_sq(T* distance, GSLMatrixType* y,
size_t L, size_t bMax,
GSLMatrixType* x, wXf wXcallback,
wXd wXDcallback) = 0;
/**
* @brief calculate modified van mises distance based on x and y
*
* @param gradient pointer to gradient to be calculated
* @param y input data points
* @param L number of elements in x
* @param bMax bMax
* @param x input data points
* @param wXcallback callback for the weight function
* @param wXDcallback callback for the gradient of the weight function
* @return true, on success, false otherwise
*/
virtual void modified_van_mises_distance_sq_derivative(
GSLMatrixType* gradient, GSLMatrixType* y, size_t L, size_t bMax,
GSLMatrixType* x, wXf wXcallback, wXd wXDcallback) = 0;
};
#endif // DIRAC_TO_DIRAC_APPROX_FUNCTION_I_H