Reference documentation for deal.II version 8.5.1
dof_accessor_set.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 OutputVector, typename number>
43 void
46  OutputVector &values,
47  const unsigned int fe_index) const
48 {
49  if (!this->has_children() && !this->is_artificial ())
50  {
52  (this->dof_handler)
53  != 0)
54  ||
55  // for hp-DoFHandlers, we need to require that on
56  // active cells, you either don't specify an fe_index,
57  // or that you specify the correct one
58  (fe_index == this->active_fe_index())
59  ||
60  (fe_index == DoFHandlerType::default_fe_index))
61  // simply set the values on this cell
62  this->set_dof_values (local_values, values);
63  else
64  {
65  Assert (local_values.size() == this->dof_handler->get_fe()[fe_index].dofs_per_cell,
66  ExcMessage ("Incorrect size of local_values vector.") );
67 
68  FullMatrix<double> interpolation (this->get_fe().dofs_per_cell, this->dof_handler->get_fe()[fe_index].dofs_per_cell);
69 
70  this->get_fe().get_interpolation_matrix (this->dof_handler->get_fe()[fe_index],
71  interpolation);
72 
73  // do the interpolation to the target space. for historical
74  // reasons, matrices are set to size 0x0 internally even
75  // we reinit as 4x0, so we have to treat this case specially
76  Vector<number> tmp (this->get_fe().dofs_per_cell);
77  if ((tmp.size() > 0) && (local_values.size() > 0))
78  interpolation.vmult (tmp, local_values);
79 
80  // now set the dof values in the global vector
81  this->set_dof_values (tmp, values);
82  }
83  }
84  else
85  // otherwise distribute them to the children
86  {
88  (this->dof_handler)
89  != 0)
90  ||
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."));
98 
99  const FiniteElement<dim,spacedim> &fe = this->get_dof_handler().get_fe()[fe_index];
100  const unsigned int dofs_per_cell = fe.dofs_per_cell;
101 
102  Assert (this->dof_handler != 0,
103  typename BaseClass::ExcInvalidObject());
104  Assert (local_values.size() == dofs_per_cell,
105  typename BaseClass::ExcVectorDoesNotMatch());
106  Assert (values.size() == this->dof_handler->n_dofs(),
107  typename BaseClass::ExcVectorDoesNotMatch());
108 
109  Vector<number> tmp(dofs_per_cell);
110 
111  for (unsigned int child=0; child<this->n_children(); ++child)
112  {
113  if (tmp.size() > 0)
114  fe.get_prolongation_matrix(child, this->refinement_case())
115  .vmult (tmp, local_values);
116  this->child(child)->set_dof_values_by_interpolation (tmp, values, fe_index);
117  }
118  }
119 }
120 
121 
122 // --------------------------------------------------------------------------
123 // explicit instantiations
124 #include "dof_accessor_set.inst"
125 
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
Definition: fe.cc:325
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:313
std::size_t size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:295