--- /dev/null
+Changed: FiniteElement::adjust_quad_dof_index_for_face_orientation() and
+FiniteElement::face_to_cell_index() now use the standardized
+<tt>combined_orientation</tt> orientation encoding as input arguments rather
+than three booleans.
+<br>
+(David Wells, 2023/12/21)
* index must be between zero and dofs_per_face.
* @param face The number of the face this degree of freedom lives on. This
* number must be between zero and GeometryInfo::faces_per_cell.
- * @param face_orientation One part of the description of the orientation of
- * the face. See
- * @ref GlossFaceOrientation.
- * @param face_flip One part of the description of the orientation of the
- * face. See
- * @ref GlossFaceOrientation.
- * @param face_rotation One part of the description of the orientation of
- * the face. See
- * @ref GlossFaceOrientation.
+ * @param combined_orientation The combined orientation flag containing the
+ * orientation, rotation, and flip of the face. See @ref GlossFaceOrientation.
* @return The index of this degree of freedom within the set of degrees of
* freedom on the entire cell. The returned value will be between zero and
* dofs_per_cell.
* freedom actually represent.
*/
virtual unsigned int
- face_to_cell_index(const unsigned int face_dof_index,
- const unsigned int face,
- const bool face_orientation = true,
- const bool face_flip = false,
- const bool face_rotation = false) const;
+ face_to_cell_index(
+ const unsigned int face_dof_index,
+ const unsigned int face,
+ const unsigned char combined_orientation =
+ ReferenceCell::default_combined_face_orientation()) const;
/**
* For lines with non-standard line_orientation in 3d, the dofs on lines
* index must be between zero and dofs_per_face.
* @param face The number of the face this degree of freedom lives on. This
* number must be between zero and GeometryInfo::faces_per_cell.
- * @param face_orientation One part of the description of the orientation of
- * the face. See
- * @ref GlossFaceOrientation.
- * @param face_flip One part of the description of the orientation of the
- * face. See
- * @ref GlossFaceOrientation.
- * @param face_rotation One part of the description of the orientation of
- * the face. See
- * @ref GlossFaceOrientation.
+ * @param combined_orientation The combined orientation flag containing the
+ * orientation, rotation, and flip of the face. See @ref GlossFaceOrientation.
* @return The index of this degree of freedom within the set of degrees of
* freedom on the entire cell. The returned value will be between zero and
* dofs_per_cell.
*/
virtual unsigned int
- face_to_cell_index(const unsigned int face_dof_index,
- const unsigned int face,
- const bool face_orientation = true,
- const bool face_flip = false,
- const bool face_rotation = false) const override;
+ face_to_cell_index(
+ const unsigned int face_dof_index,
+ const unsigned int face,
+ const unsigned char combined_orientation =
+ ReferenceCell::default_combined_face_orientation()) const override;
/**
* Return a list of constant modes of the element. For this element, the
* index must be between zero and dofs_per_face.
* @param face The number of the face this degree of freedom lives on. This
* number must be between zero and GeometryInfo::faces_per_cell.
- * @param face_orientation One part of the description of the orientation of
- * the face. See
- * @ref GlossFaceOrientation.
- * @param face_flip One part of the description of the orientation of the
- * face. See
- * @ref GlossFaceOrientation.
- * @param face_rotation One part of the description of the orientation of
- * the face. See
- * @ref GlossFaceOrientation.
+ * @param combined_orientation The combined orientation flag containing the
+ * orientation, rotation, and flip of the face. See @ref GlossFaceOrientation.
* @return The index of this degree of freedom within the set of degrees of
* freedom on the entire cell. The returned value will be between zero and
* dofs_per_cell.
*/
virtual unsigned int
- face_to_cell_index(const unsigned int face_dof_index,
- const unsigned int face,
- const bool face_orientation = true,
- const bool face_flip = false,
- const bool face_rotation = false) const override;
+ face_to_cell_index(
+ const unsigned int face_dof_index,
+ const unsigned int face,
+ const unsigned char combined_orientation =
+ ReferenceCell::default_combined_face_orientation()) const override;
/**
* Implementation of the respective function in the base class.
dof_values(i) +=
tmp * (normals[q_point] *
fe_values[vec].value(
- fe.face_to_cell_index(i,
- face,
- cell->face_orientation(face),
- cell->face_flip(face),
- cell->face_rotation(face)),
+ fe.face_to_cell_index(
+ i, face, cell->combined_face_orientation(face)),
q_point));
}
fe.get_nonzero_components(fe.face_to_cell_index(
i,
face,
- cell->face_orientation(face),
- cell->face_flip(face),
- cell->face_rotation(face)))[first_vector_component])
+ cell->combined_face_orientation(face)))[first_vector_component])
constraints.add_constraint(face_dof_indices[i],
{},
(std::abs(dof_values[i]) > 1e-14 ?
dof_values_local(i) +=
tmp * (normals[q_point] *
fe_values[vec].value(
- fe.face_to_cell_index(i,
- face,
- cell->face_orientation(face),
- cell->face_flip(face),
- cell->face_rotation(face)),
+ fe.face_to_cell_index(
+ i, face, cell->combined_face_orientation(face)),
q_point));
}
fe.get_nonzero_components(fe.face_to_cell_index(
i,
face,
- cell->face_orientation(face),
- cell->face_flip(face),
- cell->face_rotation(face)))[first_vector_component])
+ cell->combined_face_orientation(face)))[first_vector_component])
{
dof_values[face_dof_indices[i]] = dof_values_local(i);
projected_dofs[face_dof_indices[i]] = fe.degree;
static const int dim = FaceIterator::AccessorType::dimension;
static const int spacedim = FaceIterator::AccessorType::space_dimension;
+ const unsigned char combined_face_orientation =
+ ::dealii::internal::combined_face_orientation(face_orientation,
+ face_rotation,
+ face_flip);
+
const bool use_mg = (level != numbers::invalid_unsigned_int);
// If we don't use multigrid, we should be in the case where face_1 is
// Build up a cell to face index for face_2:
for (unsigned int i = 0; i < dofs_per_face; ++i)
{
- const unsigned int cell_index =
- fe.face_to_cell_index(i,
- 0, /* It doesn't really matter, just
- * assume we're on the first face...
- */
- true,
- false,
- false // default orientation
- );
+ const unsigned int cell_index = fe.face_to_cell_index(
+ i,
+ // It doesn't really matter, just assume we're on the first face...
+ 0);
cell_to_rotated_face_index[cell_index] = i;
}
// given orientation:
const unsigned int j =
cell_to_rotated_face_index[fe.face_to_cell_index(
- jj, 0, face_orientation, face_flip, face_rotation)];
+ jj, 0, combined_face_orientation)];
if (std::abs(transformation(i, jj)) > eps)
constraint_entries.emplace_back(dofs_1[j],
// orientation:
const unsigned int j =
cell_to_rotated_face_index[fe.face_to_cell_index(
- target, 0, face_orientation, face_flip, face_rotation)];
+ target, 0, combined_face_orientation)];
auto dof_left = dofs_1[j];
auto dof_right = dofs_2[i];
template <int dim, int spacedim>
unsigned int
-FiniteElement<dim, spacedim>::face_to_cell_index(const unsigned int face_index,
- const unsigned int face,
- const bool face_orientation,
- const bool face_flip,
- const bool face_rotation) const
+FiniteElement<dim, spacedim>::face_to_cell_index(
+ const unsigned int face_index,
+ const unsigned int face,
+ const unsigned char combined_orientation) const
{
AssertIndexRange(face_index, this->n_dofs_per_face(face));
AssertIndexRange(face, this->reference_cell().n_faces());
// see the function's documentation for an explanation of this
// assertion -- in essence, derived classes have to implement
// an overloaded version of this function if we are to use any
- // other than standard orientation
- if ((face_orientation != true) || (face_flip != false) ||
- (face_rotation != false))
+ // other than default (standard) orientation
+ if (combined_orientation !=
+ ReferenceCell::default_combined_face_orientation())
Assert((this->n_dofs_per_line() <= 1) && (this->n_dofs_per_quad(face) <= 1),
ExcMessage(
"The function in this base class can not handle this case. "
// then get the number of this vertex on the cell and translate
// this to a DoF number on the cell
return (this->reference_cell().face_to_cell_vertices(
- face,
- face_vertex,
- internal::combined_face_orientation(face_orientation,
- face_rotation,
- face_flip)) *
+ face, face_vertex, combined_orientation) *
this->n_dofs_per_vertex() +
dof_index_on_vertex);
}
const unsigned int dof_index_on_line = index % this->n_dofs_per_line();
return (this->get_first_line_index() +
- this->reference_cell().face_to_cell_lines(
- face,
- face_line,
- internal::combined_face_orientation(face_orientation,
- face_rotation,
- face_flip)) *
+ this->reference_cell().face_to_cell_lines(face,
+ face_line,
+ combined_orientation) *
this->n_dofs_per_line() +
dof_index_on_line);
}
template <int dim, int spacedim>
unsigned int
-FE_Q_Base<dim, spacedim>::face_to_cell_index(const unsigned int face_index,
- const unsigned int face,
- const bool face_orientation,
- const bool face_flip,
- const bool face_rotation) const
+FE_Q_Base<dim, spacedim>::face_to_cell_index(
+ const unsigned int face_index,
+ const unsigned int face,
+ const unsigned char combined_orientation) const
{
AssertIndexRange(face_index, this->n_dofs_per_face(face));
AssertIndexRange(face, GeometryInfo<dim>::faces_per_cell);
// then get the number of this vertex on the cell and translate
// this to a DoF number on the cell
- return (GeometryInfo<dim>::face_to_cell_vertices(
- face, face_vertex, face_orientation, face_flip, face_rotation) *
- this->n_dofs_per_vertex() +
- dof_index_on_vertex);
+ return this->reference_cell().face_to_cell_vertices(
+ face, face_vertex, combined_orientation) *
+ this->n_dofs_per_vertex() +
+ dof_index_on_vertex;
}
else if (face_index < this->get_first_face_quad_index(face))
// DoF is on a face
break;
case 2:
- // in 2d, only face_orientation has a meaning. if it is false (i.e.,
- // the non-default case), then consider dofs in reverse order
- if (face_orientation == true)
+ if (combined_orientation ==
+ ReferenceCell::default_combined_face_orientation())
adjusted_dof_index_on_line = dof_index_on_line;
else
adjusted_dof_index_on_line =
// that said, the Q2 case is easy enough to implement, as is the
// case where everything is in standard orientation
Assert((this->n_dofs_per_line() <= 1) ||
- ((face_orientation == true) && (face_flip == false) &&
- (face_rotation == false)),
+ combined_orientation ==
+ ReferenceCell::default_combined_face_orientation(),
ExcNotImplemented());
adjusted_dof_index_on_line = dof_index_on_line;
break;
}
return (this->get_first_line_index() +
- GeometryInfo<dim>::face_to_cell_lines(
- face, face_line, face_orientation, face_flip, face_rotation) *
+ this->reference_cell().face_to_cell_lines(face,
+ face_line,
+ combined_orientation) *
this->n_dofs_per_line() +
adjusted_dof_index_on_line);
}
// just have to draw a bunch of pictures. in the meantime,
// we can implement the Q2 case in which it is simple
Assert((this->n_dofs_per_quad(face) <= 1) ||
- ((face_orientation == true) && (face_flip == false) &&
- (face_rotation == false)),
+ combined_orientation ==
+ ReferenceCell::default_combined_face_orientation(),
ExcNotImplemented());
return (this->get_first_quad_index(face) + index);
}
template <int dim, int spacedim>
unsigned int
-FESystem<dim, spacedim>::face_to_cell_index(const unsigned int face_dof_index,
- const unsigned int face,
- const bool face_orientation,
- const bool face_flip,
- const bool face_rotation) const
+FESystem<dim, spacedim>::face_to_cell_index(
+ const unsigned int face_dof_index,
+ const unsigned int face,
+ const unsigned char combined_orientation) const
{
// we need to ask the base elements how they want to translate
// the DoFs within their own numbering. thus, translate to
const unsigned int base_face_to_cell_index =
this->base_element(face_base_index.first.first)
- .face_to_cell_index(face_base_index.second,
- face,
- face_orientation,
- face_flip,
- face_rotation);
+ .face_to_cell_index(face_base_index.second, face, combined_orientation);
// it would be nice if we had a base_to_system_index function, but
// all that exists is a component_to_system_index function. we can't do
{
deallog << "Face=" << face << std::endl;
- for (int orientation = 0; orientation < 2; ++orientation)
+ for (unsigned char orientation = 0; orientation < 2; ++orientation)
{
deallog << " orientation=" << (orientation == 0 ? "false" : "true")
<< std::endl
<< " ";
for (unsigned int i = 0; i < dofs_per_face; ++i)
- deallog << fe.face_to_cell_index(i,
- face,
- (orientation == 0 ? false : true),
- false,
- false)
- << " - ";
+ deallog << fe.face_to_cell_index(i, face, orientation) << " - ";
deallog << std::endl;
}
}
{
deallog << "Face=" << face << std::endl;
- for (int orientation = 0; orientation < 2; ++orientation)
+ for (unsigned char orientation = 0; orientation < 2; ++orientation)
{
deallog << " orientation=" << (orientation == 0 ? "false" : "true")
<< std::endl
<< " ";
for (unsigned int i = 0; i < dofs_per_face; ++i)
- deallog << fe.face_to_cell_index(i,
- face,
- (orientation == 0 ? false : true),
- false,
- false)
- << " - ";
+ deallog << fe.face_to_cell_index(i, face, orientation) << " - ";
deallog << std::endl;
}
}
{
deallog << "Face=" << face << std::endl;
- for (int orientation = 0; orientation < 2; ++orientation)
+ for (unsigned char orientation = 0; orientation < 2; ++orientation)
{
deallog << " orientation=" << (orientation == 0 ? "false" : "true")
<< std::endl
<< " ";
for (unsigned int i = 0; i < dofs_per_face; ++i)
- deallog << fe.face_to_cell_index(i,
- face,
- (orientation == 0 ? false : true),
- false,
- false)
- << " - ";
+ deallog << fe.face_to_cell_index(i, face, orientation) << " - ";
deallog << std::endl;
}
}
{
deallog << "Face=" << face << std::endl;
- for (int orientation = 0; orientation < 2; ++orientation)
+ for (unsigned char orientation = 0; orientation < 2; ++orientation)
{
deallog << " orientation=" << (orientation == 0 ? "false" : "true")
<< std::endl
<< " ";
for (unsigned int i = 0; i < dofs_per_face; ++i)
- deallog << fe.face_to_cell_index(i,
- face,
- (orientation == 0 ? false : true),
- false,
- false)
- << " - ";
+ deallog << fe.face_to_cell_index(i, face, orientation) << " - ";
deallog << std::endl;
}
}
{
deallog << "Face=" << face << std::endl;
- for (int orientation = 0; orientation < 2; ++orientation)
+ for (unsigned char orientation = 0; orientation < 2; ++orientation)
{
deallog << " orientation=" << (orientation == 0 ? "false" : "true")
<< std::endl
<< " ";
for (unsigned int i = 0; i < dofs_per_face; ++i)
- deallog << fe.face_to_cell_index(i,
- face,
- (orientation == 0 ? false : true),
- false,
- false)
- << " - ";
+ deallog << fe.face_to_cell_index(i, face, orientation) << " - ";
deallog << std::endl;
}
}
{
deallog << "Face=" << face << std::endl;
- for (int orientation = 0; orientation < 2; ++orientation)
+ for (unsigned char orientation = 0; orientation < 2; ++orientation)
{
deallog << " orientation=" << (orientation == 0 ? "false" : "true")
<< std::endl
<< " ";
for (unsigned int i = 0; i < dofs_per_face; ++i)
- deallog << fe.face_to_cell_index(i,
- face,
- (orientation == 0 ? false : true),
- false,
- false)
- << " - ";
+ deallog << fe.face_to_cell_index(i, face, orientation) << " - ";
deallog << std::endl;
}
}