Reference documentation for deal.II version 9.0.0
solution_transfer.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2009 - 2017 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 
17 #include <deal.II/base/config.h>
18 
19 #ifdef DEAL_II_WITH_P4EST
20 
21 #include <deal.II/lac/vector.h>
22 #include <deal.II/lac/block_vector.h>
23 #include <deal.II/lac/la_parallel_vector.h>
24 #include <deal.II/lac/la_parallel_block_vector.h>
25 #include <deal.II/lac/petsc_parallel_vector.h>
26 #include <deal.II/lac/petsc_parallel_block_vector.h>
27 #include <deal.II/lac/trilinos_vector.h>
28 #include <deal.II/lac/trilinos_parallel_block_vector.h>
29 
30 #include <deal.II/distributed/solution_transfer.h>
31 #include <deal.II/distributed/tria.h>
32 #include <deal.II/dofs/dof_tools.h>
33 #include <deal.II/dofs/dof_accessor.h>
34 #include <deal.II/grid/tria_accessor.h>
35 #include <deal.II/grid/tria_iterator.h>
36 
37 #include <functional>
38 
39 DEAL_II_NAMESPACE_OPEN
40 
41 namespace parallel
42 {
43  namespace distributed
44  {
45 
46  template <int dim, typename VectorType, typename DoFHandlerType>
48  :
49  dof_handler(&dof, typeid(*this).name()),
50  offset (numbers::invalid_unsigned_int)
51  {
53  (&dof_handler->get_triangulation()) != nullptr),
54  ExcMessage("parallel::distributed::SolutionTransfer requires a parallel::distributed::Triangulation object."));
55  }
56 
57 
58 
59  template <int dim, typename VectorType, typename DoFHandlerType>
60  void
62  prepare_for_coarsening_and_refinement (const std::vector<const VectorType *> &all_in)
63  {
64  input_vectors = all_in;
65  register_data_attach( get_data_size() * input_vectors.size() );
66  }
67 
68 
69 
70  template <int dim, typename VectorType, typename DoFHandlerType>
71  void
73  {
74  Assert(size > 0, ExcMessage("Please transfer at least one vector!"));
75 
76 //TODO: casting away constness is bad
80  (&dof_handler->get_triangulation())));
81  Assert (tria != nullptr, ExcInternalError());
82 
83  offset
84  = tria->register_data_attach(size,
85  std::bind(&SolutionTransfer<dim, VectorType,
86  DoFHandlerType>::pack_callback,
87  this,
88  std::placeholders::_1,
89  std::placeholders::_2,
90  std::placeholders::_3));
91 
92  }
93 
94 
95 
96  template <int dim, typename VectorType, typename DoFHandlerType>
97  void
100  {
101  std::vector<const VectorType *> all_in(1, &in);
102  prepare_for_coarsening_and_refinement(all_in);
103  }
104 
105 
106 
107  template <int dim, typename VectorType, typename DoFHandlerType>
108  void
110  {
111  std::vector<const VectorType *> all_in(1, &in);
112  prepare_serialization(all_in);
113  }
114 
115 
116 
117  template <int dim, typename VectorType, typename DoFHandlerType>
118  void
120  (const std::vector<const VectorType *> &all_in)
121  {
122  prepare_for_coarsening_and_refinement (all_in);
123  }
124 
125 
126 
127  template <int dim, typename VectorType, typename DoFHandlerType>
128  void
130  {
131  std::vector<VectorType *> all_in(1, &in);
132  deserialize(all_in);
133  }
134 
135 
136 
137  template <int dim, typename VectorType, typename DoFHandlerType>
138  void
140  {
141  register_data_attach( get_data_size() * all_in.size() );
142 
143  // this makes interpolate() happy
144  input_vectors.resize(all_in.size());
145 
146  interpolate(all_in);
147  }
148 
149 
150  template <int dim, typename VectorType, typename DoFHandlerType>
151  void
153  {
154  Assert(input_vectors.size()==all_out.size(),
155  ExcDimensionMismatch(input_vectors.size(), all_out.size()) );
156 
157 //TODO: casting away constness is bad
161  (&dof_handler->get_triangulation())));
162  Assert (tria != nullptr, ExcInternalError());
163 
164  tria->notify_ready_to_unpack(offset,
165  std::bind(&SolutionTransfer<dim, VectorType,
166  DoFHandlerType>::unpack_callback,
167  this,
168  std::placeholders::_1,
169  std::placeholders::_2,
170  std::placeholders::_3,
171  std::ref(all_out)));
172 
173 
174  for (typename std::vector<VectorType *>::iterator it=all_out.begin();
175  it !=all_out.end();
176  ++it)
177  (*it)->compress(::VectorOperation::insert);
178 
179  input_vectors.clear();
180  }
181 
182 
183 
184  template <int dim, typename VectorType, typename DoFHandlerType>
185  void
187  {
188  std::vector<VectorType *> all_out(1, &out);
189  interpolate(all_out);
190  }
191 
192 
193 
194  template <int dim, typename VectorType, typename DoFHandlerType>
195  unsigned int
197  {
198  return sizeof(typename VectorType::value_type)* DoFTools::max_dofs_per_cell(*dof_handler);
199  }
200 
201 
202  template <int dim, typename VectorType, typename DoFHandlerType>
203  void
206  const typename Triangulation<dim,DoFHandlerType::space_dimension>::CellStatus /*status*/,
207  void *data)
208  {
209  typename VectorType::value_type *data_store = reinterpret_cast<typename VectorType::value_type *>(data);
210 
211  typename DoFHandlerType::cell_iterator cell(*cell_, dof_handler);
212 
213  const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell;
214  ::Vector<typename VectorType::value_type> dofvalues(dofs_per_cell);
215  for (typename std::vector<const VectorType *>::iterator it=input_vectors.begin();
216  it !=input_vectors.end();
217  ++it)
218  {
219  cell->get_interpolated_dof_values(*(*it), dofvalues);
220  std::memcpy(data_store, &dofvalues(0), sizeof(typename VectorType::value_type)*dofs_per_cell);
221  data_store += dofs_per_cell;
222  }
223  }
224 
225 
226  template <int dim, typename VectorType, typename DoFHandlerType>
227  void
230  const typename Triangulation<dim,DoFHandlerType::space_dimension>::CellStatus /*status*/,
231  const void *data,
232  std::vector<VectorType *> &all_out)
233  {
234  typename DoFHandlerType::cell_iterator
235  cell(*cell_, dof_handler);
236 
237  const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell;
238  ::Vector<typename VectorType::value_type> dofvalues(dofs_per_cell);
239  const typename VectorType::value_type *data_store = reinterpret_cast<const typename VectorType::value_type *>(data);
240 
241  for (typename std::vector<VectorType *>::iterator it = all_out.begin();
242  it != all_out.end();
243  ++it)
244  {
245  std::memcpy(&dofvalues(0), data_store, sizeof(typename VectorType::value_type)*dofs_per_cell);
246  cell->set_dof_values_by_interpolation(dofvalues, *(*it));
247  data_store += dofs_per_cell;
248  }
249  }
250 
251 
252  }
253 }
254 
255 
256 // explicit instantiations
257 #include "solution_transfer.inst"
258 
259 DEAL_II_NAMESPACE_CLOSE
260 
261 #endif
void prepare_for_coarsening_and_refinement(const std::vector< const VectorType *> &all_in)
SolutionTransfer(const DoFHandlerType &dof)
void pack_callback(const typename Triangulation< dim, DoFHandlerType::space_dimension >::cell_iterator &cell, const typename Triangulation< dim, DoFHandlerType::space_dimension >::CellStatus status, void *data)
void interpolate(std::vector< VectorType *> &all_out)
unsigned int register_data_attach(const std::size_t size, const std::function< void(const cell_iterator &, const CellStatus, void *)> &pack_callback)
Definition: tria.cc:3224
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
unsigned int max_dofs_per_cell(const DoFHandler< dim, spacedim > &dh)
void notify_ready_to_unpack(const unsigned int handle, const std::function< void(const cell_iterator &, const CellStatus, const void *)> &unpack_callback)
Definition: tria.cc:3248
SmartPointer< const DoFHandlerType, SolutionTransfer< dim, VectorType, DoFHandlerType > > dof_handler
void prepare_serialization(const VectorType &in)
void unpack_callback(const typename Triangulation< dim, DoFHandlerType::space_dimension >::cell_iterator &cell, const typename Triangulation< dim, DoFHandlerType::space_dimension >::CellStatus status, const void *data, std::vector< VectorType *> &all_out)
::Triangulation< dim, spacedim >::cell_iterator cell_iterator
Definition: tria.h:265
static ::ExceptionBase & ExcInternalError()