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> 33 DEAL_II_NAMESPACE_OPEN
37 template<
int dim,
typename VectorType,
typename DoFHandlerType>
40 dof_handler(&dof, typeid(*this).name()),
47 template<
int dim,
typename VectorType,
typename DoFHandlerType>
55 template<
int dim,
typename VectorType,
typename DoFHandlerType>
67 template<
int dim,
typename VectorType,
typename DoFHandlerType>
72 ExcAlreadyPrepForCoarseAndRef());
76 const unsigned int n_active_cells =
dof_handler->get_triangulation().n_active_cells();
80 std::vector<std::vector<types::global_dof_index> > (n_active_cells)
83 typename DoFHandlerType::active_cell_iterator cell =
dof_handler->begin_active(),
86 for (
unsigned int i=0; cell!=endc; ++cell, ++i)
96 cell_map[std::make_pair(cell->level(),cell->index())]
104 template<
int dim,
typename VectorType,
typename DoFHandlerType>
107 (
const VectorType &in,
108 VectorType &out)
const 113 ExcDimensionMismatch(out.size(),
dof_handler->n_dofs()));
115 ExcMessage (
"Vectors cannot be used as input and output" 116 " at the same time!"));
120 typename DoFHandlerType::cell_iterator cell =
dof_handler->begin(),
123 typename std::map<std::pair<unsigned int, unsigned int>,
Pointerstruct>::const_iterator
127 for (; cell!=endc; ++cell)
129 pointerstruct=
cell_map.find(std::make_pair(cell->level(),cell->index()));
131 if (pointerstruct!=cell_map_end)
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);
146 Assert(dofs_per_cell==(*pointerstruct->second.indices_ptr).size(),
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,
173 template <
typename DoFHandlerType>
178 template <
int dim,
int spacedim>
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)
188 matrices(i,j).reinit (fe[i].dofs_per_cell, fe[j].dofs_per_cell);
199 fe[i].get_interpolation_matrix (fe[j], matrices(i,j));
203 matrices(i,j).reinit (0,0);
209 template <
int dim,
int spacedim>
211 std::vector<std::vector<bool> > &)
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)
218 restriction_is_additive.resize (fe.size());
219 for (
unsigned int f=0; f<fe.size(); ++f)
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);
230 template<
int dim,
typename VectorType,
typename DoFHandlerType>
237 ExcAlreadyPrepForCoarseAndRef());
239 const unsigned int in_size=all_in.size();
241 ExcMessage(
"The array of input vectors you pass to this " 242 "function has no elements. This is not useful."));
246 const unsigned int n_active_cells =
dof_handler->get_triangulation().n_active_cells();
247 (void)n_active_cells;
250 for (
unsigned int i=0; i<in_size; ++i)
253 ExcDimensionMismatch(all_in[i].size(),
n_dofs_old));
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();
264 if (act_cell->coarsen_flag_set())
265 ++n_cells_to_coarsen;
267 ++n_cells_to_stay_or_refine;
269 Assert((n_cells_to_coarsen+n_cells_to_stay_or_refine)==n_active_cells,
272 unsigned int n_coarsen_fathers=0;
273 for (
typename DoFHandlerType::cell_iterator cell=
dof_handler->begin();
275 if (!cell->active() && cell->child(0)->coarsen_flag_set())
277 Assert(n_cells_to_coarsen>=2*n_coarsen_fathers, ExcInternalError());
282 std::vector<std::vector<types::global_dof_index> >(n_cells_to_stay_or_refine)
285 std::vector<std::vector<Vector<typename VectorType::value_type> > >
287 std::vector<Vector<typename VectorType::value_type> > (in_size))
291 std::vector<std::vector<bool> > restriction_is_additive;
294 internal::restriction_additive (
dof_handler->get_fe(), restriction_is_additive);
299 unsigned int n_sr=0, n_cf=0;
300 for (
typename DoFHandlerType::cell_iterator cell=
dof_handler->begin();
304 if (cell->active() && !cell->coarsen_flag_set())
306 const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell;
313 cell_map[std::make_pair(cell->level(), cell->index())]
319 else if (cell->has_children() && cell->child(0)->coarsen_flag_set())
325 for (
unsigned int i=1; i<cell->n_children(); ++i)
326 Assert(cell->child(i)->coarsen_flag_set(),
328 "Triangulation::prepare_coarsening_and_refinement before " 329 "calling the current function. This can't work."));
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())
341 different_fe_on_children =
true;
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();
355 const unsigned int dofs_per_cell=cell->get_dof_handler().get_fe()[target_fe_index].dofs_per_cell;
357 std::vector<Vector<typename VectorType::value_type> >(in_size,
366 for (
unsigned int j=0; j<in_size; ++j)
367 cell->get_interpolated_dof_values(all_in[j],
370 cell_map[std::make_pair(cell->level(), cell->index())]
375 Assert(n_sr==n_cells_to_stay_or_refine, ExcInternalError());
376 Assert(n_cf==n_coarsen_fathers, ExcInternalError());
383 template<
int dim,
typename VectorType,
typename DoFHandlerType>
386 (
const VectorType &in)
388 std::vector<VectorType> all_in=std::vector<VectorType>(1, in);
394 template<
int dim,
typename VectorType,
typename DoFHandlerType>
397 std::vector<VectorType> &all_out)
const 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)
404 ExcDimensionMismatch(all_in[i].size(),
n_dofs_old));
405 for (
unsigned int i=0; i<all_out.size(); ++i)
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!"));
415 std::vector<types::global_dof_index> dofs;
417 typename std::map<std::pair<unsigned int, unsigned int>,
Pointerstruct>::const_iterator
425 typename DoFHandlerType::cell_iterator cell =
dof_handler->begin(),
427 for (; cell!=endc; ++cell)
429 pointerstruct=
cell_map.find(std::make_pair(cell->level(),cell->index()));
431 if (pointerstruct!=cell_map_end)
433 const std::vector<types::global_dof_index> *
const indexptr
434 =pointerstruct->second.indices_ptr;
436 const std::vector<Vector<typename VectorType::value_type> > *
const valuesptr
437 =pointerstruct->second.dof_values_ptr;
445 const unsigned int old_fe_index =
446 pointerstruct->second.active_fe_index;
453 unsigned int in_size = indexptr->size();
454 for (
unsigned int j=0; j<size; ++j)
456 tmp.
reinit (in_size,
true);
457 for (
unsigned int i=0; i<in_size; ++i)
458 tmp(i) = all_in[j]((*indexptr)[i]);
460 cell->set_dof_values_by_interpolation (tmp, all_out[j],
468 Assert (!cell->has_children(), ExcInternalError());
472 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
473 dofs.resize(dofs_per_cell);
476 cell->get_dof_indices(dofs);
481 for (
unsigned int j=0; j<size; ++j)
494 const unsigned int active_fe_index = cell->active_fe_index();
495 if (active_fe_index != pointerstruct->second.active_fe_index)
497 const unsigned int old_index = pointerstruct->second.active_fe_index;
498 tmp.
reinit (dofs_per_cell,
true);
500 interpolation_hp(active_fe_index,old_index).n());
502 interpolation_hp(active_fe_index,old_index).m());
503 interpolation_hp(active_fe_index,old_index).vmult (tmp, (*valuesptr)[j]);
507 data = &(*valuesptr)[j];
510 for (
unsigned int i=0; i<dofs_per_cell; ++i)
511 all_out[j](dofs[i])=(*data)(i);
516 Assert(
false, ExcInternalError());
523 template<
int dim,
typename VectorType,
typename DoFHandlerType>
525 (
const VectorType &in,
526 VectorType &out)
const 531 ExcDimensionMismatch(out.size(),
dof_handler->n_dofs()));
533 std::vector<VectorType> all_in(1);
535 std::vector<VectorType> all_out(1);
544 template<
int dim,
typename VectorType,
typename DoFHandlerType>
561 template<
int dim,
typename VectorType,
typename DoFHandlerType>
565 return sizeof(*this);
570 #define SPLIT_INSTANTIATIONS_COUNT 4 571 #ifndef SPLIT_INSTANTIATIONS_INDEX 572 #define SPLIT_INSTANTIATIONS_INDEX 0 574 #include "solution_transfer.inst" 576 DEAL_II_NAMESPACE_CLOSE
void extract_interpolation_matrices(const DoFHandlerType &, ::Table< 2, FullMatrix< double > > &)
#define AssertDimension(dim1, dim2)
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)
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::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