trilinos_epetra_matrix.C
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 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(libMesh::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 CommWorld.sum (summed_m_l); 00067 CommWorld.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 (libMesh::COMM_WORLD)); 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 { 00126 if ((m==0) || (n==0)) 00127 return; 00128 00129 { 00130 // Clear initialized matrices 00131 if (this->initialized()) 00132 this->clear(); 00133 00134 libmesh_assert (!this->_mat); 00135 libmesh_assert (!this->_map); 00136 00137 this->_is_initialized = true; 00138 } 00139 00140 // error checking 00141 #ifndef NDEBUG 00142 { 00143 libmesh_assert_equal_to (n, m); 00144 libmesh_assert_equal_to (n_l, m_l); 00145 00146 numeric_index_type 00147 summed_m_l = m_l, 00148 summed_n_l = n_l; 00149 00150 CommWorld.sum (summed_m_l); 00151 CommWorld.sum (summed_n_l); 00152 00153 libmesh_assert_equal_to (m, summed_m_l); 00154 libmesh_assert_equal_to (n, summed_n_l); 00155 } 00156 #endif 00157 00158 // build a map defining the data distribution 00159 _map = new Epetra_Map (static_cast<int>(m), 00160 m_l, 00161 0, 00162 Epetra_MpiComm (libMesh::COMM_WORLD)); 00163 00164 libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->NumGlobalPoints()), m); 00165 libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->MaxAllGID()+1), m); 00166 00167 _mat = new Epetra_FECrsMatrix (Copy, *_map, nnz + noz); 00168 } 00169 00170 00171 00172 00173 template <typename T> 00174 void EpetraMatrix<T>::init () 00175 { 00176 libmesh_assert(this->_dof_map); 00177 00178 { 00179 // Clear initialized matrices 00180 if (this->initialized()) 00181 this->clear(); 00182 00183 this->_is_initialized = true; 00184 } 00185 00186 00187 _mat = new Epetra_FECrsMatrix (Copy, *_graph); 00188 } 00189 00190 00191 00192 template <typename T> 00193 void EpetraMatrix<T>::zero () 00194 { 00195 libmesh_assert (this->initialized()); 00196 00197 _mat->Scale(0.0); 00198 } 00199 00200 00201 00202 template <typename T> 00203 void EpetraMatrix<T>::clear () 00204 { 00205 // delete _mat; 00206 // delete _map; 00207 00208 this->_is_initialized = false; 00209 00210 libmesh_assert (!this->initialized()); 00211 } 00212 00213 00214 00215 template <typename T> 00216 Real EpetraMatrix<T>::l1_norm () const 00217 { 00218 libmesh_assert (this->initialized()); 00219 00220 libmesh_assert(_mat); 00221 00222 return static_cast<Real>(_mat->NormOne()); 00223 } 00224 00225 00226 00227 template <typename T> 00228 Real EpetraMatrix<T>::linfty_norm () const 00229 { 00230 libmesh_assert (this->initialized()); 00231 00232 00233 libmesh_assert(_mat); 00234 00235 return static_cast<Real>(_mat->NormInf()); 00236 } 00237 00238 00239 00240 template <typename T> 00241 void EpetraMatrix<T>::print_matlab (const std::string) const 00242 { 00243 libmesh_assert (this->initialized()); 00244 00245 // libmesh_assert (this->closed()); 00246 this->close(); 00247 00248 libmesh_not_implemented(); 00249 } 00250 00251 00252 00253 00254 template <typename T> 00255 void EpetraMatrix<T>::add_matrix(const DenseMatrix<T>& dm, 00256 const std::vector<numeric_index_type>& rows, 00257 const std::vector<numeric_index_type>& cols) 00258 { 00259 libmesh_assert (this->initialized()); 00260 00261 const numeric_index_type m = dm.m(); 00262 const numeric_index_type n = dm.n(); 00263 00264 libmesh_assert_equal_to (rows.size(), m); 00265 libmesh_assert_equal_to (cols.size(), n); 00266 00267 _mat->SumIntoGlobalValues(m, (int *)&rows[0], n, (int *)&cols[0], &dm.get_values()[0]); 00268 } 00269 00270 00271 00272 00273 00274 00275 template <typename T> 00276 void EpetraMatrix<T>::get_diagonal (NumericVector<T>& dest) const 00277 { 00278 // Convert vector to EpetraVector. 00279 EpetraVector<T>* epetra_dest = libmesh_cast_ptr<EpetraVector<T>*>(&dest); 00280 00281 // Call Epetra function. 00282 _mat->ExtractDiagonalCopy(*(epetra_dest->vec())); 00283 } 00284 00285 00286 00287 template <typename T> 00288 void EpetraMatrix<T>::get_transpose (SparseMatrix<T>& dest) const 00289 { 00290 // Make sure the SparseMatrix passed in is really a EpetraMatrix 00291 EpetraMatrix<T>& epetra_dest = libmesh_cast_ref<EpetraMatrix<T>&>(dest); 00292 00293 if(&epetra_dest != this) 00294 epetra_dest = *this; 00295 00296 epetra_dest._use_transpose = !epetra_dest._use_transpose; 00297 epetra_dest._mat->SetUseTranspose(epetra_dest._use_transpose); 00298 } 00299 00300 00301 00302 template <typename T> 00303 EpetraMatrix<T>::EpetraMatrix() 00304 : _destroy_mat_on_exit(true), 00305 _use_transpose(false) 00306 {} 00307 00308 00309 00310 00311 template <typename T> 00312 EpetraMatrix<T>::EpetraMatrix(Epetra_FECrsMatrix * m) 00313 : _destroy_mat_on_exit(false), 00314 _use_transpose(false) // dumb guess is the best we can do... 00315 { 00316 this->_mat = m; 00317 this->_is_initialized = true; 00318 } 00319 00320 00321 00322 00323 template <typename T> 00324 EpetraMatrix<T>::~EpetraMatrix() 00325 { 00326 this->clear(); 00327 } 00328 00329 00330 00331 template <typename T> 00332 void EpetraMatrix<T>::close () const 00333 { 00334 libmesh_assert(_mat); 00335 00336 _mat->GlobalAssemble(); 00337 } 00338 00339 00340 00341 template <typename T> 00342 numeric_index_type EpetraMatrix<T>::m () const 00343 { 00344 libmesh_assert (this->initialized()); 00345 00346 return static_cast<numeric_index_type>(_mat->NumGlobalRows()); 00347 } 00348 00349 00350 00351 template <typename T> 00352 numeric_index_type EpetraMatrix<T>::n () const 00353 { 00354 libmesh_assert (this->initialized()); 00355 00356 return static_cast<numeric_index_type>(_mat->NumGlobalCols()); 00357 } 00358 00359 00360 00361 template <typename T> 00362 numeric_index_type EpetraMatrix<T>::row_start () const 00363 { 00364 libmesh_assert (this->initialized()); 00365 libmesh_assert(_map); 00366 00367 return static_cast<numeric_index_type>(_map->MinMyGID()); 00368 } 00369 00370 00371 00372 template <typename T> 00373 numeric_index_type EpetraMatrix<T>::row_stop () const 00374 { 00375 libmesh_assert (this->initialized()); 00376 libmesh_assert(_map); 00377 00378 return static_cast<numeric_index_type>(_map->MaxMyGID())+1; 00379 } 00380 00381 00382 00383 template <typename T> 00384 void EpetraMatrix<T>::set (const numeric_index_type i, 00385 const numeric_index_type j, 00386 const T value) 00387 { 00388 libmesh_assert (this->initialized()); 00389 00390 int 00391 epetra_i = static_cast<int>(i), 00392 epetra_j = static_cast<int>(j); 00393 00394 T epetra_value = value; 00395 00396 if (_mat->Filled()) 00397 _mat->ReplaceGlobalValues (epetra_i, 1, &epetra_value, &epetra_j); 00398 else 00399 _mat->InsertGlobalValues (epetra_i, 1, &epetra_value, &epetra_j); 00400 } 00401 00402 00403 00404 template <typename T> 00405 void EpetraMatrix<T>::add (const numeric_index_type i, 00406 const numeric_index_type j, 00407 const T value) 00408 { 00409 libmesh_assert (this->initialized()); 00410 00411 int 00412 epetra_i = static_cast<int>(i), 00413 epetra_j = static_cast<int>(j); 00414 00415 T epetra_value = value; 00416 00417 _mat->SumIntoGlobalValues (epetra_i, 1, &epetra_value, &epetra_j); 00418 } 00419 00420 00421 00422 template <typename T> 00423 void EpetraMatrix<T>::add_matrix(const DenseMatrix<T>& dm, 00424 const std::vector<numeric_index_type>& dof_indices) 00425 { 00426 this->add_matrix (dm, dof_indices, dof_indices); 00427 } 00428 00429 00430 00431 template <typename T> 00432 void EpetraMatrix<T>::add (const T a_in, SparseMatrix<T> &X_in) 00433 { 00434 libmesh_assert (this->initialized()); 00435 00436 // sanity check. but this cannot avoid 00437 // crash due to incompatible sparsity structure... 00438 libmesh_assert_equal_to (this->m(), X_in.m()); 00439 libmesh_assert_equal_to (this->n(), X_in.n()); 00440 00441 EpetraMatrix<T>* X = libmesh_cast_ptr<EpetraMatrix<T>*> (&X_in); 00442 00443 EpetraExt::MatrixMatrix::Add (*X->_mat, false, a_in, *_mat, 1.); 00444 } 00445 00446 00447 00448 00449 template <typename T> 00450 T EpetraMatrix<T>::operator () (const numeric_index_type i, 00451 const numeric_index_type j) const 00452 { 00453 libmesh_assert (this->initialized()); 00454 libmesh_assert(this->_mat); 00455 libmesh_assert (this->_mat->MyGlobalRow(static_cast<int>(i))); 00456 libmesh_assert_greater_equal (i, this->row_start()); 00457 libmesh_assert_less (i, this->row_stop()); 00458 00459 00460 int row_length, *row_indices; 00461 double *values; 00462 00463 _mat->ExtractMyRowView (i-this->row_start(), 00464 row_length, 00465 values, 00466 row_indices); 00467 00468 //libMesh::out << "row_length=" << row_length << std::endl; 00469 00470 int *index = std::lower_bound (row_indices, row_indices+row_length, j); 00471 00472 libmesh_assert_less (*index, row_length); 00473 libmesh_assert_equal_to (static_cast<numeric_index_type>(row_indices[*index]), j); 00474 00475 //libMesh::out << "val=" << values[*index] << std::endl; 00476 00477 return values[*index]; 00478 } 00479 00480 00481 00482 00483 template <typename T> 00484 bool EpetraMatrix<T>::closed() const 00485 { 00486 libmesh_assert (this->initialized()); 00487 libmesh_assert(this->_mat); 00488 00489 return this->_mat->Filled(); 00490 } 00491 00492 00493 template <typename T> 00494 void EpetraMatrix<T>::swap(EpetraMatrix<T> & m) 00495 { 00496 std::swap(_mat, m._mat); 00497 std::swap(_destroy_mat_on_exit, m._destroy_mat_on_exit); 00498 } 00499 00500 00501 00502 00503 00504 template <typename T> 00505 void EpetraMatrix<T>::print_personal(std::ostream& os) const 00506 { 00507 libmesh_assert (this->initialized()); 00508 libmesh_assert(_mat); 00509 00510 os << *_mat; 00511 } 00512 00513 00514 00515 //------------------------------------------------------------------ 00516 // Explicit instantiations 00517 template class EpetraMatrix<Number>; 00518 00519 } // namespace libMesh 00520 00521 00522 #endif // #ifdef LIBMESH_HAVE_TRILINOS
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:49 UTC
Hosted By: