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