deal.II version GIT relicensing-2657-g43ecde61a7 2025-02-18 02:30: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
trilinos_tpetra_block_vector.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2024 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_trilinos_tpetra_block_vector_h
16#define dealii_trilinos_tpetra_block_vector_h
17
18#include <deal.II/base/config.h>
19
20#ifdef DEAL_II_TRILINOS_WITH_TPETRA
21
25
27
28namespace LinearAlgebra
29{
34 namespace TpetraWrappers
35 {
54 template <typename Number, typename MemorySpace = ::MemorySpace::Host>
56 TpetraWrappers::Vector<Number, MemorySpace>>
57 {
58 public:
62 using BaseClass =
64
69
74 using pointer = typename BaseClass::pointer;
79 using iterator = typename BaseClass::iterator;
81
85 BlockVector() = default;
86
93 BlockVector(const std::vector<IndexSet> &parallel_partitioning,
94 const MPI_Comm communicator = MPI_COMM_WORLD);
95
101 BlockVector(const std::vector<IndexSet> &parallel_partitioning,
102 const std::vector<IndexSet> &ghost_values,
103 const MPI_Comm communicator,
104 const bool vector_writable = false);
105
111
117 explicit BlockVector(const size_type num_blocks);
118
122 ~BlockVector() override = default;
123
129 operator=(const Number s);
130
143
152 void
153 reinit(const std::vector<IndexSet> &parallel_partitioning,
154 const MPI_Comm communicator = MPI_COMM_WORLD,
155 const bool omit_zeroing_entries = false);
156
174 void
175 reinit(const std::vector<IndexSet> &partitioning,
176 const std::vector<IndexSet> &ghost_values,
177 const MPI_Comm communicator = MPI_COMM_WORLD,
178 const bool vector_writable = false);
179
194 void
196 const bool omit_zeroing_entries = false);
197
204 void
205 reinit(const size_type num_blocks);
206
213 bool
215
219 void
220 print(std::ostream &out,
221 const unsigned int precision = 3,
222 const bool scientific = true,
223 const bool across = true) const;
224 };
225 } // namespace TpetraWrappers
226
229} // namespace LinearAlgebra
230
232
233#endif // DEAL_II_TRILINOS_WITH_TPETRA
234
235#endif // dealii_trilinos_tpetra_block_vector_h
::internal::BlockVectorIterators::Iterator< BlockVectorBase, false > iterator
const value_type * const_pointer
typename BlockType::const_reference const_reference
types::global_dof_index size_type
typename BlockType::value_type value_type
::internal::BlockVectorIterators::Iterator< BlockVectorBase, true > const_iterator
typename BlockType::reference reference
BlockVector< Number, MemorySpace > & operator=(const Number s)
BlockVector< Number, MemorySpace > & operator=(const BlockVector< Number, MemorySpace > &v)
BlockVector(const BlockVector< Number, MemorySpace > &v)
BlockVector(const size_type num_blocks)
void reinit(const BlockVector< Number, MemorySpace > &V, const bool omit_zeroing_entries=false)
void reinit(const std::vector< IndexSet > &partitioning, const std::vector< IndexSet > &ghost_values, const MPI_Comm communicator=MPI_COMM_WORLD, const bool vector_writable=false)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
void reinit(const size_type num_blocks)
typename BaseClass::const_reference const_reference
BlockVector(const std::vector< IndexSet > &parallel_partitioning, const std::vector< IndexSet > &ghost_values, const MPI_Comm communicator, const bool vector_writable=false)
BlockVector(const std::vector< IndexSet > &parallel_partitioning, const MPI_Comm communicator=MPI_COMM_WORLD)
void reinit(const std::vector< IndexSet > &parallel_partitioning, const MPI_Comm communicator=MPI_COMM_WORLD, const bool omit_zeroing_entries=false)
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:518
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:519