laspack_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 00020 // C++ includes 00021 00022 // Local includes 00023 #include "libmesh/libmesh_config.h" 00024 00025 #ifdef LIBMESH_HAVE_LASPACK 00026 00027 #include "libmesh/laspack_matrix.h" 00028 #include "libmesh/dense_matrix.h" 00029 #include "libmesh/dof_map.h" 00030 #include "libmesh/sparsity_pattern.h" 00031 00032 namespace libMesh 00033 { 00034 00035 00036 //----------------------------------------------------------------------- 00037 // LaspackMatrix members 00038 template <typename T> 00039 void LaspackMatrix<T>::update_sparsity_pattern (const SparsityPattern::Graph &sparsity_pattern) 00040 { 00041 // clear data, start over 00042 this->clear (); 00043 00044 // big trouble if this fails! 00045 libmesh_assert(this->_dof_map); 00046 00047 const numeric_index_type n_rows = sparsity_pattern.size(); 00048 00049 // Initialize the _row_start data structure, 00050 // allocate storage for the _csr array 00051 { 00052 std::size_t size = 0; 00053 00054 for (numeric_index_type row=0; row<n_rows; row++) 00055 size += sparsity_pattern[row].size(); 00056 00057 _csr.resize (size); 00058 _row_start.reserve(n_rows + 1); 00059 } 00060 00061 00062 // Initize the _csr data structure. 00063 { 00064 std::vector<numeric_index_type>::iterator position = _csr.begin(); 00065 00066 _row_start.push_back (position); 00067 00068 for (numeric_index_type row=0; row<n_rows; row++) 00069 { 00070 // insert the row indices 00071 for (SparsityPattern::Row::const_iterator col = sparsity_pattern[row].begin(); 00072 col != sparsity_pattern[row].end(); ++col) 00073 { 00074 libmesh_assert (position != _csr.end()); 00075 *position = *col; 00076 ++position; 00077 } 00078 00079 _row_start.push_back (position); 00080 } 00081 } 00082 00083 00084 // Initialize the matrix 00085 libmesh_assert (!this->initialized()); 00086 this->init (); 00087 libmesh_assert (this->initialized()); 00088 //libMesh::out << "n_rows=" << n_rows << std::endl; 00089 //libMesh::out << "m()=" << m() << std::endl; 00090 libmesh_assert_equal_to (n_rows, this->m()); 00091 00092 // Tell the matrix about its structure. Initialize it 00093 // to zero. 00094 for (numeric_index_type i=0; i<n_rows; i++) 00095 { 00096 const std::vector<numeric_index_type>::const_iterator 00097 rs = _row_start[i]; 00098 00099 const numeric_index_type length = _row_start[i+1] - rs; 00100 00101 Q_SetLen (&_QMat, i+1, length); 00102 00103 for (numeric_index_type l=0; l<length; l++) 00104 { 00105 const numeric_index_type j = *(rs+l); 00106 00107 // sanity check 00108 //libMesh::out << "m()=" << m() << std::endl; 00109 //libMesh::out << "(i,j,l) = (" << i 00110 // << "," << j 00111 // << "," << l 00112 // << ")" << std::endl; 00113 //libMesh::out << "pos(i,j)=" << pos(i,j) 00114 // << std::endl; 00115 libmesh_assert_equal_to (this->pos(i,j), l); 00116 Q_SetEntry (&_QMat, i+1, l, j+1, 0.); 00117 } 00118 } 00119 00120 // That's it! 00121 //libmesh_here(); 00122 } 00123 00124 00125 00126 template <typename T> 00127 void LaspackMatrix<T>::init (const numeric_index_type libmesh_dbg_var(m_in), 00128 const numeric_index_type libmesh_dbg_var(n_in), 00129 const numeric_index_type libmesh_dbg_var(m_l), 00130 const numeric_index_type libmesh_dbg_var(n_l), 00131 const numeric_index_type libmesh_dbg_var(nnz), 00132 const numeric_index_type) 00133 { 00134 // noz ignored... only used for multiple processors! 00135 libmesh_assert_equal_to (m_in, m_l); 00136 libmesh_assert_equal_to (n_in, n_l); 00137 libmesh_assert_equal_to (m_in, n_in); 00138 libmesh_assert_greater (nnz, 0); 00139 00140 00141 libMesh::err << "ERROR: Only the init() member that uses the" << std::endl 00142 << "DofMap is implemented for Laspack matrices!" << std::endl; 00143 libmesh_error(); 00144 00145 this->_is_initialized = true; 00146 } 00147 00148 00149 00150 template <typename T> 00151 void LaspackMatrix<T>::init () 00152 { 00153 // Ignore calls on initialized objects 00154 if (this->initialized()) 00155 return; 00156 00157 // We need the DofMap for this! 00158 libmesh_assert(this->_dof_map); 00159 00160 // Clear intialized matrices 00161 if (this->initialized()) 00162 this->clear(); 00163 00164 const numeric_index_type n_rows = this->_dof_map->n_dofs(); 00165 #ifndef NDEBUG 00166 // The following variables are only used for assertions, 00167 // so avoid declaring them when asserts are inactive. 00168 const numeric_index_type n_cols = n_rows; 00169 const numeric_index_type n_l = this->_dof_map->n_dofs_on_processor(0); 00170 const numeric_index_type m_l = n_l; 00171 #endif 00172 00173 // Laspack Matrices only work for uniprocessor cases 00174 libmesh_assert_equal_to (n_rows, n_cols); 00175 libmesh_assert_equal_to (m_l, n_rows); 00176 libmesh_assert_equal_to (n_l, n_cols); 00177 00178 #ifndef NDEBUG 00179 // The following variables are only used for assertions, 00180 // so avoid declaring them when asserts are inactive. 00181 const std::vector<numeric_index_type>& n_nz = this->_dof_map->get_n_nz(); 00182 const std::vector<numeric_index_type>& n_oz = this->_dof_map->get_n_oz(); 00183 #endif 00184 00185 // Make sure the sparsity pattern isn't empty 00186 libmesh_assert_equal_to (n_nz.size(), n_l); 00187 libmesh_assert_equal_to (n_oz.size(), n_l); 00188 00189 if (n_rows==0) 00190 return; 00191 00192 Q_Constr(&_QMat, const_cast<char*>("Mat"), n_rows, _LPFalse, Rowws, Normal, _LPTrue); 00193 00194 this->_is_initialized = true; 00195 00196 libmesh_assert_equal_to (n_rows, this->m()); 00197 } 00198 00199 00200 00201 template <typename T> 00202 void LaspackMatrix<T>::add_matrix(const DenseMatrix<T>& dm, 00203 const std::vector<numeric_index_type>& rows, 00204 const std::vector<numeric_index_type>& cols) 00205 00206 { 00207 libmesh_assert (this->initialized()); 00208 unsigned int n_rows = libmesh_cast_int<unsigned int>(rows.size()); 00209 unsigned int n_cols = libmesh_cast_int<unsigned int>(cols.size()); 00210 libmesh_assert_equal_to (dm.m(), n_rows); 00211 libmesh_assert_equal_to (dm.n(), n_cols); 00212 00213 00214 for (unsigned int i=0; i<n_rows; i++) 00215 for (unsigned int j=0; j<n_cols; j++) 00216 this->add(rows[i],cols[j],dm(i,j)); 00217 } 00218 00219 00220 00221 template <typename T> 00222 void LaspackMatrix<T>::get_diagonal (NumericVector<T>& /*dest*/) const 00223 { 00224 libmesh_not_implemented(); 00225 } 00226 00227 00228 00229 template <typename T> 00230 void LaspackMatrix<T>::get_transpose (SparseMatrix<T>& /*dest*/) const 00231 { 00232 libmesh_not_implemented(); 00233 } 00234 00235 00236 00237 template <typename T> 00238 LaspackMatrix<T>::LaspackMatrix () : 00239 _closed (false) 00240 { 00241 } 00242 00243 00244 00245 template <typename T> 00246 LaspackMatrix<T>::~LaspackMatrix () 00247 { 00248 this->clear (); 00249 } 00250 00251 00252 00253 template <typename T> 00254 void LaspackMatrix<T>::clear () 00255 { 00256 if (this->initialized()) 00257 { 00258 Q_Destr(&_QMat); 00259 } 00260 00261 _csr.clear(); 00262 _row_start.clear(); 00263 _closed = false; 00264 this->_is_initialized = false; 00265 } 00266 00267 00268 00269 template <typename T> 00270 void LaspackMatrix<T>::zero () 00271 { 00272 const numeric_index_type n_rows = this->m(); 00273 00274 for (numeric_index_type row=0; row<n_rows; row++) 00275 { 00276 const std::vector<numeric_index_type>::const_iterator 00277 r_start = _row_start[row]; 00278 00279 const numeric_index_type len = (_row_start[row+1] - _row_start[row]); 00280 00281 // Make sure we agree on the row length 00282 libmesh_assert_equal_to (len, Q_GetLen(&_QMat, row+1)); 00283 00284 for (numeric_index_type l=0; l<len; l++) 00285 { 00286 const numeric_index_type j = *(r_start + l); 00287 00288 // Make sure the data structures are working 00289 libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, row+1, l)); 00290 00291 Q_SetEntry (&_QMat, row+1, l, j+1, 0.); 00292 } 00293 } 00294 } 00295 00296 00297 00298 template <typename T> 00299 numeric_index_type LaspackMatrix<T>::m () const 00300 { 00301 libmesh_assert (this->initialized()); 00302 00303 return static_cast<numeric_index_type>(Q_GetDim(const_cast<QMatrix*>(&_QMat))); 00304 } 00305 00306 00307 00308 template <typename T> 00309 numeric_index_type LaspackMatrix<T>::n () const 00310 { 00311 libmesh_assert (this->initialized()); 00312 00313 return static_cast<numeric_index_type>(Q_GetDim(const_cast<QMatrix*>(&_QMat))); 00314 } 00315 00316 00317 00318 template <typename T> 00319 numeric_index_type LaspackMatrix<T>::row_start () const 00320 { 00321 return 0; 00322 } 00323 00324 00325 00326 template <typename T> 00327 numeric_index_type LaspackMatrix<T>::row_stop () const 00328 { 00329 return this->m(); 00330 } 00331 00332 00333 00334 template <typename T> 00335 void LaspackMatrix<T>::set (const numeric_index_type i, 00336 const numeric_index_type j, 00337 const T value) 00338 { 00339 libmesh_assert (this->initialized()); 00340 libmesh_assert_less (i, this->m()); 00341 libmesh_assert_less (j, this->n()); 00342 00343 const numeric_index_type position = this->pos(i,j); 00344 00345 // Sanity check 00346 libmesh_assert_equal_to (*(_row_start[i]+position), j); 00347 libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, i+1, position)); 00348 00349 Q_SetEntry (&_QMat, i+1, position, j+1, value); 00350 } 00351 00352 00353 00354 template <typename T> 00355 void LaspackMatrix<T>::add (const numeric_index_type i, 00356 const numeric_index_type j, 00357 const T value) 00358 { 00359 libmesh_assert (this->initialized()); 00360 libmesh_assert_less (i, this->m()); 00361 libmesh_assert_less (j, this->n()); 00362 00363 const numeric_index_type position = this->pos(i,j); 00364 00365 // Sanity check 00366 libmesh_assert_equal_to (*(_row_start[i]+position), j); 00367 00368 Q_AddVal (&_QMat, i+1, position, value); 00369 } 00370 00371 00372 00373 template <typename T> 00374 void LaspackMatrix<T>::add_matrix(const DenseMatrix<T>& dm, 00375 const std::vector<numeric_index_type>& dof_indices) 00376 { 00377 this->add_matrix (dm, dof_indices, dof_indices); 00378 } 00379 00380 00381 00382 template <typename T> 00383 void LaspackMatrix<T>::add (const T a_in, SparseMatrix<T> &X_in) 00384 { 00385 libmesh_assert (this->initialized()); 00386 libmesh_assert_equal_to (this->m(), X_in.m()); 00387 libmesh_assert_equal_to (this->n(), X_in.n()); 00388 00389 LaspackMatrix<T>* X = libmesh_cast_ptr<LaspackMatrix<T>*> (&X_in); 00390 _LPNumber a = static_cast<_LPNumber> (a_in); 00391 00392 libmesh_assert(X); 00393 00394 // loops taken from LaspackMatrix<T>::zero () 00395 00396 const numeric_index_type n_rows = this->m(); 00397 00398 for (numeric_index_type row=0; row<n_rows; row++) 00399 { 00400 const std::vector<numeric_index_type>::const_iterator 00401 r_start = _row_start[row]; 00402 00403 const numeric_index_type len = (_row_start[row+1] - _row_start[row]); 00404 00405 // Make sure we agree on the row length 00406 libmesh_assert_equal_to (len, Q_GetLen(&_QMat, row+1)); 00407 // compare matrix sparsity structures 00408 libmesh_assert_equal_to (len, Q_GetLen(&(X->_QMat), row+1)); 00409 00410 00411 for (numeric_index_type l=0; l<len; l++) 00412 { 00413 const numeric_index_type j = *(r_start + l); 00414 00415 // Make sure the data structures are working 00416 libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, row+1, l)); 00417 00418 const _LPNumber value = a * Q_GetEl(const_cast<QMatrix*>(&(X->_QMat)), row+1, j+1); 00419 Q_AddVal (&_QMat, row+1, l, value); 00420 } 00421 } 00422 } 00423 00424 00425 00426 00427 template <typename T> 00428 T LaspackMatrix<T>::operator () (const numeric_index_type i, 00429 const numeric_index_type j) const 00430 { 00431 libmesh_assert (this->initialized()); 00432 libmesh_assert_less (i, this->m()); 00433 libmesh_assert_less (j, this->n()); 00434 00435 return Q_GetEl (const_cast<QMatrix*>(&_QMat), i+1, j+1); 00436 } 00437 00438 00439 00440 template <typename T> 00441 numeric_index_type LaspackMatrix<T>::pos (const numeric_index_type i, 00442 const numeric_index_type j) const 00443 { 00444 libmesh_assert_less (i, this->m()); 00445 libmesh_assert_less (j, this->n()); 00446 libmesh_assert_less (i+1, _row_start.size()); 00447 libmesh_assert (_row_start.back() == _csr.end()); 00448 00449 // note this requires the _csr to be 00450 std::pair<std::vector<numeric_index_type>::const_iterator, 00451 std::vector<numeric_index_type>::const_iterator> p = 00452 std::equal_range (_row_start[i], 00453 _row_start[i+1], 00454 j); 00455 00456 // Make sure the row contains the element j 00457 libmesh_assert (p.first != p.second); 00458 00459 // Make sure the values match 00460 libmesh_assert (*p.first == j); 00461 00462 // Return the position in the compressed row 00463 return std::distance (_row_start[i], p.first); 00464 } 00465 00466 00467 00468 //------------------------------------------------------------------ 00469 // Explicit instantiations 00470 template class LaspackMatrix<Number>; 00471 00472 } // namespace libMesh 00473 00474 00475 #endif // #ifdef LIBMESH_HAVE_LASPACK
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:47 UTC
Hosted By: