laspack_matrix.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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 
19 
20 // C++ includes
21 
22 // Local includes
23 #include "libmesh/libmesh_config.h"
24 
25 #ifdef LIBMESH_HAVE_LASPACK
26 
27 #include "libmesh/laspack_matrix.h"
28 #include "libmesh/dense_matrix.h"
29 #include "libmesh/dof_map.h"
31 
32 namespace libMesh
33 {
34 
35 
36 //-----------------------------------------------------------------------
37 // LaspackMatrix members
38 template <typename T>
40 {
41  // clear data, start over
42  this->clear ();
43 
44  // big trouble if this fails!
45  libmesh_assert(this->_dof_map);
46 
47  const numeric_index_type n_rows = sparsity_pattern.size();
48 
49  // Initialize the _row_start data structure,
50  // allocate storage for the _csr array
51  {
52  std::size_t size = 0;
53 
54  for (numeric_index_type row=0; row<n_rows; row++)
55  size += sparsity_pattern[row].size();
56 
57  _csr.resize (size);
58  _row_start.reserve(n_rows + 1);
59  }
60 
61 
62  // Initialize the _csr data structure.
63  {
64  std::vector<numeric_index_type>::iterator position = _csr.begin();
65 
66  _row_start.push_back (position);
67 
68  for (numeric_index_type row=0; row<n_rows; row++)
69  {
70  // insert the row indices
71  for (const auto & col : sparsity_pattern[row])
72  {
73  libmesh_assert (position != _csr.end());
74  *position = col;
75  ++position;
76  }
77 
78  _row_start.push_back (position);
79  }
80  }
81 
82 
83  // Initialize the matrix
84  libmesh_assert (!this->initialized());
85  this->init ();
86  libmesh_assert (this->initialized());
87  //libMesh::out << "n_rows=" << n_rows << std::endl;
88  //libMesh::out << "m()=" << m() << std::endl;
89  libmesh_assert_equal_to (n_rows, this->m());
90 
91  // Tell the matrix about its structure. Initialize it
92  // to zero.
93  for (numeric_index_type i=0; i<n_rows; i++)
94  {
95  auto rs = _row_start[i];
96 
97  const numeric_index_type length = _row_start[i+1] - rs;
98 
99  Q_SetLen (&_QMat, i+1, length);
100 
101  for (numeric_index_type l=0; l<length; l++)
102  {
103  const numeric_index_type j = *(rs+l);
104 
105  // sanity check
106  //libMesh::out << "m()=" << m() << std::endl;
107  //libMesh::out << "(i,j,l) = (" << i
108  // << "," << j
109  // << "," << l
110  // << ")" << std::endl;
111  //libMesh::out << "pos(i,j)=" << pos(i,j)
112  // << std::endl;
113  libmesh_assert_equal_to (this->pos(i,j), l);
114  Q_SetEntry (&_QMat, i+1, l, j+1, 0.);
115  }
116  }
117 
118  // That's it!
119  //libmesh_here();
120 }
121 
122 
123 
124 template <typename T>
125 void LaspackMatrix<T>::init (const numeric_index_type libmesh_dbg_var(m_in),
126  const numeric_index_type libmesh_dbg_var(n_in),
127  const numeric_index_type libmesh_dbg_var(m_l),
128  const numeric_index_type libmesh_dbg_var(n_l),
129  const numeric_index_type libmesh_dbg_var(nnz),
130  const numeric_index_type,
131  const numeric_index_type)
132 {
133  // noz ignored... only used for multiple processors!
134  libmesh_assert_equal_to (m_in, m_l);
135  libmesh_assert_equal_to (n_in, n_l);
136  libmesh_assert_equal_to (m_in, n_in);
137  libmesh_assert_greater (nnz, 0);
138 
139  libmesh_error_msg("ERROR: Only the init() member that uses the DofMap is implemented for Laspack matrices!");
140 
141  this->_is_initialized = true;
142 }
143 
144 
145 
146 template <typename T>
148 {
149  // Ignore calls on initialized objects
150  if (this->initialized())
151  return;
152 
153  // We need the DofMap for this!
154  libmesh_assert(this->_dof_map);
155 
156  // Clear initialized matrices
157  if (this->initialized())
158  this->clear();
159 
160  const numeric_index_type n_rows = this->_dof_map->n_dofs();
161 #ifndef NDEBUG
162  // The following variables are only used for assertions,
163  // so avoid declaring them when asserts are inactive.
164  const numeric_index_type n_cols = n_rows;
165  const numeric_index_type n_l = this->_dof_map->n_dofs_on_processor(0);
166  const numeric_index_type m_l = n_l;
167 #endif
168 
169  // Laspack Matrices only work for uniprocessor cases
170  libmesh_assert_equal_to (n_rows, n_cols);
171  libmesh_assert_equal_to (m_l, n_rows);
172  libmesh_assert_equal_to (n_l, n_cols);
173 
174 #ifndef NDEBUG
175  // The following variables are only used for assertions,
176  // so avoid declaring them when asserts are inactive.
177  const std::vector<numeric_index_type> & n_nz = this->_dof_map->get_n_nz();
178  const std::vector<numeric_index_type> & n_oz = this->_dof_map->get_n_oz();
179 #endif
180 
181  // Make sure the sparsity pattern isn't empty
182  libmesh_assert_equal_to (n_nz.size(), n_l);
183  libmesh_assert_equal_to (n_oz.size(), n_l);
184 
185  if (n_rows==0)
186  return;
187 
188  Q_Constr(&_QMat, const_cast<char *>("Mat"), n_rows, _LPFalse, Rowws, Normal, _LPTrue);
189 
190  this->_is_initialized = true;
191 
192  libmesh_assert_equal_to (n_rows, this->m());
193 }
194 
195 
196 
197 template <typename T>
199  const std::vector<numeric_index_type> & rows,
200  const std::vector<numeric_index_type> & cols)
201 
202 {
203  libmesh_assert (this->initialized());
204  unsigned int n_rows = cast_int<unsigned int>(rows.size());
205  unsigned int n_cols = cast_int<unsigned int>(cols.size());
206  libmesh_assert_equal_to (dm.m(), n_rows);
207  libmesh_assert_equal_to (dm.n(), n_cols);
208 
209 
210  for (unsigned int i=0; i<n_rows; i++)
211  for (unsigned int j=0; j<n_cols; j++)
212  this->add(rows[i],cols[j],dm(i,j));
213 }
214 
215 
216 
217 template <typename T>
219 {
220  libmesh_not_implemented();
221 }
222 
223 
224 
225 template <typename T>
227 {
228  libmesh_not_implemented();
229 }
230 
231 
232 
233 template <typename T>
235  SparseMatrix<T>(comm),
236  _closed (false)
237 {
238 }
239 
240 
241 
242 template <typename T>
244 {
245  this->clear ();
246 }
247 
248 
249 
250 template <typename T>
252 {
253  if (this->initialized())
254  {
255  Q_Destr(&_QMat);
256  }
257 
258  _csr.clear();
259  _row_start.clear();
260  _closed = false;
261  this->_is_initialized = false;
262 }
263 
264 
265 
266 template <typename T>
268 {
269  const numeric_index_type n_rows = this->m();
270 
271  for (numeric_index_type row=0; row<n_rows; row++)
272  {
273  auto r_start = _row_start[row];
274 
275  const numeric_index_type len = (_row_start[row+1] - _row_start[row]);
276 
277  // Make sure we agree on the row length
278  libmesh_assert_equal_to (len, Q_GetLen(&_QMat, row+1));
279 
280  for (numeric_index_type l=0; l<len; l++)
281  {
282  const numeric_index_type j = *(r_start + l);
283 
284  // Make sure the data structures are working
285  libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, row+1, l));
286 
287  Q_SetEntry (&_QMat, row+1, l, j+1, 0.);
288  }
289  }
290 
291  this->close();
292 }
293 
294 
295 
296 template <typename T>
298 {
299  libmesh_assert (this->initialized());
300 
301  return static_cast<numeric_index_type>(Q_GetDim(const_cast<QMatrix*>(&_QMat)));
302 }
303 
304 
305 
306 template <typename T>
308 {
309  libmesh_assert (this->initialized());
310 
311  return static_cast<numeric_index_type>(Q_GetDim(const_cast<QMatrix*>(&_QMat)));
312 }
313 
314 
315 
316 template <typename T>
318 {
319  return 0;
320 }
321 
322 
323 
324 template <typename T>
326 {
327  return this->m();
328 }
329 
330 
331 
332 template <typename T>
334  const numeric_index_type j,
335  const T value)
336 {
337  libmesh_assert (this->initialized());
338  libmesh_assert_less (i, this->m());
339  libmesh_assert_less (j, this->n());
340 
341  const numeric_index_type position = this->pos(i,j);
342 
343  // Sanity check
344  libmesh_assert_equal_to (*(_row_start[i]+position), j);
345  libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, i+1, position));
346 
347  Q_SetEntry (&_QMat, i+1, position, j+1, value);
348 }
349 
350 
351 
352 template <typename T>
354  const numeric_index_type j,
355  const T value)
356 {
357  libmesh_assert (this->initialized());
358  libmesh_assert_less (i, this->m());
359  libmesh_assert_less (j, this->n());
360 
361  const numeric_index_type position = this->pos(i,j);
362 
363  // Sanity check
364  libmesh_assert_equal_to (*(_row_start[i]+position), j);
365 
366  Q_AddVal (&_QMat, i+1, position, value);
367 }
368 
369 
370 
371 template <typename T>
373  const std::vector<numeric_index_type> & dof_indices)
374 {
375  this->add_matrix (dm, dof_indices, dof_indices);
376 }
377 
378 
379 
380 template <typename T>
381 void LaspackMatrix<T>::add (const T a_in, const SparseMatrix<T> & X_in)
382 {
383  libmesh_assert (this->initialized());
384  libmesh_assert_equal_to (this->m(), X_in.m());
385  libmesh_assert_equal_to (this->n(), X_in.n());
386 
387  const LaspackMatrix<T> * X =
388  cast_ptr<const LaspackMatrix<T> *> (&X_in);
389 
390  _LPNumber a = static_cast<_LPNumber> (a_in);
391 
392  libmesh_assert(X);
393 
394  // loops taken from LaspackMatrix<T>::zero ()
395 
396  const numeric_index_type n_rows = this->m();
397 
398  for (numeric_index_type row=0; row<n_rows; row++)
399  {
400  auto r_start = _row_start[row];
401 
402  const numeric_index_type len = (_row_start[row+1] - _row_start[row]);
403 
404  // Make sure we agree on the row length
405  libmesh_assert_equal_to (len, Q_GetLen(&_QMat, row+1));
406  // compare matrix sparsity structures
407  libmesh_assert_equal_to (len, Q_GetLen(&(X->_QMat), row+1));
408 
409 
410  for (numeric_index_type l=0; l<len; l++)
411  {
412  const numeric_index_type j = *(r_start + l);
413 
414  // Make sure the data structures are working
415  libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, row+1, l));
416 
417  const _LPNumber value = a * Q_GetEl(const_cast<QMatrix*>(&(X->_QMat)), row+1, j+1);
418  Q_AddVal (&_QMat, row+1, l, value);
419  }
420  }
421 }
422 
423 
424 
425 
426 template <typename T>
428  const numeric_index_type j) const
429 {
430  libmesh_assert (this->initialized());
431  libmesh_assert_less (i, this->m());
432  libmesh_assert_less (j, this->n());
433 
434  return Q_GetEl (const_cast<QMatrix*>(&_QMat), i+1, j+1);
435 }
436 
437 
438 
439 template <typename T>
441  const numeric_index_type j) const
442 {
443  libmesh_assert_less (i, this->m());
444  libmesh_assert_less (j, this->n());
445  libmesh_assert_less (i+1, _row_start.size());
446  libmesh_assert (_row_start.back() == _csr.end());
447 
448  // note this requires the _csr to be sorted
449  auto p = std::equal_range (_row_start[i], _row_start[i+1], j);
450 
451  // Make sure the row contains the element j
452  libmesh_assert (p.first != p.second);
453 
454  // Make sure the values match
455  libmesh_assert (*p.first == j);
456 
457  // Return the position in the compressed row
458  return std::distance (_row_start[i], p.first);
459 }
460 
461 
462 
463 template <typename T>
465 {
466  libmesh_assert(this->initialized());
467 
468  this->_closed = true;
469 
470  // We've probably changed some entries so we need to tell LASPACK
471  // that cached data is now invalid.
472  *_QMat.DiagElAlloc = _LPFalse;
473  *_QMat.ElSorted = _LPFalse;
474  if (*_QMat.ILUExists)
475  {
476  *_QMat.ILUExists = _LPFalse;
477  Q_Destr(_QMat.ILU);
478  }
479 }
480 
481 
482 
483 //------------------------------------------------------------------
484 // Explicit instantiations
485 template class LaspackMatrix<Number>;
486 
487 } // namespace libMesh
488 
489 
490 #endif // #ifdef LIBMESH_HAVE_LASPACK
virtual void get_diagonal(NumericVector< T > &dest) const override
virtual void init() override
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
Definition: diff_context.h:40
unsigned int m() const
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value) override
numeric_index_type pos(const numeric_index_type i, const numeric_index_type j) const
virtual void get_transpose(SparseMatrix< T > &dest) const override
virtual void close() override
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type m() const =0
virtual T operator()(const numeric_index_type i, const numeric_index_type j) const override
void init(triangulateio &t)
virtual numeric_index_type row_start() const override
virtual void zero() override
virtual void clear() override
virtual numeric_index_type m() const override
virtual numeric_index_type n() const override
static const bool value
Definition: xdr_io.C:109
bool initialized()
Definition: libmesh.C:258
virtual void update_sparsity_pattern(const SparsityPattern::Graph &) override
virtual numeric_index_type row_stop() const override
unsigned int n() const
virtual void set(const numeric_index_type i, const numeric_index_type j, const T value) override
A matrix object used for finite element assembly and numerics.
Definition: dense_matrix.h:54
virtual numeric_index_type n() const =0
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) override
LaspackMatrix(const Parallel::Communicator &comm)