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