$extrastylesheet
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 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 const numeric_index_type) 00134 { 00135 // noz ignored... only used for multiple processors! 00136 libmesh_assert_equal_to (m_in, m_l); 00137 libmesh_assert_equal_to (n_in, n_l); 00138 libmesh_assert_equal_to (m_in, n_in); 00139 libmesh_assert_greater (nnz, 0); 00140 00141 libmesh_error_msg("ERROR: Only the init() member that uses the DofMap is implemented for Laspack matrices!"); 00142 00143 this->_is_initialized = true; 00144 } 00145 00146 00147 00148 template <typename T> 00149 void LaspackMatrix<T>::init () 00150 { 00151 // Ignore calls on initialized objects 00152 if (this->initialized()) 00153 return; 00154 00155 // We need the DofMap for this! 00156 libmesh_assert(this->_dof_map); 00157 00158 // Clear intialized matrices 00159 if (this->initialized()) 00160 this->clear(); 00161 00162 const numeric_index_type n_rows = this->_dof_map->n_dofs(); 00163 #ifndef NDEBUG 00164 // The following variables are only used for assertions, 00165 // so avoid declaring them when asserts are inactive. 00166 const numeric_index_type n_cols = n_rows; 00167 const numeric_index_type n_l = this->_dof_map->n_dofs_on_processor(0); 00168 const numeric_index_type m_l = n_l; 00169 #endif 00170 00171 // Laspack Matrices only work for uniprocessor cases 00172 libmesh_assert_equal_to (n_rows, n_cols); 00173 libmesh_assert_equal_to (m_l, n_rows); 00174 libmesh_assert_equal_to (n_l, n_cols); 00175 00176 #ifndef NDEBUG 00177 // The following variables are only used for assertions, 00178 // so avoid declaring them when asserts are inactive. 00179 const std::vector<numeric_index_type>& n_nz = this->_dof_map->get_n_nz(); 00180 const std::vector<numeric_index_type>& n_oz = this->_dof_map->get_n_oz(); 00181 #endif 00182 00183 // Make sure the sparsity pattern isn't empty 00184 libmesh_assert_equal_to (n_nz.size(), n_l); 00185 libmesh_assert_equal_to (n_oz.size(), n_l); 00186 00187 if (n_rows==0) 00188 return; 00189 00190 Q_Constr(&_QMat, const_cast<char*>("Mat"), n_rows, _LPFalse, Rowws, Normal, _LPTrue); 00191 00192 this->_is_initialized = true; 00193 00194 libmesh_assert_equal_to (n_rows, this->m()); 00195 } 00196 00197 00198 00199 template <typename T> 00200 void LaspackMatrix<T>::add_matrix(const DenseMatrix<T>& dm, 00201 const std::vector<numeric_index_type>& rows, 00202 const std::vector<numeric_index_type>& cols) 00203 00204 { 00205 libmesh_assert (this->initialized()); 00206 unsigned int n_rows = cast_int<unsigned int>(rows.size()); 00207 unsigned int n_cols = cast_int<unsigned int>(cols.size()); 00208 libmesh_assert_equal_to (dm.m(), n_rows); 00209 libmesh_assert_equal_to (dm.n(), n_cols); 00210 00211 00212 for (unsigned int i=0; i<n_rows; i++) 00213 for (unsigned int j=0; j<n_cols; j++) 00214 this->add(rows[i],cols[j],dm(i,j)); 00215 } 00216 00217 00218 00219 template <typename T> 00220 void LaspackMatrix<T>::get_diagonal (NumericVector<T>& /*dest*/) const 00221 { 00222 libmesh_not_implemented(); 00223 } 00224 00225 00226 00227 template <typename T> 00228 void LaspackMatrix<T>::get_transpose (SparseMatrix<T>& /*dest*/) const 00229 { 00230 libmesh_not_implemented(); 00231 } 00232 00233 00234 00235 template <typename T> 00236 LaspackMatrix<T>::LaspackMatrix (const Parallel::Communicator &comm) : 00237 SparseMatrix<T>(comm), 00238 _closed (false) 00239 { 00240 } 00241 00242 00243 00244 template <typename T> 00245 LaspackMatrix<T>::~LaspackMatrix () 00246 { 00247 this->clear (); 00248 } 00249 00250 00251 00252 template <typename T> 00253 void LaspackMatrix<T>::clear () 00254 { 00255 if (this->initialized()) 00256 { 00257 Q_Destr(&_QMat); 00258 } 00259 00260 _csr.clear(); 00261 _row_start.clear(); 00262 _closed = false; 00263 this->_is_initialized = false; 00264 } 00265 00266 00267 00268 template <typename T> 00269 void LaspackMatrix<T>::zero () 00270 { 00271 const numeric_index_type n_rows = this->m(); 00272 00273 for (numeric_index_type row=0; row<n_rows; row++) 00274 { 00275 const std::vector<numeric_index_type>::const_iterator 00276 r_start = _row_start[row]; 00277 00278 const numeric_index_type len = (_row_start[row+1] - _row_start[row]); 00279 00280 // Make sure we agree on the row length 00281 libmesh_assert_equal_to (len, Q_GetLen(&_QMat, row+1)); 00282 00283 for (numeric_index_type l=0; l<len; l++) 00284 { 00285 const numeric_index_type j = *(r_start + l); 00286 00287 // Make sure the data structures are working 00288 libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, row+1, l)); 00289 00290 Q_SetEntry (&_QMat, row+1, l, j+1, 0.); 00291 } 00292 } 00293 } 00294 00295 00296 00297 template <typename T> 00298 numeric_index_type LaspackMatrix<T>::m () const 00299 { 00300 libmesh_assert (this->initialized()); 00301 00302 return static_cast<numeric_index_type>(Q_GetDim(const_cast<QMatrix*>(&_QMat))); 00303 } 00304 00305 00306 00307 template <typename T> 00308 numeric_index_type LaspackMatrix<T>::n () const 00309 { 00310 libmesh_assert (this->initialized()); 00311 00312 return static_cast<numeric_index_type>(Q_GetDim(const_cast<QMatrix*>(&_QMat))); 00313 } 00314 00315 00316 00317 template <typename T> 00318 numeric_index_type LaspackMatrix<T>::row_start () const 00319 { 00320 return 0; 00321 } 00322 00323 00324 00325 template <typename T> 00326 numeric_index_type LaspackMatrix<T>::row_stop () const 00327 { 00328 return this->m(); 00329 } 00330 00331 00332 00333 template <typename T> 00334 void LaspackMatrix<T>::set (const numeric_index_type i, 00335 const numeric_index_type j, 00336 const T value) 00337 { 00338 libmesh_assert (this->initialized()); 00339 libmesh_assert_less (i, this->m()); 00340 libmesh_assert_less (j, this->n()); 00341 00342 const numeric_index_type position = this->pos(i,j); 00343 00344 // Sanity check 00345 libmesh_assert_equal_to (*(_row_start[i]+position), j); 00346 libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, i+1, position)); 00347 00348 Q_SetEntry (&_QMat, i+1, position, j+1, value); 00349 } 00350 00351 00352 00353 template <typename T> 00354 void LaspackMatrix<T>::add (const numeric_index_type i, 00355 const numeric_index_type j, 00356 const T value) 00357 { 00358 libmesh_assert (this->initialized()); 00359 libmesh_assert_less (i, this->m()); 00360 libmesh_assert_less (j, this->n()); 00361 00362 const numeric_index_type position = this->pos(i,j); 00363 00364 // Sanity check 00365 libmesh_assert_equal_to (*(_row_start[i]+position), j); 00366 00367 Q_AddVal (&_QMat, i+1, position, value); 00368 } 00369 00370 00371 00372 template <typename T> 00373 void LaspackMatrix<T>::add_matrix(const DenseMatrix<T>& dm, 00374 const std::vector<numeric_index_type>& dof_indices) 00375 { 00376 this->add_matrix (dm, dof_indices, dof_indices); 00377 } 00378 00379 00380 00381 template <typename T> 00382 void LaspackMatrix<T>::add (const T a_in, SparseMatrix<T> &X_in) 00383 { 00384 libmesh_assert (this->initialized()); 00385 libmesh_assert_equal_to (this->m(), X_in.m()); 00386 libmesh_assert_equal_to (this->n(), X_in.n()); 00387 00388 LaspackMatrix<T>* X = cast_ptr<LaspackMatrix<T>*> (&X_in); 00389 _LPNumber a = static_cast<_LPNumber> (a_in); 00390 00391 libmesh_assert(X); 00392 00393 // loops taken from LaspackMatrix<T>::zero () 00394 00395 const numeric_index_type n_rows = this->m(); 00396 00397 for (numeric_index_type row=0; row<n_rows; row++) 00398 { 00399 const std::vector<numeric_index_type>::const_iterator 00400 r_start = _row_start[row]; 00401 00402 const numeric_index_type len = (_row_start[row+1] - _row_start[row]); 00403 00404 // Make sure we agree on the row length 00405 libmesh_assert_equal_to (len, Q_GetLen(&_QMat, row+1)); 00406 // compare matrix sparsity structures 00407 libmesh_assert_equal_to (len, Q_GetLen(&(X->_QMat), row+1)); 00408 00409 00410 for (numeric_index_type l=0; l<len; l++) 00411 { 00412 const numeric_index_type j = *(r_start + l); 00413 00414 // Make sure the data structures are working 00415 libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, row+1, l)); 00416 00417 const _LPNumber value = a * Q_GetEl(const_cast<QMatrix*>(&(X->_QMat)), row+1, j+1); 00418 Q_AddVal (&_QMat, row+1, l, value); 00419 } 00420 } 00421 } 00422 00423 00424 00425 00426 template <typename T> 00427 T LaspackMatrix<T>::operator () (const numeric_index_type i, 00428 const numeric_index_type j) const 00429 { 00430 libmesh_assert (this->initialized()); 00431 libmesh_assert_less (i, this->m()); 00432 libmesh_assert_less (j, this->n()); 00433 00434 return Q_GetEl (const_cast<QMatrix*>(&_QMat), i+1, j+1); 00435 } 00436 00437 00438 00439 template <typename T> 00440 numeric_index_type LaspackMatrix<T>::pos (const numeric_index_type i, 00441 const numeric_index_type j) const 00442 { 00443 libmesh_assert_less (i, this->m()); 00444 libmesh_assert_less (j, this->n()); 00445 libmesh_assert_less (i+1, _row_start.size()); 00446 libmesh_assert (_row_start.back() == _csr.end()); 00447 00448 // note this requires the _csr to be 00449 std::pair<std::vector<numeric_index_type>::const_iterator, 00450 std::vector<numeric_index_type>::const_iterator> p = 00451 std::equal_range (_row_start[i], 00452 _row_start[i+1], 00453 j); 00454 00455 // Make sure the row contains the element j 00456 libmesh_assert (p.first != p.second); 00457 00458 // Make sure the values match 00459 libmesh_assert (*p.first == j); 00460 00461 // Return the position in the compressed row 00462 return std::distance (_row_start[i], p.first); 00463 } 00464 00465 00466 00467 //------------------------------------------------------------------ 00468 // Explicit instantiations 00469 template class LaspackMatrix<Number>; 00470 00471 } // namespace libMesh 00472 00473 00474 #endif // #ifdef LIBMESH_HAVE_LASPACK