std_cxx20::ranges::iota_view<unsigned int, unsigned int>
dof_indices() const;
+ /**
+ * Return whether the face associated to this FEFaceEvaluation object
+ * is at the boundary.
+ */
+ bool
+ at_boundary() const;
+
+ /**
+ * Return the boundary indicator of the face associated to this
+ * FEFaceEvaluation object.
+ *
+ * If the return value is the special value
+ * numbers::internal_face_boundary_id, then the face is in the interior of
+ * the domain.
+ *
+ * @note Alternatively to this function, you can use
+ * MatrixFree::get_boundary_id() to get the same information if no
+ * FEFaceEvaluation object has been set up.
+ */
+ types::boundary_id
+ boundary_id() const;
+
/**
* The number of degrees of freedom of a single component on the cell for
* the underlying evaluation object. Usually close to
+template <int dim,
+ int fe_degree,
+ int n_q_points_1d,
+ int n_components_,
+ typename Number,
+ typename VectorizedArrayType>
+bool
+FEFaceEvaluation<dim,
+ fe_degree,
+ n_q_points_1d,
+ n_components_,
+ Number,
+ VectorizedArrayType>::at_boundary() const
+{
+ Assert(this->dof_access_index !=
+ internal::MatrixFreeFunctions::DoFInfo::dof_access_cell,
+ ExcNotImplemented());
+
+ if (this->is_interior_face() == false)
+ return false;
+ else if (this->cell < this->matrix_free->n_inner_face_batches())
+ return false;
+ else if (this->cell < (this->matrix_free->n_inner_face_batches() +
+ this->matrix_free->n_boundary_face_batches()))
+ return true;
+ else
+ return false;
+}
+
+
+
+template <int dim,
+ int fe_degree,
+ int n_q_points_1d,
+ int n_components_,
+ typename Number,
+ typename VectorizedArrayType>
+types::boundary_id
+FEFaceEvaluation<dim,
+ fe_degree,
+ n_q_points_1d,
+ n_components_,
+ Number,
+ VectorizedArrayType>::boundary_id() const
+{
+ Assert(this->dof_access_index !=
+ internal::MatrixFreeFunctions::DoFInfo::dof_access_cell,
+ ExcNotImplemented());
+
+ if (at_boundary())
+ return this->matrix_free->get_boundary_id(this->cell);
+ else
+ return numbers::internal_face_boundary_id;
+}
+
+
+
/*------------------------- end FEFaceEvaluation ------------------------- */
* this method can be used to query the boundary id of a given face in the
* faces' own sorting by lanes in a VectorizedArray. Only valid for an index
* indicating a boundary face.
+ *
+ * @note Alternatively to this function, you can use
+ * FEFaceEvaluation::boundary_id() to get the same information if a
+ * FEFaceEvaluation object has been set up already.
*/
types::boundary_id
get_boundary_id(const unsigned int face_batch_index) const;