Reference documentation for deal.II version 8.5.1
solution_transfer.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2009 - 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 
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_vector.h>
26 #include <deal.II/lac/petsc_block_vector.h>
27 #include <deal.II/lac/trilinos_vector.h>
28 #include <deal.II/lac/trilinos_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 <deal.II/base/std_cxx11/bind.h>
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()) != 0,
54  ExcMessage("parallel::distributed::SolutionTransfer requires a parallel::distributed::Triangulation object."));
55  }
56 
57 
58 
59  template<int dim, typename VectorType, typename DoFHandlerType>
61  {}
62 
63 
64 
65  template<int dim, typename VectorType, typename DoFHandlerType>
66  void
68  prepare_for_coarsening_and_refinement (const std::vector<const VectorType *> &all_in)
69  {
70  input_vectors = all_in;
71  register_data_attach( get_data_size() * input_vectors.size() );
72  }
73 
74 
75 
76  template<int dim, typename VectorType, typename DoFHandlerType>
77  void
79  {
80  Assert(size > 0, ExcMessage("Please transfer at least one vector!"));
81 
82 //TODO: casting away constness is bad
86  (&dof_handler->get_triangulation())));
87  Assert (tria != 0, ExcInternalError());
88 
89  offset
90  = tria->register_data_attach(size,
91  std_cxx11::bind(&SolutionTransfer<dim, VectorType,
92  DoFHandlerType>::pack_callback,
93  this,
94  std_cxx11::_1,
95  std_cxx11::_2,
96  std_cxx11::_3));
97 
98  }
99 
100 
101 
102  template<int dim, typename VectorType, typename DoFHandlerType>
103  void
106  {
107  std::vector<const VectorType *> all_in(1, &in);
108  prepare_for_coarsening_and_refinement(all_in);
109  }
110 
111 
112 
113  template<int dim, typename VectorType, typename DoFHandlerType>
114  void
116  {
117  std::vector<const VectorType *> all_in(1, &in);
118  prepare_serialization(all_in);
119  }
120 
121 
122 
123  template<int dim, typename VectorType, typename DoFHandlerType>
124  void
126  (const std::vector<const VectorType *> &all_in)
127  {
128  prepare_for_coarsening_and_refinement (all_in);
129  }
130 
131 
132 
133  template<int dim, typename VectorType, typename DoFHandlerType>
134  void
136  {
137  std::vector<VectorType *> all_in(1, &in);
138  deserialize(all_in);
139  }
140 
141 
142 
143  template<int dim, typename VectorType, typename DoFHandlerType>
144  void
146  {
147  register_data_attach( get_data_size() * all_in.size() );
148 
149  // this makes interpolate() happy
150  input_vectors.resize(all_in.size());
151 
152  interpolate(all_in);
153  }
154 
155 
156  template<int dim, typename VectorType, typename DoFHandlerType>
157  void
159  {
160  Assert(input_vectors.size()==all_out.size(),
161  ExcDimensionMismatch(input_vectors.size(), all_out.size()) );
162 
163 //TODO: casting away constness is bad
167  (&dof_handler->get_triangulation())));
168  Assert (tria != 0, ExcInternalError());
169 
170  tria->notify_ready_to_unpack(offset,
171  std_cxx11::bind(&SolutionTransfer<dim, VectorType,
172  DoFHandlerType>::unpack_callback,
173  this,
174  std_cxx11::_1,
175  std_cxx11::_2,
176  std_cxx11::_3,
177  std_cxx11::ref(all_out)));
178 
179 
180  for (typename std::vector<VectorType *>::iterator it=all_out.begin();
181  it !=all_out.end();
182  ++it)
183  (*it)->compress(::VectorOperation::insert);
184 
185  input_vectors.clear();
186  }
187 
188 
189 
190  template<int dim, typename VectorType, typename DoFHandlerType>
191  void
193  {
194  std::vector<VectorType *> all_out(1, &out);
195  interpolate(all_out);
196  }
197 
198 
199 
200  template<int dim, typename VectorType, typename DoFHandlerType>
201  unsigned int
203  {
204  return sizeof(typename VectorType::value_type)* DoFTools::max_dofs_per_cell(*dof_handler);
205  }
206 
207 
208  template<int dim, typename VectorType, typename DoFHandlerType>
209  void
212  const typename Triangulation<dim,DoFHandlerType::space_dimension>::CellStatus /*status*/,
213  void *data)
214  {
215  typename VectorType::value_type *data_store = reinterpret_cast<typename VectorType::value_type *>(data);
216 
217  typename DoFHandlerType::cell_iterator cell(*cell_, dof_handler);
218 
219  const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell;
220  ::Vector<typename VectorType::value_type> dofvalues(dofs_per_cell);
221  for (typename std::vector<const VectorType *>::iterator it=input_vectors.begin();
222  it !=input_vectors.end();
223  ++it)
224  {
225  cell->get_interpolated_dof_values(*(*it), dofvalues);
226  std::memcpy(data_store, &dofvalues(0), sizeof(typename VectorType::value_type)*dofs_per_cell);
227  data_store += dofs_per_cell;
228  }
229  }
230 
231 
232  template<int dim, typename VectorType, typename DoFHandlerType>
233  void
236  const typename Triangulation<dim,DoFHandlerType::space_dimension>::CellStatus /*status*/,
237  const void *data,
238  std::vector<VectorType *> &all_out)
239  {
240  typename DoFHandlerType::cell_iterator
241  cell(*cell_, dof_handler);
242 
243  const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell;
244  ::Vector<typename VectorType::value_type> dofvalues(dofs_per_cell);
245  const typename VectorType::value_type *data_store = reinterpret_cast<const typename VectorType::value_type *>(data);
246 
247  for (typename std::vector<VectorType *>::iterator it = all_out.begin();
248  it != all_out.end();
249  ++it)
250  {
251  std::memcpy(&dofvalues(0), data_store, sizeof(typename VectorType::value_type)*dofs_per_cell);
252  cell->set_dof_values_by_interpolation(dofvalues, *(*it));
253  data_store += dofs_per_cell;
254  }
255  }
256 
257 
258  }
259 }
260 
261 
262 // explicit instantiations
263 #include "solution_transfer.inst"
264 
265 DEAL_II_NAMESPACE_CLOSE
266 
267 #endif
void prepare_for_coarsening_and_refinement(const std::vector< const VectorType *> &all_in)
unsigned int register_data_attach(const std::size_t size, const std_cxx11::function< void(const cell_iterator &, const CellStatus, void *)> &pack_callback)
Definition: tria.cc:3277
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)
void notify_ready_to_unpack(const unsigned int offset, const std_cxx11::function< void(const cell_iterator &, const CellStatus, const void *)> &unpack_callback)
Definition: tria.cc:3300
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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()