Medusa  1.1
Coordinate Free Mehless Method implementation
Gaussian_fwd.hpp
Go to the documentation of this file.
1 #ifndef MEDUSA_BITS_APPROXIMATIONS_GAUSSIAN_FWD_HPP_
2 #define MEDUSA_BITS_APPROXIMATIONS_GAUSSIAN_FWD_HPP_
3 
11 #include <medusa/Config.hpp>
12 #include <iosfwd>
13 #include <cmath>
14 #include "Operators_fwd.hpp"
15 namespace mm {
16 
28 template <class scal_t>
29 class Gaussian {
30  public:
31  typedef scal_t scalar_t;
32 
33  private:
35 
36  public:
38  Gaussian(scalar_t shape = 1.0);
39 
41  scalar_t shape() const { return shape_; }
42 
45 
52  inline scalar_t operator()(scalar_t r2, int derivative) const;
53 
59  template <int dimension>
60  inline scalar_t operator()(scalar_t r2, Lap<dimension> lap) const;
61 
63  inline scalar_t operator()(scalar_t r2) const { return std::exp(-r2/shape_/shape_); }
64 
66  template <typename V>
67  friend std::ostream& operator<<(std::ostream& os, const Gaussian<V>& m);
68 };
69 
70 } // namespace mm
71 
72 #endif // MEDUSA_BITS_APPROXIMATIONS_GAUSSIAN_FWD_HPP_
mm
Root namespace for the whole library.
Definition: Gaussian.hpp:14
mm::Gaussian::shape_
scalar_t shape_
Shape parameter.
Definition: Gaussian_fwd.hpp:34
scalar_t
Scalar scalar_t
Type of the elements, alias of Scalar.
Definition: MatrixBaseAddons.hpp:16
mm::Gaussian::operator()
scalar_t operator()(scalar_t r2) const
Evaluate this RBF given squared radial distance.
Definition: Gaussian_fwd.hpp:63
mm::sh::lap
static const shape_flags lap
Indicates to calculate laplace shapes.
Definition: shape_flags.hpp:24
mm::Gaussian::operator<<
friend std::ostream & operator<<(std::ostream &os, const Gaussian< V > &m)
Output basic information about given Gaussian RBF.
mm::Gaussian::setShape
void setShape(scalar_t shape)
Sets shape parameter to a new value.
Definition: Gaussian_fwd.hpp:44
mm::Gaussian::shape
scalar_t shape() const
Returns shape parameter.
Definition: Gaussian_fwd.hpp:41
mm::Gaussian::Gaussian
Gaussian(scalar_t shape=1.0)
Creates a Gaussian RBF with shape parameter shape. The shape should be positive.
Definition: Gaussian.hpp:17
mm::Gaussian::scalar_t
scal_t scalar_t
Scalar type used for computations.
Definition: Gaussian_fwd.hpp:31
Config.hpp
mm::Gaussian::operator()
scalar_t operator()(scalar_t r2, int derivative) const
Evaluate derivative of this RBF wrt.
Definition: Gaussian.hpp:22
Operators_fwd.hpp
mm::Gaussian
Gaussian Radial Basis Function.
Definition: Gaussian_fwd.hpp:29
mm::Lap
Represents the Laplacian operator.
Definition: Monomials_fwd.hpp:20