deal.II version GIT relicensing-2206-gaa53ff9447 2024-12-02 09: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
trilinos_epetra_vector.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2016 - 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_epetra_vector_h
16#define dealii_trilinos_epetra_vector_h
17
18
19#include <deal.II/base/config.h>
20
21#ifdef DEAL_II_WITH_TRILINOS
22
25
30
31# include <Epetra_FEVector.h>
32
33# include <memory>
34
36
37namespace LinearAlgebra
38{
39 // Forward declaration
40 template <typename Number>
41 class ReadWriteVector;
42
50 namespace EpetraWrappers
51 {
57 {
58 public:
59 using value_type = double;
61 };
62
73 namespace internal
74 {
84 class VectorReference
85 {
86 private:
87 using value_type = VectorTraits::value_type;
88 using size_type = VectorTraits::size_type;
89
94 VectorReference(Vector &vector, const size_type index);
95
96 public:
100 VectorReference(const VectorReference &) = default;
101
113 const VectorReference &
114 operator=(const VectorReference &r) const;
115
119 VectorReference &
120 operator=(const VectorReference &r);
121
125 const VectorReference &
126 operator=(const value_type &s) const;
127
131 const VectorReference &
132 operator+=(const value_type &s) const;
133
137 const VectorReference &
138 operator-=(const value_type &s) const;
139
143 const VectorReference &
144 operator*=(const value_type &s) const;
145
149 const VectorReference &
150 operator/=(const value_type &s) const;
151
156 operator value_type() const;
157
161 DeclException1(ExcTrilinosError,
162 int,
163 << "An error with error number " << arg1
164 << " occurred while calling a Trilinos function");
165
166 /*
167 * Access to a an element that is not (locally-)owned.
168 *
169 * @ingroup Exceptions
170 */
172 ExcAccessToNonLocalElement,
173 size_type,
174 size_type,
175 size_type,
176 size_type,
177 << "You are trying to access element " << arg1
178 << " of a distributed vector, but this element is not stored "
179 << "on the current processor. Note: There are " << arg2
180 << " elements stored "
181 << "on the current processor from within the range [" << arg3 << ','
182 << arg4 << "] but Trilinos vectors need not store contiguous "
183 << "ranges on each processor, and not every element in "
184 << "this range may in fact be stored locally."
185 << "\n\n"
186 << "A common source for this kind of problem is that you "
187 << "are passing a 'fully distributed' vector into a function "
188 << "that needs read access to vector elements that correspond "
189 << "to degrees of freedom on ghost cells (or at least to "
190 << "'locally active' degrees of freedom that are not also "
191 << "'locally owned'). You need to pass a vector that has these "
192 << "elements as ghost entries.");
193
194 private:
198 Vector &vector;
199
203 const size_type index;
204
205 // Make the vector class a friend, so that it can create objects of the
206 // present type.
207 friend class ::LinearAlgebra::EpetraWrappers::Vector;
208 }; // class VectorReference
209
210 } // namespace internal
224 class Vector : public ReadVector<VectorTraits::value_type>
225 {
226 public:
232
241 Vector();
242
247 Vector(const Vector &V);
248
255 explicit Vector(const IndexSet &parallel_partitioner,
256 const MPI_Comm communicator);
257
264 void
266 const MPI_Comm communicator,
267 const bool omit_zeroing_entries = false);
268
273 void
274 reinit(const Vector &V, const bool omit_zeroing_entries = false);
275
279 virtual void
282 const ArrayView<double> &elements) const override;
283
313 Vector &
314 operator=(const Vector &V);
315
320 Vector &
321 operator=(const double s);
322
331 void
335 const std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
337
342 void
343 import(
346 const std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
348 {
350 }
351
368 operator()(const size_type index);
369
378 operator()(const size_type index) const;
379
386 operator[](const size_type index);
387
394 operator[](const size_type index) const;
395
406 Vector &
407 operator*=(const double factor);
408
412 Vector &
413 operator/=(const double factor);
414
418 Vector &
419 operator+=(const Vector &V);
420
424 Vector &
425 operator-=(const Vector &V);
426
431 double
432 operator*(const Vector &V) const;
433
437 void
438 add(const double a);
439
444 void
445 add(const double a, const Vector &V);
446
451 void
452 add(const double a, const Vector &V, const double b, const Vector &W);
453
458 void
459 sadd(const double s, const double a, const Vector &V);
460
467 void
469
473 void
474 equ(const double a, const Vector &V);
475
479 bool
480 all_zero() const;
481
493 double
494 mean_value() const;
495
500 double
501 l1_norm() const;
502
507 double
508 l2_norm() const;
509
514 double
515 linfty_norm() const;
516
539 double
540 add_and_dot(const double a, const Vector &V, const Vector &W);
541
553 bool
554 has_ghost_elements() const;
555
560 virtual size_type
561 size() const override;
562
568 locally_owned_size() const;
569
574 get_mpi_communicator() const;
575
589
608 void
610
615 const Epetra_FEVector &
616 trilinos_vector() const;
617
622 Epetra_FEVector &
624
628 void
629 print(std::ostream &out,
630 const unsigned int precision = 3,
631 const bool scientific = true,
632 const bool across = true) const;
633
637 std::size_t
638 memory_consumption() const;
639
647
654
661 int,
662 << "An error with error number " << arg1
663 << " occurred while calling a Trilinos function");
664
665 private:
671 void
673 const MPI_Comm mpi_comm);
674
678 std::unique_ptr<Epetra_FEVector> vector;
679
684
689 std::shared_ptr<const CommunicationPattern> epetra_comm_pattern;
690
691 // Make the reference class a friend.
693 };
694
695# ifndef DOXYGEN
696
697 // VectorReference
698 namespace internal
699 {
700 inline VectorReference::VectorReference(Vector &vector,
701 const size_type index)
702 : vector(vector)
703 , index(index)
704 {}
705
706
707
708 inline const VectorReference &
709 VectorReference::operator=(const VectorReference &r) const
710 {
711 // as explained in the class
712 // documentation, this is not the copy
713 // operator. so simply pass on to the
714 // "correct" assignment operator
715 *this = static_cast<value_type>(r);
716
717 return *this;
718 }
719
720
721
722 inline VectorReference &
723 VectorReference::operator=(const VectorReference &r)
724 {
725 // as above
726 *this = static_cast<value_type>(r);
727
728 return *this;
729 }
730
731
732
733 inline const VectorReference &
734 VectorReference::operator=(const value_type &value) const
735 {
736 (*vector.vector)[0][index] = value;
737
738 return *this;
739 }
740
741
742
743 inline const VectorReference &
744 VectorReference::operator+=(const value_type &value) const
745 {
746 value_type new_value = static_cast<value_type>(*this) + value;
747 (*vector.vector)[0][index] = new_value;
748
749 return *this;
750 }
751
752
753
754 inline const VectorReference &
755 VectorReference::operator-=(const value_type &value) const
756 {
757 value_type new_value = static_cast<value_type>(*this) - value;
758 (*vector.vector)[0][index] = new_value;
759
760 return *this;
761 }
762
763
764
765 inline const VectorReference &
766 VectorReference::operator*=(const value_type &value) const
767 {
768 value_type new_value = static_cast<value_type>(*this) * value;
769 (*vector.vector)[0][index] = new_value;
770
771 return *this;
772 }
773
774
775
776 inline const VectorReference &
777 VectorReference::operator/=(const value_type &value) const
778 {
779 value_type new_value = static_cast<value_type>(*this) / value;
780 (*vector.vector)[0][index] = new_value;
781
782 return *this;
783 }
784 } // namespace internal
785
786# endif /* DOXYGEN */
787
788
789 inline internal::VectorReference
791 {
792 return internal::VectorReference(*this, index);
793 }
794
795 inline internal::VectorReference
797 {
798 return operator()(index);
799 }
800
801 inline Vector::value_type
802 Vector::operator[](const size_type index) const
803 {
804 return operator()(index);
805 }
806
807
808 inline bool
810 {
811 return false;
812 }
813 } // namespace EpetraWrappers
814} // namespace LinearAlgebra
815
816
820template <>
821struct is_serial_vector<LinearAlgebra::EpetraWrappers::Vector> : std::false_type
822{};
823
825
826#endif
827
828#endif
void equ(const double a, const Vector &V)
void compress(const VectorOperation::values operation)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
std::unique_ptr< Epetra_FEVector > vector
const Epetra_FEVector & trilinos_vector() const
void import_elements(const ReadWriteVector< double > &V, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
void sadd(const double s, const double a, const Vector &V)
value_type operator()(const size_type index) const
virtual size_type size() const override
reference operator[](const size_type index)
void scale(const Vector &scaling_factors)
void create_epetra_comm_pattern(const IndexSet &source_index_set, const MPI_Comm mpi_comm)
std::shared_ptr< const CommunicationPattern > epetra_comm_pattern
void reinit(const IndexSet &parallel_partitioner, const MPI_Comm communicator, const bool omit_zeroing_entries=false)
reference operator()(const size_type index)
virtual void extract_subvector_to(const ArrayView< const types::global_dof_index > &indices, const ArrayView< double > &elements) const override
const internal::VectorReference const_reference
double add_and_dot(const double a, const Vector &V, const Vector &W)
bool has_ghost_elements() const
Number operator[](const size_type i) const
Number operator()(const size_type i) const
#define DEAL_II_DEPRECATED
Definition config.h:205
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
static ::ExceptionBase & ExcVectorTypeNotCompatible()
#define DeclException0(Exception0)
Definition exceptions.h:466
static ::ExceptionBase & ExcTrilinosError(int arg1)
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
Definition exceptions.h:580
static ::ExceptionBase & ExcDifferentParallelPartitioning()
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:511
unsigned int global_dof_index
Definition types.h:81