(this->dof_handler)
!= 0)
||
+ // for hp-DoFHandlers, we need to require that on
+ // active cells, you either don't specify an fe_index,
+ // or that you specify the correct one
(fe_index == this->active_fe_index())
||
(fe_index == DH::default_fe_index))
FullMatrix<double> interpolation (this->dof_handler->get_fe()[fe_index].dofs_per_cell,
this->get_fe().dofs_per_cell);
this->dof_handler->get_fe()[fe_index].get_interpolation_matrix (this->get_fe(),
- interpolation);
+ interpolation);
interpolation.vmult (interpolated_values, tmp);
}
}
// space to this cell's (unknown) FE space unless an explicit
// fe_index is given
Assert ((dynamic_cast<DoFHandler<DH::dimension,DH::space_dimension>*>
- (this->dof_handler)
- != 0)
+ (this->dof_handler)
+ != 0)
||
- (fe_index != DH::default_fe_index),
- ExcMessage ("You cannot call this function on non-active cells "
- "of hp::DoFHandler objects unless you provide an explicit "
- "finite element index because they do not have naturally "
- "associated finite element spaces associated: degrees "
- "of freedom are only distributed on active cells for which "
- "the active_fe_index has been set."));
+ (fe_index != DH::default_fe_index),
+ ExcMessage ("You cannot call this function on non-active cells "
+ "of hp::DoFHandler objects unless you provide an explicit "
+ "finite element index because they do not have naturally "
+ "associated finite element spaces associated: degrees "
+ "of freedom are only distributed on active cells for which "
+ "the active_fe_index has been set."));
const FiniteElement<dim,spacedim> &fe = this->get_dof_handler().get_fe()[fe_index];
const unsigned int dofs_per_cell = fe.dofs_per_cell;
Assert (this->dof_handler != 0,
- typename BaseClass::ExcInvalidObject());
+ typename BaseClass::ExcInvalidObject());
Assert (&fe != 0,
- typename BaseClass::ExcInvalidObject());
+ typename BaseClass::ExcInvalidObject());
Assert (interpolated_values.size() == dofs_per_cell,
- typename BaseClass::ExcVectorDoesNotMatch());
+ typename BaseClass::ExcVectorDoesNotMatch());
Assert (values.size() == this->dof_handler->n_dofs(),
- typename BaseClass::ExcVectorDoesNotMatch());
+ typename BaseClass::ExcVectorDoesNotMatch());
Vector<number> tmp1(dofs_per_cell);
void
DoFCellAccessor<DH,lda>::
set_dof_values_by_interpolation (const Vector<number> &local_values,
- OutputVector &values,
- const unsigned int fe_index) const
+ OutputVector &values,
+ const unsigned int fe_index) const
+{
+ if (!this->has_children())
{
- if (!this->has_children())
- {
- if ((dynamic_cast<DoFHandler<DH::dimension,DH::space_dimension>*>
- (this->dof_handler)
- != 0)
- ||
- (fe_index == this->active_fe_index())
- ||
- (fe_index == DH::default_fe_index))
- // simply set the values on this cell
- this->set_dof_values (local_values, values);
- else
- {
- Assert (local_values.size() == this->dof_handler->get_fe()[fe_index].dofs_per_cell,
- ExcMessage ("Incorrect size of local_values vector.") );
-
- FullMatrix<double> interpolation (this->get_fe().dofs_per_cell, this->dof_handler->get_fe()[fe_index].dofs_per_cell);
-
- this->get_fe().get_interpolation_matrix (this->dof_handler->get_fe()[fe_index],
- interpolation);
-
- Vector<number> tmp (this->get_fe().dofs_per_cell);
- interpolation.vmult (tmp, local_values);
-
- //now set the dof values in the global vector
- this->set_dof_values (tmp, values);
- }
- }
- else
- // otherwise distribute them to the children
- {
- Assert ((dynamic_cast<DoFHandler<DH::dimension,DH::space_dimension>*>
- (this->dof_handler)
- != 0)
- ||
- (fe_index != DH::default_fe_index),
- ExcMessage ("You cannot call this function on non-active cells "
- "of hp::DoFHandler objects unless you provide an explicit "
- "finite element index because they do not have naturally "
- "associated finite element spaces associated: degrees "
- "of freedom are only distributed on active cells for which "
- "the active_fe_index has been set."));
-
- const FiniteElement<dim,spacedim> &fe = this->get_dof_handler().get_fe()[fe_index];
- const unsigned int dofs_per_cell = fe.dofs_per_cell;
-
- //const unsigned int dofs_per_cell = this->get_fe().dofs_per_cell;
-
- Assert (this->dof_handler != 0,
- typename BaseClass::ExcInvalidObject());
- Assert (&this->get_fe() != 0,
- typename BaseClass::ExcInvalidObject());
- Assert (local_values.size() == dofs_per_cell,
- typename BaseClass::ExcVectorDoesNotMatch());
- Assert (values.size() == this->dof_handler->n_dofs(),
- typename BaseClass::ExcVectorDoesNotMatch());
-
- Vector<number> tmp(dofs_per_cell);
-
- for (unsigned int child=0; child<this->n_children(); ++child)
- {
- fe.get_prolongation_matrix(child, this->refinement_case()).vmult (tmp, local_values);
- this->child(child)->set_dof_values_by_interpolation (tmp, values, fe_index);
- }
- }
+ if ((dynamic_cast<DoFHandler<DH::dimension,DH::space_dimension>*>
+ (this->dof_handler)
+ != 0)
+ ||
+ // for hp-DoFHandlers, we need to require that on
+ // active cells, you either don't specify an fe_index,
+ // or that you specify the correct one
+ (fe_index == this->active_fe_index())
+ ||
+ (fe_index == DH::default_fe_index))
+ // simply set the values on this cell
+ this->set_dof_values (local_values, values);
+ else
+ {
+ Assert (local_values.size() == this->dof_handler->get_fe()[fe_index].dofs_per_cell,
+ ExcMessage ("Incorrect size of local_values vector.") );
+
+ FullMatrix<double> interpolation (this->get_fe().dofs_per_cell, this->dof_handler->get_fe()[fe_index].dofs_per_cell);
+
+ this->get_fe().get_interpolation_matrix (this->dof_handler->get_fe()[fe_index],
+ interpolation);
+
+ Vector<number> tmp (this->get_fe().dofs_per_cell);
+ interpolation.vmult (tmp, local_values);
+
+ //now set the dof values in the global vector
+ this->set_dof_values (tmp, values);
+ }
}
+ else
+ // otherwise distribute them to the children
+ {
+ Assert ((dynamic_cast<DoFHandler<DH::dimension,DH::space_dimension>*>
+ (this->dof_handler)
+ != 0)
+ ||
+ (fe_index != DH::default_fe_index),
+ ExcMessage ("You cannot call this function on non-active cells "
+ "of hp::DoFHandler objects unless you provide an explicit "
+ "finite element index because they do not have naturally "
+ "associated finite element spaces associated: degrees "
+ "of freedom are only distributed on active cells for which "
+ "the active_fe_index has been set."));
+
+ const FiniteElement<dim,spacedim> &fe = this->get_dof_handler().get_fe()[fe_index];
+ const unsigned int dofs_per_cell = fe.dofs_per_cell;
+
+ Assert (this->dof_handler != 0,
+ typename BaseClass::ExcInvalidObject());
+ Assert (&this->get_fe() != 0,
+ typename BaseClass::ExcInvalidObject());
+ Assert (local_values.size() == dofs_per_cell,
+ typename BaseClass::ExcVectorDoesNotMatch());
+ Assert (values.size() == this->dof_handler->n_dofs(),
+ typename BaseClass::ExcVectorDoesNotMatch());
+
+ Vector<number> tmp(dofs_per_cell);
+
+ for (unsigned int child=0; child<this->n_children(); ++child)
+ {
+ fe.get_prolongation_matrix(child, this->refinement_case()).vmult (tmp, local_values);
+ this->child(child)->set_dof_values_by_interpolation (tmp, values, fe_index);
+ }
+ }
+}
// --------------------------------------------------------------------------