/**
* Constructor taking vectors on the multigrid levels rather than the active
* cells only. The vector of vectors is expected to have as many entries as
- * there are levels in the triangulation and provide valid data on each
- * level, i.e., be of compatible length DoFHandler::n_dofs(level). Apart
- * from the source of the vectors, the same arguments as in the other
- * constructor need to be provided.
+ * there are global levels in the triangulation and provide valid data on
+ * each level, i.e., be of compatible length DoFHandler::n_dofs(level). A
+ * prerequisite of this constructor is that DoFHandler::distribute_mg_dofs()
+ * has been called. Apart from the level vectors, the same arguments as in
+ * the other constructor need to be provided.
*/
MappingFEField(const DoFHandlerType & euler_dof_handler,
const std::vector<VectorType> &euler_vector,
const ComponentMask & mask = ComponentMask());
/**
- * Constructor taking vectors on the multigrid levels rather than the active
- * cells only, variant with MGLevelObject instead of std::vector. The vector
- * of vectors is expected to have as many entries as there are levels in the
- * triangulation and provide valid data on each level, i.e., be of
- * compatible length DoFHandler::n_dofs(level). Apart from the source of the
- * vectors, the same arguments as in the other constructor need to be
- * provided.
+ * Constructor with MGLevelObject instead of std::vector, otherwise the same
+ * as above. It is required that `euler_vector.max_level()+1` equals the
+ * global number of levels in the triangulation. The minimum level may be
+ * zero or more — it only needs to be consistent between what is set
+ * here and later used for evaluation of the mapping.
*/
MappingFEField(const DoFHandlerType & euler_dof_handler,
const MGLevelObject<VectorType> &euler_vector,
AssertDimension(size, spacedim);
Assert(euler_dof_handler.has_level_dofs(),
- ExcMessage("The underlying did not call distribute_mg_dofs(). "
- "In this case, the construction via level vectors does not "
- "make sense."));
+ ExcMessage("The underlying DoFHandler object did not call "
+ "distribute_mg_dofs(). In this case, the construction via "
+ "level vectors does not make sense."));
AssertDimension(euler_vector.size(),
euler_dof_handler.get_triangulation().n_global_levels());
this->euler_vector.clear();
AssertDimension(size, spacedim);
Assert(euler_dof_handler.has_level_dofs(),
- ExcMessage("The underlying did not call distribute_mg_dofs(). "
- "In this case, the construction via level vectors does not "
- "make sense."));
- AssertDimension(euler_vector.max_level() + 1 - euler_vector.min_level(),
+ ExcMessage("The underlying DoFHandler object did not call "
+ "distribute_mg_dofs(). In this case, the construction via "
+ "level vectors does not make sense."));
+ AssertDimension(euler_vector.max_level() + 1,
euler_dof_handler.get_triangulation().n_global_levels());
this->euler_vector.clear();
this->euler_vector.resize(
euler_dof_handler.get_triangulation().n_global_levels());
for (unsigned int i = euler_vector.min_level(); i <= euler_vector.max_level();
++i)
- this->euler_vector[i - euler_vector.min_level()] = &euler_vector[i];
+ this->euler_vector[i] = &euler_vector[i];
}
fe_values.n_quadrature_points);
AssertDimension(fe_to_real.size(),
euler_dof_handler->get_fe().n_components());
+ if (uses_level_dofs)
+ {
+ AssertIndexRange(cell->level(), euler_vector.size());
+ AssertDimension(euler_vector[cell->level()]->size(),
+ euler_dof_handler->n_dofs(cell->level()));
+ }
+ else
+ AssertDimension(euler_vector[0]->size(), euler_dof_handler->n_dofs());
{
std::lock_guard<std::mutex> lock(fe_values_mutex);
const VectorType &vector =
uses_level_dofs ? *euler_vector[cell->level()] : *euler_vector[0];
+
std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell> vertices;
for (unsigned int i = 0; i < dofs_per_cell; ++i)
{
typename DoFHandlerType::cell_iterator dof_cell(*cell, euler_dof_handler);
Assert(uses_level_dofs || dof_cell->active() == true, ExcInactiveCell());
+ if (uses_level_dofs)
+ {
+ AssertIndexRange(cell->level(), euler_vector.size());
+ AssertDimension(euler_vector[cell->level()]->size(),
+ euler_dof_handler->n_dofs(cell->level()));
+ }
+ else
+ AssertDimension(euler_vector[0]->size(), euler_dof_handler->n_dofs());
if (uses_level_dofs)
dof_cell->get_mg_dof_indices(data.local_dof_indices);
const VectorType &vector =
uses_level_dofs ? *euler_vector[cell->level()] : *euler_vector[0];
+
for (unsigned int i = 0; i < data.local_dof_values.size(); ++i)
data.local_dof_values[i] =
internal::ElementAccess<VectorType>::get(vector,