OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
TriangleQuadratureRule.hpp
1// Note that this code is mostly from the Drake project: https://drake.mit.edu/
2// and is governed by their BSD 3-Clause license:
3// https://github.com/RobotLocomotion/drake/blob/master/LICENSE.TXT
4// Only minor modifications have been made to make it work in OpenMD.
5
6#ifndef MATH_INTEGRATION_TRIANGLEQUADRATURERULE_HPP
7#define MATH_INTEGRATION_TRIANGLEQUADRATURERULE_HPP
8
9#include <vector>
10
11#include "math/Vector2.hpp"
12#include "math/Vector3.hpp"
13
14namespace OpenMD {
15
16 /// A "rule" (weights and quadrature points) for computing quadrature over
17 /// triangular domains.
19 public:
21 TriangleQuadratureRule& operator=(const TriangleQuadratureRule&) = default;
23 TriangleQuadratureRule& operator=(TriangleQuadratureRule&&) = default;
24 TriangleQuadratureRule() = default;
25 virtual ~TriangleQuadratureRule() {}
26
27 /// Returns the order of this rule.
28 int order() const {
29 int rule_order = do_order();
30 assert(rule_order >= 1);
31 return rule_order;
32 }
33
34 /// Returns the vector of quadrature points. These are returned as
35 /// the first two barycentric coordinates b0 b1; the third is just
36 /// b2 = 1 - b0 - b1. Each of these has a corresponding weight
37 /// returned by weights().
38 const std::vector<Vector2d>& quadrature_points() const {
39 return do_quadrature_points();
40 }
41
42 /// Returns the vector of weights. These sum to 1 and there is one
43 /// weight for each point returned by quadrature_points().
44 const std::vector<RealType>& weights() const { return do_weights(); }
45
46 protected:
47 /// Derived classes shall return the order (>= 1) of this rule.
48 virtual int do_order() const = 0;
49
50 /// Derived classes shall return the vector of quadrature
51 /// points. Each of these Vector2d objects represents
52 /// the barycentric coordinates of a triangle (the third
53 /// barycentric coordinate is implicit: it is the difference
54 /// between unity and the sum of the other two coordinates).
55 virtual const std::vector<Vector2d>& do_quadrature_points() const = 0;
56
57 /// Derived classes shall return the vector of weights. The sum of
58 /// all weights must equal 1.0.
59 virtual const std::vector<RealType>& do_weights() const = 0;
60 };
61} // namespace OpenMD
62#endif
A "rule" (weights and quadrature points) for computing quadrature over triangular domains.
int order() const
Returns the order of this rule.
const std::vector< RealType > & weights() const
Returns the vector of weights. These sum to 1 and there is one weight for each point returned by quad...
virtual int do_order() const =0
Derived classes shall return the order (>= 1) of this rule.
virtual const std::vector< Vector2d > & do_quadrature_points() const =0
Derived classes shall return the vector of quadrature points. Each of these Vector2d objects represen...
virtual const std::vector< RealType > & do_weights() const =0
Derived classes shall return the vector of weights. The sum of all weights must equal 1....
const std::vector< Vector2d > & quadrature_points() const
Returns the vector of quadrature points. These are returned as the first two barycentric coordinates ...
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.