Reference documentation for deal.II version GIT acb749f4f5 2024-02-21 16:10:02+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
mg_transfer_matrix_free.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2016 - 2023 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
55template <int dim, typename Number>
57 : public MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number>>
58{
59public:
65
71
75 virtual ~MGTransferMatrixFree() override = default;
76
80 void
82
86 void
87 clear();
88
104 void
105 build(const DoFHandler<dim> &dof_handler,
106 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
107 &external_partitioners =
108 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>());
109
114 void
115 build(const DoFHandler<dim> &dof_handler,
116 const std::function<void(const unsigned int,
119
134 virtual void
136 const unsigned int to_level,
138 const LinearAlgebra::distributed::Vector<Number> &src) const override;
139
140 virtual void
142 const unsigned int to_level,
144 const LinearAlgebra::distributed::Vector<Number> &src) const override;
145
164 virtual void
166 const unsigned int from_level,
168 const LinearAlgebra::distributed::Vector<Number> &src) const override;
169
184 template <typename BlockVectorType2>
185 void
187 const DoFHandler<dim> &dof_handler,
189 const BlockVectorType2 &src) const;
190
194 std::size_t
195 memory_consumption() const;
196
197private:
203 unsigned int fe_degree;
204
210
215 unsigned int n_components;
216
222 unsigned int n_child_cell_dofs;
223
234 std::vector<std::vector<unsigned int>> level_dof_indices;
235
240 std::vector<std::vector<std::pair<unsigned int, unsigned int>>>
242
247 std::vector<unsigned int> n_owned_level_cells;
248
254
259
271 std::vector<AlignedVector<VectorizedArray<Number>>> weights_on_refined;
272
278 std::vector<std::vector<std::vector<unsigned short>>> dirichlet_indices;
279
288
292 template <int degree>
293 void
295 const unsigned int to_level,
298
302 template <int degree>
303 void
304 do_restrict_add(const unsigned int from_level,
307};
308
309
310
316template <int dim, typename Number, typename TransferType>
318 : public MGTransferBase<LinearAlgebra::distributed::BlockVector<Number>>
319{
320public:
324
339 virtual void
341 const unsigned int to_level,
343 const LinearAlgebra::distributed::BlockVector<Number> &src) const override;
344
345 virtual void
347 const unsigned int to_level,
349 const LinearAlgebra::distributed::BlockVector<Number> &src) const override;
350
369 virtual void
371 const unsigned int from_level,
373 const LinearAlgebra::distributed::BlockVector<Number> &src) const override;
374
385 template <typename BlockVectorType2>
386 void
388 const DoFHandler<dim> &dof_handler,
390 const BlockVectorType2 &src) const;
391
395 template <typename BlockVectorType2>
396 void
398 const std::vector<const DoFHandler<dim> *> &dof_handler,
400 const BlockVectorType2 &src) const;
401
405 template <typename BlockVectorType2>
406 void
408 const DoFHandler<dim> &dof_handler,
409 BlockVectorType2 &dst,
411 const;
412
416 template <typename BlockVectorType2>
417 void
419 const std::vector<const DoFHandler<dim> *> &dof_handler,
420 BlockVectorType2 &dst,
422 const;
423
428 static const bool supports_dof_handler_vector = true;
429
430protected:
435 virtual const TransferType &
436 get_matrix_free_transfer(const unsigned int b) const = 0;
437
442 const bool same_for_all;
443};
444
445
446
460template <int dim, typename Number>
463 Number,
464 MGTransferMatrixFree<dim, Number>>
465{
466public:
472
477 MGTransferBlockMatrixFree(const MGConstrainedDoFs &mg_constrained_dofs);
478
483 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
484
488 virtual ~MGTransferBlockMatrixFree() override = default;
489
493 void
494 initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs);
495
499 void
501 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
502
506 void
507 clear();
508
512 void
513 build(const DoFHandler<dim> &dof_handler);
514
518 void
519 build(const std::vector<const DoFHandler<dim> *> &dof_handler);
520
524 std::size_t
525 memory_consumption() const;
526
527protected:
529 get_matrix_free_transfer(const unsigned int b) const override;
530
531private:
535 std::vector<MGTransferMatrixFree<dim, Number>> matrix_free_transfer_vector;
536};
537
538
542//------------------------ templated functions -------------------------
543#ifndef DOXYGEN
544
545
546template <int dim, typename Number>
547template <typename BlockVectorType2>
548void
550 const DoFHandler<dim> &dof_handler,
552 const BlockVectorType2 &src) const
553{
554 const unsigned int min_level = dst.min_level();
555 const unsigned int max_level = dst.max_level();
556
557 Assert(max_level == dof_handler.get_triangulation().n_global_levels() - 1,
559 max_level, dof_handler.get_triangulation().n_global_levels() - 1));
560
561 const auto &fe = dof_handler.get_fe();
562
563 for (unsigned int level = min_level; level <= max_level; ++level)
564 if (dst[level].size() != dof_handler.n_dofs(level) ||
565 dst[level].locally_owned_size() !=
567 dst[level].reinit(this->vector_partitioners[level]);
568
569 // copy fine level vector to active cells in MG hierarchy
570 this->copy_to_mg(dof_handler, dst, src, true);
571
572 // FIXME: maybe need to store hanging nodes constraints per level?
573 // MGConstrainedDoFs does NOT keep this info right now, only periodicity
574 // constraints...
575
576 // do the transfer from level to level-1:
577 dst[max_level].update_ghost_values();
578
579 for (unsigned int level = max_level; level > min_level; --level)
580 {
581 // auxiliary vector which always has ghost elements
582 const LinearAlgebra::distributed::Vector<Number> *input = nullptr;
584 if (dst[level].get_partitioner().get() ==
585 this->vector_partitioners[level].get())
586 input = &dst[level];
587 else
588 {
589 ghosted_fine.reinit(this->vector_partitioners[level]);
590 ghosted_fine.copy_locally_owned_data_from(dst[level]);
591 ghosted_fine.update_ghost_values();
592 input = &ghosted_fine;
593 }
594
595 std::vector<Number> dof_values_coarse(fe.n_dofs_per_cell());
596 Vector<Number> dof_values_fine(fe.n_dofs_per_cell());
597 Vector<Number> tmp(fe.n_dofs_per_cell());
598 std::vector<types::global_dof_index> dof_indices(fe.n_dofs_per_cell());
599 for (const auto &cell : dof_handler.cell_iterators_on_level(level - 1))
600 if (cell->is_locally_owned_on_level())
601 {
602 // if we get to a cell without children (== active), we can
603 // skip it as there values should be already set by the
604 // equivalent of copy_to_mg()
605 if (cell->is_active())
606 continue;
607
608 std::fill(dof_values_coarse.begin(), dof_values_coarse.end(), 0.);
609 for (unsigned int child = 0; child < cell->n_children(); ++child)
610 {
611 cell->child(child)->get_mg_dof_indices(dof_indices);
612
614 this->mg_constrained_dofs, level, dof_indices);
615
616 for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
617 dof_values_fine(i) = (*input)(dof_indices[i]);
618 fe.get_restriction_matrix(child, cell->refinement_case())
619 .vmult(tmp, dof_values_fine);
620 for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
621 if (fe.restriction_is_additive(i))
622 dof_values_coarse[i] += tmp[i];
623 else if (tmp(i) != 0.)
624 dof_values_coarse[i] = tmp[i];
625 }
626 cell->get_mg_dof_indices(dof_indices);
627 for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
628 if (dof_handler.locally_owned_mg_dofs(level - 1).is_element(
629 dof_indices[i]))
630 dst[level - 1](dof_indices[i]) = dof_values_coarse[i];
631 }
632
633 dst[level - 1].update_ghost_values();
634 }
635
636 for (unsigned int level = min_level; level <= max_level; ++level)
637 dst[level].zero_out_ghost_values();
638}
639
640
641
642template <int dim, typename Number, typename TransferType>
643template <typename BlockVectorType2>
644void
646 const DoFHandler<dim> &dof_handler,
648 const BlockVectorType2 &src) const
649{
650 Assert(same_for_all,
652 "This object was initialized with support for usage with one "
653 "DoFHandler for each block, but this method assumes that "
654 "the same DoFHandler is used for all the blocks!"));
655 const std::vector<const DoFHandler<dim> *> mg_dofs(src.n_blocks(),
656 &dof_handler);
657
658 copy_to_mg(mg_dofs, dst, src);
659}
660
661
662
663template <int dim, typename Number, typename TransferType>
664template <typename BlockVectorType2>
665void
667 const std::vector<const DoFHandler<dim> *> &dof_handler,
669 const BlockVectorType2 &src) const
670{
671 const unsigned int n_blocks = src.n_blocks();
672 AssertDimension(dof_handler.size(), n_blocks);
673
674 if (n_blocks == 0)
675 return;
676
677 const unsigned int min_level = dst.min_level();
678 const unsigned int max_level = dst.max_level();
679
680 for (unsigned int level = min_level; level <= max_level; ++level)
681 if (dst[level].n_blocks() != n_blocks)
682 dst[level].reinit(n_blocks);
683
684 // FIXME: this a quite ugly as we need a temporary object:
686 min_level, max_level);
687
688 for (unsigned int b = 0; b < n_blocks; ++b)
689 {
690 const unsigned int data_block = same_for_all ? 0 : b;
691 get_matrix_free_transfer(data_block)
692 .copy_to_mg(*dof_handler[b], dst_non_block, src.block(b));
693
694 for (unsigned int l = min_level; l <= max_level; ++l)
695 dst[l].block(b) = dst_non_block[l];
696 }
697
698 for (unsigned int level = min_level; level <= max_level; ++level)
699 dst[level].collect_sizes();
700}
701
702template <int dim, typename Number, typename TransferType>
703template <typename BlockVectorType2>
704void
706 const DoFHandler<dim> &dof_handler,
707 BlockVectorType2 &dst,
709 const
710{
711 Assert(same_for_all,
713 "This object was initialized with support for usage with one "
714 "DoFHandler for each block, but this method assumes that "
715 "the same DoFHandler is used for all the blocks!"));
716 const std::vector<const DoFHandler<dim> *> mg_dofs(dst.n_blocks(),
717 &dof_handler);
718
719 copy_from_mg(mg_dofs, dst, src);
720}
721
722template <int dim, typename Number, typename TransferType>
723template <typename BlockVectorType2>
724void
726 const std::vector<const DoFHandler<dim> *> &dof_handler,
727 BlockVectorType2 &dst,
729 const
730{
731 const unsigned int n_blocks = dst.n_blocks();
732 AssertDimension(dof_handler.size(), n_blocks);
733
734 if (n_blocks == 0)
735 return;
736
737 const unsigned int min_level = src.min_level();
738 const unsigned int max_level = src.max_level();
739
740 for (unsigned int l = min_level; l <= max_level; ++l)
741 AssertDimension(src[l].n_blocks(), dst.n_blocks());
742
743 // FIXME: this a quite ugly as we need a temporary object:
745 min_level, max_level);
746
747 for (unsigned int b = 0; b < n_blocks; ++b)
748 {
749 for (unsigned int l = min_level; l <= max_level; ++l)
750 {
751 src_non_block[l].reinit(src[l].block(b));
752 src_non_block[l] = src[l].block(b);
753 }
754 const unsigned int data_block = same_for_all ? 0 : b;
755 get_matrix_free_transfer(data_block)
756 .copy_from_mg(*dof_handler[b], dst.block(b), src_non_block);
757 }
758}
759
760
761
762#endif // DOXYGEN
763
764
766
767#endif
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index 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
size_type n_elements() const
Definition index_set.h:1924
bool is_element(const size_type index) const
Definition index_set.h:1884
void copy_locally_owned_data_from(const Vector< Number2, MemorySpace > &src)
void reinit(const size_type size, const bool omit_zeroing_entries=false)
std::function< void(const unsigned int, LinearAlgebra::distributed::Vector< Number > &)> initialize_dof_vector
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
virtual void prolongate_and_add(const unsigned int to_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
virtual const TransferType & get_matrix_free_transfer(const unsigned int b) const =0
void copy_from_mg(const std::vector< const DoFHandler< dim > * > &dof_handler, BlockVectorType2 &dst, const MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &src) const
void copy_to_mg(const DoFHandler< dim > &dof_handler, MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &dst, const BlockVectorType2 &src) const
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
MGTransferBlockMatrixFreeBase(const bool same_for_all)
void copy_from_mg(const DoFHandler< dim > &dof_handler, BlockVectorType2 &dst, const MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &src) const
void copy_to_mg(const std::vector< const DoFHandler< dim > * > &dof_handler, MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &dst, const BlockVectorType2 &src) const
std::vector< MGTransferMatrixFree< dim, Number > > matrix_free_transfer_vector
void build(const DoFHandler< dim > &dof_handler)
MGTransferBlockMatrixFree()=default
virtual ~MGTransferBlockMatrixFree() override=default
const MGTransferMatrixFree< dim, Number > & get_matrix_free_transfer(const unsigned int b) const override
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
void interpolate_to_mg(const DoFHandler< dim > &dof_handler, MGLevelObject< LinearAlgebra::distributed::Vector< Number > > &dst, const BlockVectorType2 &src) const
MGLevelObject< std::shared_ptr< const Utilities::MPI::Partitioner > > vector_partitioners
AlignedVector< VectorizedArray< Number > > evaluation_data
void build(const DoFHandler< dim > &dof_handler, const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &external_partitioners=std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > >())
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
virtual unsigned int n_global_levels() const
void update_ghost_values() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
unsigned int level
Definition grid_out.cc:4617
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
std::enable_if_t< IsBlockVector< VectorType >::value, void > collect_sizes(VectorType &vector)
Definition operators.h:96
std::enable_if_t< IsBlockVector< VectorType >::value, unsigned int > n_blocks(const VectorType &vector)
Definition operators.h:50
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 resolve_identity_constraints(const MGConstrainedDoFs *mg_constrained_dofs, const unsigned int level, std::vector< types::global_dof_index > &dof_indices)