17 #include <deal.II/base/quadrature_lib.h> 18 #include <deal.II/lac/vector.h> 19 #include <deal.II/fe/fe_q_iso_q1.h> 20 #include <deal.II/fe/fe_nothing.h> 24 #include <deal.II/base/std_cxx14/memory.h> 26 DEAL_II_NAMESPACE_OPEN
32 template <
int dim,
int spacedim>
37 (
Polynomials::generate_complete_Lagrange_basis_on_subdivisions(subdivisions, 1)),
41 std::vector<bool> (1, false))
44 ExcMessage (
"This element can only be used with a positive number of " 55 template <
int dim,
int spacedim>
63 std::ostringstream namebuf;
64 namebuf <<
"FE_Q_iso_Q1<" 66 <<
">(" << this->degree <<
")";
72 template <
int dim,
int spacedim>
76 std::vector<double> &nodal_values)
const 79 this->get_unit_support_points().size());
85 for (
unsigned int i=0; i<this->dofs_per_cell; ++i)
89 nodal_values[i] = support_point_values[i](0);
95 template <
int dim,
int spacedim>
96 std::unique_ptr<FiniteElement<dim,spacedim> >
99 return std_cxx14::make_unique<FE_Q_iso_Q1<dim,spacedim>>(*this);
104 template <
int dim,
int spacedim>
115 if (this->degree < fe_q_iso_q1_other->degree &&
116 fe_q_iso_q1_other->degree % this->degree == 0)
118 else if (this->degree == fe_q_iso_q1_other->degree)
120 else if (this->degree > fe_q_iso_q1_other->degree &&
121 this->degree % fe_q_iso_q1_other->degree == 0)
128 if (fe_nothing->is_dominating())
146 #include "fe_q_iso_q1.inst" 148 DEAL_II_NAMESPACE_CLOSE
#define AssertDimension(dim1, dim2)
const std::vector< Point< dim > > & get_points() const
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::unique_ptr< FiniteElement< dim, spacedim > > clone() const
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
std::string dim_string(const int dim, const int spacedim)
void initialize(const std::vector< Point< 1 > > &support_points_1d)
FE_Q_iso_Q1(const unsigned int n_subdivisions)
virtual std::string get_name() const
static ::ExceptionBase & ExcNotImplemented()
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const