Mercurial > hg > kwantix
diff include/interpolator.hpp @ 15:5144dd3c5468
Almost completed adding all Doxygen docstrings.
author | Jordi GutiƩrrez Hermoso <jordigh@gmail.com> |
---|---|
date | Fri, 08 Aug 2008 00:08:57 -0500 (2008-08-08) |
parents | e6d2cd3b6e77 |
children | 29a7b95c2805 |
line wrap: on
line diff
--- a/include/interpolator.hpp +++ b/include/interpolator.hpp @@ -1,3 +1,10 @@ +/*!\file interpolator.hpp + * + * \brief An interpolator \f$ \displaystyle u(x) = \sum_{\xi \in \Xi} + * \lambda_\xi \varphi_\xi(x) \f$, templated by RBFs, and its helper + * classes are declared here. + * + */ #ifndef __INTERPOLATOR_HPP__ #define __INTERPOLATOR_HPP__ @@ -15,6 +22,26 @@ class interp_values; + /*! \brief An interpolator \f$ \displaystyle u(x) = \sum_{\xi \in + * \Xi} \lambda_\xi \varphi_\xi(x) \f$, where the \f$ \varphi_\xi + * \f$ are RBFs. + * + * An interpolator \f$ \displaystyle u(x) = \sum_{\xi \in \Xi} + * \lambda_\xi \varphi_\xi(x) \f$, where the \f$\varphi_\xi\f$ are + * RBFs centred at \f$\xi\f$, and \f$\Xi\f$ is a set of points on + * the interior and boundary of some domain \f$\Omega\f$. + * + * The interpolator can be used in various ways, initialised by + * either a BVP::linear_bvp2 or by interpolation data + * (i.e. std::map< linalg::point, double>). Once an interpolator is + * initialised, it becomes a bona fide bvp::realfunc, so that its + * evaluations and derivatives become available. For certain + * problems, it is convenient to be able to evaluate the + * interpolator or its derivatives at all points of the domain. It + * is also useful to be able to interpolate again on the same domain + * or with a slightly modified BVP. Interfaces to both of these + * functions are provided. + */ template<typename RBF> class interpolator : public realfunc{ public: @@ -39,7 +66,11 @@ * at those points * * Must provide domain information. The values of Xi must match - * points on the given domain Omega. + * points on the given domain Omega. The reason why the domain + * might be relevant is in case it is important to specify which + * points are on the interior and which are on the boundary, so + * that it's possible to interpolate again by changing the values + * on either the interior or the boundary. */ interpolator(shared_ptr<domain> Omega, const map<point, double>& Xi); //@} @@ -81,7 +112,7 @@ interp_values operator()() const; ///Full domain evaluation interp_values at() const; - /// Full domain first erivative + /// Full domain first derivative interp_values d(size_t k) const; /// Full domain second derivatives interp_values d2(size_t k1, size_t k2) const; @@ -154,11 +185,11 @@ ///Is the interpolator ready for use? bool initted; - /** @name Exception throwers + /** @name Exception builders */ //@{ - void not_initted(int line, string file) const; - void not_precomputed(int line, string file) const; + error_handling::badArgument not_initted(int line, string file) const; + error_handling::badArgument not_precomputed(int line, string file) const; //@} ///Coefficients of the RBFs @@ -235,17 +266,23 @@ linalg::vector v; ///Two interp_values can be added iff these two values matches. size_t rbfs_hash; - ///Exception thrower - void different_rbfs(int line, string file) const; + ///Exception builder + badArgument different_rbfs(int line, string file) const; template<typename RBF> friend class interpolator; - friend interp_values operator*(const interp_values& v, const - interp_values& w); + + /// For comfortable syntax + friend interp_values operator*(double a, const interp_values& w); }; ///For comfortable syntax + //@{ template <typename RBF> - interpolator<RBF> operator*(double a, const interpolator<RBF>& u); + interpolator<RBF> operator*(double a, const interpolator<RBF>& u) + { + return u*a; + } + //@} } #endif