$extrastylesheet
trilinos_epetra_matrix.C
Go to the documentation of this file.
00001 // The libMesh Finite Element Library.
00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
00003 
00004 // This library is free software; you can redistribute it and/or
00005 // modify it under the terms of the GNU Lesser General Public
00006 // License as published by the Free Software Foundation; either
00007 // version 2.1 of the License, or (at your option) any later version.
00008 
00009 // This library is distributed in the hope that it will be useful,
00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00012 // Lesser General Public License for more details.
00013 
00014 // You should have received a copy of the GNU Lesser General Public
00015 // License along with this library; if not, write to the Free Software
00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00017 
00018 
00019 // C++ includes
00020 #include "libmesh/libmesh_config.h"
00021 
00022 #ifdef LIBMESH_HAVE_TRILINOS
00023 
00024 // Local includes
00025 #include "libmesh/trilinos_epetra_matrix.h"
00026 #include "libmesh/trilinos_epetra_vector.h"
00027 #include "libmesh/dof_map.h"
00028 #include "libmesh/dense_matrix.h"
00029 #include "libmesh/parallel.h"
00030 #include "libmesh/sparsity_pattern.h"
00031 
00032 namespace libMesh
00033 {
00034 
00035 
00036 
00037 //-----------------------------------------------------------------------
00038 //EpetraMatrix members
00039 
00040 template <typename T>
00041 void EpetraMatrix<T>::update_sparsity_pattern (const SparsityPattern::Graph &sparsity_pattern)
00042 {
00043   // clear data, start over
00044   this->clear ();
00045 
00046   // big trouble if this fails!
00047   libmesh_assert(this->_dof_map);
00048 
00049   const numeric_index_type n_rows = sparsity_pattern.size();
00050 
00051   const numeric_index_type m   = this->_dof_map->n_dofs();
00052   const numeric_index_type n   = m;
00053   const numeric_index_type n_l = this->_dof_map->n_dofs_on_processor(this->processor_id());
00054   const numeric_index_type m_l = n_l;
00055 
00056   // error checking
00057 #ifndef NDEBUG
00058   {
00059     libmesh_assert_equal_to (n, m);
00060     libmesh_assert_equal_to (n_l, m_l);
00061 
00062     numeric_index_type
00063       summed_m_l = m_l,
00064       summed_n_l = n_l;
00065 
00066     this->comm().sum (summed_m_l);
00067     this->comm().sum (summed_n_l);
00068 
00069     libmesh_assert_equal_to (m, summed_m_l);
00070     libmesh_assert_equal_to (n, summed_n_l);
00071   }
00072 #endif
00073 
00074   // build a map defining the data distribution
00075   _map = new Epetra_Map (static_cast<int>(m),
00076                          m_l,
00077                          0,
00078                          Epetra_MpiComm (this->comm().get()));
00079 
00080   libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->NumGlobalPoints()), m);
00081   libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->MaxAllGID()+1), m);
00082 
00083   const std::vector<numeric_index_type>& n_nz = this->_dof_map->get_n_nz();
00084   const std::vector<numeric_index_type>& n_oz = this->_dof_map->get_n_oz();
00085 
00086   // Make sure the sparsity pattern isn't empty
00087   libmesh_assert_equal_to (n_nz.size(), n_l);
00088   libmesh_assert_equal_to (n_oz.size(), n_l);
00089 
00090   // Epetra wants the total number of nonzeros, both local and remote.
00091   std::vector<int> n_nz_tot;  n_nz_tot.reserve(n_nz.size());
00092 
00093   for (numeric_index_type i=0; i<n_nz.size(); i++)
00094     n_nz_tot.push_back(std::min(n_nz[i] + n_oz[i], n));
00095 
00096   if (m==0)
00097     return;
00098 
00099   _graph = new Epetra_CrsGraph(Copy, *_map, &n_nz_tot[0]);
00100 
00101   // Tell the matrix about its structure.  Initialize it
00102   // to zero.
00103   for (numeric_index_type i=0; i<n_rows; i++)
00104     _graph->InsertGlobalIndices(_graph->GRID(i),
00105                                 sparsity_pattern[i].size(),
00106                                 const_cast<int *>((const int *)&sparsity_pattern[i][0]));
00107 
00108   _graph->FillComplete();
00109 
00110   //Initialize the matrix
00111   libmesh_assert (!this->initialized());
00112   this->init ();
00113   libmesh_assert (this->initialized());
00114 }
00115 
00116 
00117 
00118 template <typename T>
00119 void EpetraMatrix<T>::init (const numeric_index_type m,
00120                             const numeric_index_type n,
00121                             const numeric_index_type m_l,
00122                             const numeric_index_type n_l,
00123                             const numeric_index_type nnz,
00124                             const numeric_index_type noz,
00125                             const numeric_index_type /* blocksize */)
00126 {
00127   if ((m==0) || (n==0))
00128     return;
00129 
00130   {
00131     // Clear initialized matrices
00132     if (this->initialized())
00133       this->clear();
00134 
00135     libmesh_assert (!this->_mat);
00136     libmesh_assert (!this->_map);
00137 
00138     this->_is_initialized = true;
00139   }
00140 
00141   // error checking
00142 #ifndef NDEBUG
00143   {
00144     libmesh_assert_equal_to (n, m);
00145     libmesh_assert_equal_to (n_l, m_l);
00146 
00147     numeric_index_type
00148       summed_m_l = m_l,
00149       summed_n_l = n_l;
00150 
00151     this->comm().sum (summed_m_l);
00152     this->comm().sum (summed_n_l);
00153 
00154     libmesh_assert_equal_to (m, summed_m_l);
00155     libmesh_assert_equal_to (n, summed_n_l);
00156   }
00157 #endif
00158 
00159   // build a map defining the data distribution
00160   _map = new Epetra_Map (static_cast<int>(m),
00161                          m_l,
00162                          0,
00163                          Epetra_MpiComm (this->comm().get()));
00164 
00165   libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->NumGlobalPoints()), m);
00166   libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->MaxAllGID()+1), m);
00167 
00168   _mat = new Epetra_FECrsMatrix (Copy, *_map, nnz + noz);
00169 }
00170 
00171 
00172 
00173 
00174 template <typename T>
00175 void EpetraMatrix<T>::init ()
00176 {
00177   libmesh_assert(this->_dof_map);
00178 
00179   {
00180     // Clear initialized matrices
00181     if (this->initialized())
00182       this->clear();
00183 
00184     this->_is_initialized = true;
00185   }
00186 
00187 
00188   _mat = new Epetra_FECrsMatrix (Copy, *_graph);
00189 }
00190 
00191 
00192 
00193 template <typename T>
00194 void EpetraMatrix<T>::zero ()
00195 {
00196   libmesh_assert (this->initialized());
00197 
00198   _mat->Scale(0.0);
00199 }
00200 
00201 
00202 
00203 template <typename T>
00204 void EpetraMatrix<T>::clear ()
00205 {
00206   //  delete _mat;
00207   //  delete _map;
00208 
00209   this->_is_initialized = false;
00210 
00211   libmesh_assert (!this->initialized());
00212 }
00213 
00214 
00215 
00216 template <typename T>
00217 Real EpetraMatrix<T>::l1_norm () const
00218 {
00219   libmesh_assert (this->initialized());
00220 
00221   libmesh_assert(_mat);
00222 
00223   return static_cast<Real>(_mat->NormOne());
00224 }
00225 
00226 
00227 
00228 template <typename T>
00229 Real EpetraMatrix<T>::linfty_norm () const
00230 {
00231   libmesh_assert (this->initialized());
00232 
00233 
00234   libmesh_assert(_mat);
00235 
00236   return static_cast<Real>(_mat->NormInf());
00237 }
00238 
00239 
00240 
00241 template <typename T>
00242 void EpetraMatrix<T>::add_matrix(const DenseMatrix<T>& dm,
00243                                  const std::vector<numeric_index_type>& rows,
00244                                  const std::vector<numeric_index_type>& cols)
00245 {
00246   libmesh_assert (this->initialized());
00247 
00248   const numeric_index_type m = dm.m();
00249   const numeric_index_type n = dm.n();
00250 
00251   libmesh_assert_equal_to (rows.size(), m);
00252   libmesh_assert_equal_to (cols.size(), n);
00253 
00254   _mat->SumIntoGlobalValues(m, (int *)&rows[0], n, (int *)&cols[0], &dm.get_values()[0]);
00255 }
00256 
00257 
00258 
00259 
00260 
00261 
00262 template <typename T>
00263 void EpetraMatrix<T>::get_diagonal (NumericVector<T>& dest) const
00264 {
00265   // Convert vector to EpetraVector.
00266   EpetraVector<T>* epetra_dest = cast_ptr<EpetraVector<T>*>(&dest);
00267 
00268   // Call Epetra function.
00269   _mat->ExtractDiagonalCopy(*(epetra_dest->vec()));
00270 }
00271 
00272 
00273 
00274 template <typename T>
00275 void EpetraMatrix<T>::get_transpose (SparseMatrix<T>& dest) const
00276 {
00277   // Make sure the SparseMatrix passed in is really a EpetraMatrix
00278   EpetraMatrix<T>& epetra_dest = cast_ref<EpetraMatrix<T>&>(dest);
00279 
00280   if(&epetra_dest != this)
00281     epetra_dest = *this;
00282 
00283   epetra_dest._use_transpose = !epetra_dest._use_transpose;
00284   epetra_dest._mat->SetUseTranspose(epetra_dest._use_transpose);
00285 }
00286 
00287 
00288 
00289 template <typename T>
00290 EpetraMatrix<T>::EpetraMatrix(const Parallel::Communicator &comm)
00291   : SparseMatrix<T>(comm),
00292     _destroy_mat_on_exit(true),
00293     _use_transpose(false)
00294 {}
00295 
00296 
00297 
00298 
00299 template <typename T>
00300 EpetraMatrix<T>::EpetraMatrix(Epetra_FECrsMatrix * m,
00301                               const Parallel::Communicator &comm)
00302   : SparseMatrix<T>(comm),
00303     _destroy_mat_on_exit(false),
00304     _use_transpose(false) // dumb guess is the best we can do...
00305 {
00306   this->_mat = m;
00307   this->_is_initialized = true;
00308 }
00309 
00310 
00311 
00312 
00313 template <typename T>
00314 EpetraMatrix<T>::~EpetraMatrix()
00315 {
00316   this->clear();
00317 }
00318 
00319 
00320 
00321 template <typename T>
00322 void EpetraMatrix<T>::close () const
00323 {
00324   libmesh_assert(_mat);
00325 
00326   _mat->GlobalAssemble();
00327 }
00328 
00329 
00330 
00331 template <typename T>
00332 numeric_index_type EpetraMatrix<T>::m () const
00333 {
00334   libmesh_assert (this->initialized());
00335 
00336   return static_cast<numeric_index_type>(_mat->NumGlobalRows());
00337 }
00338 
00339 
00340 
00341 template <typename T>
00342 numeric_index_type EpetraMatrix<T>::n () const
00343 {
00344   libmesh_assert (this->initialized());
00345 
00346   return static_cast<numeric_index_type>(_mat->NumGlobalCols());
00347 }
00348 
00349 
00350 
00351 template <typename T>
00352 numeric_index_type EpetraMatrix<T>::row_start () const
00353 {
00354   libmesh_assert (this->initialized());
00355   libmesh_assert(_map);
00356 
00357   return static_cast<numeric_index_type>(_map->MinMyGID());
00358 }
00359 
00360 
00361 
00362 template <typename T>
00363 numeric_index_type EpetraMatrix<T>::row_stop () const
00364 {
00365   libmesh_assert (this->initialized());
00366   libmesh_assert(_map);
00367 
00368   return static_cast<numeric_index_type>(_map->MaxMyGID())+1;
00369 }
00370 
00371 
00372 
00373 template <typename T>
00374 void EpetraMatrix<T>::set (const numeric_index_type i,
00375                            const numeric_index_type j,
00376                            const T value)
00377 {
00378   libmesh_assert (this->initialized());
00379 
00380   int
00381     epetra_i = static_cast<int>(i),
00382     epetra_j = static_cast<int>(j);
00383 
00384   T epetra_value = value;
00385 
00386   if (_mat->Filled())
00387     _mat->ReplaceGlobalValues (epetra_i, 1, &epetra_value, &epetra_j);
00388   else
00389     _mat->InsertGlobalValues (epetra_i, 1, &epetra_value, &epetra_j);
00390 }
00391 
00392 
00393 
00394 template <typename T>
00395 void EpetraMatrix<T>::add (const numeric_index_type i,
00396                            const numeric_index_type j,
00397                            const T value)
00398 {
00399   libmesh_assert (this->initialized());
00400 
00401   int
00402     epetra_i = static_cast<int>(i),
00403     epetra_j = static_cast<int>(j);
00404 
00405   T epetra_value = value;
00406 
00407   _mat->SumIntoGlobalValues (epetra_i, 1, &epetra_value, &epetra_j);
00408 }
00409 
00410 
00411 
00412 template <typename T>
00413 void EpetraMatrix<T>::add_matrix(const DenseMatrix<T>& dm,
00414                                  const std::vector<numeric_index_type>& dof_indices)
00415 {
00416   this->add_matrix (dm, dof_indices, dof_indices);
00417 }
00418 
00419 
00420 
00421 template <typename T>
00422 void EpetraMatrix<T>::add (const T a_in, SparseMatrix<T> &X_in)
00423 {
00424   libmesh_assert (this->initialized());
00425 
00426   // sanity check. but this cannot avoid
00427   // crash due to incompatible sparsity structure...
00428   libmesh_assert_equal_to (this->m(), X_in.m());
00429   libmesh_assert_equal_to (this->n(), X_in.n());
00430 
00431   EpetraMatrix<T>* X = cast_ptr<EpetraMatrix<T>*> (&X_in);
00432 
00433   EpetraExt::MatrixMatrix::Add (*X->_mat, false, a_in, *_mat, 1.);
00434 }
00435 
00436 
00437 
00438 
00439 template <typename T>
00440 T EpetraMatrix<T>::operator () (const numeric_index_type i,
00441                                 const numeric_index_type j) const
00442 {
00443   libmesh_assert (this->initialized());
00444   libmesh_assert(this->_mat);
00445   libmesh_assert (this->_mat->MyGlobalRow(static_cast<int>(i)));
00446   libmesh_assert_greater_equal (i, this->row_start());
00447   libmesh_assert_less (i, this->row_stop());
00448 
00449 
00450   int row_length, *row_indices;
00451   double *values;
00452 
00453   _mat->ExtractMyRowView (i-this->row_start(),
00454                           row_length,
00455                           values,
00456                           row_indices);
00457 
00458   //libMesh::out << "row_length=" << row_length << std::endl;
00459 
00460   int *index = std::lower_bound (row_indices, row_indices+row_length, j);
00461 
00462   libmesh_assert_less (*index, row_length);
00463   libmesh_assert_equal_to (static_cast<numeric_index_type>(row_indices[*index]), j);
00464 
00465   //libMesh::out << "val=" << values[*index] << std::endl;
00466 
00467   return values[*index];
00468 }
00469 
00470 
00471 
00472 
00473 template <typename T>
00474 bool EpetraMatrix<T>::closed() const
00475 {
00476   libmesh_assert (this->initialized());
00477   libmesh_assert(this->_mat);
00478 
00479   return this->_mat->Filled();
00480 }
00481 
00482 
00483 template <typename T>
00484 void EpetraMatrix<T>::swap(EpetraMatrix<T> & m)
00485 {
00486   std::swap(_mat, m._mat);
00487   std::swap(_destroy_mat_on_exit, m._destroy_mat_on_exit);
00488 }
00489 
00490 
00491 
00492 
00493 
00494 template <typename T>
00495 void EpetraMatrix<T>::print_personal(std::ostream& os) const
00496 {
00497   libmesh_assert (this->initialized());
00498   libmesh_assert(_mat);
00499 
00500   os << *_mat;
00501 }
00502 
00503 
00504 
00505 //------------------------------------------------------------------
00506 // Explicit instantiations
00507 template class EpetraMatrix<Number>;
00508 
00509 } // namespace libMesh
00510 
00511 
00512 #endif // #ifdef LIBMESH_HAVE_TRILINOS