deal.II version GIT relicensing-2169-gec1b43f35b 2024-11-22 07:10:00+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// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2016 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_mg_transfer_matrix_free_h
16#define dealii_mg_transfer_matrix_free_h
17
18#include <deal.II/base/config.h>
19
22
24
27
32
33
35
36
54template <int dim, typename Number>
56 : public MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number>>
57{
58public:
64
70
74 virtual ~MGTransferMatrixFree() override = default;
75
79 void
81
85 void
86 clear();
87
103 void
104 build(const DoFHandler<dim> &dof_handler,
105 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
106 &external_partitioners =
107 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>());
108
113 void
114 build(const DoFHandler<dim> &dof_handler,
115 const std::function<void(const unsigned int,
118
133 virtual void
135 const unsigned int to_level,
137 const LinearAlgebra::distributed::Vector<Number> &src) const override;
138
139 virtual void
141 const unsigned int to_level,
143 const LinearAlgebra::distributed::Vector<Number> &src) const override;
144
163 virtual void
165 const unsigned int from_level,
167 const LinearAlgebra::distributed::Vector<Number> &src) const override;
168
183 template <typename BlockVectorType2>
184 void
186 const DoFHandler<dim> &dof_handler,
188 const BlockVectorType2 &src) const;
189
193 std::size_t
194 memory_consumption() const;
195
196private:
202 unsigned int fe_degree;
203
209
214 unsigned int n_components;
215
221 unsigned int n_child_cell_dofs;
222
233 std::vector<std::vector<unsigned int>> level_dof_indices;
234
239 std::vector<std::vector<std::pair<unsigned int, unsigned int>>>
241
246 std::vector<unsigned int> n_owned_level_cells;
247
253
258
270 std::vector<AlignedVector<VectorizedArray<Number>>> weights_on_refined;
271
277 std::vector<std::vector<std::vector<unsigned short>>> dirichlet_indices;
278
287
291 template <int degree>
292 void
294 const unsigned int to_level,
297
301 template <int degree>
302 void
303 do_restrict_add(const unsigned int from_level,
306};
307
308
309
315template <int dim, typename Number, typename TransferType>
317 : public MGTransferBase<LinearAlgebra::distributed::BlockVector<Number>>
318{
319public:
323
338 virtual void
340 const unsigned int to_level,
342 const LinearAlgebra::distributed::BlockVector<Number> &src) const override;
343
344 virtual void
346 const unsigned int to_level,
348 const LinearAlgebra::distributed::BlockVector<Number> &src) const override;
349
368 virtual void
370 const unsigned int from_level,
372 const LinearAlgebra::distributed::BlockVector<Number> &src) const override;
373
384 template <typename BlockVectorType2>
385 void
387 const DoFHandler<dim> &dof_handler,
389 const BlockVectorType2 &src) const;
390
394 template <typename BlockVectorType2>
395 void
397 const std::vector<const DoFHandler<dim> *> &dof_handler,
399 const BlockVectorType2 &src) const;
400
404 template <typename BlockVectorType2>
405 void
407 const DoFHandler<dim> &dof_handler,
408 BlockVectorType2 &dst,
410 const;
411
415 template <typename BlockVectorType2>
416 void
418 const std::vector<const DoFHandler<dim> *> &dof_handler,
419 BlockVectorType2 &dst,
421 const;
422
427 static const bool supports_dof_handler_vector = true;
428
429protected:
434 virtual const TransferType &
435 get_matrix_free_transfer(const unsigned int b) const = 0;
436
441 const bool same_for_all;
442};
443
444
445
459template <int dim, typename Number>
462 Number,
463 MGTransferMatrixFree<dim, Number>>
464{
465public:
471
476 MGTransferBlockMatrixFree(const MGConstrainedDoFs &mg_constrained_dofs);
477
482 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
483
487 virtual ~MGTransferBlockMatrixFree() override = default;
488
492 void
493 initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs);
494
498 void
500 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
501
505 void
506 clear();
507
511 void
512 build(const DoFHandler<dim> &dof_handler);
513
517 void
518 build(const std::vector<const DoFHandler<dim> *> &dof_handler);
519
523 std::size_t
524 memory_consumption() const;
525
526protected:
528 get_matrix_free_transfer(const unsigned int b) const override;
529
530private:
534 std::vector<MGTransferMatrixFree<dim, Number>> matrix_free_transfer_vector;
535};
536
537
541//------------------------ templated functions -------------------------
542#ifndef DOXYGEN
543
544
545template <int dim, typename Number>
546template <typename BlockVectorType2>
547void
549 const DoFHandler<dim> &dof_handler,
551 const BlockVectorType2 &src) const
552{
553 const unsigned int min_level = dst.min_level();
554 const unsigned int max_level = dst.max_level();
555
556 Assert(max_level == dof_handler.get_triangulation().n_global_levels() - 1,
558 max_level, dof_handler.get_triangulation().n_global_levels() - 1));
559
560 const auto &fe = dof_handler.get_fe();
561
562 for (unsigned int level = min_level; level <= max_level; ++level)
563 if (dst[level].size() != dof_handler.n_dofs(level) ||
564 dst[level].locally_owned_size() !=
566 dst[level].reinit(this->vector_partitioners[level]);
567
568 // copy fine level vector to active cells in MG hierarchy
569 this->copy_to_mg(dof_handler, dst, src, true);
570
571 // FIXME: maybe need to store hanging nodes constraints per level?
572 // MGConstrainedDoFs does NOT keep this info right now, only periodicity
573 // constraints...
574
575 // do the transfer from level to level-1:
576 dst[max_level].update_ghost_values();
577
578 for (unsigned int level = max_level; level > min_level; --level)
579 {
580 // auxiliary vector which always has ghost elements
581 const LinearAlgebra::distributed::Vector<Number> *input = nullptr;
583 if (dst[level].get_partitioner().get() ==
584 this->vector_partitioners[level].get())
585 input = &dst[level];
586 else
587 {
588 ghosted_fine.reinit(this->vector_partitioners[level]);
589 ghosted_fine.copy_locally_owned_data_from(dst[level]);
590 ghosted_fine.update_ghost_values();
591 input = &ghosted_fine;
592 }
593
594 std::vector<Number> dof_values_coarse(fe.n_dofs_per_cell());
595 Vector<Number> dof_values_fine(fe.n_dofs_per_cell());
596 Vector<Number> tmp(fe.n_dofs_per_cell());
597 std::vector<types::global_dof_index> dof_indices(fe.n_dofs_per_cell());
598 for (const auto &cell : dof_handler.cell_iterators_on_level(level - 1))
599 if (cell->is_locally_owned_on_level())
600 {
601 // if we get to a cell without children (== active), we can
602 // skip it as there values should be already set by the
603 // equivalent of copy_to_mg()
604 if (cell->is_active())
605 continue;
606
607 std::fill(dof_values_coarse.begin(), dof_values_coarse.end(), 0.);
608 for (unsigned int child = 0; child < cell->n_children(); ++child)
609 {
610 cell->child(child)->get_mg_dof_indices(dof_indices);
611
613 this->mg_constrained_dofs, level, dof_indices);
614
615 for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
616 dof_values_fine(i) = (*input)(dof_indices[i]);
617 fe.get_restriction_matrix(child, cell->refinement_case())
618 .vmult(tmp, dof_values_fine);
619 for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
620 if (fe.restriction_is_additive(i))
621 dof_values_coarse[i] += tmp[i];
622 else if (tmp(i) != 0.)
623 dof_values_coarse[i] = tmp[i];
624 }
625 cell->get_mg_dof_indices(dof_indices);
626 for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
627 if (dof_handler.locally_owned_mg_dofs(level - 1).is_element(
628 dof_indices[i]))
629 dst[level - 1](dof_indices[i]) = dof_values_coarse[i];
630 }
631
632 dst[level - 1].update_ghost_values();
633 }
634
635 for (unsigned int level = min_level; level <= max_level; ++level)
636 dst[level].zero_out_ghost_values();
637}
638
639
640
641template <int dim, typename Number, typename TransferType>
642template <typename BlockVectorType2>
643void
645 const DoFHandler<dim> &dof_handler,
647 const BlockVectorType2 &src) const
648{
649 Assert(same_for_all,
651 "This object was initialized with support for usage with one "
652 "DoFHandler for each block, but this method assumes that "
653 "the same DoFHandler is used for all the blocks!"));
654 const std::vector<const DoFHandler<dim> *> mg_dofs(src.n_blocks(),
655 &dof_handler);
656
657 copy_to_mg(mg_dofs, dst, src);
658}
659
660
661
662template <int dim, typename Number, typename TransferType>
663template <typename BlockVectorType2>
664void
666 const std::vector<const DoFHandler<dim> *> &dof_handler,
668 const BlockVectorType2 &src) const
669{
670 const unsigned int n_blocks = src.n_blocks();
671 AssertDimension(dof_handler.size(), n_blocks);
672
673 if (n_blocks == 0)
674 return;
675
676 const unsigned int min_level = dst.min_level();
677 const unsigned int max_level = dst.max_level();
678
679 for (unsigned int level = min_level; level <= max_level; ++level)
680 if (dst[level].n_blocks() != n_blocks)
681 dst[level].reinit(n_blocks);
682
683 // FIXME: this a quite ugly as we need a temporary object:
685 min_level, max_level);
686
687 for (unsigned int b = 0; b < n_blocks; ++b)
688 {
689 const unsigned int data_block = same_for_all ? 0 : b;
690 get_matrix_free_transfer(data_block)
691 .copy_to_mg(*dof_handler[b], dst_non_block, src.block(b));
692
693 for (unsigned int l = min_level; l <= max_level; ++l)
694 dst[l].block(b) = dst_non_block[l];
695 }
696
697 for (unsigned int level = min_level; level <= max_level; ++level)
698 dst[level].collect_sizes();
699}
700
701template <int dim, typename Number, typename TransferType>
702template <typename BlockVectorType2>
703void
705 const DoFHandler<dim> &dof_handler,
706 BlockVectorType2 &dst,
708 const
709{
710 Assert(same_for_all,
712 "This object was initialized with support for usage with one "
713 "DoFHandler for each block, but this method assumes that "
714 "the same DoFHandler is used for all the blocks!"));
715 const std::vector<const DoFHandler<dim> *> mg_dofs(dst.n_blocks(),
716 &dof_handler);
717
718 copy_from_mg(mg_dofs, dst, src);
719}
720
721template <int dim, typename Number, typename TransferType>
722template <typename BlockVectorType2>
723void
725 const std::vector<const DoFHandler<dim> *> &dof_handler,
726 BlockVectorType2 &dst,
728 const
729{
730 const unsigned int n_blocks = dst.n_blocks();
731 AssertDimension(dof_handler.size(), n_blocks);
732
733 if (n_blocks == 0)
734 return;
735
736 const unsigned int min_level = src.min_level();
737 const unsigned int max_level = src.max_level();
738
739 for (unsigned int l = min_level; l <= max_level; ++l)
740 AssertDimension(src[l].n_blocks(), dst.n_blocks());
741
742 // FIXME: this a quite ugly as we need a temporary object:
744 min_level, max_level);
745
746 for (unsigned int b = 0; b < n_blocks; ++b)
747 {
748 for (unsigned int l = min_level; l <= max_level; ++l)
749 {
750 src_non_block[l].reinit(src[l].block(b));
751 src_non_block[l] = src[l].block(b);
752 }
753 const unsigned int data_block = same_for_all ? 0 : b;
754 get_matrix_free_transfer(data_block)
755 .copy_from_mg(*dof_handler[b], dst.block(b), src_non_block);
756 }
757}
758
759
760
761#endif // DOXYGEN
762
763
765
766#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:1934
bool is_element(const size_type index) const
Definition index_set.h:1894
void copy_locally_owned_data_from(const Vector< Number2, MemorySpace > &src)
void reinit(const size_type size, const bool omit_zeroing_entries=false)
ObserverPointer< const MGConstrainedDoFs > mg_constrained_dofs
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:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
unsigned int level
Definition grid_out.cc:4626
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::size_t size
Definition mpi.cc:734
std::enable_if_t< IsBlockVector< VectorType >::value, void > collect_sizes(VectorType &vector)
Definition operators.h:95
std::enable_if_t< IsBlockVector< VectorType >::value, unsigned int > 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 resolve_identity_constraints(const MGConstrainedDoFs *mg_constrained_dofs, const unsigned int level, std::vector< types::global_dof_index > &dof_indices)