trilinos_epetra_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 // C++ includes
20 #include "libmesh/libmesh_config.h"
21 
22 #ifdef LIBMESH_TRILINOS_HAVE_EPETRA
23 
24 // Local includes
27 #include "libmesh/dof_map.h"
28 #include "libmesh/dense_matrix.h"
29 #include "libmesh/parallel.h"
31 
32 namespace libMesh
33 {
34 
35 
36 
37 //-----------------------------------------------------------------------
38 //EpetraMatrix members
39 
40 template <typename T>
42 {
43  // clear data, start over
44  this->clear ();
45 
46  // big trouble if this fails!
47  libmesh_assert(this->_dof_map);
48 
49  const numeric_index_type n_rows = cast_int<numeric_index_type>
50  (sparsity_pattern.size());
51 
52  const numeric_index_type m = this->_dof_map->n_dofs();
53  const numeric_index_type n = m;
54  const numeric_index_type n_l = this->_dof_map->n_dofs_on_processor(this->processor_id());
55  const numeric_index_type m_l = n_l;
56 
57  // error checking
58 #ifndef NDEBUG
59  {
60  libmesh_assert_equal_to (n, m);
61  libmesh_assert_equal_to (n_l, m_l);
62 
64  summed_m_l = m_l,
65  summed_n_l = n_l;
66 
67  this->comm().sum (summed_m_l);
68  this->comm().sum (summed_n_l);
69 
70  libmesh_assert_equal_to (m, summed_m_l);
71  libmesh_assert_equal_to (n, summed_n_l);
72  }
73 #endif
74 
75  // build a map defining the data distribution
76  _map = new Epetra_Map (static_cast<int>(m),
77  m_l,
78  0,
79  Epetra_MpiComm (this->comm().get()));
80 
81  libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->NumGlobalPoints()), m);
82  libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->MaxAllGID()+1), m);
83 
84  const std::vector<numeric_index_type> & n_nz = this->_dof_map->get_n_nz();
85  const std::vector<numeric_index_type> & n_oz = this->_dof_map->get_n_oz();
86 
87  // Make sure the sparsity pattern isn't empty
88  libmesh_assert_equal_to (n_nz.size(), n_l);
89  libmesh_assert_equal_to (n_oz.size(), n_l);
90 
91  // Epetra wants the total number of nonzeros, both local and remote.
92  std::vector<int> n_nz_tot; n_nz_tot.reserve(n_nz.size());
93 
94  for (numeric_index_type i=0; i<n_nz.size(); i++)
95  n_nz_tot.push_back(std::min(n_nz[i] + n_oz[i], n));
96 
97  if (m==0)
98  return;
99 
100  _graph = new Epetra_CrsGraph(Copy, *_map, n_nz_tot.data());
101 
102  // Tell the matrix about its structure. Initialize it
103  // to zero.
104  for (numeric_index_type i=0; i<n_rows; i++)
105  _graph->InsertGlobalIndices(_graph->GRID(i),
106  cast_int<numeric_index_type>(sparsity_pattern[i].size()),
107  const_cast<int *>(reinterpret_cast<const int *>(sparsity_pattern[i].data())));
108 
109  _graph->FillComplete();
110 
111  //Initialize the matrix
112  libmesh_assert (!this->initialized());
113  this->init ();
114  libmesh_assert (this->initialized());
115 }
116 
117 
118 
119 template <typename T>
121  const numeric_index_type n,
122  const numeric_index_type m_l,
123  const numeric_index_type libmesh_dbg_var(n_l),
124  const numeric_index_type nnz,
125  const numeric_index_type noz,
126  const numeric_index_type /* blocksize */)
127 {
128  if ((m==0) || (n==0))
129  return;
130 
131  {
132  // Clear initialized matrices
133  if (this->initialized())
134  this->clear();
135 
136  libmesh_assert (!this->_mat);
137  libmesh_assert (!this->_map);
138 
139  this->_is_initialized = true;
140  }
141 
142  // error checking
143 #ifndef NDEBUG
144  {
145  libmesh_assert_equal_to (n, m);
146  libmesh_assert_equal_to (n_l, m_l);
147 
149  summed_m_l = m_l,
150  summed_n_l = n_l;
151 
152  this->comm().sum (summed_m_l);
153  this->comm().sum (summed_n_l);
154 
155  libmesh_assert_equal_to (m, summed_m_l);
156  libmesh_assert_equal_to (n, summed_n_l);
157  }
158 #endif
159 
160  // build a map defining the data distribution
161  _map = new Epetra_Map (static_cast<int>(m),
162  m_l,
163  0,
164  Epetra_MpiComm (this->comm().get()));
165 
166  libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->NumGlobalPoints()), m);
167  libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->MaxAllGID()+1), m);
168 
169  _mat = new Epetra_FECrsMatrix (Copy, *_map, nnz + noz);
170 }
171 
172 
173 
174 
175 template <typename T>
177 {
178  libmesh_assert(this->_dof_map);
179 
180  {
181  // Clear initialized matrices
182  if (this->initialized())
183  this->clear();
184 
185  this->_is_initialized = true;
186  }
187 
188 
189  _mat = new Epetra_FECrsMatrix (Copy, *_graph);
190 }
191 
192 
193 
194 template <typename T>
196 {
197  libmesh_assert (this->initialized());
198 
199  _mat->Scale(0.0);
200 }
201 
202 
203 
204 template <typename T>
206 {
207  // delete _mat;
208  // delete _map;
209 
210  this->_is_initialized = false;
211 
212  libmesh_assert (!this->initialized());
213 }
214 
215 
216 
217 template <typename T>
219 {
220  libmesh_assert (this->initialized());
221 
222  libmesh_assert(_mat);
223 
224  return static_cast<Real>(_mat->NormOne());
225 }
226 
227 
228 
229 template <typename T>
231 {
232  libmesh_assert (this->initialized());
233 
234 
235  libmesh_assert(_mat);
236 
237  return static_cast<Real>(_mat->NormInf());
238 }
239 
240 
241 
242 template <typename T>
244  const std::vector<numeric_index_type> & rows,
245  const std::vector<numeric_index_type> & cols)
246 {
247  libmesh_assert (this->initialized());
248 
249  const numeric_index_type m = dm.m();
250  const numeric_index_type n = dm.n();
251 
252  libmesh_assert_equal_to (rows.size(), m);
253  libmesh_assert_equal_to (cols.size(), n);
254 
255  _mat->SumIntoGlobalValues(m, numeric_trilinos_cast(rows.data()),
256  n, numeric_trilinos_cast(cols.data()),
257  dm.get_values().data());
258 }
259 
260 
261 
262 
263 
264 
265 template <typename T>
267 {
268  // Convert vector to EpetraVector.
269  EpetraVector<T> * epetra_dest = cast_ptr<EpetraVector<T> *>(&dest);
270 
271  // Call Epetra function.
272  _mat->ExtractDiagonalCopy(*(epetra_dest->vec()));
273 }
274 
275 
276 
277 template <typename T>
279 {
280  // Make sure the SparseMatrix passed in is really a EpetraMatrix
281  EpetraMatrix<T> & epetra_dest = cast_ref<EpetraMatrix<T> &>(dest);
282 
283  // We currently only support calling get_transpose() with ourself
284  // as the destination. Previously, this called the default copy
285  // constructor which was not safe because this class manually
286  // manages memory.
287  if (&epetra_dest != this)
288  libmesh_not_implemented();
289 
290  epetra_dest._use_transpose = !epetra_dest._use_transpose;
291  epetra_dest._mat->SetUseTranspose(epetra_dest._use_transpose);
292 }
293 
294 
295 
296 template <typename T>
298  SparseMatrix<T>(comm),
299  _destroy_mat_on_exit(true),
300  _use_transpose(false)
301 {}
302 
303 
304 
305 
306 template <typename T>
307 EpetraMatrix<T>::EpetraMatrix(Epetra_FECrsMatrix * m,
308  const Parallel::Communicator & comm) :
309  SparseMatrix<T>(comm),
310  _destroy_mat_on_exit(false),
311  _use_transpose(false) // dumb guess is the best we can do...
312 {
313  this->_mat = m;
314  this->_is_initialized = true;
315 }
316 
317 
318 
319 
320 template <typename T>
322 {
323  this->clear();
324 }
325 
326 
327 
328 template <typename T>
330 {
331  libmesh_assert(_mat);
332 
333  _mat->GlobalAssemble();
334 }
335 
336 
337 
338 template <typename T>
340 {
341  libmesh_assert (this->initialized());
342 
343  return static_cast<numeric_index_type>(_mat->NumGlobalRows());
344 }
345 
346 
347 
348 template <typename T>
350 {
351  libmesh_assert (this->initialized());
352 
353  return static_cast<numeric_index_type>(_mat->NumGlobalCols());
354 }
355 
356 
357 
358 template <typename T>
360 {
361  libmesh_assert (this->initialized());
362  libmesh_assert(_map);
363 
364  return static_cast<numeric_index_type>(_map->MinMyGID());
365 }
366 
367 
368 
369 template <typename T>
371 {
372  libmesh_assert (this->initialized());
373  libmesh_assert(_map);
374 
375  return static_cast<numeric_index_type>(_map->MaxMyGID())+1;
376 }
377 
378 
379 
380 template <typename T>
382  const numeric_index_type j,
383  const T value)
384 {
385  libmesh_assert (this->initialized());
386 
387  int
388  epetra_i = static_cast<int>(i),
389  epetra_j = static_cast<int>(j);
390 
391  T epetra_value = value;
392 
393  if (_mat->Filled())
394  _mat->ReplaceGlobalValues (epetra_i, 1, &epetra_value, &epetra_j);
395  else
396  _mat->InsertGlobalValues (epetra_i, 1, &epetra_value, &epetra_j);
397 }
398 
399 
400 
401 template <typename T>
403  const numeric_index_type j,
404  const T value)
405 {
406  libmesh_assert (this->initialized());
407 
408  int
409  epetra_i = static_cast<int>(i),
410  epetra_j = static_cast<int>(j);
411 
412  T epetra_value = value;
413 
414  _mat->SumIntoGlobalValues (epetra_i, 1, &epetra_value, &epetra_j);
415 }
416 
417 
418 
419 template <typename T>
421  const std::vector<numeric_index_type> & dof_indices)
422 {
423  this->add_matrix (dm, dof_indices, dof_indices);
424 }
425 
426 
427 
428 template <typename T>
429 void EpetraMatrix<T>::add (const T a_in, const SparseMatrix<T> & X_in)
430 {
431 #ifdef LIBMESH_TRILINOS_HAVE_EPETRAEXT
432  libmesh_assert (this->initialized());
433 
434  // sanity check. but this cannot avoid
435  // crash due to incompatible sparsity structure...
436  libmesh_assert_equal_to (this->m(), X_in.m());
437  libmesh_assert_equal_to (this->n(), X_in.n());
438 
439  const EpetraMatrix<T> * X =
440  cast_ptr<const EpetraMatrix<T> *> (&X_in);
441 
442  EpetraExt::MatrixMatrix::Add (*X->_mat, false, a_in, *_mat, 1.);
443 #else
444  libmesh_error_msg("ERROR: EpetraExt is required for EpetraMatrix::add()!");
445 #endif
446 }
447 
448 
449 
450 
451 template <typename T>
453  const numeric_index_type j) const
454 {
455  libmesh_assert (this->initialized());
456  libmesh_assert(this->_mat);
457  libmesh_assert (this->_mat->MyGlobalRow(static_cast<int>(i)));
458  libmesh_assert_greater_equal (i, this->row_start());
459  libmesh_assert_less (i, this->row_stop());
460 
461 
462  int row_length;
463  int * row_indices;
464  double * values;
465 
466  _mat->ExtractMyRowView (i-this->row_start(),
467  row_length,
468  values,
469  row_indices);
470 
471  //libMesh::out << "row_length=" << row_length << std::endl;
472 
473  int * index = std::lower_bound (row_indices, row_indices+row_length, j);
474 
475  libmesh_assert_less (*index, row_length);
476  libmesh_assert_equal_to (static_cast<numeric_index_type>(row_indices[*index]), j);
477 
478  //libMesh::out << "val=" << values[*index] << std::endl;
479 
480  return values[*index];
481 }
482 
483 
484 
485 
486 template <typename T>
488 {
489  libmesh_assert (this->initialized());
490  libmesh_assert(this->_mat);
491 
492  return this->_mat->Filled();
493 }
494 
495 
496 template <typename T>
498 {
499  std::swap(_mat, m._mat);
500  std::swap(_destroy_mat_on_exit, m._destroy_mat_on_exit);
501 }
502 
503 
504 
505 
506 
507 template <typename T>
508 void EpetraMatrix<T>::print_personal(std::ostream & os) const
509 {
510  libmesh_assert (this->initialized());
511  libmesh_assert(_mat);
512 
513  os << *_mat;
514 }
515 
516 
517 
518 //------------------------------------------------------------------
519 // Explicit instantiations
520 template class EpetraMatrix<Number>;
521 
522 } // namespace libMesh
523 
524 
525 #endif // LIBMESH_TRILINOS_HAVE_EPETRA
virtual numeric_index_type row_stop() const override
virtual void get_transpose(SparseMatrix< T > &dest) const 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_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) override
virtual numeric_index_type n() const override
virtual void set(const numeric_index_type i, const numeric_index_type j, const T value) override
virtual void init() override
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type m() const override
virtual numeric_index_type m() const =0
void init(triangulateio &t)
std::vector< T > & get_values()
Definition: dense_matrix.h:341
virtual void clear() override
virtual numeric_index_type row_start() const override
virtual void close() override
virtual void print_personal(std::ostream &os=libMesh::out) const override
EpetraMatrix(const Parallel::Communicator &comm)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real l1_norm() const override
virtual Real linfty_norm() const override
virtual void get_diagonal(NumericVector< T > &dest) const override
void swap(Iterator &lhs, Iterator &rhs)
virtual bool closed() const override
virtual void zero() override
static const bool value
Definition: xdr_io.C:109
bool initialized()
Definition: libmesh.C:258
int * numeric_trilinos_cast(const numeric_index_type *p)
Epetra_FECrsMatrix * _mat
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value) override
IterBase * data
unsigned int n() const
A matrix object used for finite element assembly and numerics.
Definition: dense_matrix.h:54
void swap(EpetraMatrix< T > &)
long double min(long double a, double b)
virtual T operator()(const numeric_index_type i, const numeric_index_type j) const override
virtual void update_sparsity_pattern(const SparsityPattern::Graph &) override
virtual numeric_index_type n() const =0