23#include <deal.II/grid/tria_iterator.templates.h>
45template <
int dim,
int spacedim,
bool lda>
46template <
class OutputVector,
typename number>
51 const unsigned int fe_index_)
const
53 const unsigned int fe_index =
54 (this->dof_handler->hp_capability_enabled ==
false &&
59 if (this->is_active() && !this->is_artificial())
61 if ((this->dof_handler->hp_capability_enabled ==
false) ||
65 (fe_index == this->active_fe_index()) ||
68 this->set_dof_values(local_values,
values);
72 this->dof_handler->get_fe(fe_index).n_dofs_per_cell(),
73 ExcMessage(
"Incorrect size of local_values vector."));
76 this->get_fe().n_dofs_per_cell(),
77 this->dof_handler->get_fe(fe_index).n_dofs_per_cell());
79 this->get_fe().get_interpolation_matrix(
80 this->dof_handler->get_fe(fe_index), interpolation);
86 if ((tmp.
size() > 0) && (local_values.
size() > 0))
87 interpolation.
vmult(tmp, local_values);
90 this->set_dof_values(tmp,
values);
96 Assert((this->dof_handler->hp_capability_enabled ==
false) ||
99 "You cannot call this function on non-active cells "
100 "of DoFHandler objects unless you provide an explicit "
101 "finite element index because they do not have naturally "
102 "associated finite element spaces associated: degrees "
103 "of freedom are only distributed on active cells for which "
104 "the active FE index has been set."));
107 this->get_dof_handler().get_fe(fe_index);
110 Assert(this->dof_handler !=
nullptr,
111 typename BaseClass::ExcInvalidObject());
113 typename BaseClass::ExcVectorDoesNotMatch());
115 typename BaseClass::ExcVectorDoesNotMatch());
119 for (
unsigned int child = 0; child < this->n_children(); ++child)
123 .vmult(tmp, local_values);
124 this->child(child)->set_dof_values_by_interpolation(tmp,
134#include "dof_accessor_set.inst"
void set_dof_values_by_interpolation(const Vector< number > &local_values, OutputVector &values, const unsigned int fe_index=DoFHandler< dimension_, space_dimension_ >::invalid_fe_index) const
unsigned int n_dofs_per_cell() const
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
void vmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcMessage(std::string arg1)