#include <deal.II/base/polynomial.h>
#include <deal.II/base/polynomial_space.h>
#include <deal.II/base/tensor.h>
+#include <deal.II/base/tensor_polynomials_base.h>
#include <deal.II/base/tensor_product_polynomials.h>
#include <deal.II/base/thread_management.h>
* @date 2006
*/
template <int dim>
-class PolynomialsABF
+class PolynomialsABF : public TensorPolynomialsBase<dim>
{
public:
/**
std::vector<Tensor<2, dim>> &grads,
std::vector<Tensor<3, dim>> &grad_grads,
std::vector<Tensor<4, dim>> &third_derivatives,
- std::vector<Tensor<5, dim>> &fourth_derivatives) const;
-
- /**
- * Return the number of ABF polynomials.
- */
- unsigned int
- n() const;
-
- /**
- * Return the degree of the ABF space, which is two less than the highest
- * polynomial degree.
- */
- unsigned int
- degree() const;
+ std::vector<Tensor<5, dim>> &fourth_derivatives) const override;
/**
* Return the name of the space, which is <tt>ABF</tt>.
*/
std::string
- name() const;
+ name() const override;
/**
* Return the number of polynomials in the space <tt>RT(degree)</tt> without
static unsigned int
compute_n_pols(unsigned int degree);
-private:
/**
- * The degree of this object as given to the constructor.
+ * @copydoc TensorPolynomialsBase<dim>::clone()
*/
- const unsigned int my_degree;
+ virtual std::unique_ptr<TensorPolynomialsBase<dim>>
+ clone() const override;
+private:
/**
* An object representing the polynomial space for a single component. We
* can re-use it for the other vector components by rotating the
*/
const AnisotropicPolynomials<dim> polynomial_space;
- /**
- * Number of Raviart-Thomas polynomials.
- */
- unsigned int n_pols;
-
/**
* A mutex that guards the following scratch arrays.
*/
};
-template <int dim>
-inline unsigned int
-PolynomialsABF<dim>::n() const
-{
- return n_pols;
-}
-
-
-template <int dim>
-inline unsigned int
-PolynomialsABF<dim>::degree() const
-{
- return my_degree;
-}
-
-
template <int dim>
inline std::string
PolynomialsABF<dim>::name() const
#include <deal.II/base/polynomial.h>
#include <deal.II/base/polynomial_space.h>
#include <deal.II/base/tensor.h>
+#include <deal.II/base/tensor_polynomials_base.h>
#include <deal.II/base/thread_management.h>
#include <vector>
* @date 2003, 2005, 2009
*/
template <int dim>
-class PolynomialsBDM
+class PolynomialsBDM : public TensorPolynomialsBase<dim>
{
public:
/**
std::vector<Tensor<2, dim>> &grads,
std::vector<Tensor<3, dim>> &grad_grads,
std::vector<Tensor<4, dim>> &third_derivatives,
- std::vector<Tensor<5, dim>> &fourth_derivatives) const;
-
- /**
- * Return the number of BDM polynomials.
- */
- unsigned int
- n() const;
+ std::vector<Tensor<5, dim>> &fourth_derivatives) const override;
/**
* Return the degree of the BDM space, which is one less than the highest
* Return the name of the space, which is <tt>BDM</tt>.
*/
std::string
- name() const;
+ name() const override;
/**
* Return the number of polynomials in the space <tt>BDM(degree)</tt>
static unsigned int
compute_n_pols(unsigned int degree);
+ /**
+ * @copydoc TensorPolynomialsBase<dim>::clone()
+ */
+ virtual std::unique_ptr<TensorPolynomialsBase<dim>>
+ clone() const override;
+
private:
/**
* An object representing the polynomial space used here. The constructor
*/
std::vector<Polynomials::Polynomial<double>> monomials;
- /**
- * Number of BDM polynomials.
- */
- unsigned int n_pols;
-
/**
* A mutex that guards the following scratch arrays.
*/
};
-template <int dim>
-inline unsigned int
-PolynomialsBDM<dim>::n() const
-{
- return n_pols;
-}
-
-
template <int dim>
inline unsigned int
PolynomialsBDM<dim>::degree() const
#include <deal.II/base/polynomial.h>
#include <deal.II/base/polynomial_space.h>
#include <deal.II/base/tensor.h>
+#include <deal.II/base/tensor_polynomials_base.h>
#include <deal.II/base/tensor_product_polynomials.h>
#include <vector>
* @date 2018
*/
template <int dim>
-class PolynomialsBernardiRaugel
+class PolynomialsBernardiRaugel : public TensorPolynomialsBase<dim>
{
public:
/**
*/
PolynomialsBernardiRaugel(const unsigned int k);
- /**
- * Return the number of Bernardi-Raugel polynomials.
- */
- unsigned int
- n() const;
-
-
- /**
- * Return the degree of Bernardi-Raugel polynomials.
- * Since the bubble functions are quadratic in at least one variable,
- * the degree of the Bernardi-Raugel polynomials is two.
- */
- unsigned int
- degree() const;
-
/**
* Return the name of the space, which is <tt>BernardiRaugel</tt>.
*/
std::string
- name() const;
+ name() const override;
/**
* Compute the value and derivatives of each Bernardi-Raugel
std::vector<Tensor<2, dim>> &grads,
std::vector<Tensor<3, dim>> &grad_grads,
std::vector<Tensor<4, dim>> &third_derivatives,
- std::vector<Tensor<5, dim>> &fourth_derivatives) const;
+ std::vector<Tensor<5, dim>> &fourth_derivatives) const override;
/**
* Return the number of polynomials in the space <tt>BR(degree)</tt> without
static unsigned int
compute_n_pols(const unsigned int k);
-private:
- /**
- * The degree of this object given to the constructor (must be 1).
- */
- const unsigned int my_degree;
-
/**
- * The number of Bernardi-Raugel polynomials.
+ * @copydoc TensorPolynomialsBase<dim>::clone()
*/
- const unsigned int n_pols;
+ virtual std::unique_ptr<TensorPolynomialsBase<dim>>
+ clone() const override;
+private:
/**
* An object representing the polynomial space of Q
* functions which forms the <tt>BR</tt> polynomials through
};
-template <int dim>
-inline unsigned int
-PolynomialsBernardiRaugel<dim>::n() const
-{
- return n_pols;
-}
-
-
-
-template <int dim>
-inline unsigned int
-PolynomialsBernardiRaugel<dim>::degree() const
-{
- return 2;
-}
-
-
-
template <int dim>
inline std::string
PolynomialsBernardiRaugel<dim>::name() const
#include <deal.II/base/polynomial.h>
#include <deal.II/base/polynomial_space.h>
#include <deal.II/base/tensor.h>
+#include <deal.II/base/tensor_polynomials_base.h>
#include <deal.II/base/tensor_product_polynomials.h>
#include <vector>
* @date 2009, 2010
*/
template <int dim>
-class PolynomialsNedelec
+class PolynomialsNedelec : public TensorPolynomialsBase<dim>
{
public:
/**
std::vector<Tensor<2, dim>> &grads,
std::vector<Tensor<3, dim>> &grad_grads,
std::vector<Tensor<4, dim>> &third_derivatives,
- std::vector<Tensor<5, dim>> &fourth_derivatives) const;
-
- /**
- * Return the number of Nédélec polynomials.
- */
- unsigned int
- n() const;
-
- /**
- * Return the degree of the Nédélec space, which is one less than the
- * highest polynomial degree.
- */
- unsigned int
- degree() const;
+ std::vector<Tensor<5, dim>> &fourth_derivatives) const override;
/**
* Return the name of the space, which is <tt>Nedelec</tt>.
*/
std::string
- name() const;
+ name() const override;
/**
* Return the number of polynomials in the space <tt>N(degree)</tt> without
static unsigned int
compute_n_pols(unsigned int degree);
-private:
/**
- * The degree of this object as given to the constructor.
+ * @copydoc TensorPolynomialsBase<dim>::clone()
*/
- const unsigned int my_degree;
+ virtual std::unique_ptr<TensorPolynomialsBase<dim>>
+ clone() const override;
+private:
/**
* An object representing the polynomial space for a single component. We
* can re-use it by rotating the coordinates of the evaluation point.
*/
const AnisotropicPolynomials<dim> polynomial_space;
- /**
- * Number of Nédélec polynomials.
- */
- const unsigned int n_pols;
-
/**
* A static member function that creates the polynomial space we use to
* initialize the #polynomial_space member variable.
};
-template <int dim>
-inline unsigned int
-PolynomialsNedelec<dim>::n() const
-{
- return n_pols;
-}
-
-
-template <int dim>
-inline unsigned int
-PolynomialsNedelec<dim>::degree() const
-{
- return my_degree;
-}
-
-
template <int dim>
inline std::string
PolynomialsNedelec<dim>::name() const
#include <deal.II/base/polynomial.h>
#include <deal.II/base/polynomial_space.h>
#include <deal.II/base/tensor.h>
+#include <deal.II/base/tensor_polynomials_base.h>
#include <deal.II/base/tensor_product_polynomials.h>
#include <vector>
* @date 2005
*/
template <int dim>
-class PolynomialsRaviartThomas
+class PolynomialsRaviartThomas : public TensorPolynomialsBase<dim>
{
public:
/**
std::vector<Tensor<2, dim>> &grads,
std::vector<Tensor<3, dim>> &grad_grads,
std::vector<Tensor<4, dim>> &third_derivatives,
- std::vector<Tensor<5, dim>> &fourth_derivatives) const;
-
- /**
- * Return the number of Raviart-Thomas polynomials.
- */
- unsigned int
- n() const;
-
- /**
- * Return the degree of the Raviart-Thomas space, which is one less than
- * the highest polynomial degree.
- */
- unsigned int
- degree() const;
+ std::vector<Tensor<5, dim>> &fourth_derivatives) const override;
/**
* Return the name of the space, which is <tt>RaviartThomas</tt>.
*/
std::string
- name() const;
+ name() const override;
/**
* Return the number of polynomials in the space <tt>RT(degree)</tt> without
static unsigned int
compute_n_pols(unsigned int degree);
-private:
/**
- * The degree of this object as given to the constructor.
+ * @copydoc TensorPolynomialsBase<dim>::clone()
*/
- const unsigned int my_degree;
+ virtual std::unique_ptr<TensorPolynomialsBase<dim>>
+ clone() const override;
+private:
/**
* An object representing the polynomial space for a single component. We
* can re-use it by rotating the coordinates of the evaluation point.
*/
const AnisotropicPolynomials<dim> polynomial_space;
- /**
- * Number of Raviart-Thomas polynomials.
- */
- const unsigned int n_pols;
-
/**
* A static member function that creates the polynomial space we use to
* initialize the #polynomial_space member variable.
};
-
-template <int dim>
-inline unsigned int
-PolynomialsRaviartThomas<dim>::n() const
-{
- return n_pols;
-}
-
-
-
-template <int dim>
-inline unsigned int
-PolynomialsRaviartThomas<dim>::degree() const
-{
- return my_degree;
-}
-
-
-
template <int dim>
inline std::string
PolynomialsRaviartThomas<dim>::name() const
#include <deal.II/base/polynomial.h>
#include <deal.II/base/polynomials_raviart_thomas.h>
#include <deal.II/base/tensor.h>
+#include <deal.II/base/tensor_polynomials_base.h>
#include <vector>
*/
template <int dim>
-class PolynomialsRT_Bubbles
+class PolynomialsRT_Bubbles : public TensorPolynomialsBase<dim>
{
public:
/**
std::vector<Tensor<2, dim>> &grads,
std::vector<Tensor<3, dim>> &grad_grads,
std::vector<Tensor<4, dim>> &third_derivatives,
- std::vector<Tensor<5, dim>> &fourth_derivatives) const;
-
- /**
- * Returns the number of enhanced Raviart-Thomas polynomials.
- */
- unsigned int
- n() const;
-
- /**
- * Returns the degree of the RT_bubble space, which is one less than the
- * highest polynomial degree.
- */
- unsigned int
- degree() const;
+ std::vector<Tensor<5, dim>> &fourth_derivatives) const override;
/**
* Return the name of the space, which is <tt>RT_Bubbles</tt>.
*/
std::string
- name() const;
+ name() const override;
/**
* Return the number of polynomials in the space <tt>RT_Bubbles(degree)</tt>
static unsigned int
compute_n_pols(const unsigned int degree);
-private:
/**
- * The degree of this object as given to the constructor.
+ * @copydoc TensorPolynomialsBase<dim>::clone()
*/
- const unsigned int my_degree;
+ virtual std::unique_ptr<TensorPolynomialsBase<dim>>
+ clone() const override;
+private:
/**
* An object representing the Raviart-Thomas part of the space
*/
* to <i>k+1</i>.
*/
std::vector<Polynomials::Polynomial<double>> monomials;
-
- /**
- * Number of RT_Bubbles polynomials.
- */
- unsigned int n_pols;
};
-
-template <int dim>
-inline unsigned int
-PolynomialsRT_Bubbles<dim>::n() const
-{
- return n_pols;
-}
-
-
-
-template <int dim>
-inline unsigned int
-PolynomialsRT_Bubbles<dim>::degree() const
-{
- return my_degree;
-}
-
-
-
template <int dim>
inline std::string
PolynomialsRT_Bubbles<dim>::name() const
#include <deal.II/base/polynomials_abf.h>
#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/base/std_cxx14/memory.h>
#include <iomanip>
#include <iostream>
template <int dim>
PolynomialsABF<dim>::PolynomialsABF(const unsigned int k)
- : my_degree(k)
+ : TensorPolynomialsBase<dim>(k, compute_n_pols(k))
, polynomial_space(get_abf_polynomials<dim>(k))
- , n_pols(compute_n_pols(k))
{
// check that the dimensions match. we only store one of the 'dim'
// anisotropic polynomials that make up the vector-valued space, so
std::vector<Tensor<4, dim>> &third_derivatives,
std::vector<Tensor<5, dim>> &fourth_derivatives) const
{
- Assert(values.size() == n_pols || values.size() == 0,
- ExcDimensionMismatch(values.size(), n_pols));
- Assert(grads.size() == n_pols || grads.size() == 0,
- ExcDimensionMismatch(grads.size(), n_pols));
- Assert(grad_grads.size() == n_pols || grad_grads.size() == 0,
- ExcDimensionMismatch(grad_grads.size(), n_pols));
- Assert(third_derivatives.size() == n_pols || third_derivatives.size() == 0,
- ExcDimensionMismatch(third_derivatives.size(), n_pols));
- Assert(fourth_derivatives.size() == n_pols || fourth_derivatives.size() == 0,
- ExcDimensionMismatch(fourth_derivatives.size(), n_pols));
+ Assert(values.size() == this->n() || values.size() == 0,
+ ExcDimensionMismatch(values.size(), this->n()));
+ Assert(grads.size() == this->n() || grads.size() == 0,
+ ExcDimensionMismatch(grads.size(), this->n()));
+ Assert(grad_grads.size() == this->n() || grad_grads.size() == 0,
+ ExcDimensionMismatch(grad_grads.size(), this->n()));
+ Assert(third_derivatives.size() == this->n() || third_derivatives.size() == 0,
+ ExcDimensionMismatch(third_derivatives.size(), this->n()));
+ Assert(fourth_derivatives.size() == this->n() ||
+ fourth_derivatives.size() == 0,
+ ExcDimensionMismatch(fourth_derivatives.size(), this->n()));
const unsigned int n_sub = polynomial_space.n();
// guard access to the scratch
}
+template <int dim>
+std::unique_ptr<TensorPolynomialsBase<dim>>
+PolynomialsABF<dim>::clone() const
+{
+ return std_cxx14::make_unique<PolynomialsABF<dim>>(*this);
+}
+
+
template class PolynomialsABF<1>;
template class PolynomialsABF<2>;
template class PolynomialsABF<3>;
#include <deal.II/base/polynomial_space.h>
#include <deal.II/base/polynomials_bdm.h>
#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/base/std_cxx14/memory.h>
#include <iomanip>
#include <iostream>
template <int dim>
PolynomialsBDM<dim>::PolynomialsBDM(const unsigned int k)
- : polynomial_space(Polynomials::Legendre::generate_complete_basis(k))
+ : TensorPolynomialsBase<dim>(k, compute_n_pols(k))
+ , polynomial_space(Polynomials::Legendre::generate_complete_basis(k))
, monomials((dim == 2) ? (1) : (k + 2))
- , n_pols(compute_n_pols(k))
, p_values(polynomial_space.n())
, p_grads(polynomial_space.n())
, p_grad_grads(polynomial_space.n())
std::vector<Tensor<4, dim>> &third_derivatives,
std::vector<Tensor<5, dim>> &fourth_derivatives) const
{
- Assert(values.size() == n_pols || values.size() == 0,
- ExcDimensionMismatch(values.size(), n_pols));
- Assert(grads.size() == n_pols || grads.size() == 0,
- ExcDimensionMismatch(grads.size(), n_pols));
- Assert(grad_grads.size() == n_pols || grad_grads.size() == 0,
- ExcDimensionMismatch(grad_grads.size(), n_pols));
- Assert(third_derivatives.size() == n_pols || third_derivatives.size() == 0,
- ExcDimensionMismatch(third_derivatives.size(), n_pols));
- Assert(fourth_derivatives.size() == n_pols || fourth_derivatives.size() == 0,
- ExcDimensionMismatch(fourth_derivatives.size(), n_pols));
+ Assert(values.size() == this->n() || values.size() == 0,
+ ExcDimensionMismatch(values.size(), this->n()));
+ Assert(grads.size() == this->n() || grads.size() == 0,
+ ExcDimensionMismatch(grads.size(), this->n()));
+ Assert(grad_grads.size() == this->n() || grad_grads.size() == 0,
+ ExcDimensionMismatch(grad_grads.size(), this->n()));
+ Assert(third_derivatives.size() == this->n() || third_derivatives.size() == 0,
+ ExcDimensionMismatch(third_derivatives.size(), this->n()));
+ Assert(fourth_derivatives.size() == this->n() ||
+ fourth_derivatives.size() == 0,
+ ExcDimensionMismatch(fourth_derivatives.size(), this->n()));
// third and fourth derivatives not implemented
(void)third_derivatives;
grad_grads[start + 2][2][2][2] = 0.;
}
}
- Assert(start == n_pols, ExcInternalError());
+ Assert(start == this->n(), ExcInternalError());
}
}
}
+template <int dim>
+std::unique_ptr<TensorPolynomialsBase<dim>>
+PolynomialsBDM<dim>::clone() const
+{
+ return std_cxx14::make_unique<PolynomialsBDM<dim>>(*this);
+}
+
+
template class PolynomialsBDM<1>;
template class PolynomialsBDM<2>;
template class PolynomialsBDM<3>;
#include <deal.II/base/polynomials_bernardi_raugel.h>
+#include <deal.II/base/std_cxx14/memory.h>
DEAL_II_NAMESPACE_OPEN
template <int dim>
PolynomialsBernardiRaugel<dim>::PolynomialsBernardiRaugel(const unsigned int k)
- : my_degree(k)
- , n_pols(compute_n_pols(k))
+ : TensorPolynomialsBase<dim>(k + 1, compute_n_pols(k))
, polynomial_space_Q(create_polynomials_Q())
, polynomial_space_bubble(create_polynomials_bubble())
{}
std::vector<Tensor<4, dim>> &third_derivatives,
std::vector<Tensor<5, dim>> &fourth_derivatives) const
{
- Assert(values.size() == n_pols || values.size() == 0,
- ExcDimensionMismatch(values.size(), n_pols));
- Assert(grads.size() == n_pols || grads.size() == 0,
- ExcDimensionMismatch(grads.size(), n_pols));
- Assert(grad_grads.size() == n_pols || grad_grads.size() == 0,
- ExcDimensionMismatch(grad_grads.size(), n_pols));
- Assert(third_derivatives.size() == n_pols || third_derivatives.size() == 0,
- ExcDimensionMismatch(third_derivatives.size(), n_pols));
- Assert(fourth_derivatives.size() == n_pols || fourth_derivatives.size() == 0,
- ExcDimensionMismatch(fourth_derivatives.size(), n_pols));
+ Assert(values.size() == this->n() || values.size() == 0,
+ ExcDimensionMismatch(values.size(), this->n()));
+ Assert(grads.size() == this->n() || grads.size() == 0,
+ ExcDimensionMismatch(grads.size(), this->n()));
+ Assert(grad_grads.size() == this->n() || grad_grads.size() == 0,
+ ExcDimensionMismatch(grad_grads.size(), this->n()));
+ Assert(third_derivatives.size() == this->n() || third_derivatives.size() == 0,
+ ExcDimensionMismatch(third_derivatives.size(), this->n()));
+ Assert(fourth_derivatives.size() == this->n() ||
+ fourth_derivatives.size() == 0,
+ ExcDimensionMismatch(fourth_derivatives.size(), this->n()));
std::vector<double> Q_values;
std::vector<Tensor<1, dim>> Q_grads;
return 0;
}
+
+template <int dim>
+std::unique_ptr<TensorPolynomialsBase<dim>>
+PolynomialsBernardiRaugel<dim>::clone() const
+{
+ return std_cxx14::make_unique<PolynomialsBernardiRaugel<dim>>(*this);
+}
+
template class PolynomialsBernardiRaugel<1>; // to prevent errors
template class PolynomialsBernardiRaugel<2>;
template class PolynomialsBernardiRaugel<3>;
#include <deal.II/base/polynomial.h>
#include <deal.II/base/polynomials_nedelec.h>
#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/base/std_cxx14/memory.h>
#include <iomanip>
#include <iostream>
template <int dim>
PolynomialsNedelec<dim>::PolynomialsNedelec(const unsigned int k)
- : my_degree(k)
+ : TensorPolynomialsBase<dim>(k, compute_n_pols(k))
, polynomial_space(create_polynomials(k))
- , n_pols(compute_n_pols(k))
{}
template <int dim>
std::vector<Tensor<4, dim>> &third_derivatives,
std::vector<Tensor<5, dim>> &fourth_derivatives) const
{
- Assert(values.size() == n_pols || values.size() == 0,
- ExcDimensionMismatch(values.size(), n_pols));
- Assert(grads.size() == n_pols || grads.size() == 0,
- ExcDimensionMismatch(grads.size(), n_pols));
- Assert(grad_grads.size() == n_pols || grad_grads.size() == 0,
- ExcDimensionMismatch(grad_grads.size(), n_pols));
- Assert(third_derivatives.size() == n_pols || third_derivatives.size() == 0,
- ExcDimensionMismatch(third_derivatives.size(), n_pols));
- Assert(fourth_derivatives.size() == n_pols || fourth_derivatives.size() == 0,
- ExcDimensionMismatch(fourth_derivatives.size(), n_pols));
+ Assert(values.size() == this->n() || values.size() == 0,
+ ExcDimensionMismatch(values.size(), this->n()));
+ Assert(grads.size() == this->n() || grads.size() == 0,
+ ExcDimensionMismatch(grads.size(), this->n()));
+ Assert(grad_grads.size() == this->n() || grad_grads.size() == 0,
+ ExcDimensionMismatch(grad_grads.size(), this->n()));
+ Assert(third_derivatives.size() == this->n() || third_derivatives.size() == 0,
+ ExcDimensionMismatch(third_derivatives.size(), this->n()));
+ Assert(fourth_derivatives.size() == this->n() ||
+ fourth_derivatives.size() == 0,
+ ExcDimensionMismatch(fourth_derivatives.size(), this->n()));
// third and fourth derivatives not implemented
(void)third_derivatives;
// and second derivatives vectors of
// <tt>polynomial_space</tt> at
// <tt>unit_point</tt>
- const unsigned int n_basis = polynomial_space.n();
+ const unsigned int n_basis = polynomial_space.n();
+ const unsigned int my_degree = this->degree();
std::vector<double> unit_point_values((values.size() == 0) ? 0 : n_basis);
std::vector<Tensor<1, dim>> unit_point_grads((grads.size() == 0) ? 0 :
n_basis);
}
+template <int dim>
+std::unique_ptr<TensorPolynomialsBase<dim>>
+PolynomialsNedelec<dim>::clone() const
+{
+ return std_cxx14::make_unique<PolynomialsNedelec<dim>>(*this);
+}
+
+
template class PolynomialsNedelec<1>;
template class PolynomialsNedelec<2>;
template class PolynomialsNedelec<3>;
#include <deal.II/base/polynomials_raviart_thomas.h>
#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/base/std_cxx14/memory.h>
#include <deal.II/base/thread_management.h>
#include <iomanip>
template <int dim>
PolynomialsRaviartThomas<dim>::PolynomialsRaviartThomas(const unsigned int k)
- : my_degree(k)
+ : TensorPolynomialsBase<dim>(k, compute_n_pols(k))
, polynomial_space(create_polynomials(k))
- , n_pols(compute_n_pols(k))
{}
std::vector<Tensor<4, dim>> &third_derivatives,
std::vector<Tensor<5, dim>> &fourth_derivatives) const
{
- Assert(values.size() == n_pols || values.size() == 0,
- ExcDimensionMismatch(values.size(), n_pols));
- Assert(grads.size() == n_pols || grads.size() == 0,
- ExcDimensionMismatch(grads.size(), n_pols));
- Assert(grad_grads.size() == n_pols || grad_grads.size() == 0,
- ExcDimensionMismatch(grad_grads.size(), n_pols));
- Assert(third_derivatives.size() == n_pols || third_derivatives.size() == 0,
- ExcDimensionMismatch(third_derivatives.size(), n_pols));
- Assert(fourth_derivatives.size() == n_pols || fourth_derivatives.size() == 0,
- ExcDimensionMismatch(fourth_derivatives.size(), n_pols));
+ Assert(values.size() == this->n() || values.size() == 0,
+ ExcDimensionMismatch(values.size(), this->n()));
+ Assert(grads.size() == this->n() || grads.size() == 0,
+ ExcDimensionMismatch(grads.size(), this->n()));
+ Assert(grad_grads.size() == this->n() || grad_grads.size() == 0,
+ ExcDimensionMismatch(grad_grads.size(), this->n()));
+ Assert(third_derivatives.size() == this->n() || third_derivatives.size() == 0,
+ ExcDimensionMismatch(third_derivatives.size(), this->n()));
+ Assert(fourth_derivatives.size() == this->n() ||
+ fourth_derivatives.size() == 0,
+ ExcDimensionMismatch(fourth_derivatives.size(), this->n()));
// have a few scratch
// arrays. because we don't want to
}
+template <int dim>
+std::unique_ptr<TensorPolynomialsBase<dim>>
+PolynomialsRaviartThomas<dim>::clone() const
+{
+ return std_cxx14::make_unique<PolynomialsRaviartThomas<dim>>(*this);
+}
+
+
template class PolynomialsRaviartThomas<1>;
template class PolynomialsRaviartThomas<2>;
template class PolynomialsRaviartThomas<3>;
#include "deal.II/base/polynomials_rt_bubbles.h"
#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/base/std_cxx14/memory.h>
#include <deal.II/base/thread_management.h>
#include <iomanip>
template <int dim>
PolynomialsRT_Bubbles<dim>::PolynomialsRT_Bubbles(const unsigned int k)
- : my_degree(k)
+ : TensorPolynomialsBase<dim>(k, compute_n_pols(k))
, raviart_thomas_space(k - 1)
, monomials(k + 2)
- , n_pols(compute_n_pols(k))
{
Assert(dim >= 2, ExcImpossibleInDim(dim));
std::vector<Tensor<4, dim>> &third_derivatives,
std::vector<Tensor<5, dim>> &fourth_derivatives) const
{
- Assert(values.size() == n_pols || values.size() == 0,
- ExcDimensionMismatch(values.size(), n_pols));
- Assert(grads.size() == n_pols || grads.size() == 0,
- ExcDimensionMismatch(grads.size(), n_pols));
- Assert(grad_grads.size() == n_pols || grad_grads.size() == 0,
- ExcDimensionMismatch(grad_grads.size(), n_pols));
- Assert(third_derivatives.size() == n_pols || third_derivatives.size() == 0,
- ExcDimensionMismatch(third_derivatives.size(), n_pols));
- Assert(fourth_derivatives.size() == n_pols || fourth_derivatives.size() == 0,
- ExcDimensionMismatch(fourth_derivatives.size(), n_pols));
+ Assert(values.size() == this->n() || values.size() == 0,
+ ExcDimensionMismatch(values.size(), this->n()));
+ Assert(grads.size() == this->n() || grads.size() == 0,
+ ExcDimensionMismatch(grads.size(), this->n()));
+ Assert(grad_grads.size() == this->n() || grad_grads.size() == 0,
+ ExcDimensionMismatch(grad_grads.size(), this->n()));
+ Assert(third_derivatives.size() == this->n() || third_derivatives.size() == 0,
+ ExcDimensionMismatch(third_derivatives.size(), this->n()));
+ Assert(fourth_derivatives.size() == this->n() ||
+ fourth_derivatives.size() == 0,
+ ExcDimensionMismatch(fourth_derivatives.size(), this->n()));
// Third and fourth derivatives are not implemented
(void)third_derivatives;
(void)fourth_derivatives;
Assert(fourth_derivatives.size() == 0, ExcNotImplemented());
- const unsigned int n_sub = raviart_thomas_space.n();
+ const unsigned int n_sub = raviart_thomas_space.n();
+ const unsigned int my_degree = this->degree();
// Guard access to the scratch arrays in the following block
// using a mutex to make sure they are not used by multiple threads
monoval_plus[0][1] * monoval_i[1][2];
}
}
- Assert(start == n_pols - my_degree - 1, ExcInternalError());
+ Assert(start == this->n() - my_degree - 1, ExcInternalError());
}
else if (dim == 3)
{
start += 2;
}
}
- Assert(start == n_pols - 2 * n_curls, ExcInternalError());
+ Assert(start == this->n() - 2 * n_curls, ExcInternalError());
}
}
}
+template <int dim>
+std::unique_ptr<TensorPolynomialsBase<dim>>
+PolynomialsRT_Bubbles<dim>::clone() const
+{
+ return std_cxx14::make_unique<PolynomialsRT_Bubbles<dim>>(*this);
+}
+
+
template class PolynomialsRT_Bubbles<1>;
template class PolynomialsRT_Bubbles<2>;
template class PolynomialsRT_Bubbles<3>;