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> 34 DEAL_II_NAMESPACE_OPEN
45 std::vector<bool>(dim,true)))
48 Assert (deg > 0,
ExcMessage(
"Lowest order BDM element are degree 1, but you asked for degree 0"));
77 for (
unsigned int i=0; i<GeometryInfo<dim>::max_children_per_face; ++i)
82 unsigned int target_row=0;
83 for (
unsigned int d=0; d<GeometryInfo<dim>::max_children_per_face; ++d)
84 for (
unsigned int i=0; i<face_embeddings[d].
m(); ++i)
86 for (
unsigned int j=0; j<face_embeddings[d].
n(); ++j)
109 std::ostringstream namebuf;
110 namebuf <<
"FE_BDM<" << dim <<
">(" << this->degree-1 <<
")";
112 return namebuf.str();
129 std::vector<double> &nodal_values)
const 131 Assert (support_point_values.size() == this->generalized_support_points.size(),
134 Assert (nodal_values.size() == this->dofs_per_cell,
141 unsigned int dbase = 0;
143 unsigned int pbase = 0;
144 for (
unsigned int f = 0; f<GeometryInfo<dim>::faces_per_cell; ++f)
148 if (test_values_face.size() == 0)
150 for (
unsigned int i=0; i<this->dofs_per_face; ++i)
152 pbase += this->dofs_per_face;
156 for (
unsigned int i=0; i<this->dofs_per_face; ++i)
159 for (
unsigned int k=0; k<test_values_face.size(); ++k)
161 nodal_values[dbase+i] = s;
163 pbase += test_values_face.size();
165 dbase += this->dofs_per_face;
169 AssertDimension (pbase, this->generalized_support_points.size() - test_values_cell.size());
172 if (dbase == this->dofs_per_cell)
return;
180 for (
unsigned int d=0; d<dim; ++d, dbase += test_values_cell[0].size())
182 for (
unsigned int i=0; i<test_values_cell[0].size(); ++i)
185 for (
unsigned int k=0; k<test_values_cell.size(); ++k)
186 s += support_point_values[pbase+k][d] * test_values_cell[k][i];
187 nodal_values[dbase+i] = s;
199 std::vector<double> &,
200 const std::vector<double> &)
const 209 std::vector<double> &,
221 std::vector<double> &local_dofs,
222 const VectorSlice<
const std::vector<std::vector<double> > > &values)
const 225 Assert (values[0].size() == this->generalized_support_points.size(),
227 Assert (local_dofs.size() == this->dofs_per_cell,
234 unsigned int dbase = 0;
236 unsigned int pbase = 0;
237 for (
unsigned int f = 0; f<GeometryInfo<dim>::faces_per_cell; ++f)
241 if (test_values_face.size() == 0)
243 for (
unsigned int i=0; i<this->dofs_per_face; ++i)
245 pbase += this->dofs_per_face;
249 for (
unsigned int i=0; i<this->dofs_per_face; ++i)
252 for (
unsigned int k=0; k<test_values_face.size(); ++k)
254 local_dofs[dbase+i] = s;
256 pbase += test_values_face.size();
258 dbase += this->dofs_per_face;
262 AssertDimension (pbase, this->generalized_support_points.size() - test_values_cell.size());
265 if (dbase == this->dofs_per_cell)
return;
273 for (
unsigned int d=0; d<dim; ++d, dbase += test_values_cell[0].size())
275 for (
unsigned int i=0; i<test_values_cell[0].size(); ++i)
278 for (
unsigned int k=0; k<test_values_cell.size(); ++k)
279 s += values[d][pbase+k] * test_values_cell[k][i];
280 local_dofs[dbase+i] = s;
291 std::vector<unsigned int>
304 std::vector<unsigned int> dpo(dim+1, 0u);
305 dpo[dim] = (deg > 1 ?
322 return std::vector<bool>();
336 std::vector<bool> ret_val(dofs_per_cell,
false);
338 i < dofs_per_cell; ++i)
354 initialize_test_values (std::vector<std::vector<double> > &test_values,
356 const unsigned int deg)
359 std::vector<Tensor<1,dim> > dummy1;
360 std::vector<Tensor<2,dim> > dummy2;
361 std::vector<Tensor<3,dim> > dummy3;
362 std::vector<Tensor<4,dim> > dummy4;
364 test_values.resize(quadrature.
size());
366 for (
unsigned int k=0; k<quadrature.
size(); ++k)
368 test_values[k].resize(poly.n());
369 poly.compute(quadrature.
point(k), test_values[k], dummy1, dummy2,
371 for (
unsigned int i=0; i < poly.n(); ++i)
373 test_values[k][i] *= quadrature.
weight(k);
383 initialize_test_values (std::vector<std::vector<double> > &,
399 QGauss<dim-1> face_points (deg+1);
402 this->generalized_face_support_points.resize (face_points.size());
403 for (
unsigned int k=0; k<face_points.size(); ++k)
404 this->generalized_face_support_points[k] = face_points.point(k);
413 const unsigned int npoints
416 this->generalized_support_points.resize (npoints);
420 this->generalized_support_points[k]
422 ::DataSetDescriptor::face(0,
true,
false,
false,
423 this->dofs_per_face));
429 initialize_test_values(test_values_face, face_points, deg);
435 for (
unsigned int k=0; k<cell_points.
size(); ++k)
437 this->generalized_support_points[ibase+k] = cell_points.
point(k);
443 initialize_test_values(test_values_cell, cell_points, deg-2);
449 #include "fe_bdm.inst" 451 DEAL_II_NAMESPACE_CLOSE
void initialize_support_points(const unsigned int bdm_degree)
#define AssertDimension(dim1, dim2)
FullMatrix< double > interface_constraints
virtual void convert_generalized_support_point_values_to_nodal_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) 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)
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
virtual FiniteElement< dim > * clone() const
static unsigned int compute_n_pols(const unsigned int n)
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const 1
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)
static ::ExceptionBase & ExcNotImplemented()
FullMatrix< double > inverse_node_matrix
double weight(const unsigned int i) const
static ::ExceptionBase & ExcInternalError()