Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
dof_accessor_set.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1998 - 2020 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
19
20#include <deal.II/fe/fe.h>
21
23#include <deal.II/grid/tria_iterator.templates.h>
24
26
38#include <deal.II/lac/vector.h>
39
40#include <vector>
41
43
44
45template <int dim, int spacedim, bool lda>
46template <class OutputVector, typename number>
47void
49 const Vector<number> &local_values,
50 OutputVector & values,
51 const unsigned int fe_index_) const
52{
53 const unsigned int fe_index =
54 (this->dof_handler->hp_capability_enabled == false &&
57 fe_index_;
58
59 if (this->is_active() && !this->is_artificial())
60 {
61 if ((this->dof_handler->hp_capability_enabled == false) ||
62 // for hp-DoFHandlers, we need to require that on
63 // active cells, you either don't specify an fe_index,
64 // or that you specify the correct one
65 (fe_index == this->active_fe_index()) ||
67 // simply set the values on this cell
68 this->set_dof_values(local_values, values);
69 else
70 {
71 Assert(local_values.size() ==
72 this->dof_handler->get_fe(fe_index).n_dofs_per_cell(),
73 ExcMessage("Incorrect size of local_values vector."));
74
75 FullMatrix<double> interpolation(
76 this->get_fe().n_dofs_per_cell(),
77 this->dof_handler->get_fe(fe_index).n_dofs_per_cell());
78
79 this->get_fe().get_interpolation_matrix(
80 this->dof_handler->get_fe(fe_index), interpolation);
81
82 // do the interpolation to the target space. for historical
83 // reasons, matrices are set to size 0x0 internally even
84 // we reinit as 4x0, so we have to treat this case specially
85 Vector<number> tmp(this->get_fe().n_dofs_per_cell());
86 if ((tmp.size() > 0) && (local_values.size() > 0))
87 interpolation.vmult(tmp, local_values);
88
89 // now set the dof values in the global vector
90 this->set_dof_values(tmp, values);
91 }
92 }
93 else
94 // otherwise distribute them to the children
95 {
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."));
105
107 this->get_dof_handler().get_fe(fe_index);
108 const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
109
110 Assert(this->dof_handler != nullptr,
111 typename BaseClass::ExcInvalidObject());
112 Assert(local_values.size() == dofs_per_cell,
113 typename BaseClass::ExcVectorDoesNotMatch());
114 Assert(values.size() == this->dof_handler->n_dofs(),
115 typename BaseClass::ExcVectorDoesNotMatch());
116
117 Vector<number> tmp(dofs_per_cell);
118
119 for (unsigned int child = 0; child < this->n_children(); ++child)
120 {
121 if (tmp.size() > 0)
122 fe.get_prolongation_matrix(child, this->refinement_case())
123 .vmult(tmp, local_values);
124 this->child(child)->set_dof_values_by_interpolation(tmp,
125 values,
126 fe_index);
127 }
128 }
129}
130
131
132// --------------------------------------------------------------------------
133// explicit instantiations
134#include "dof_accessor_set.inst"
135
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
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcMessage(std::string arg1)
size_type size() const