Flow123d
release_3.0.0-973-g92f55e826
|
Go to the documentation of this file.
19 #ifndef QUADRATURE_HH_
20 #define QUADRATURE_HH_
45 template<
unsigned int dim>
52 Quadrature(
const unsigned int n_quadrature_points = 0);
71 void resize(
const unsigned int n_q_points);
74 const unsigned int size()
const;
77 const arma::vec::fixed<dim> &
point(
const unsigned int i)
const;
87 void set_point(
const unsigned int i,
const arma::vec::fixed<dim> &p);
90 double weight(
const unsigned int i)
const;
96 void set_weight(
const unsigned int i,
const double w);
117 template<
unsigned int dim>
123 template<
unsigned int dim>
125 quadrature_points(q.quadrature_points),
129 template<
unsigned int dim>
132 arma::vec::fixed<dim> v;
134 quadrature_points.resize(n_q, v);
135 weights.resize(n_q, 0);
138 template<
unsigned int dim>
140 return weights.size();
143 template<
unsigned int dim>
145 const unsigned int i)
const {
146 return quadrature_points[i];
149 template<
unsigned int dim>
151 return quadrature_points;
154 template<
unsigned int dim>
157 quadrature_points[i] = p;
160 template<
unsigned int dim>
165 template<
unsigned int dim>
170 template<
unsigned int dim>
176 template<
unsigned int dim>
180 template<
unsigned int dim>
inline
217 arma::vec::fixed<dim+1> el_bar_coords, final_bar;
219 for (
unsigned int k=0; k<subq.
size(); k++)
223 arma::vec::fixed<dim> pp;
226 for (
unsigned int i=0; i<RefElement<dim>::n_nodes_per_side; i++) {
234 weights[k] = subq.
weight(k);
242 arma::vec::fixed<1> p({(double)sid});
std::vector< arma::vec::fixed< dim > > quadrature_points
List of quadrature points.
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc.
const unsigned int size() const
Returns number of quadrature points.
void set_weight(const unsigned int i, const double w)
Sets individual quadrature weight.
double weight(const unsigned int i) const
Returns the ith weight.
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
static LocalPoint bary_to_local(const BaryPoint &bp)
Converts from barycentric to local coordinates.
std::vector< double > weights
List of weights to the quadrature points.
Quadrature(const unsigned int n_quadrature_points=0)
Constructor.
void set_point(const unsigned int i, const arma::vec::fixed< dim > &p)
Sets individual quadrature point coordinates.
const std::vector< arma::vec::fixed< dim > > & get_points() const
Return a reference to the whole array of quadrature points.
virtual ~Quadrature()
Virtual destructor.
static BaryPoint local_to_bary(const LocalPoint &lp)
Converts from local to barycentric coordinates.
void resize(const unsigned int n_q_points)
Modify the number of quadrature points.
Base class for quadrature rules on simplices in arbitrary dimensions.
const std::vector< double > & get_weights() const
Return a reference to the whole array of weights.