17 #include <deal.II/base/quadrature_lib.h> 18 #include <deal.II/base/qprojector.h> 19 #include <deal.II/base/polynomials_p.h> 20 #include <deal.II/base/table.h> 21 #include <deal.II/grid/tria.h> 22 #include <deal.II/grid/tria_iterator.h> 23 #include <deal.II/dofs/dof_accessor.h> 24 #include <deal.II/fe/fe.h> 25 #include <deal.II/fe/mapping.h> 26 #include <deal.II/fe/fe_bdm.h> 27 #include <deal.II/fe/fe_values.h> 28 #include <deal.II/fe/fe_tools.h> 32 #include <deal.II/base/std_cxx14/memory.h> 35 DEAL_II_NAMESPACE_OPEN
46 std::vector<bool>(dim,true)))
49 Assert (deg > 0,
ExcMessage(
"Lowest order BDM element are degree 1, but you asked for degree 0"));
78 for (
unsigned int i=0; i<GeometryInfo<dim>::max_children_per_face; ++i)
83 unsigned int target_row=0;
84 for (
unsigned int d=0; d<GeometryInfo<dim>::max_children_per_face; ++d)
85 for (
unsigned int i=0; i<face_embeddings[d].
m(); ++i)
87 for (
unsigned int j=0; j<face_embeddings[d].
n(); ++j)
110 std::ostringstream namebuf;
111 namebuf <<
"FE_BDM<" << dim <<
">(" << this->degree-1 <<
")";
113 return namebuf.str();
118 std::unique_ptr<FiniteElement<dim,dim> >
121 return std_cxx14::make_unique<FE_BDM<dim>>(*this);
130 std::vector<double> &nodal_values)
const 132 Assert (support_point_values.size() == this->generalized_support_points.size(),
135 Assert (nodal_values.size() == this->dofs_per_cell,
142 unsigned int dbase = 0;
144 unsigned int pbase = 0;
145 for (
unsigned int f = 0; f<GeometryInfo<dim>::faces_per_cell; ++f)
149 if (test_values_face.size() == 0)
151 for (
unsigned int i=0; i<this->dofs_per_face; ++i)
153 pbase += this->dofs_per_face;
157 for (
unsigned int i=0; i<this->dofs_per_face; ++i)
160 for (
unsigned int k=0; k<test_values_face.size(); ++k)
162 nodal_values[dbase+i] = s;
164 pbase += test_values_face.size();
166 dbase += this->dofs_per_face;
170 AssertDimension (pbase, this->generalized_support_points.size() - test_values_cell.size());
173 if (dbase == this->dofs_per_cell)
return;
181 for (
unsigned int d=0; d<dim; ++d, dbase += test_values_cell[0].size())
183 for (
unsigned int i=0; i<test_values_cell[0].size(); ++i)
186 for (
unsigned int k=0; k<test_values_cell.size(); ++k)
187 s += support_point_values[pbase+k][d] * test_values_cell[k][i];
188 nodal_values[dbase+i] = s;
198 std::vector<unsigned int>
211 std::vector<unsigned int> dpo(dim+1, 0u);
212 dpo[dim] = (deg > 1 ?
229 return std::vector<bool>();
243 std::vector<bool> ret_val(dofs_per_cell,
false);
245 i < dofs_per_cell; ++i)
265 initialize_test_values (std::vector<std::vector<double> > &test_values,
267 const unsigned int deg)
270 std::vector<Tensor<1,dim> > dummy1;
271 std::vector<Tensor<2,dim> > dummy2;
272 std::vector<Tensor<3,dim> > dummy3;
273 std::vector<Tensor<4,dim> > dummy4;
275 test_values.resize(quadrature.
size());
277 for (
unsigned int k=0; k<quadrature.
size(); ++k)
279 test_values[k].resize(poly.n());
280 poly.compute(quadrature.
point(k), test_values[k], dummy1, dummy2,
282 for (
unsigned int i=0; i < poly.n(); ++i)
284 test_values[k][i] *= quadrature.
weight(k);
294 initialize_test_values (std::vector<std::vector<double> > &,
312 QGauss<dim-1> face_points (deg+1);
315 this->generalized_face_support_points.resize (face_points.size());
316 for (
unsigned int k=0; k<face_points.size(); ++k)
317 this->generalized_face_support_points[k] = face_points.point(k);
326 const unsigned int npoints
329 this->generalized_support_points.resize (npoints);
333 this->generalized_support_points[k]
335 ::DataSetDescriptor::face(0,
true,
false,
false,
336 this->dofs_per_face));
342 internal::FE_BDM::initialize_test_values(test_values_face, face_points, deg);
348 for (
unsigned int k=0; k<cell_points.
size(); ++k)
350 this->generalized_support_points[ibase+k] = cell_points.
point(k);
356 internal::FE_BDM::initialize_test_values(test_values_cell, cell_points, deg-2);
362 #include "fe_bdm.inst" 364 DEAL_II_NAMESPACE_CLOSE
void initialize_support_points(const unsigned int bdm_degree)
#define AssertDimension(dim1, dim2)
FullMatrix< double > interface_constraints
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const
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)
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
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)
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::string get_name() const
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()