From: Wolfgang Bangerth Date: Thu, 6 Jul 2023 13:02:54 +0000 (-0600) Subject: Remove derivation of LinearAlgebra::distributed::Vector from VectorSpaceVector. X-Git-Tag: relicensing~729^2~1 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=a91b3f4d9d10ae6d8934e1e41dbfd0476fc2ca14;p=dealii.git Remove derivation of LinearAlgebra::distributed::Vector from VectorSpaceVector. --- diff --git a/include/deal.II/lac/la_parallel_vector.h b/include/deal.II/lac/la_parallel_vector.h index 0d6dc96e31..6b117fd4a7 100644 --- a/include/deal.II/lac/la_parallel_vector.h +++ b/include/deal.II/lac/la_parallel_vector.h @@ -248,9 +248,7 @@ namespace LinearAlgebra * @endcode */ template - class Vector : public ::dealii::LinearAlgebra::VectorSpaceVector, - public ::dealii::ReadVector, - public Subscriptor + class Vector : public ::dealii::ReadVector, public Subscriptor { public: using memory_space = MemorySpace; @@ -338,7 +336,7 @@ namespace LinearAlgebra /** * Destructor. */ - virtual ~Vector() override; + ~Vector(); /** * Set the global size of the vector to @p size without any actual @@ -537,8 +535,8 @@ namespace LinearAlgebra * with VectorOperation::max two times subsequently, the maximal value * after the second calculation will be zero. */ - virtual void - compress(VectorOperation::values operation) override; + void + compress(VectorOperation::values operation); /** * Fills the data field for ghost indices with the values stored in the @@ -710,41 +708,33 @@ namespace LinearAlgebra /** @} */ /** - * @name 3: Implementation of VectorSpaceVector + * @name 3: Implementation of vector space operations */ /** @{ */ - /** - * Change the dimension to that of the vector V. The elements of V are not - * copied. - */ - virtual void - reinit(const VectorSpaceVector &V, - const bool omit_zeroing_entries = false) override; - /** * Multiply the entire vector by a fixed factor. */ - virtual Vector & - operator*=(const Number factor) override; + Vector & + operator*=(const Number factor); /** * Divide the entire vector by a fixed factor. */ - virtual Vector & - operator/=(const Number factor) override; + Vector & + operator/=(const Number factor); /** * Add the vector @p V to the present one. */ - virtual Vector & - operator+=(const VectorSpaceVector &V) override; + Vector & + operator+=(const Vector &V); /** * Subtract the vector @p V from the present one. */ - virtual Vector & - operator-=(const VectorSpaceVector &V) override; + Vector & + operator-=(const Vector &V); /** * Import all the elements present in the vector's IndexSet from the input @@ -757,21 +747,21 @@ namespace LinearAlgebra * @note If the MemorySpace is Default, the data in the ReadWriteVector will * be moved to the @ref GlossDevice "device". */ - virtual void + void import_elements( const LinearAlgebra::ReadWriteVector &V, VectorOperation::values operation, std::shared_ptr - communication_pattern = {}) override; + communication_pattern = {}); /** * @deprecated Use import_elements() instead. */ - DEAL_II_DEPRECATED virtual void + DEAL_II_DEPRECATED void import(const LinearAlgebra::ReadWriteVector &V, VectorOperation::values operation, std::shared_ptr - communication_pattern = {}) override + communication_pattern = {}) { import_elements(V, operation, communication_pattern); } @@ -779,35 +769,35 @@ namespace LinearAlgebra /** * Return the scalar product of two vectors. */ - virtual Number - operator*(const VectorSpaceVector &V) const override; + Number + operator*(const Vector &V) const; /** * Add @p a to all components. Note that @p a is a scalar not a vector. */ - virtual void - add(const Number a) override; + void + add(const Number a); /** * Simple addition of a multiple of a vector, i.e. *this += a*V. */ - virtual void - add(const Number a, const VectorSpaceVector &V) override; + void + add(const Number a, const Vector &V); /** * Multiple addition of scaled vectors, i.e. *this += a*V+b*W. */ - virtual void - add(const Number a, - const VectorSpaceVector &V, - const Number b, - const VectorSpaceVector &W) override; + void + add(const Number a, + const Vector &V, + const Number b, + const Vector &W); /** * A collective add operation: This function adds a whole set of values * stored in @p values to the vector components specified by @p indices. */ - virtual void + void add(const std::vector &indices, const std::vector & values); @@ -815,38 +805,38 @@ namespace LinearAlgebra * Scaling and simple addition of a multiple of a vector, i.e. *this = * s*(*this)+a*V. */ - virtual void - sadd(const Number s, - const Number a, - const VectorSpaceVector &V) override; + void + sadd(const Number s, + const Number a, + const Vector &V); /** * Scale each element of this vector by the corresponding element in the * argument. This function is mostly meant to simulate multiplication (and * immediate re-assignment) by a diagonal scaling matrix. */ - virtual void - scale(const VectorSpaceVector &scaling_factors) override; + void + scale(const Vector &scaling_factors); /** * Assignment *this = a*V. */ - virtual void - equ(const Number a, const VectorSpaceVector &V) override; + void + equ(const Number a, const Vector &V); /** * Return the l1 norm of the vector (i.e., the sum of the * absolute values of all entries among all processors). */ - virtual real_type - l1_norm() const override; + real_type + l1_norm() const; /** * Return the $l_2$ norm of the vector (i.e., the square root of * the sum of the square of all entries among all processors). */ - virtual real_type - l2_norm() const override; + real_type + l2_norm() const; /** * Return the square of the $l_2$ norm of the vector. @@ -858,8 +848,8 @@ namespace LinearAlgebra * Return the maximum norm of the vector (i.e., the maximum absolute value * among all entries and among all processors). */ - virtual real_type - linfty_norm() const override; + real_type + linfty_norm() const; /** * Perform a combined operation of a vector addition and a subsequent @@ -881,10 +871,10 @@ namespace LinearAlgebra * implemented as * $\left=\sum_i v_i \bar{w_i}$. */ - virtual Number - add_and_dot(const Number a, - const VectorSpaceVector &V, - const VectorSpaceVector &W) override; + Number + add_and_dot(const Number a, + const Vector &V, + const Vector &W); /** * Return the global size of the vector, equal to the sum of the number of @@ -904,27 +894,27 @@ namespace LinearAlgebra * vec.locally_owned_elements() == complete_index_set(vec.size()) * @endcode */ - virtual dealii::IndexSet - locally_owned_elements() const override; + dealii::IndexSet + locally_owned_elements() const; /** * Print the vector to the output stream @p out. */ - virtual void + void print(std::ostream & out, const unsigned int precision = 3, const bool scientific = true, - const bool across = true) const override; + const bool across = true) const; /** * Return the memory consumption of this class in bytes. */ - virtual std::size_t - memory_consumption() const override; + std::size_t + memory_consumption() const; /** @} */ /** - * @name 4: Other vector operations not included in VectorSpaceVector + * @name 4: Other vector operations */ /** @{ */ @@ -933,8 +923,8 @@ namespace LinearAlgebra * zero, also ghost elements are set to zero, otherwise they remain * unchanged. */ - virtual Vector & - operator=(const Number s) override; + Vector & + operator=(const Number s); /** * This is a collective add operation that adds a whole set of values @@ -1168,14 +1158,14 @@ namespace LinearAlgebra * This is a collective operation. This function is expensive, because * potentially all elements have to be checked. */ - virtual bool - all_zero() const override; + bool + all_zero() const; /** * Compute the mean value of all the entries in the vector. */ - virtual Number - mean_value() const override; + Number + mean_value() const; /** * $l_p$-norm of the vector. The pth root of the sum of the pth powers @@ -1301,16 +1291,16 @@ namespace LinearAlgebra * without MPI communication. */ void - add_local(const Number a, const VectorSpaceVector &V); + add_local(const Number a, const Vector &V); /** * Scaling and simple addition of a multiple of a vector, i.e. *this = * s*(*this)+a*V without MPI communication. */ void - sadd_local(const Number s, - const Number a, - const VectorSpaceVector &V); + sadd_local(const Number s, + const Number a, + const Vector &V); /** * Local part of the inner product of two vectors. diff --git a/include/deal.II/lac/la_parallel_vector.templates.h b/include/deal.II/lac/la_parallel_vector.templates.h index 967a6bbb25..84bbd3b1f3 100644 --- a/include/deal.II/lac/la_parallel_vector.templates.h +++ b/include/deal.II/lac/la_parallel_vector.templates.h @@ -1423,33 +1423,11 @@ namespace LinearAlgebra - template - void - Vector::reinit(const VectorSpaceVector &V, - const bool omit_zeroing_entries) - { - // Downcast. Throws an exception if invalid. - using VectorType = Vector; - Assert(dynamic_cast(&V) != nullptr, - ExcVectorTypeNotCompatible()); - const VectorType &down_V = dynamic_cast(V); - - reinit(down_V, omit_zeroing_entries); - } - - - template Vector & Vector::operator+=( - const VectorSpaceVector &vv) + const Vector &v) { - // Downcast. Throws an exception if invalid. - using VectorType = Vector; - Assert(dynamic_cast(&vv) != nullptr, - ExcVectorTypeNotCompatible()); - const VectorType &v = dynamic_cast(vv); - AssertDimension(locally_owned_size(), v.locally_owned_size()); dealii::internal::VectorOperations:: @@ -1470,14 +1448,8 @@ namespace LinearAlgebra template Vector & Vector::operator-=( - const VectorSpaceVector &vv) + const Vector &v) { - // Downcast. Throws an exception if invalid. - using VectorType = Vector; - Assert(dynamic_cast(&vv) != nullptr, - ExcVectorTypeNotCompatible()); - const VectorType &v = dynamic_cast(vv); - AssertDimension(locally_owned_size(), v.locally_owned_size()); dealii::internal::VectorOperations:: @@ -1514,15 +1486,9 @@ namespace LinearAlgebra template void Vector::add_local( - const Number a, - const VectorSpaceVector &vv) + const Number a, + const Vector &v) { - // Downcast. Throws an exception if invalid. - using VectorType = Vector; - Assert(dynamic_cast(&vv) != nullptr, - ExcVectorTypeNotCompatible()); - const VectorType &v = dynamic_cast(vv); - AssertIsFinite(a); AssertDimension(locally_owned_size(), v.locally_owned_size()); @@ -1543,8 +1509,9 @@ namespace LinearAlgebra template void - Vector::add(const Number a, - const VectorSpaceVector &vv) + Vector::add( + const Number a, + const Vector &vv) { add_local(a, vv); @@ -1556,20 +1523,12 @@ namespace LinearAlgebra template void - Vector::add(const Number a, - const VectorSpaceVector &vv, - const Number b, - const VectorSpaceVector &ww) - { - // Downcast. Throws an exception if invalid. - using VectorType = Vector; - Assert(dynamic_cast(&vv) != nullptr, - ExcVectorTypeNotCompatible()); - const VectorType &v = dynamic_cast(vv); - Assert(dynamic_cast(&ww) != nullptr, - ExcVectorTypeNotCompatible()); - const VectorType &w = dynamic_cast(ww); - + Vector::add( + const Number a, + const Vector &v, + const Number b, + const Vector &w) + { AssertIsFinite(a); AssertIsFinite(b); @@ -1631,16 +1590,10 @@ namespace LinearAlgebra template void Vector::sadd_local( - const Number x, - const Number a, - const VectorSpaceVector &vv) + const Number x, + const Number a, + const Vector &v) { - // Downcast. Throws an exception if invalid. - using VectorType = Vector; - Assert((dynamic_cast(&vv) != nullptr), - ExcVectorTypeNotCompatible()); - const VectorType &v = dynamic_cast(vv); - AssertIsFinite(x); AssertIsFinite(a); AssertDimension(locally_owned_size(), v.locally_owned_size()); @@ -1659,11 +1612,12 @@ namespace LinearAlgebra template void - Vector::sadd(const Number x, - const Number a, - const VectorSpaceVector &vv) + Vector::sadd( + const Number x, + const Number a, + const Vector &v) { - sadd_local(x, a, vv); + sadd_local(x, a, v); if (vector_is_ghosted) update_ghost_values(); @@ -1704,14 +1658,9 @@ namespace LinearAlgebra template void - Vector::scale(const VectorSpaceVector &vv) + Vector::scale( + const Vector &v) { - // Downcast. Throws an exception if invalid. - using VectorType = Vector; - Assert(dynamic_cast(&vv) != nullptr, - ExcVectorTypeNotCompatible()); - const VectorType &v = dynamic_cast(vv); - AssertDimension(locally_owned_size(), v.locally_owned_size()); dealii::internal::VectorOperations:: @@ -1726,15 +1675,10 @@ namespace LinearAlgebra template void - Vector::equ(const Number a, - const VectorSpaceVector &vv) + Vector::equ( + const Number a, + const Vector &v) { - // Downcast. Throws an exception if invalid. - using VectorType = Vector; - Assert(dynamic_cast(&vv) != nullptr, - ExcVectorTypeNotCompatible()); - const VectorType &v = dynamic_cast(vv); - AssertIsFinite(a); AssertDimension(locally_owned_size(), v.locally_owned_size()); @@ -1803,14 +1747,8 @@ namespace LinearAlgebra template Number Vector::operator*( - const VectorSpaceVector &vv) const + const Vector &v) const { - // Downcast. Throws an exception if invalid. - using VectorType = Vector; - Assert((dynamic_cast(&vv) != nullptr), - ExcVectorTypeNotCompatible()); - const VectorType &v = dynamic_cast(vv); - Number local_result = inner_product_local(v); if (partitioner->n_mpi_processes() > 1) return Utilities::MPI::sum(local_result, @@ -2021,19 +1959,10 @@ namespace LinearAlgebra template Number Vector::add_and_dot( - const Number a, - const VectorSpaceVector &vv, - const VectorSpaceVector &ww) - { - // Downcast. Throws an exception if invalid. - using VectorType = Vector; - Assert((dynamic_cast(&vv) != nullptr), - ExcVectorTypeNotCompatible()); - const VectorType &v = dynamic_cast(vv); - Assert((dynamic_cast(&ww) != nullptr), - ExcVectorTypeNotCompatible()); - const VectorType &w = dynamic_cast(ww); - + const Number a, + const Vector &v, + const Vector &w) + { Number local_result = add_and_dot_local(a, v, w); if (partitioner->n_mpi_processes() > 1) return Utilities::MPI::sum(local_result,