Reference documentation for deal.II version 8.5.1
dof_accessor_get.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2016 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
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>
26 
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>
35 
36 #include <vector>
37 
38 DEAL_II_NAMESPACE_OPEN
39 
40 
41 template <typename DoFHandlerType, bool lda>
42 template <class InputVector, typename number>
43 void
45 get_interpolated_dof_values (const InputVector &values,
46  Vector<number> &interpolated_values,
47  const unsigned int fe_index) const
48 {
49  if (!this->has_children())
50  // if this cell has no children: simply return the exact values on this
51  // cell unless the finite element we need to interpolate to is different than
52  // the one we have on the current cell
53  {
55  (this->dof_handler)
56  != 0)
57  ||
58  // for hp-DoFHandlers, we need to require that on
59  // active cells, you either don't specify an fe_index,
60  // or that you specify the correct one
61  (fe_index == this->active_fe_index())
62  ||
63  (fe_index == DoFHandlerType::default_fe_index))
64  this->get_dof_values (values, interpolated_values);
65  else
66  {
67  // well, here we need to first get the values from the current
68  // cell and then interpolate it to the element requested. this
69  // can clearly only happen for hp::DoFHandler objects
70  Vector<number> tmp (this->get_fe().dofs_per_cell);
71  this->get_dof_values (values, tmp);
72 
73  FullMatrix<double> interpolation (this->dof_handler->get_fe()[fe_index].dofs_per_cell,
74  this->get_fe().dofs_per_cell);
75  this->dof_handler->get_fe()[fe_index].get_interpolation_matrix (this->get_fe(),
76  interpolation);
77  interpolation.vmult (interpolated_values, tmp);
78  }
79  }
80  else
81  // otherwise obtain them from the children
82  {
83  // we are on a non-active cell. these do not have any finite
84  // element associated with them in the hp context (in the non-hp
85  // context, we can simply assume that the FE space to which we
86  // want to interpolate is the same as for all elements in the
87  // mesh). consequently, we cannot interpolate from children's FE
88  // space to this cell's (unknown) FE space unless an explicit
89  // fe_index is given
91  (this->dof_handler)
92  != 0)
93  ||
94  (fe_index != DoFHandlerType::default_fe_index),
95  ExcMessage ("You cannot call this function on non-active cells "
96  "of hp::DoFHandler objects unless you provide an explicit "
97  "finite element index because they do not have naturally "
98  "associated finite element spaces associated: degrees "
99  "of freedom are only distributed on active cells for which "
100  "the active_fe_index has been set."));
101 
102  const FiniteElement<dim,spacedim> &fe = this->get_dof_handler().get_fe()[fe_index];
103  const unsigned int dofs_per_cell = fe.dofs_per_cell;
104 
105  Assert (this->dof_handler != 0,
106  typename BaseClass::ExcInvalidObject());
107  Assert (interpolated_values.size() == dofs_per_cell,
108  typename BaseClass::ExcVectorDoesNotMatch());
109  Assert (values.size() == this->dof_handler->n_dofs(),
110  typename BaseClass::ExcVectorDoesNotMatch());
111 
112 
113  // see if the finite element we have on the current cell has any
114  // degrees of freedom to begin with; if not (e.g., when
115  // interpolating FE_Nothing), then simply skip all of the
116  // following since the output vector would be of size zero
117  // anyway (and in fact is of size zero, see the assertion above)
118  if (fe.dofs_per_cell > 0)
119  {
120  Vector<number> tmp1(dofs_per_cell);
121  Vector<number> tmp2(dofs_per_cell);
122 
123  interpolated_values = 0;
124 
125  // later on we will have to push the values interpolated from the child
126  // to the mother cell into the output vector. unfortunately, there are
127  // two types of elements: ones where you add up the contributions from
128  // the different child cells, and ones where you overwrite.
129  //
130  // an example for the first is piecewise constant (and discontinuous)
131  // elements, where we build the value on the coarse cell by averaging
132  // the values from the cell (i.e. by adding up a fraction of the values
133  // of their values)
134  //
135  // an example for the latter are the usual continuous elements. the
136  // value on a vertex of a coarse cell must there be the same,
137  // irrespective of the adjacent cell we are presently on. so we always
138  // overwrite. in fact, we must, since we cannot know in advance how many
139  // neighbors there will be, so there is no way to compute the average
140  // with fixed factors
141  //
142  // so we have to find out to which type this element belongs. the
143  // difficulty is: the finite element may be a composed one, so we can
144  // only hope to do this for each shape function individually. in fact,
145  // there are even weird finite elements (for example the Raviart-Thomas
146  // element) which have shape functions that are additive (interior ones)
147  // and others that are overwriting (face degrees of freedom that need to
148  // be continuous across the face).
149  for (unsigned int child=0; child<this->n_children(); ++child)
150  {
151  // get the values from the present child, if necessary by
152  // interpolation itself either from its own children or
153  // by interpolating from the finite element on an active
154  // child to the finite element space requested here
155  this->child(child)->get_interpolated_dof_values (values,
156  tmp1,
157  fe_index);
158  // interpolate these to the mother cell
159  fe.get_restriction_matrix(child, this->refinement_case()).vmult (tmp2, tmp1);
160 
161  // and add up or set them in the output vector
162  for (unsigned int i=0; i<dofs_per_cell; ++i)
163  if (fe.restriction_is_additive(i))
164  interpolated_values(i) += tmp2(i);
165  else if (tmp2(i) != number())
166  interpolated_values(i) = tmp2(i);
167  }
168  }
169  }
170 }
171 
172 
173 // --------------------------------------------------------------------------
174 // explicit instantiations
175 #include "dof_accessor_get.inst"
176 
177 DEAL_II_NAMESPACE_CLOSE
bool restriction_is_additive(const unsigned int index) const
Definition: fe.h:3000
void vmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
void get_interpolated_dof_values(const InputVector &values, Vector< number > &interpolated_values, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:313
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe.cc:305
std::size_t size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:295