17 #include <deal.II/base/polynomials_p.h> 18 #include <deal.II/base/qprojector.h> 19 #include <deal.II/base/quadrature_lib.h> 20 #include <deal.II/base/std_cxx14/memory.h> 22 #include <deal.II/dofs/dof_accessor.h> 24 #include <deal.II/fe/fe.h> 25 #include <deal.II/fe/fe_bdm.h> 26 #include <deal.II/fe/fe_tools.h> 27 #include <deal.II/fe/fe_values.h> 28 #include <deal.II/fe/mapping.h> 30 #include <deal.II/grid/tria.h> 31 #include <deal.II/grid/tria_iterator.h> 37 DEAL_II_NAMESPACE_OPEN
49 std::vector<bool>(dim, true)))
55 "Lowest order BDM element are degree 1, but you asked for degree 0"));
84 for (
unsigned int i = 0; i < GeometryInfo<dim>::max_children_per_face; ++i)
89 unsigned int target_row = 0;
90 for (
unsigned int d = 0; d < GeometryInfo<dim>::max_children_per_face; ++d)
91 for (
unsigned int i = 0; i < face_embeddings[d].
m(); ++i)
93 for (
unsigned int j = 0; j < face_embeddings[d].
n(); ++j)
116 std::ostringstream namebuf;
117 namebuf <<
"FE_BDM<" << dim <<
">(" << this->degree - 1 <<
")";
119 return namebuf.str();
124 std::unique_ptr<FiniteElement<dim, dim>>
127 return std_cxx14::make_unique<FE_BDM<dim>>(*this);
136 std::vector<double> & nodal_values)
const 138 Assert(support_point_values.size() == this->generalized_support_points.size(),
140 this->generalized_support_points.size()));
142 Assert(nodal_values.size() == this->dofs_per_cell,
149 unsigned int dbase = 0;
151 unsigned int pbase = 0;
152 for (
unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
156 if (test_values_face.size() == 0)
158 for (
unsigned int i = 0; i < this->dofs_per_face; ++i)
159 nodal_values[dbase + i] =
160 support_point_values[pbase + i]
162 pbase += this->dofs_per_face;
166 for (
unsigned int i = 0; i < this->dofs_per_face; ++i)
169 for (
unsigned int k = 0; k < test_values_face.size(); ++k)
173 test_values_face[k][i];
174 nodal_values[dbase + i] = s;
176 pbase += test_values_face.size();
178 dbase += this->dofs_per_face;
184 this->generalized_support_points.size() -
185 test_values_cell.size());
188 if (dbase == this->dofs_per_cell)
197 for (
unsigned int d = 0; d < dim; ++d, dbase += test_values_cell[0].size())
199 for (
unsigned int i = 0; i < test_values_cell[0].size(); ++i)
202 for (
unsigned int k = 0; k < test_values_cell.size(); ++k)
203 s += support_point_values[pbase + k][d] * test_values_cell[k][i];
204 nodal_values[dbase + i] = s;
214 std::vector<unsigned int>
227 std::vector<unsigned int> dpo(dim + 1, 0u);
244 return std::vector<bool>();
248 const unsigned int dofs_per_face =
259 std::vector<bool> ret_val(dofs_per_cell,
false);
282 initialize_test_values(std::vector<std::vector<double>> &test_values,
284 const unsigned int deg)
287 std::vector<Tensor<1, dim>> dummy1;
288 std::vector<Tensor<2, dim>> dummy2;
289 std::vector<Tensor<3, dim>> dummy3;
290 std::vector<Tensor<4, dim>> dummy4;
292 test_values.resize(quadrature.
size());
294 for (
unsigned int k = 0; k < quadrature.
size(); ++k)
296 test_values[k].resize(poly.n());
297 poly.compute(quadrature.
point(k),
303 for (
unsigned int i = 0; i < poly.n(); ++i)
305 test_values[k][i] *= quadrature.
weight(k);
315 initialize_test_values(std::vector<std::vector<double>> &,
333 QGauss<dim - 1> face_points(deg + 1);
336 this->generalized_face_support_points.resize(face_points.size());
337 for (
unsigned int k = 0; k < face_points.size(); ++k)
338 this->generalized_face_support_points[k] = face_points.point(k);
347 const unsigned int npoints =
350 this->generalized_support_points.resize(npoints);
353 for (
unsigned int k = 0;
356 this->generalized_support_points[k] =
358 0,
true,
false,
false, this->dofs_per_face));
364 internal::FE_BDM::initialize_test_values(test_values_face,
372 const unsigned int ibase =
374 for (
unsigned int k = 0; k < cell_points.
size(); ++k)
376 this->generalized_support_points[ibase + k] = cell_points.
point(k);
382 internal::FE_BDM::initialize_test_values(test_values_cell,
390 #include "fe_bdm.inst" 392 DEAL_II_NAMESPACE_CLOSE
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_support_points(const unsigned int bdm_degree)
#define AssertDimension(dim1, dim2)
FullMatrix< double > interface_constraints
const Point< dim > & point(const unsigned int i) const
static std::vector< bool > get_ria_vector(const unsigned int degree)
FE_BDM(const unsigned int p)
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
std::vector< std::vector< FullMatrix< double > > > prolongation
void invert(const FullMatrix< number2 > &M)
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
virtual std::string get_name() const override
unsigned int size() const
const unsigned int dofs_per_cell
static unsigned int compute_n_pols(const unsigned int n)
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
static unsigned int compute_n_pols(unsigned int degree)
const unsigned int dofs_per_face
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
FullMatrix< double > inverse_node_matrix
double weight(const unsigned int i) const
static ::ExceptionBase & ExcInternalError()