17 #include <deal.II/fe/fe_rannacher_turek.h> 18 #include <deal.II/base/quadrature_lib.h> 19 #include <deal.II/lac/vector.h> 23 #include <deal.II/base/std_cxx14/memory.h> 26 DEAL_II_NAMESPACE_OPEN
31 const unsigned int n_face_support_points) :
38 std::vector<bool>(4, false),
41 n_face_support_points(n_face_support_points)
53 std::vector<unsigned int> dpo(dim + 1, 0);
64 std::ostringstream namebuf;
65 namebuf <<
"FE_RannacherTurek" 67 <<
"(" << this->order <<
", " << this->n_face_support_points <<
")";
74 std::unique_ptr<FiniteElement<dim,dim> >
77 return std_cxx14::make_unique<FE_RannacherTurek<dim>>(this->order, this->n_face_support_points);
86 ::QGauss<dim-1> face_quadrature(this->n_face_support_points);
87 this->weights = face_quadrature.get_weights();
88 this->generalized_support_points.resize(4*face_quadrature.size());
89 for (
unsigned int q = 0;
90 q < face_quadrature.size();
93 this->generalized_support_points[0*face_quadrature.size() + q] =
95 this->generalized_support_points[1*face_quadrature.size() + q] =
97 this->generalized_support_points[2*face_quadrature.size() + q] =
99 this->generalized_support_points[3*face_quadrature.size() + q] =
110 std::vector<double> &nodal_values)
const 112 AssertDimension(support_point_values.size(), this->generalized_support_points.size());
115 const unsigned int q_points_per_face = this->weights.size();
116 std::fill(nodal_values.begin(), nodal_values.end(), 0.0);
118 std::vector<Vector<double> >::const_iterator value = support_point_values.begin();
119 for (
unsigned int face = 0;
120 face < ::GeometryInfo<dim>::faces_per_cell;
123 for (
unsigned int q = 0;
124 q < q_points_per_face;
127 nodal_values[face] += (*value)[0] * this->weights[q];
136 #include "fe_rannacher_turek.inst" 138 DEAL_II_NAMESPACE_CLOSE
#define AssertDimension(dim1, dim2)
FE_RannacherTurek(const unsigned int order=0, const unsigned int n_face_support_points=2)
void initialize_support_points()
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const
#define Assert(cond, exc)
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
virtual std::string get_name() const
static ::ExceptionBase & ExcNotImplemented()
std::vector< unsigned int > get_dpo_vector()