Reference documentation for deal.II version 9.3.3
\(\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.h
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
16#ifndef dealii_solution_transfer_h
17# define dealii_solution_transfer_h
18
19
20/*---------------------------- solutiontransfer.h ----------------------*/
21
22
23# include <deal.II/base/config.h>
24
27
29
30# include <deal.II/lac/vector.h>
31
32# include <vector>
33
35
335template <int dim,
336 typename VectorType = Vector<double>,
337 typename DoFHandlerType = DoFHandler<dim>>
339{
340# ifndef DEAL_II_MSVC
341 static_assert(dim == DoFHandlerType::dimension,
342 "The dimension explicitly provided as a template "
343 "argument, and the dimension of the DoFHandlerType "
344 "template argument must match.");
345# endif
346public:
350 SolutionTransfer(const DoFHandlerType &dof);
351
356
361 void
362 clear();
363
369 void
371
379 void
380 prepare_for_coarsening_and_refinement(const std::vector<VectorType> &all_in);
381
386 void
387 prepare_for_coarsening_and_refinement(const VectorType &in);
388
400 void
401 refine_interpolate(const VectorType &in, VectorType &out) const;
402
422 void
423 interpolate(const std::vector<VectorType> &all_in,
424 std::vector<VectorType> & all_out) const;
425
435 void
436 interpolate(const VectorType &in, VectorType &out) const;
437
442 std::size_t
443 memory_consumption() const;
444
449 "You are attempting an operation for which this object is "
450 "not prepared. This may be because you either did not call "
451 "one of the prepare_*() functions at all, or because you "
452 "called the wrong one for the operation you are currently "
453 "attempting.");
454
460 "You are attempting to call one of the prepare_*() functions "
461 "of this object to prepare it for an operation for which it "
462 "is already prepared. Specifically, the object was "
463 "previously prepared for pure refinement.");
464
470 "You are attempting to call one of the prepare_*() functions "
471 "of this object to prepare it for an operation for which it "
472 "is already prepared. Specifically, the object was "
473 "previously prepared for both coarsening and refinement.");
474
475private:
479 SmartPointer<const DoFHandlerType,
482
487
494 {
507 };
508
513
514
520 std::vector<std::vector<types::global_dof_index>> indices_on_cell;
521
534 {
536 : indices_ptr(nullptr)
537 , dof_values_ptr(nullptr)
538 , active_fe_index(0)
539 {}
540 Pointerstruct(std::vector<types::global_dof_index> *indices_ptr_in,
541 const unsigned int active_fe_index_in = 0)
542 : indices_ptr(indices_ptr_in)
543 , dof_values_ptr(nullptr)
544 , active_fe_index(active_fe_index_in)
545 {}
547 std::vector<Vector<typename VectorType::value_type>> *dof_values_ptr_in,
548 const unsigned int active_fe_index_in = 0)
549 : indices_ptr(nullptr)
550 , dof_values_ptr(dof_values_ptr_in)
551 , active_fe_index(active_fe_index_in)
552 {}
553 std::size_t
554 memory_consumption() const;
555
556 std::vector<types::global_dof_index> * indices_ptr;
557 std::vector<Vector<typename VectorType::value_type>> *dof_values_ptr;
558 unsigned int active_fe_index;
559 };
560
567 std::map<std::pair<unsigned int, unsigned int>, Pointerstruct> cell_map;
568
573 std::vector<std::vector<Vector<typename VectorType::value_type>>>
575};
576
577namespace Legacy
578{
585 template <int dim,
586 typename VectorType = Vector<double>,
587 typename DoFHandlerType = DoFHandler<dim>>
590} // namespace Legacy
591
592
594
595#endif // dealii_solutiontransfer_h
596/*---------------------------- solutiontransfer.h ---------------------------*/
types::global_dof_index n_dofs_old
std::vector< std::vector< Vector< typename VectorType::value_type > > > dof_values_on_cell
std::vector< std::vector< types::global_dof_index > > indices_on_cell
void prepare_for_pure_refinement()
SolutionTransfer(const DoFHandlerType &dof)
std::size_t memory_consumption() const
std::map< std::pair< unsigned int, unsigned int >, Pointerstruct > cell_map
SmartPointer< const DoFHandlerType, SolutionTransfer< dim, VectorType, DoFHandlerType > > dof_handler
void refine_interpolate(const VectorType &in, VectorType &out) const
PreparationState prepared_for
void interpolate(const std::vector< VectorType > &all_in, std::vector< VectorType > &all_out) const
void prepare_for_coarsening_and_refinement(const std::vector< VectorType > &all_in)
Definition: vector.h:110
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
static ::ExceptionBase & ExcNotPrepared()
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:493
static ::ExceptionBase & ExcAlreadyPrepForCoarseAndRef()
static ::ExceptionBase & ExcAlreadyPrepForRef()
std::vector< types::global_dof_index > * indices_ptr
Pointerstruct(std::vector< Vector< typename VectorType::value_type > > *dof_values_ptr_in, const unsigned int active_fe_index_in=0)
std::size_t memory_consumption() const
Pointerstruct(std::vector< types::global_dof_index > *indices_ptr_in, const unsigned int active_fe_index_in=0)
std::vector< Vector< typename VectorType::value_type > > * dof_values_ptr