Deterministic Gaussian Sampling
Loading...
Searching...
No Matches
gm_to_dirac_approx_standard_normal_distribution_i.h
1#ifndef GM_TO_DIRAC_APPROX_STANDARD_NORMAL_DISTRIBUTION_I_H
2#define GM_TO_DIRAC_APPROX_STANDARD_NORMAL_DISTRIBUTION_I_H
3
4#include <gsl/gsl_matrix.h>
5#include <gsl/gsl_matrix_float.h>
6#include <gsl/gsl_matrix_long_double.h>
7
8#include <type_traits>
9
10#include "approximate_options.h"
11#include "gsl_minimizer.h"
12#include "gsl_vector_matrix_types.h"
13
19template <typename T>
21 public:
22 using GSLVectorType = typename GSLTemplateTypeAlias<T>::VectorType;
23 using GSLVectorViewType = typename GSLTemplateTypeAlias<T>::VectorViewType;
24 using GSLMatrixType = typename GSLTemplateTypeAlias<T>::MatrixType;
25
27
39 virtual bool approximate(size_t L, size_t N, size_t bMax, T* x, const T* wX,
41 const ApproximateOptions& options) = 0;
42
56 virtual void modified_van_mises_distance_sq(T* distance, size_t L, size_t N,
57 size_t bMax, T* x,
58 const T* wX) = 0;
59
74 size_t N, size_t bMax,
75 T* x, const T* wX) = 0;
76
89 virtual bool approximate(size_t L, size_t N, size_t bMax, GSLVectorType* x,
90 const GSLVectorType* wX, GslminimizerResult* result,
91 const ApproximateOptions& options) = 0;
92
106 virtual void modified_van_mises_distance_sq(T* distance, size_t L, size_t N,
107 size_t bMax, GSLVectorType* x,
108 const GSLVectorType* wX) = 0;
109
123 virtual void modified_van_mises_distance_sq_derivative(GSLVectorType* gradient, size_t L,
124 size_t N, size_t bMax,
125 GSLVectorType* x, const GSLVectorType* wX) = 0;
126
139 virtual bool approximate(size_t L, size_t N, size_t bMax, GSLMatrixType* x,
140 const GSLVectorType* wX, GslminimizerResult* result,
141 const ApproximateOptions& options) = 0;
142
156 virtual void modified_van_mises_distance_sq(T* distance, size_t L, size_t N,
157 size_t bMax, GSLMatrixType* x,
158 const GSLVectorType* wX) = 0;
159
173 virtual void modified_van_mises_distance_sq_derivative(GSLMatrixType* gradient, size_t L,
174 size_t N, size_t bMax,
175 GSLMatrixType* x, const GSLVectorType* wX) = 0;
176};
177
178#endif // GM_TO_DIRAC_APPROX_STANDARD_NORMAL_DISTRIBUTION_I_H
Definition dirac_to_dirac_approx_short.h:9
interface for the gausian mixture to dirac approximation
Definition gm_to_dirac_approx_standard_normal_distribution_i.h:20
virtual bool approximate(size_t L, size_t N, size_t bMax, GSLMatrixType *x, const GSLVectorType *wX, GslminimizerResult *result, const ApproximateOptions &options)=0
approximate using gsl matricies where possible
virtual bool approximate(size_t L, size_t N, size_t bMax, GSLVectorType *x, const GSLVectorType *wX, GslminimizerResult *result, const ApproximateOptions &options)=0
approximate using gsl vectors
virtual void modified_van_mises_distance_sq(T *distance, size_t L, size_t N, size_t bMax, GSLVectorType *x, const GSLVectorType *wX)=0
calculate modified van mises distance based on standard normal deviation and x
virtual void modified_van_mises_distance_sq(T *distance, size_t L, size_t N, size_t bMax, T *x, const T *wX)=0
calculate modified van mises distance based on standard normal deviation and x
virtual void modified_van_mises_distance_sq(T *distance, size_t L, size_t N, size_t bMax, GSLMatrixType *x, const GSLVectorType *wX)=0
calculate modified van mises distance based on standard normal deviation and x
virtual void modified_van_mises_distance_sq_derivative(GSLMatrixType *gradient, size_t L, size_t N, size_t bMax, GSLMatrixType *x, const GSLVectorType *wX)=0
calculate modified van mises distance based on standard normal deviation and x
virtual void modified_van_mises_distance_sq_derivative(T *gradient, size_t L, size_t N, size_t bMax, T *x, const T *wX)=0
calculate modified van mises distance based on standard normal deviation and x
virtual bool approximate(size_t L, size_t N, size_t bMax, T *x, const T *wX, GslminimizerResult *result, const ApproximateOptions &options)=0
approximate using raw pointers
virtual void modified_van_mises_distance_sq_derivative(GSLVectorType *gradient, size_t L, size_t N, size_t bMax, GSLVectorType *x, const GSLVectorType *wX)=0
calculate modified van mises distance based on standard normal deviation and x
Definition approximate_options.h:6
struct to hold the result of the minimization
Definition gsl_minimizer.h:32