Reference documentation for deal.II version 8.4.2
dof_accessor_set.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2015 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/parallel_vector.h>
19 #include <deal.II/lac/parallel_block_vector.h>
20 #include <deal.II/lac/petsc_vector.h>
21 #include <deal.II/lac/petsc_block_vector.h>
22 #include <deal.II/lac/trilinos_vector.h>
23 #include <deal.II/lac/trilinos_block_vector.h>
24 #include <deal.II/lac/sparse_matrix.h>
25 
26 #include <deal.II/dofs/dof_accessor.h>
27 #include <deal.II/dofs/dof_handler.h>
28 #include <deal.II/dofs/dof_levels.h>
29 #include <deal.II/hp/dof_handler.h>
30 #include <deal.II/grid/tria_boundary.h>
31 #include <deal.II/grid/tria_iterator.h>
32 #include <deal.II/grid/tria_iterator.templates.h>
33 #include <deal.II/fe/fe.h>
34 
35 #include <vector>
36 
37 DEAL_II_NAMESPACE_OPEN
38 
39 
40 template <typename DoFHandlerType, bool lda>
41 template <class OutputVector, typename number>
42 void
45  OutputVector &values,
46  const unsigned int fe_index) const
47 {
48  if (!this->has_children() && !this->is_artificial ())
49  {
51  (this->dof_handler)
52  != 0)
53  ||
54  // for hp-DoFHandlers, we need to require that on
55  // active cells, you either don't specify an fe_index,
56  // or that you specify the correct one
57  (fe_index == this->active_fe_index())
58  ||
59  (fe_index == DoFHandlerType::default_fe_index))
60  // simply set the values on this cell
61  this->set_dof_values (local_values, values);
62  else
63  {
64  Assert (local_values.size() == this->dof_handler->get_fe()[fe_index].dofs_per_cell,
65  ExcMessage ("Incorrect size of local_values vector.") );
66 
67  FullMatrix<double> interpolation (this->get_fe().dofs_per_cell, this->dof_handler->get_fe()[fe_index].dofs_per_cell);
68 
69  this->get_fe().get_interpolation_matrix (this->dof_handler->get_fe()[fe_index],
70  interpolation);
71 
72  // do the interpolation to the target space. for historical
73  // reasons, matrices are set to size 0x0 internally even
74  // we reinit as 4x0, so we have to treat this case specially
75  Vector<number> tmp (this->get_fe().dofs_per_cell);
76  if ((tmp.size() > 0) && (local_values.size() > 0))
77  interpolation.vmult (tmp, local_values);
78 
79  // now set the dof values in the global vector
80  this->set_dof_values (tmp, values);
81  }
82  }
83  else
84  // otherwise distribute them to the children
85  {
87  (this->dof_handler)
88  != 0)
89  ||
90  (fe_index != DoFHandlerType::default_fe_index),
91  ExcMessage ("You cannot call this function on non-active cells "
92  "of hp::DoFHandler objects unless you provide an explicit "
93  "finite element index because they do not have naturally "
94  "associated finite element spaces associated: degrees "
95  "of freedom are only distributed on active cells for which "
96  "the active_fe_index has been set."));
97 
98  const FiniteElement<dim,spacedim> &fe = this->get_dof_handler().get_fe()[fe_index];
99  const unsigned int dofs_per_cell = fe.dofs_per_cell;
100 
101  Assert (this->dof_handler != 0,
102  typename BaseClass::ExcInvalidObject());
103  Assert (local_values.size() == dofs_per_cell,
104  typename BaseClass::ExcVectorDoesNotMatch());
105  Assert (values.size() == this->dof_handler->n_dofs(),
106  typename BaseClass::ExcVectorDoesNotMatch());
107 
108  Vector<number> tmp(dofs_per_cell);
109 
110  for (unsigned int child=0; child<this->n_children(); ++child)
111  {
112  if (tmp.size() > 0)
113  fe.get_prolongation_matrix(child, this->refinement_case())
114  .vmult (tmp, local_values);
115  this->child(child)->set_dof_values_by_interpolation (tmp, values, fe_index);
116  }
117  }
118 }
119 
120 
121 // --------------------------------------------------------------------------
122 // explicit instantiations
123 #include "dof_accessor_set.inst"
124 
125 DEAL_II_NAMESPACE_CLOSE
::ExceptionBase & ExcMessage(std::string arg1)
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:326
#define Assert(cond, exc)
Definition: exceptions.h:294
std::size_t size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:283