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\}}\)
mg_transfer_matrix_free.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2016 - 2021 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_mg_transfer_matrix_free_h
17#define dealii_mg_transfer_matrix_free_h
18
19#include <deal.II/base/config.h>
20
23
25
28
33
34
36
37
40
53template <int dim, typename Number>
55 : public MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number>>
56{
57public:
63
69
73 virtual ~MGTransferMatrixFree() override = default;
74
78 void
80
84 void
85 clear();
86
102 void
103 build(const DoFHandler<dim, dim> &dof_handler,
104 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
105 &external_partitioners =
106 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>());
107
122 virtual void
124 const unsigned int to_level,
126 const LinearAlgebra::distributed::Vector<Number> &src) const override;
127
128 virtual void
130 const unsigned int to_level,
132 const LinearAlgebra::distributed::Vector<Number> &src) const override;
133
152 virtual void
154 const unsigned int from_level,
156 const LinearAlgebra::distributed::Vector<Number> &src) const override;
157
172 template <typename Number2, int spacedim>
173 void
175 const DoFHandler<dim, spacedim> & dof_handler,
178
183
187 std::size_t
188 memory_consumption() const;
189
190private:
196 unsigned int fe_degree;
197
203
208 unsigned int n_components;
209
215 unsigned int n_child_cell_dofs;
216
227 std::vector<std::vector<unsigned int>> level_dof_indices;
228
233 std::vector<std::vector<std::pair<unsigned int, unsigned int>>>
235
240 std::vector<unsigned int> n_owned_level_cells;
241
247
252
264 std::vector<AlignedVector<VectorizedArray<Number>>> weights_on_refined;
265
271 std::vector<std::vector<std::vector<unsigned short>>> dirichlet_indices;
272
281
285 template <int degree>
286 void
288 const unsigned int to_level,
291
295 template <int degree>
296 void
297 do_restrict_add(const unsigned int from_level,
300};
301
302
316template <int dim, typename Number>
318 : public MGTransferBase<LinearAlgebra::distributed::BlockVector<Number>>
319{
320public:
326
331 MGTransferBlockMatrixFree(const MGConstrainedDoFs &mg_constrained_dofs);
332
337 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
338
342 virtual ~MGTransferBlockMatrixFree() override = default;
343
347 void
348 initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs);
349
353 void
355 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
356
360 void
361 clear();
362
366 void
367 build(const DoFHandler<dim, dim> &dof_handler);
368
372 void
373 build(const std::vector<const DoFHandler<dim, dim> *> &dof_handler);
374
389 virtual void
391 const unsigned int to_level,
393 const LinearAlgebra::distributed::BlockVector<Number> &src) const override;
394
395 virtual void
397 const unsigned int to_level,
399 const LinearAlgebra::distributed::BlockVector<Number> &src) const override;
400
419 virtual void
421 const unsigned int from_level,
423 const LinearAlgebra::distributed::BlockVector<Number> &src) const override;
424
435 template <typename Number2, int spacedim>
436 void
438 const DoFHandler<dim, spacedim> & dof_handler,
441
445 template <typename Number2, int spacedim>
446 void
448 const std::vector<const DoFHandler<dim, spacedim> *> & dof_handler,
451
455 template <typename Number2, int spacedim>
456 void
458 const DoFHandler<dim, spacedim> & dof_handler,
461 const;
462
466 template <typename Number2, int spacedim>
467 void
469 const std::vector<const DoFHandler<dim, spacedim> *> &dof_handler,
472 const;
473
477 std::size_t
478 memory_consumption() const;
479
484 static const bool supports_dof_handler_vector = true;
485
486private:
490 std::vector<MGTransferMatrixFree<dim, Number>> matrix_free_transfer_vector;
491
496 const bool same_for_all;
497};
498
499
503//------------------------ templated functions -------------------------
504#ifndef DOXYGEN
505
506
507template <int dim, typename Number>
508template <typename Number2, int spacedim>
509void
511 const DoFHandler<dim, spacedim> & dof_handler,
514{
515 const unsigned int min_level = dst.min_level();
516 const unsigned int max_level = dst.max_level();
517
518 Assert(max_level == dof_handler.get_triangulation().n_global_levels() - 1,
520 max_level, dof_handler.get_triangulation().n_global_levels() - 1));
521
522 const FiniteElement<dim, spacedim> &fe = dof_handler.get_fe();
523
524 for (unsigned int level = min_level; level <= max_level; ++level)
525 if (dst[level].size() != dof_handler.n_dofs(level) ||
526 dst[level].locally_owned_size() !=
528 dst[level].reinit(this->vector_partitioners[level]);
529
530 // copy fine level vector to active cells in MG hierarchy
531 this->copy_to_mg(dof_handler, dst, src, true);
532
533 // FIXME: maybe need to store hanging nodes constraints per level?
534 // MGConstrainedDoFs does NOT keep this info right now, only periodicity
535 // constraints...
536
537 // do the transfer from level to level-1:
538 dst[max_level].update_ghost_values();
539 for (unsigned int level = max_level; level > min_level; --level)
540 {
541 // auxiliary vector which always has ghost elements
542 const LinearAlgebra::distributed::Vector<Number> *input = nullptr;
544 if (dst[level].get_partitioner().get() ==
545 this->vector_partitioners[level].get())
546 input = &dst[level];
547 else
548 {
549 ghosted_fine.reinit(this->vector_partitioners[level]);
550 ghosted_fine.copy_locally_owned_data_from(dst[level]);
551 ghosted_fine.update_ghost_values();
552 input = &ghosted_fine;
553 }
554
555 std::vector<Number> dof_values_coarse(fe.n_dofs_per_cell());
556 Vector<Number> dof_values_fine(fe.n_dofs_per_cell());
558 std::vector<types::global_dof_index> dof_indices(fe.n_dofs_per_cell());
559 for (const auto &cell : dof_handler.cell_iterators_on_level(level - 1))
560 if (cell->is_locally_owned_on_level())
561 {
562 // if we get to a cell without children (== active), we can
563 // skip it as there values should be already set by the
564 // equivalent of copy_to_mg()
565 if (cell->is_active())
566 continue;
567
568 std::fill(dof_values_coarse.begin(), dof_values_coarse.end(), 0.);
569 for (unsigned int child = 0; child < cell->n_children(); ++child)
570 {
571 cell->child(child)->get_mg_dof_indices(dof_indices);
572 for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
573 dof_values_fine(i) = (*input)(dof_indices[i]);
574 fe.get_restriction_matrix(child, cell->refinement_case())
575 .vmult(tmp, dof_values_fine);
576 for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
577 if (fe.restriction_is_additive(i))
578 dof_values_coarse[i] += tmp[i];
579 else if (tmp(i) != 0.)
580 dof_values_coarse[i] = tmp[i];
581 }
582 cell->get_mg_dof_indices(dof_indices);
583 for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
584 if (dof_handler.locally_owned_mg_dofs(level - 1).is_element(
585 dof_indices[i]))
586 dst[level - 1](dof_indices[i]) = dof_values_coarse[i];
587 }
588
589 dst[level - 1].update_ghost_values();
590 }
591}
592
593
594
595template <int dim, typename Number>
596template <typename Number2, int spacedim>
597void
599 const DoFHandler<dim, spacedim> & dof_handler,
602{
603 AssertDimension(matrix_free_transfer_vector.size(), 1);
604 Assert(same_for_all,
606 "This object was initialized with support for usage with one "
607 "DoFHandler for each block, but this method assumes that "
608 "the same DoFHandler is used for all the blocks!"));
609 const std::vector<const DoFHandler<dim, spacedim> *> mg_dofs(src.n_blocks(),
610 &dof_handler);
611
612 copy_to_mg(mg_dofs, dst, src);
613}
614
615
616
617template <int dim, typename Number>
618template <typename Number2, int spacedim>
619void
621 const std::vector<const DoFHandler<dim, spacedim> *> & dof_handler,
624{
625 const unsigned int n_blocks = src.n_blocks();
626 AssertDimension(dof_handler.size(), n_blocks);
627
628 if (n_blocks == 0)
629 return;
630
631 const unsigned int min_level = dst.min_level();
632 const unsigned int max_level = dst.max_level();
633
634 // this function is normally called within the Multigrid class with
635 // dst == defect level block vector. At first run this vector is not
636 // initialized. Do this below:
637 {
639 (dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
640 &(dof_handler[0]->get_triangulation())));
641 for (unsigned int i = 1; i < n_blocks; ++i)
643 (dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
644 &(dof_handler[0]->get_triangulation())) == tria),
645 ExcMessage("The DoFHandler use different Triangulations!"));
646
647 MGLevelObject<bool> do_reinit;
648 do_reinit.resize(min_level, max_level);
649 for (unsigned int level = min_level; level <= max_level; ++level)
650 {
651 do_reinit[level] = false;
652 if (dst[level].n_blocks() != n_blocks)
653 {
654 do_reinit[level] = true;
655 continue; // level
656 }
657 for (unsigned int b = 0; b < n_blocks; ++b)
658 {
660 if (v.size() !=
661 dof_handler[b]->locally_owned_mg_dofs(level).size() ||
662 v.locally_owned_size() !=
663 dof_handler[b]->locally_owned_mg_dofs(level).n_elements())
664 {
665 do_reinit[level] = true;
666 break; // b
667 }
668 }
669 }
670
671 for (unsigned int level = min_level; level <= max_level; ++level)
672 {
673 if (do_reinit[level])
674 {
675 dst[level].reinit(n_blocks);
676 for (unsigned int b = 0; b < n_blocks; ++b)
677 {
679 dst[level].block(b);
680 v.reinit(dof_handler[b]->locally_owned_mg_dofs(level),
681 dof_handler[b]->get_communicator());
682 }
683 dst[level].collect_sizes();
684 }
685 else
686 dst[level] = 0;
687 }
688 }
689
690 // FIXME: this a quite ugly as we need a temporary object:
692 min_level, max_level);
693
694 for (unsigned int b = 0; b < n_blocks; ++b)
695 {
696 for (unsigned int l = min_level; l <= max_level; ++l)
697 dst_non_block[l].reinit(dst[l].block(b));
698 const unsigned int data_block = same_for_all ? 0 : b;
699 matrix_free_transfer_vector[data_block].copy_to_mg(*dof_handler[b],
700 dst_non_block,
701 src.block(b));
702
703 for (unsigned int l = min_level; l <= max_level; ++l)
704 dst[l].block(b) = dst_non_block[l];
705 }
706}
707
708template <int dim, typename Number>
709template <typename Number2, int spacedim>
710void
712 const DoFHandler<dim, spacedim> & dof_handler,
715 const
716{
717 AssertDimension(matrix_free_transfer_vector.size(), 1);
718 const std::vector<const DoFHandler<dim, spacedim> *> mg_dofs(dst.n_blocks(),
719 &dof_handler);
720
721 copy_from_mg(mg_dofs, dst, src);
722}
723
724template <int dim, typename Number>
725template <typename Number2, int spacedim>
726void
728 const std::vector<const DoFHandler<dim, spacedim> *> &dof_handler,
731 const
732{
733 const unsigned int n_blocks = dst.n_blocks();
734 AssertDimension(dof_handler.size(), n_blocks);
735
736 if (n_blocks == 0)
737 return;
738
739 const unsigned int min_level = src.min_level();
740 const unsigned int max_level = src.max_level();
741
742 for (unsigned int l = min_level; l <= max_level; ++l)
743 AssertDimension(src[l].n_blocks(), dst.n_blocks());
744
745 // FIXME: this a quite ugly as we need a temporary object:
747 min_level, max_level);
748
749 for (unsigned int b = 0; b < n_blocks; ++b)
750 {
751 for (unsigned int l = min_level; l <= max_level; ++l)
752 {
753 src_non_block[l].reinit(src[l].block(b));
754 src_non_block[l] = src[l].block(b);
755 }
756 const unsigned int data_block = same_for_all ? 0 : b;
757 matrix_free_transfer_vector[data_block].copy_from_mg(*dof_handler[b],
758 dst.block(b),
759 src_non_block);
760 }
761}
762
763
764
765#endif // DOXYGEN
766
767
769
770#endif
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
const Triangulation< dim, spacedim > & get_triangulation() const
types::global_dof_index n_dofs() const
MPI_Comm get_communicator() const
unsigned int n_dofs_per_cell() const
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
bool restriction_is_additive(const unsigned int index) const
size_type size() const
Definition: index_set.h:1634
size_type n_elements() const
Definition: index_set.h:1832
bool is_element(const size_type index) const
Definition: index_set.h:1765
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > > > mg_constrained_dofs
Definition: mg_transfer.h:586
void resize(const unsigned int new_minlevel, const unsigned int new_maxlevel, Args &&... args)
void copy_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &dst, const LinearAlgebra::distributed::BlockVector< Number2 > &src) const
std::vector< MGTransferMatrixFree< dim, Number > > matrix_free_transfer_vector
virtual void prolongate_and_add(const unsigned int to_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
static const bool supports_dof_handler_vector
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
std::size_t memory_consumption() const
MGTransferBlockMatrixFree()=default
void copy_to_mg(const std::vector< const DoFHandler< dim, spacedim > * > &dof_handler, MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &dst, const LinearAlgebra::distributed::BlockVector< Number2 > &src) const
void copy_from_mg(const DoFHandler< dim, spacedim > &dof_handler, LinearAlgebra::distributed::BlockVector< Number2 > &dst, const MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &src) const
virtual ~MGTransferBlockMatrixFree() override=default
void build(const DoFHandler< dim, dim > &dof_handler)
void copy_from_mg(const std::vector< const DoFHandler< dim, spacedim > * > &dof_handler, LinearAlgebra::distributed::BlockVector< Number2 > &dst, const MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &src) const
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
void do_restrict_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
void do_prolongate_add(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
virtual ~MGTransferMatrixFree() override=default
std::vector< std::vector< unsigned int > > level_dof_indices
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const override
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
std::vector< std::vector< std::vector< unsigned short > > > dirichlet_indices
MGLevelObject< std::shared_ptr< const Utilities::MPI::Partitioner > > vector_partitioners
void build(const DoFHandler< dim, dim > &dof_handler, const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &external_partitioners=std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > >())
AlignedVector< VectorizedArray< Number > > evaluation_data
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > parent_child_connect
virtual void prolongate_and_add(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const override
std::size_t memory_consumption() const
std::vector< AlignedVector< VectorizedArray< Number > > > weights_on_refined
std::vector< unsigned int > n_owned_level_cells
AlignedVector< VectorizedArray< Number > > prolongation_matrix_1d
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const override
void interpolate_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< LinearAlgebra::distributed::Vector< Number > > &dst, const LinearAlgebra::distributed::Vector< Number2 > &src) const
virtual unsigned int n_global_levels() const
Definition: vector.h:110
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
unsigned int level
Definition: grid_out.cc:4590
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
#define DeclException0(Exception0)
Definition: exceptions.h:470
static ::ExceptionBase & ExcNoProlongation()
unsigned int n_blocks() const
#define Assert(cond, exc)
Definition: exceptions.h:1465
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
BlockType & block(const unsigned int i)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
size_type locally_owned_size() const
virtual size_type size() const override
void copy_locally_owned_data_from(const Vector< Number2, MemorySpace > &src)
void reinit(const size_type size, const bool omit_zeroing_entries=false)
void update_ghost_values() const
std::enable_if< IsBlockVector< VectorType >::value, unsignedint >::type n_blocks(const VectorType &vector)
Definition: operators.h:49
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:618