Reference documentation for deal.II version 8.5.1
solution_transfer.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 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/base/memory_consumption.h>
17 #include <deal.II/grid/tria.h>
18 #include <deal.II/grid/tria_accessor.h>
19 #include <deal.II/grid/tria_iterator.h>
20 #include <deal.II/distributed/tria.h>
21 #include <deal.II/dofs/dof_handler.h>
22 #include <deal.II/dofs/dof_accessor.h>
23 #include <deal.II/fe/fe.h>
24 #include <deal.II/lac/vector.h>
25 #include <deal.II/lac/la_vector.h>
26 #include <deal.II/lac/la_parallel_vector.h>
27 #include <deal.II/lac/petsc_vector.h>
28 #include <deal.II/lac/trilinos_vector.h>
29 #include <deal.II/lac/block_vector.h>
30 #include <deal.II/lac/la_parallel_block_vector.h>
31 #include <deal.II/lac/petsc_block_vector.h>
32 #include <deal.II/lac/trilinos_block_vector.h>
33 #include <deal.II/numerics/solution_transfer.h>
34 
35 DEAL_II_NAMESPACE_OPEN
36 
37 
38 
39 template<int dim, typename VectorType, typename DoFHandlerType>
41 SolutionTransfer(const DoFHandlerType &dof)
42  :
43  dof_handler(&dof, typeid(*this).name()),
44  n_dofs_old(0),
45  prepared_for(none)
46 {
48  (&dof_handler->get_tria())
49  == 0),
50  ExcMessage ("You are calling the ::SolutionTransfer class "
51  "with a DoF handler that is built on a "
52  "parallel::distributed::Triangulation. This will not "
53  "work for parallel computations. You probably want to "
54  "use the parallel::distributed::SolutionTransfer class."));
55 }
56 
57 
58 
59 template<int dim, typename VectorType, typename DoFHandlerType>
61 {
62  clear ();
63 }
64 
65 
66 
67 template<int dim, typename VectorType, typename DoFHandlerType>
69 {
70  indices_on_cell.clear();
71  dof_values_on_cell.clear();
72  cell_map.clear();
73 
74  prepared_for=none;
75 }
76 
77 
78 
79 template<int dim, typename VectorType, typename DoFHandlerType>
81 {
82  Assert(prepared_for!=pure_refinement, ExcAlreadyPrepForRef());
83  Assert(prepared_for!=coarsening_and_refinement,
84  ExcAlreadyPrepForCoarseAndRef());
85 
86  clear();
87 
88  const unsigned int n_active_cells = dof_handler->get_triangulation().n_active_cells();
89  n_dofs_old=dof_handler->n_dofs();
90 
91  // efficient reallocation of indices_on_cell
92  std::vector<std::vector<types::global_dof_index> > (n_active_cells)
93  .swap(indices_on_cell);
94 
95  typename DoFHandlerType::active_cell_iterator cell = dof_handler->begin_active(),
96  endc = dof_handler->end();
97 
98  for (unsigned int i=0; cell!=endc; ++cell, ++i)
99  {
100  indices_on_cell[i].resize(cell->get_fe().dofs_per_cell);
101  // on each cell store the indices of the
102  // dofs. after refining we get the values
103  // on the children by taking these
104  // indices, getting the respective values
105  // out of the data vectors and prolonging
106  // them to the children
107  cell->get_dof_indices(indices_on_cell[i]);
108  cell_map[std::make_pair(cell->level(),cell->index())]
109  = Pointerstruct(&indices_on_cell[i], cell->active_fe_index());
110  }
111  prepared_for=pure_refinement;
112 }
113 
114 
115 
116 template<int dim, typename VectorType, typename DoFHandlerType>
117 void
119 (const VectorType &in,
120  VectorType &out) const
121 {
122  Assert(prepared_for==pure_refinement, ExcNotPrepared());
123  Assert(in.size()==n_dofs_old, ExcDimensionMismatch(in.size(),n_dofs_old));
124  Assert(out.size()==dof_handler->n_dofs(),
125  ExcDimensionMismatch(out.size(),dof_handler->n_dofs()));
126  Assert(&in != &out,
127  ExcMessage ("Vectors cannot be used as input and output"
128  " at the same time!"));
129 
131 
132  typename DoFHandlerType::cell_iterator cell = dof_handler->begin(),
133  endc = dof_handler->end();
134 
135  typename std::map<std::pair<unsigned int, unsigned int>, Pointerstruct>::const_iterator
136  pointerstruct,
137  cell_map_end=cell_map.end();
138 
139  for (; cell!=endc; ++cell)
140  {
141  pointerstruct=cell_map.find(std::make_pair(cell->level(),cell->index()));
142 
143  if (pointerstruct!=cell_map_end)
144  // this cell was refined or not
145  // touched at all, so we can get
146  // the new values by just setting
147  // or interpolating to the children,
148  // which is both done by one
149  // function
150  {
151  const unsigned int this_fe_index = pointerstruct->second.active_fe_index;
152  const unsigned int dofs_per_cell=cell->get_dof_handler().get_fe()[this_fe_index].dofs_per_cell;
153  local_values.reinit(dofs_per_cell, true);
154 
155  // make sure that the size of the stored indices is the same as
156  // dofs_per_cell. since we store the desired fe_index, we know
157  // what this size should be
158  Assert(dofs_per_cell==(*pointerstruct->second.indices_ptr).size(),
159  ExcInternalError());
160  for (unsigned int i=0; i<dofs_per_cell; ++i)
161  local_values(i)=in((*pointerstruct->second.indices_ptr)[i]);
162  cell->set_dof_values_by_interpolation(local_values, out,
163  this_fe_index);
164  }
165  }
166 }
167 
168 
169 
170 namespace internal
171 {
185  template <typename DoFHandlerType>
186  void extract_interpolation_matrices (const DoFHandlerType &,
187  ::Table<2,FullMatrix<double> > &)
188  {}
189 
190  template <int dim, int spacedim>
191  void extract_interpolation_matrices (const ::hp::DoFHandler<dim,spacedim> &dof,
192  ::Table<2,FullMatrix<double> > &matrices)
193  {
194  const ::hp::FECollection<dim,spacedim> &fe = dof.get_fe();
195  matrices.reinit (fe.size(), fe.size());
196  for (unsigned int i=0; i<fe.size(); ++i)
197  for (unsigned int j=0; j<fe.size(); ++j)
198  if (i != j)
199  {
200  matrices(i,j).reinit (fe[i].dofs_per_cell, fe[j].dofs_per_cell);
201 
202  // see if we can get the interpolation matrices for this
203  // combination of elements. if not, reset the matrix sizes to zero
204  // to indicate that this particular combination isn't
205  // supported. this isn't an outright error right away since we may
206  // never need to actually interpolate between these two elements
207  // on actual cells; we simply have to trigger an error if someone
208  // actually tries
209  try
210  {
211  fe[i].get_interpolation_matrix (fe[j], matrices(i,j));
212  }
214  {
215  matrices(i,j).reinit (0,0);
216  }
217  }
218  }
219 
220 
221  template <int dim, int spacedim>
222  void restriction_additive (const FiniteElement<dim,spacedim> &,
223  std::vector<std::vector<bool> > &)
224  {}
225 
226  template <int dim, int spacedim>
227  void restriction_additive (const ::hp::FECollection<dim,spacedim> &fe,
228  std::vector<std::vector<bool> > &restriction_is_additive)
229  {
230  restriction_is_additive.resize (fe.size());
231  for (unsigned int f=0; f<fe.size(); ++f)
232  {
233  restriction_is_additive[f].resize (fe[f].dofs_per_cell);
234  for (unsigned int i=0; i<fe[f].dofs_per_cell; ++i)
235  restriction_is_additive[f][i] = fe[f].restriction_is_additive(i);
236  }
237  }
238 }
239 
240 
241 
242 template<int dim, typename VectorType, typename DoFHandlerType>
243 void
245 prepare_for_coarsening_and_refinement(const std::vector<VectorType> &all_in)
246 {
247  Assert (prepared_for!=pure_refinement, ExcAlreadyPrepForRef());
248  Assert (prepared_for!=coarsening_and_refinement,
249  ExcAlreadyPrepForCoarseAndRef());
250 
251  const unsigned int in_size=all_in.size();
252  Assert(in_size!=0,
253  ExcMessage("The array of input vectors you pass to this "
254  "function has no elements. This is not useful."));
255 
256  clear();
257 
258  const unsigned int n_active_cells = dof_handler->get_triangulation().n_active_cells();
259  (void)n_active_cells;
260  n_dofs_old = dof_handler->n_dofs();
261 
262  for (unsigned int i=0; i<in_size; ++i)
263  {
264  Assert(all_in[i].size()==n_dofs_old,
265  ExcDimensionMismatch(all_in[i].size(),n_dofs_old));
266  }
267 
268  // first count the number
269  // of cells that will be coarsened
270  // and that'll stay or be refined
271  unsigned int n_cells_to_coarsen=0;
272  unsigned int n_cells_to_stay_or_refine=0;
273  for (typename DoFHandlerType::active_cell_iterator act_cell = dof_handler->begin_active();
274  act_cell!=dof_handler->end(); ++act_cell)
275  {
276  if (act_cell->coarsen_flag_set())
277  ++n_cells_to_coarsen;
278  else
279  ++n_cells_to_stay_or_refine;
280  }
281  Assert((n_cells_to_coarsen+n_cells_to_stay_or_refine)==n_active_cells,
282  ExcInternalError());
283 
284  unsigned int n_coarsen_fathers=0;
285  for (typename DoFHandlerType::cell_iterator cell=dof_handler->begin();
286  cell!=dof_handler->end(); ++cell)
287  if (!cell->active() && cell->child(0)->coarsen_flag_set())
288  ++n_coarsen_fathers;
289  Assert(n_cells_to_coarsen>=2*n_coarsen_fathers, ExcInternalError());
290 
291  // allocate the needed memory. initialize
292  // the following arrays in an efficient
293  // way, without copying much
294  std::vector<std::vector<types::global_dof_index> >(n_cells_to_stay_or_refine)
295  .swap(indices_on_cell);
296 
297  std::vector<std::vector<Vector<typename VectorType::value_type> > >
298  (n_coarsen_fathers,
299  std::vector<Vector<typename VectorType::value_type> > (in_size))
300  .swap(dof_values_on_cell);
301 
302  Table<2,FullMatrix<double> > interpolation_hp;
303  std::vector<std::vector<bool> > restriction_is_additive;
304 
305  internal::extract_interpolation_matrices (*dof_handler, interpolation_hp);
306  internal::restriction_additive (dof_handler->get_fe(), restriction_is_additive);
307 
308  // we need counters for
309  // the 'to_stay_or_refine' cells 'n_sr' and
310  // the 'coarsen_fathers' cells 'n_cf',
311  unsigned int n_sr=0, n_cf=0;
312  for (typename DoFHandlerType::cell_iterator cell=dof_handler->begin();
313  cell!=dof_handler->end(); ++cell)
314  {
315  // CASE 1: active cell that remains as it is
316  if (cell->active() && !cell->coarsen_flag_set())
317  {
318  const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell;
319  indices_on_cell[n_sr].resize(dofs_per_cell);
320  // cell will not be coarsened,
321  // so we get away by storing the
322  // dof indices and later
323  // interpolating to the children
324  cell->get_dof_indices(indices_on_cell[n_sr]);
325  cell_map[std::make_pair(cell->level(), cell->index())]
326  = Pointerstruct(&indices_on_cell[n_sr], cell->active_fe_index());
327  ++n_sr;
328  }
329 
330  // CASE 2: cell is inactive but will become active
331  else if (cell->has_children() && cell->child(0)->coarsen_flag_set())
332  {
333  // note that if one child has the
334  // coarsen flag, then all should
335  // have if Tria::prepare_* has
336  // worked correctly
337  for (unsigned int i=1; i<cell->n_children(); ++i)
338  Assert(cell->child(i)->coarsen_flag_set(),
339  ExcMessage("It looks like you didn't call "
340  "Triangulation::prepare_coarsening_and_refinement before "
341  "calling the current function. This can't work."));
342 
343  // we will need to interpolate from the children of this cell
344  // to the current one. in the hp context, this also means
345  // we need to figure out which finite element space to interpolate
346  // to since that is not implied by the global FE as in the non-hp
347  // case.
348  bool different_fe_on_children = false;
349  for (unsigned int child=1; child<cell->n_children(); ++child)
350  if (cell->child(child)->active_fe_index()
351  != cell->child(0)->active_fe_index())
352  {
353  different_fe_on_children = true;
354  break;
355  }
356 
357  // take FE index from the child with most
358  // degrees of freedom locally
359  unsigned int most_general_child = 0;
360  if (different_fe_on_children == true)
361  for (unsigned int child=1; child<cell->n_children(); ++child)
362  if (cell->child(child)->get_fe().dofs_per_cell >
363  cell->child(most_general_child)->get_fe().dofs_per_cell)
364  most_general_child = child;
365  const unsigned int target_fe_index = cell->child(most_general_child)->active_fe_index();
366 
367  const unsigned int dofs_per_cell=cell->get_dof_handler().get_fe()[target_fe_index].dofs_per_cell;
368 
369  std::vector<Vector<typename VectorType::value_type> >(in_size,
371  .swap(dof_values_on_cell[n_cf]);
372 
373 
374  // store the data of each of the input vectors. get this data
375  // as interpolated onto a finite element space that encompasses
376  // that of all the children. note that cell->get_interpolated_dof_values
377  // already does all of the interpolations between spaces
378  for (unsigned int j=0; j<in_size; ++j)
379  cell->get_interpolated_dof_values(all_in[j],
380  dof_values_on_cell[n_cf][j],
381  target_fe_index);
382  cell_map[std::make_pair(cell->level(), cell->index())]
383  = Pointerstruct(&dof_values_on_cell[n_cf], target_fe_index);
384  ++n_cf;
385  }
386  }
387  Assert(n_sr==n_cells_to_stay_or_refine, ExcInternalError());
388  Assert(n_cf==n_coarsen_fathers, ExcInternalError());
389 
390  prepared_for=coarsening_and_refinement;
391 }
392 
393 
394 
395 template<int dim, typename VectorType, typename DoFHandlerType>
396 void
398 (const VectorType &in)
399 {
400  std::vector<VectorType> all_in=std::vector<VectorType>(1, in);
401  prepare_for_coarsening_and_refinement(all_in);
402 }
403 
404 
405 
406 template<int dim, typename VectorType, typename DoFHandlerType>
408 interpolate (const std::vector<VectorType> &all_in,
409  std::vector<VectorType> &all_out) const
410 {
411  Assert(prepared_for==coarsening_and_refinement, ExcNotPrepared());
412  const unsigned int size=all_in.size();
413  Assert(all_out.size()==size, ExcDimensionMismatch(all_out.size(), size));
414  for (unsigned int i=0; i<size; ++i)
415  Assert (all_in[i].size() == n_dofs_old,
416  ExcDimensionMismatch(all_in[i].size(), n_dofs_old));
417  for (unsigned int i=0; i<all_out.size(); ++i)
418  Assert (all_out[i].size() == dof_handler->n_dofs(),
419  ExcDimensionMismatch(all_out[i].size(), dof_handler->n_dofs()));
420  for (unsigned int i=0; i<size; ++i)
421  for (unsigned int j=0; j<size; ++j)
422  Assert(&all_in[i] != &all_out[j],
423  ExcMessage ("Vectors cannot be used as input and output"
424  " at the same time!"));
425 
427  std::vector<types::global_dof_index> dofs;
428 
429  typename std::map<std::pair<unsigned int, unsigned int>, Pointerstruct>::const_iterator
430  pointerstruct,
431  cell_map_end=cell_map.end();
432 
433  Table<2,FullMatrix<double> > interpolation_hp;
434  internal::extract_interpolation_matrices (*dof_handler, interpolation_hp);
436 
437  typename DoFHandlerType::cell_iterator cell = dof_handler->begin(),
438  endc = dof_handler->end();
439  for (; cell!=endc; ++cell)
440  {
441  pointerstruct=cell_map.find(std::make_pair(cell->level(),cell->index()));
442 
443  if (pointerstruct!=cell_map_end)
444  {
445  const std::vector<types::global_dof_index> *const indexptr
446  =pointerstruct->second.indices_ptr;
447 
448  const std::vector<Vector<typename VectorType::value_type> > *const valuesptr
449  =pointerstruct->second.dof_values_ptr;
450 
451  // cell stayed as it was or was refined
452  if (indexptr)
453  {
454  Assert (valuesptr == 0,
455  ExcInternalError());
456 
457  const unsigned int old_fe_index =
458  pointerstruct->second.active_fe_index;
459 
460  // get the values of
461  // each of the input
462  // data vectors on this
463  // cell and prolong it
464  // to its children
465  unsigned int in_size = indexptr->size();
466  for (unsigned int j=0; j<size; ++j)
467  {
468  tmp.reinit (in_size, true);
469  for (unsigned int i=0; i<in_size; ++i)
470  tmp(i) = all_in[j]((*indexptr)[i]);
471 
472  cell->set_dof_values_by_interpolation (tmp, all_out[j],
473  old_fe_index);
474  }
475  }
476  else if (valuesptr)
477  // the children of this cell were
478  // deleted
479  {
480  Assert (!cell->has_children(), ExcInternalError());
481  Assert (indexptr == 0,
482  ExcInternalError());
483 
484  const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
485  dofs.resize(dofs_per_cell);
486  // get the local
487  // indices
488  cell->get_dof_indices(dofs);
489 
490  // distribute the
491  // stored data to the
492  // new vectors
493  for (unsigned int j=0; j<size; ++j)
494  {
495  // make sure that the size of
496  // the stored indices is the
497  // same as
498  // dofs_per_cell. this is
499  // kind of a test if we use
500  // the same fe in the hp
501  // case. to really do that
502  // test we would have to
503  // store the fe_index of all
504  // cells
506  const unsigned int active_fe_index = cell->active_fe_index();
507  if (active_fe_index != pointerstruct->second.active_fe_index)
508  {
509  const unsigned int old_index = pointerstruct->second.active_fe_index;
510  tmp.reinit (dofs_per_cell, true);
511  AssertDimension ((*valuesptr)[j].size(),
512  interpolation_hp(active_fe_index,old_index).n());
513  AssertDimension (tmp.size(),
514  interpolation_hp(active_fe_index,old_index).m());
515  interpolation_hp(active_fe_index,old_index).vmult (tmp, (*valuesptr)[j]);
516  data = &tmp;
517  }
518  else
519  data = &(*valuesptr)[j];
520 
521 
522  for (unsigned int i=0; i<dofs_per_cell; ++i)
523  all_out[j](dofs[i])=(*data)(i);
524  }
525  }
526  // undefined status
527  else
528  Assert(false, ExcInternalError());
529  }
530  }
531 }
532 
533 
534 
535 template<int dim, typename VectorType, typename DoFHandlerType>
537 (const VectorType &in,
538  VectorType &out) const
539 {
540  Assert (in.size()==n_dofs_old,
541  ExcDimensionMismatch(in.size(), n_dofs_old));
542  Assert (out.size()==dof_handler->n_dofs(),
543  ExcDimensionMismatch(out.size(), dof_handler->n_dofs()));
544 
545  std::vector<VectorType> all_in(1);
546  all_in[0] = in;
547  std::vector<VectorType> all_out(1);
548  all_out[0] = out;
549  interpolate(all_in,
550  all_out);
551  out=all_out[0];
552 }
553 
554 
555 
556 template<int dim, typename VectorType, typename DoFHandlerType>
557 std::size_t
559 {
560  // at the moment we do not include the memory
561  // consumption of the cell_map as we have no
562  // real idea about memory consumption of a
563  // std::map
564  return (MemoryConsumption::memory_consumption (dof_handler) +
566  sizeof (prepared_for) +
567  MemoryConsumption::memory_consumption (indices_on_cell) +
568  MemoryConsumption::memory_consumption (dof_values_on_cell));
569 }
570 
571 
572 
573 template<int dim, typename VectorType, typename DoFHandlerType>
574 std::size_t
576 {
577  return sizeof(*this);
578 }
579 
580 
581 /*-------------- Explicit Instantiations -------------------------------*/
582 #define SPLIT_INSTANTIATIONS_COUNT 4
583 #ifndef SPLIT_INSTANTIATIONS_INDEX
584 #define SPLIT_INSTANTIATIONS_INDEX 0
585 #endif
586 #include "solution_transfer.inst"
587 
588 DEAL_II_NAMESPACE_CLOSE
void extract_interpolation_matrices(const DoFHandlerType &, ::Table< 2, FullMatrix< double > > &)
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
static ::ExceptionBase & ExcMessage(std::string arg1)
void swap(BlockIndices &u, BlockIndices &v)
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void prepare_for_pure_refinement()
std::size_t size() const
iterator begin()
SolutionTransfer(const DoFHandlerType &dof)
SmartPointer< const DoFHandlerType, SolutionTransfer< dim, VectorType, DoFHandlerType > > dof_handler
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
void prepare_for_coarsening_and_refinement(const std::vector< VectorType > &all_in)
void refine_interpolate(const VectorType &in, VectorType &out) const
virtual void reinit(const size_type N, const bool omit_zeroing_entries=false)
Definition: table.h:33
std::size_t memory_consumption() const
void interpolate(const std::vector< VectorType > &all_in, std::vector< VectorType > &all_out) const
static ::ExceptionBase & ExcInternalError()