Line data Source code
1 : // The libMesh Finite Element Library. 2 : // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 3 : 4 : // This library is free software; you can redistribute it and/or 5 : // modify it under the terms of the GNU Lesser General Public 6 : // License as published by the Free Software Foundation; either 7 : // version 2.1 of the License, or (at your option) any later version. 8 : 9 : // This library is distributed in the hope that it will be useful, 10 : // but WITHOUT ANY WARRANTY; without even the implied warranty of 11 : // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 12 : // Lesser General Public License for more details. 13 : 14 : // You should have received a copy of the GNU Lesser General Public 15 : // License along with this library; if not, write to the Free Software 16 : // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 17 : 18 : #include "libmesh/lumped_mass_matrix.h" 19 : #include "libmesh/numeric_vector.h" 20 : 21 : // C++ Includes 22 : #include <cmath> 23 : #include <memory> 24 : 25 : namespace libMesh 26 : { 27 : 28 : template <typename T> 29 213 : LumpedMassMatrix<T>::LumpedMassMatrix(const Parallel::Communicator & comm_in) 30 213 : : DiagonalMatrix<T>(comm_in) 31 : { 32 79 : } 33 : 34 : template <typename T> 35 : std::unique_ptr<SparseMatrix<T>> 36 71 : LumpedMassMatrix<T>::zero_clone() const 37 : { 38 : // Make empty copy with matching comm 39 73 : auto mat_copy = std::make_unique<LumpedMassMatrix<T>>(this->comm()); 40 : 41 : // Initialize copy with our same nonzero structure, and explicitly 42 : // zero values using fast == false. 43 71 : mat_copy->init(*this, /*fast=*/false); 44 : 45 73 : return mat_copy; 46 67 : } 47 : 48 : template <typename T> 49 : std::unique_ptr<SparseMatrix<T>> 50 71 : LumpedMassMatrix<T>::clone() const 51 : { 52 : // Make empty copy with matching comm 53 73 : auto mat_copy = std::make_unique<LumpedMassMatrix<T>>(this->comm()); 54 : 55 : // Make copy of our diagonal 56 73 : auto diag_copy = this->_diagonal->clone(); 57 : 58 : // Swap diag_copy with diagonal in mat_copy 59 4 : *mat_copy = std::move(*diag_copy); 60 : 61 73 : return mat_copy; 62 67 : } 63 : 64 : template <typename T> 65 : void 66 453 : LumpedMassMatrix<T>::set(const numeric_index_type i, 67 : const numeric_index_type libmesh_dbg_var(j), 68 : const T value) 69 : { 70 5 : libmesh_assert_msg(i == j, "Set in a lumped mass matrix really only makes sense for i == j"); 71 463 : this->_diagonal->set(i, std::abs(value)); 72 453 : } 73 : 74 : template <typename T> 75 : void 76 3911 : LumpedMassMatrix<T>::add(const numeric_index_type i, const numeric_index_type /*j*/, const T value) 77 : { 78 3924 : this->_diagonal->add(i, std::abs(value)); 79 3911 : } 80 : 81 : template <typename T> 82 : void 83 0 : LumpedMassMatrix<T>::add(const T a, const SparseMatrix<T> & X) 84 : { 85 0 : if (dynamic_cast<const LumpedMassMatrix<T> *>(&X)) 86 : { 87 0 : auto x_diagonal = this->_diagonal->zero_clone(); 88 0 : X.get_diagonal(*x_diagonal); 89 0 : this->_diagonal->add(a, *x_diagonal); 90 0 : } 91 : else 92 0 : libmesh_error_msg("Unsupported matrix type passed to LumpedMassMatrix::add"); 93 0 : } 94 : 95 : template <typename T> 96 : LumpedMassMatrix<T> & 97 0 : LumpedMassMatrix<T>::operator=(const NumericVector<T> & vec) 98 : { 99 0 : *this->_diagonal = vec; 100 0 : return *this; 101 : } 102 : 103 : template <typename T> 104 : LumpedMassMatrix<T> & 105 4 : LumpedMassMatrix<T>::operator=(NumericVector<T> && vec) 106 : { 107 71 : this->_diagonal->swap(vec); 108 71 : return *this; 109 : } 110 : 111 : template class LIBMESH_EXPORT LumpedMassMatrix<Number>; 112 : 113 : } // namespace libMesh