16 #include <deal.II/lac/vector.h> 17 #include <deal.II/lac/block_vector.h> 18 #include <deal.II/lac/la_vector.h> 19 #include <deal.II/lac/la_parallel_vector.h> 20 #include <deal.II/lac/la_parallel_block_vector.h> 21 #include <deal.II/lac/petsc_vector.h> 22 #include <deal.II/lac/petsc_block_vector.h> 23 #include <deal.II/lac/trilinos_vector.h> 24 #include <deal.II/lac/trilinos_block_vector.h> 25 #include <deal.II/lac/sparse_matrix.h> 27 #include <deal.II/dofs/dof_accessor.h> 28 #include <deal.II/dofs/dof_handler.h> 29 #include <deal.II/dofs/dof_levels.h> 30 #include <deal.II/hp/dof_handler.h> 31 #include <deal.II/grid/tria_boundary.h> 32 #include <deal.II/grid/tria_iterator.h> 33 #include <deal.II/grid/tria_iterator.templates.h> 34 #include <deal.II/fe/fe.h> 38 DEAL_II_NAMESPACE_OPEN
41 template <
typename DoFHandlerType,
bool lda>
42 template <
class OutputVector,
typename number>
47 const unsigned int fe_index)
const 49 if (!this->has_children() && !this->is_artificial ())
58 (fe_index == this->active_fe_index())
60 (fe_index == DoFHandlerType::default_fe_index))
62 this->set_dof_values (local_values, values);
65 Assert (local_values.
size() == this->dof_handler->get_fe()[fe_index].dofs_per_cell,
66 ExcMessage (
"Incorrect size of local_values vector.") );
68 FullMatrix<double> interpolation (this->get_fe().dofs_per_cell, this->dof_handler->get_fe()[fe_index].dofs_per_cell);
70 this->get_fe().get_interpolation_matrix (this->dof_handler->get_fe()[fe_index],
77 if ((tmp.
size() > 0) && (local_values.
size() > 0))
78 interpolation.vmult (tmp, local_values);
81 this->set_dof_values (tmp, values);
91 (fe_index != DoFHandlerType::default_fe_index),
92 ExcMessage (
"You cannot call this function on non-active cells " 93 "of hp::DoFHandler objects unless you provide an explicit " 94 "finite element index because they do not have naturally " 95 "associated finite element spaces associated: degrees " 96 "of freedom are only distributed on active cells for which " 97 "the active_fe_index has been set."));
102 Assert (this->dof_handler != 0,
103 typename BaseClass::ExcInvalidObject());
105 typename BaseClass::ExcVectorDoesNotMatch());
106 Assert (values.size() == this->dof_handler->n_dofs(),
107 typename BaseClass::ExcVectorDoesNotMatch());
111 for (
unsigned int child=0; child<this->n_children(); ++child)
115 .vmult (tmp, local_values);
116 this->child(child)->set_dof_values_by_interpolation (tmp, values, fe_index);
124 #include "dof_accessor_set.inst" 126 DEAL_II_NAMESPACE_CLOSE
void set_dof_values_by_interpolation(const Vector< number > &local_values, OutputVector &values, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
const unsigned int dofs_per_cell