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> 35 DEAL_II_NAMESPACE_OPEN
39 template<
int dim,
typename VectorType,
typename DoFHandlerType>
43 dof_handler(&dof, typeid(*this).name()),
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."));
59 template<
int dim,
typename VectorType,
typename DoFHandlerType>
67 template<
int dim,
typename VectorType,
typename DoFHandlerType>
70 indices_on_cell.clear();
71 dof_values_on_cell.clear();
79 template<
int dim,
typename VectorType,
typename DoFHandlerType>
82 Assert(prepared_for!=pure_refinement, ExcAlreadyPrepForRef());
83 Assert(prepared_for!=coarsening_and_refinement,
84 ExcAlreadyPrepForCoarseAndRef());
88 const unsigned int n_active_cells = dof_handler->get_triangulation().n_active_cells();
89 n_dofs_old=dof_handler->n_dofs();
92 std::vector<std::vector<types::global_dof_index> > (n_active_cells)
93 .
swap(indices_on_cell);
95 typename DoFHandlerType::active_cell_iterator cell = dof_handler->begin_active(),
96 endc = dof_handler->end();
98 for (
unsigned int i=0; cell!=endc; ++cell, ++i)
100 indices_on_cell[i].resize(cell->get_fe().dofs_per_cell);
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());
111 prepared_for=pure_refinement;
116 template<
int dim,
typename VectorType,
typename DoFHandlerType>
119 (
const VectorType &in,
120 VectorType &out)
const 122 Assert(prepared_for==pure_refinement, ExcNotPrepared());
124 Assert(out.size()==dof_handler->n_dofs(),
127 ExcMessage (
"Vectors cannot be used as input and output" 128 " at the same time!"));
132 typename DoFHandlerType::cell_iterator cell = dof_handler->begin(),
133 endc = dof_handler->end();
135 typename std::map<std::pair<unsigned int, unsigned int>,
Pointerstruct>::const_iterator
137 cell_map_end=cell_map.end();
139 for (; cell!=endc; ++cell)
141 pointerstruct=cell_map.find(std::make_pair(cell->level(),cell->index()));
143 if (pointerstruct!=cell_map_end)
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);
158 Assert(dofs_per_cell==(*pointerstruct->second.indices_ptr).size(),
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,
185 template <
typename DoFHandlerType>
190 template <
int dim,
int spacedim>
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)
200 matrices(i,j).reinit (fe[i].dofs_per_cell, fe[j].dofs_per_cell);
211 fe[i].get_interpolation_matrix (fe[j], matrices(i,j));
215 matrices(i,j).reinit (0,0);
221 template <
int dim,
int spacedim>
223 std::vector<std::vector<bool> > &)
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)
230 restriction_is_additive.resize (fe.size());
231 for (
unsigned int f=0; f<fe.size(); ++f)
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);
242 template<
int dim,
typename VectorType,
typename DoFHandlerType>
247 Assert (prepared_for!=pure_refinement, ExcAlreadyPrepForRef());
248 Assert (prepared_for!=coarsening_and_refinement,
249 ExcAlreadyPrepForCoarseAndRef());
251 const unsigned int in_size=all_in.size();
253 ExcMessage(
"The array of input vectors you pass to this " 254 "function has no elements. This is not useful."));
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();
262 for (
unsigned int i=0; i<in_size; ++i)
264 Assert(all_in[i].size()==n_dofs_old,
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)
276 if (act_cell->coarsen_flag_set())
277 ++n_cells_to_coarsen;
279 ++n_cells_to_stay_or_refine;
281 Assert((n_cells_to_coarsen+n_cells_to_stay_or_refine)==n_active_cells,
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())
294 std::vector<std::vector<types::global_dof_index> >(n_cells_to_stay_or_refine)
295 .
swap(indices_on_cell);
297 std::vector<std::vector<Vector<typename VectorType::value_type> > >
299 std::vector<Vector<typename VectorType::value_type> > (in_size))
300 .swap(dof_values_on_cell);
303 std::vector<std::vector<bool> > restriction_is_additive;
306 internal::restriction_additive (dof_handler->get_fe(), restriction_is_additive);
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)
316 if (cell->active() && !cell->coarsen_flag_set())
318 const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell;
319 indices_on_cell[n_sr].resize(dofs_per_cell);
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());
331 else if (cell->has_children() && cell->child(0)->coarsen_flag_set())
337 for (
unsigned int i=1; i<cell->n_children(); ++i)
338 Assert(cell->child(i)->coarsen_flag_set(),
340 "Triangulation::prepare_coarsening_and_refinement before " 341 "calling the current function. This can't work."));
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())
353 different_fe_on_children =
true;
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();
367 const unsigned int dofs_per_cell=cell->get_dof_handler().get_fe()[target_fe_index].dofs_per_cell;
369 std::vector<Vector<typename VectorType::value_type> >(in_size,
371 .swap(dof_values_on_cell[n_cf]);
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],
382 cell_map[std::make_pair(cell->level(), cell->index())]
390 prepared_for=coarsening_and_refinement;
395 template<
int dim,
typename VectorType,
typename DoFHandlerType>
398 (
const VectorType &in)
400 std::vector<VectorType> all_in=std::vector<VectorType>(1, in);
401 prepare_for_coarsening_and_refinement(all_in);
406 template<
int dim,
typename VectorType,
typename DoFHandlerType>
409 std::vector<VectorType> &all_out)
const 411 Assert(prepared_for==coarsening_and_refinement, ExcNotPrepared());
412 const unsigned int size=all_in.size();
414 for (
unsigned int i=0; i<size; ++i)
415 Assert (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(),
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!"));
427 std::vector<types::global_dof_index> dofs;
429 typename std::map<std::pair<unsigned int, unsigned int>,
Pointerstruct>::const_iterator
431 cell_map_end=cell_map.end();
437 typename DoFHandlerType::cell_iterator cell = dof_handler->
begin(),
438 endc = dof_handler->end();
439 for (; cell!=endc; ++cell)
441 pointerstruct=cell_map.find(std::make_pair(cell->level(),cell->index()));
443 if (pointerstruct!=cell_map_end)
445 const std::vector<types::global_dof_index> *
const indexptr
446 =pointerstruct->second.indices_ptr;
448 const std::vector<Vector<typename VectorType::value_type> > *
const valuesptr
449 =pointerstruct->second.dof_values_ptr;
457 const unsigned int old_fe_index =
458 pointerstruct->second.active_fe_index;
465 unsigned int in_size = indexptr->size();
466 for (
unsigned int j=0; j<size; ++j)
468 tmp.
reinit (in_size,
true);
469 for (
unsigned int i=0; i<in_size; ++i)
470 tmp(i) = all_in[j]((*indexptr)[i]);
472 cell->set_dof_values_by_interpolation (tmp, all_out[j],
484 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
485 dofs.resize(dofs_per_cell);
488 cell->get_dof_indices(dofs);
493 for (
unsigned int j=0; j<size; ++j)
506 const unsigned int active_fe_index = cell->active_fe_index();
507 if (active_fe_index != pointerstruct->second.active_fe_index)
509 const unsigned int old_index = pointerstruct->second.active_fe_index;
510 tmp.
reinit (dofs_per_cell,
true);
512 interpolation_hp(active_fe_index,old_index).n());
514 interpolation_hp(active_fe_index,old_index).m());
515 interpolation_hp(active_fe_index,old_index).vmult (tmp, (*valuesptr)[j]);
519 data = &(*valuesptr)[j];
522 for (
unsigned int i=0; i<dofs_per_cell; ++i)
523 all_out[j](dofs[i])=(*data)(i);
535 template<
int dim,
typename VectorType,
typename DoFHandlerType>
537 (
const VectorType &in,
538 VectorType &out)
const 540 Assert (in.size()==n_dofs_old,
542 Assert (out.size()==dof_handler->n_dofs(),
545 std::vector<VectorType> all_in(1);
547 std::vector<VectorType> all_out(1);
556 template<
int dim,
typename VectorType,
typename DoFHandlerType>
566 sizeof (prepared_for) +
573 template<
int dim,
typename VectorType,
typename DoFHandlerType>
577 return sizeof(*this);
582 #define SPLIT_INSTANTIATIONS_COUNT 4 583 #ifndef SPLIT_INSTANTIATIONS_INDEX 584 #define SPLIT_INSTANTIATIONS_INDEX 0 586 #include "solution_transfer.inst" 588 DEAL_II_NAMESPACE_CLOSE
void extract_interpolation_matrices(const DoFHandlerType &, ::Table< 2, FullMatrix< double > > &)
#define AssertDimension(dim1, dim2)
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()
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)
std::size_t memory_consumption() const
void interpolate(const std::vector< VectorType > &all_in, std::vector< VectorType > &all_out) const
static ::ExceptionBase & ExcInternalError()