16 #include <deal.II/base/memory_consumption.h> 18 #include <deal.II/distributed/tria.h> 20 #include <deal.II/dofs/dof_accessor.h> 21 #include <deal.II/dofs/dof_handler.h> 23 #include <deal.II/fe/fe.h> 25 #include <deal.II/grid/tria.h> 26 #include <deal.II/grid/tria_accessor.h> 27 #include <deal.II/grid/tria_iterator.h> 29 #include <deal.II/lac/block_vector.h> 30 #include <deal.II/lac/la_parallel_block_vector.h> 31 #include <deal.II/lac/la_parallel_vector.h> 32 #include <deal.II/lac/la_vector.h> 33 #include <deal.II/lac/petsc_block_vector.h> 34 #include <deal.II/lac/petsc_vector.h> 35 #include <deal.II/lac/trilinos_parallel_block_vector.h> 36 #include <deal.II/lac/trilinos_vector.h> 37 #include <deal.II/lac/vector.h> 38 #include <deal.II/lac/vector_element_access.h> 40 #include <deal.II/numerics/solution_transfer.h> 42 DEAL_II_NAMESPACE_OPEN
44 template <
int dim,
typename VectorType,
typename DoFHandlerType>
46 const DoFHandlerType &dof)
47 : dof_handler(&dof, typeid(*this).name())
52 DoFHandlerType::dimension,
53 DoFHandlerType::space_dimension
> *>(
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."));
64 template <
int dim,
typename VectorType,
typename DoFHandlerType>
72 template <
int dim,
typename VectorType,
typename DoFHandlerType>
76 indices_on_cell.clear();
77 dof_values_on_cell.clear();
85 template <
int dim,
typename VectorType,
typename DoFHandlerType>
89 Assert(prepared_for != pure_refinement, ExcAlreadyPrepForRef());
90 Assert(prepared_for != coarsening_and_refinement,
91 ExcAlreadyPrepForCoarseAndRef());
95 const unsigned int n_active_cells =
96 dof_handler->get_triangulation().n_active_cells();
97 n_dofs_old = dof_handler->n_dofs();
100 std::vector<std::vector<types::global_dof_index>>(n_active_cells)
101 .
swap(indices_on_cell);
103 typename DoFHandlerType::active_cell_iterator cell =
104 dof_handler->begin_active(),
105 endc = dof_handler->end();
107 for (
unsigned int i = 0; cell != endc; ++cell, ++i)
109 indices_on_cell[i].resize(cell->get_fe().dofs_per_cell);
116 cell->get_dof_indices(indices_on_cell[i]);
117 cell_map[std::make_pair(cell->level(), cell->index())] =
120 prepared_for = pure_refinement;
125 template <
int dim,
typename VectorType,
typename DoFHandlerType>
128 const VectorType &in,
129 VectorType & out)
const 131 Assert(prepared_for == pure_refinement, ExcNotPrepared());
133 Assert(out.size() == dof_handler->n_dofs(),
136 ExcMessage(
"Vectors cannot be used as input and output" 137 " at the same time!"));
141 typename DoFHandlerType::cell_iterator cell = dof_handler->begin(),
142 endc = dof_handler->end();
144 typename std::map<std::pair<unsigned int, unsigned int>,
146 cell_map_end = cell_map.end();
148 for (; cell != endc; ++cell)
151 cell_map.find(std::make_pair(cell->level(), cell->index()));
153 if (pointerstruct != cell_map_end)
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);
170 Assert(dofs_per_cell == (*pointerstruct->second.indices_ptr).size(),
172 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
173 local_values(i) = internal::ElementAccess<VectorType>::get(
174 in, (*pointerstruct->second.indices_ptr)[i]);
175 cell->set_dof_values_by_interpolation(local_values,
199 template <
typename DoFHandlerType>
205 template <
int dim,
int spacedim>
208 const ::hp::DoFHandler<dim, spacedim> &dof,
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)
217 matrices(i, j).reinit(fe[i].dofs_per_cell, fe[j].dofs_per_cell);
228 fe[i].get_interpolation_matrix(fe[j], matrices(i, j));
231 ExcInterpolationNotImplemented &)
233 matrices(i, j).reinit(0, 0);
239 template <
int dim,
int spacedim>
242 std::vector<std::vector<bool>> &)
245 template <
int dim,
int spacedim>
247 restriction_additive(const ::hp::FECollection<dim, spacedim> &fe,
248 std::vector<std::vector<bool>> &restriction_is_additive)
250 restriction_is_additive.resize(fe.size());
251 for (
unsigned int f = 0; f < fe.size(); ++f)
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);
262 template <
int dim,
typename VectorType,
typename DoFHandlerType>
267 Assert(prepared_for != pure_refinement, ExcAlreadyPrepForRef());
268 Assert(prepared_for != coarsening_and_refinement,
269 ExcAlreadyPrepForCoarseAndRef());
271 const unsigned int in_size = all_in.size();
273 ExcMessage(
"The array of input vectors you pass to this " 274 "function has no elements. This is not useful."));
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();
283 for (
unsigned int i = 0; i < in_size; ++i)
285 Assert(all_in[i].size() == n_dofs_old,
292 unsigned int n_cells_to_coarsen = 0;
293 unsigned int n_cells_to_stay_or_refine = 0;
294 for (
typename DoFHandlerType::active_cell_iterator act_cell =
295 dof_handler->begin_active();
296 act_cell != dof_handler->end();
299 if (act_cell->coarsen_flag_set())
300 ++n_cells_to_coarsen;
302 ++n_cells_to_stay_or_refine;
304 Assert((n_cells_to_coarsen + n_cells_to_stay_or_refine) == n_active_cells,
307 unsigned int n_coarsen_fathers = 0;
308 for (
typename DoFHandlerType::cell_iterator cell = dof_handler->begin();
309 cell != dof_handler->end();
311 if (!cell->active() && cell->child(0)->coarsen_flag_set())
318 std::vector<std::vector<types::global_dof_index>>(n_cells_to_stay_or_refine)
319 .
swap(indices_on_cell);
321 std::vector<std::vector<Vector<typename VectorType::value_type>>>(
323 std::vector<Vector<typename VectorType::value_type>>(in_size))
324 .swap(dof_values_on_cell);
327 std::vector<std::vector<bool>> restriction_is_additive;
330 internal::restriction_additive(dof_handler->get_fe_collection(),
331 restriction_is_additive);
336 unsigned int n_sr = 0, n_cf = 0;
337 for (
typename DoFHandlerType::cell_iterator cell = dof_handler->begin();
338 cell != dof_handler->end();
342 if (cell->active() && !cell->coarsen_flag_set())
344 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
345 indices_on_cell[n_sr].resize(dofs_per_cell);
350 cell->get_dof_indices(indices_on_cell[n_sr]);
351 cell_map[std::make_pair(cell->level(), cell->index())] =
352 Pointerstruct(&indices_on_cell[n_sr], cell->active_fe_index());
357 else if (cell->has_children() && cell->child(0)->coarsen_flag_set())
363 for (
unsigned int i = 1; i < cell->n_children(); ++i)
364 Assert(cell->child(i)->coarsen_flag_set(),
366 "It looks like you didn't call " 367 "Triangulation::prepare_coarsening_and_refinement before " 368 "calling the current function. This can't work."));
375 bool different_fe_on_children =
false;
376 for (
unsigned int child = 1; child < cell->n_children(); ++child)
377 if (cell->child(child)->active_fe_index() !=
378 cell->child(0)->active_fe_index())
380 different_fe_on_children =
true;
386 unsigned int most_general_child = 0;
387 if (different_fe_on_children ==
true)
388 for (
unsigned int child = 1; child < cell->n_children(); ++child)
389 if (cell->child(child)->get_fe().dofs_per_cell >
390 cell->child(most_general_child)->get_fe().dofs_per_cell)
391 most_general_child = child;
392 const unsigned int target_fe_index =
393 cell->child(most_general_child)->active_fe_index();
395 const unsigned int dofs_per_cell =
396 dof_handler->get_fe(target_fe_index).dofs_per_cell;
398 std::vector<Vector<typename VectorType::value_type>>(
400 .swap(dof_values_on_cell[n_cf]);
408 for (
unsigned int j = 0; j < in_size; ++j)
409 cell->get_interpolated_dof_values(all_in[j],
410 dof_values_on_cell[n_cf][j],
412 cell_map[std::make_pair(cell->level(), cell->index())] =
420 prepared_for = coarsening_and_refinement;
425 template <
int dim,
typename VectorType,
typename DoFHandlerType>
430 std::vector<VectorType> all_in(1, in);
431 prepare_for_coarsening_and_refinement(all_in);
436 template <
int dim,
typename VectorType,
typename DoFHandlerType>
439 const std::vector<VectorType> &all_in,
440 std::vector<VectorType> & all_out)
const 442 Assert(prepared_for == coarsening_and_refinement, ExcNotPrepared());
443 const unsigned int size = all_in.size();
445 for (
unsigned int i = 0; i < size; ++i)
446 Assert(all_in[i].size() == n_dofs_old,
448 for (
unsigned int i = 0; i < all_out.size(); ++i)
449 Assert(all_out[i].size() == dof_handler->n_dofs(),
451 for (
unsigned int i = 0; i < size; ++i)
452 for (
unsigned int j = 0; j < size; ++j)
453 Assert(&all_in[i] != &all_out[j],
454 ExcMessage(
"Vectors cannot be used as input and output" 455 " at the same time!"));
458 std::vector<types::global_dof_index> dofs;
460 typename std::map<std::pair<unsigned int, unsigned int>,
462 cell_map_end = cell_map.end();
468 typename DoFHandlerType::cell_iterator cell = dof_handler->
begin(),
469 endc = dof_handler->end();
470 for (; cell != endc; ++cell)
473 cell_map.find(std::make_pair(cell->level(), cell->index()));
475 if (pointerstruct != cell_map_end)
477 const std::vector<types::global_dof_index> *
const indexptr =
478 pointerstruct->second.indices_ptr;
480 const std::vector<Vector<typename VectorType::value_type>>
481 *
const valuesptr = pointerstruct->second.dof_values_ptr;
488 const unsigned int old_fe_index =
489 pointerstruct->second.active_fe_index;
496 unsigned int in_size = indexptr->size();
497 for (
unsigned int j = 0; j < size; ++j)
499 tmp.
reinit(in_size,
true);
500 for (
unsigned int i = 0; i < in_size; ++i)
502 internal::ElementAccess<VectorType>::get(all_in[j],
505 cell->set_dof_values_by_interpolation(tmp,
517 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
518 dofs.resize(dofs_per_cell);
521 cell->get_dof_indices(dofs);
526 for (
unsigned int j = 0; j < size; ++j)
539 const unsigned int active_fe_index = cell->active_fe_index();
540 if (active_fe_index != pointerstruct->second.active_fe_index)
542 const unsigned int old_index =
543 pointerstruct->second.active_fe_index;
544 tmp.
reinit(dofs_per_cell,
true);
546 (*valuesptr)[j].size(),
547 interpolation_hp(active_fe_index, old_index).n());
550 interpolation_hp(active_fe_index, old_index).m());
551 interpolation_hp(active_fe_index, old_index)
552 .vmult(tmp, (*valuesptr)[j]);
556 data = &(*valuesptr)[j];
559 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
574 template <
int dim,
typename VectorType,
typename DoFHandlerType>
577 const VectorType &in,
578 VectorType & out)
const 581 Assert(out.size() == dof_handler->n_dofs(),
584 std::vector<VectorType> all_in(1);
586 std::vector<VectorType> all_out(1);
588 interpolate(all_in, all_out);
594 template <
int dim,
typename VectorType,
typename DoFHandlerType>
604 sizeof(prepared_for) +
611 template <
int dim,
typename VectorType,
typename DoFHandlerType>
616 return sizeof(*this);
621 #define SPLIT_INSTANTIATIONS_COUNT 4 622 #ifndef SPLIT_INSTANTIATIONS_INDEX 623 # define SPLIT_INSTANTIATIONS_INDEX 0 625 #include "solution_transfer.inst" 627 DEAL_II_NAMESPACE_CLOSE
#define AssertDimension(dim1, dim2)
SmartPointer< const DoFHandlerType, SolutionTransfer< dim, VectorType, DoFHandlerType > > dof_handler
static ::ExceptionBase & ExcMessage(std::string arg1)
void swap(BlockIndices &u, BlockIndices &v)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void prepare_for_pure_refinement()
void extract_interpolation_matrices(const DoFHandlerType &, ::Table< 2, FullMatrix< double >> &)
SolutionTransfer(const DoFHandlerType &dof)
__global__ void set(Number *val, const Number s, const size_type N)
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)
std::size_t memory_consumption() const
void interpolate(const std::vector< VectorType > &all_in, std::vector< VectorType > &all_out) const
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
static ::ExceptionBase & ExcInternalError()