Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30: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
petsc_parallel_block_vector.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2004 - 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
16
17#ifdef DEAL_II_WITH_PETSC
18
20
22
23namespace PETScWrappers
24{
25 namespace MPI
26 {
27 // Check that the class we declare here satisfies the
28 // vector-space-vector concept. If we catch it here,
29 // any mistake in the vector class declaration would
30 // show up in uses of this class later on as well.
31# ifdef DEAL_II_HAVE_CXX20
33# endif
34
36
38 {
39 PetscErrorCode ierr = VecDestroy(&petsc_nest_vector);
40 AssertThrow(ierr == 0, ExcPETScError(ierr));
41 }
42
43 void
44 BlockVector::reinit(const unsigned int num_blocks)
45 {
46 std::vector<size_type> block_sizes(num_blocks, 0);
47 this->block_indices.reinit(block_sizes);
48 if (this->components.size() != this->n_blocks())
49 this->components.resize(this->n_blocks());
50
51 for (unsigned int i = 0; i < this->n_blocks(); ++i)
52 components[i].reinit(MPI_COMM_SELF, 0, 0);
53
55 }
56
57 void
59 {
60 PetscBool isnest;
61
62 PetscErrorCode ierr =
63 PetscObjectTypeCompare(reinterpret_cast<PetscObject>(v),
64 VECNEST,
65 &isnest);
66 AssertThrow(ierr == 0, ExcPETScError(ierr));
67 std::vector<Vec> sv;
68 if (isnest)
69 {
70 PetscInt nb;
71 ierr = VecNestGetSize(v, &nb);
72 AssertThrow(ierr == 0, ExcPETScError(ierr));
73 for (PetscInt i = 0; i < nb; ++i)
74 {
75 Vec vv;
76 ierr = VecNestGetSubVec(v, i, &vv);
77 sv.push_back(vv);
78 }
79 }
80 else
81 {
82 sv.push_back(v);
83 }
84
85 auto nb = sv.size();
86
87 std::vector<size_type> block_sizes(nb, 0);
88 this->block_indices.reinit(block_sizes);
89
90 this->components.resize(nb);
91 for (unsigned int i = 0; i < nb; ++i)
92 {
93 this->components[i].reinit(sv[i]);
94 }
95
97 if (!isnest)
99 else
100 {
101 ierr = PetscObjectReference(reinterpret_cast<PetscObject>(v));
102 AssertThrow(ierr == 0, ExcPETScError(ierr));
103 PetscErrorCode ierr = VecDestroy(&petsc_nest_vector);
104 AssertThrow(ierr == 0, ExcPETScError(ierr));
106 }
107 }
108
109 Vec &
114
115 BlockVector::operator const Vec &() const
116 {
117 return petsc_nest_vector;
118 }
119
120 void
126
127 void
133
134 void
136 {
137 PetscErrorCode ierr = VecDestroy(&petsc_nest_vector);
138 AssertThrow(ierr == 0, ExcPETScError(ierr));
139
140 auto n = this->n_blocks();
141
142 std::vector<Vec> pcomponents(n);
143 for (unsigned int i = 0; i < n; i++)
144 pcomponents[i] = this->components[i].petsc_vector();
145
146 MPI_Comm comm =
147 pcomponents.size() > 0 ?
148 PetscObjectComm(reinterpret_cast<PetscObject>(pcomponents[0])) :
149 PETSC_COMM_SELF;
150
151 ierr =
152 VecCreateNest(comm, n, nullptr, pcomponents.data(), &petsc_nest_vector);
153 AssertThrow(ierr == 0, ExcPETScError(ierr));
154 }
155 } // namespace MPI
156
157} // namespace PETScWrappers
158
160
161#endif // DEAL_II_WITH_PETSC
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
unsigned int n_blocks() const
void compress(VectorOperation::values operation)
void collect_sizes()
std::vector< Vector > components
void compress(VectorOperation::values operation)
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
#define AssertThrow(cond, exc)
void petsc_increment_state_counter(Vec v)
unsigned int global_dof_index
Definition types.h:81
*braid_SplitCommworld & comm