59 "Lowest order BDM element are degree 1, but you asked for degree 0"));
88 const unsigned int face_no = 0;
91 for (
unsigned int i = 0; i < GeometryInfo<dim>::max_children_per_face; ++i)
98 unsigned int target_row = 0;
99 for (
unsigned int d = 0; d < GeometryInfo<dim>::max_children_per_face; ++d)
100 for (
unsigned int i = 0; i < face_embeddings[d].
m(); ++i)
102 for (
unsigned int j = 0; j < face_embeddings[d].
n(); ++j)
141 std::ostringstream namebuf;
142 namebuf <<
"FE_BDM<" << dim <<
">(" << this->degree - 1 <<
")";
144 return namebuf.str();
149std::unique_ptr<FiniteElement<dim, dim>>
152 return std::make_unique<FE_BDM<dim>>(*this);
161 std::vector<double> & nodal_values)
const
163 Assert(support_point_values.size() == this->generalized_support_points.size(),
165 this->generalized_support_points.size()));
167 Assert(nodal_values.size() == this->n_dofs_per_cell(),
174 unsigned int dbase = 0;
176 unsigned int pbase = 0;
181 if (test_values_face.size() == 0)
183 for (
unsigned int i = 0; i < this->n_dofs_per_face(f); ++i)
184 nodal_values[dbase + i] =
185 support_point_values[pbase + i]
187 pbase += this->n_dofs_per_face(f);
191 for (
unsigned int i = 0; i < this->n_dofs_per_face(f); ++i)
194 for (
unsigned int k = 0; k < test_values_face.size(); ++k)
198 test_values_face[k][i];
199 nodal_values[dbase + i] = s;
201 pbase += test_values_face.size();
203 dbase += this->n_dofs_per_face(f);
207 const unsigned int face_no = 0;
211 this->n_dofs_per_face(face_no) *
214 this->generalized_support_points.size() -
215 test_values_cell.size());
218 if (dbase == this->n_dofs_per_cell())
227 for (
unsigned int d = 0; d < dim; ++d, dbase += test_values_cell[0].size())
229 for (
unsigned int i = 0; i < test_values_cell[0].size(); ++i)
232 for (
unsigned int k = 0; k < test_values_cell.size(); ++k)
233 s += support_point_values[pbase + k][d] * test_values_cell[k][i];
234 nodal_values[dbase + i] = s;
244std::vector<unsigned int>
257 std::vector<unsigned int> dpo(dim + 1, 0u);
274 return std::vector<bool>();
278 const unsigned int dofs_per_face =
289 std::vector<bool> ret_val(dofs_per_cell,
false);
312 initialize_test_values(std::vector<std::vector<double>> &test_values,
314 const unsigned int deg)
317 std::vector<Tensor<1, dim>> dummy1;
318 std::vector<Tensor<2, dim>> dummy2;
319 std::vector<Tensor<3, dim>> dummy3;
320 std::vector<Tensor<4, dim>> dummy4;
322 test_values.resize(quadrature.
size());
324 for (
unsigned int k = 0; k < quadrature.
size(); ++k)
326 test_values[k].resize(poly.
n());
333 for (
unsigned int i = 0; i < poly.
n(); ++i)
335 test_values[k][i] *= quadrature.
weight(k);
345 initialize_test_values(std::vector<std::vector<double>> &,
363 QGauss<dim - 1> face_points(deg + 1);
368 const unsigned int face_no = 0;
371 this->generalized_face_support_points[face_no].resize(face_points.size());
372 for (
unsigned int k = 0; k < face_points.size(); ++k)
373 this->generalized_face_support_points[face_no][k] = face_points.point(k);
382 const unsigned int npoints =
385 this->generalized_support_points.resize(npoints);
389 for (
unsigned int k = 0;
392 this->generalized_support_points[k] = faces.
point(
399 this->n_dofs_per_face(face_no)));
405 internal::FE_BDM::initialize_test_values(test_values_face,
413 const unsigned int ibase =
415 for (
unsigned int k = 0; k < cell_points.
size(); ++k)
417 this->generalized_support_points[ibase + k] = cell_points.
point(k);
423 internal::FE_BDM::initialize_test_values(test_values_cell,
431#include "fe_bdm.inst"
void initialize_support_points(const unsigned int bdm_degree)
FE_BDM(const unsigned int p)
static std::vector< bool > get_ria_vector(const unsigned int degree)
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const override
void initialize_quad_dof_index_permutation_and_sign_change()
virtual std::string get_name() const override
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
FullMatrix< double > inverse_node_matrix
std::vector< MappingKind > mapping_kind
unsigned int n_dofs_per_cell() const
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
unsigned int n_unique_faces() const
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
FullMatrix< double > interface_constraints
std::vector< std::vector< FullMatrix< double > > > prolongation
void invert(const FullMatrix< number2 > &M)
static unsigned int n_polynomials(const unsigned int n)
void evaluate(const Point< dim > &unit_point, std::vector< double > &values, std::vector< Tensor< 1, dim > > &grads, std::vector< Tensor< 2, dim > > &grad_grads, std::vector< Tensor< 3, dim > > &third_derivatives, std::vector< Tensor< 4, dim > > &fourth_derivatives) const override
static unsigned int n_polynomials(const unsigned int degree)
static Quadrature< dim > project_to_all_faces(const ReferenceCell &reference_cell, const hp::QCollection< dim - 1 > &quadrature)
const Point< dim > & point(const unsigned int i) const
double weight(const unsigned int i) const
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices()