Reference documentation for deal.II version 8.4.2
petsc_matrix_base.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2015 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/lac/petsc_matrix_base.h>
17 
18 #ifdef DEAL_II_WITH_PETSC
19 
20 # include <deal.II/lac/petsc_full_matrix.h>
21 # include <deal.II/lac/petsc_sparse_matrix.h>
22 # include <deal.II/lac/petsc_parallel_sparse_matrix.h>
23 # include <deal.II/lac/petsc_vector.h>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 namespace PETScWrappers
28 {
29  namespace MatrixIterators
30  {
31  void
34  {
35  // if we are asked to visit the
36  // past-the-end line, then simply
37  // release all our caches and go on
38  // with life
39  if (this->a_row == matrix->m())
40  {
41  colnum_cache.reset ();
42  value_cache.reset ();
43 
44  return;
45  }
46 
47  // get a representation of the present row
48  PetscInt ncols;
49  const PetscInt *colnums;
50  const PetscScalar *values;
51 
52  int ierr;
53  (void)ierr;
54  ierr = MatGetRow(*matrix, this->a_row, &ncols, &colnums, &values);
55  AssertThrow (ierr == 0, MatrixBase::ExcPETScError(ierr));
56 
57  // copy it into our caches if the line
58  // isn't empty. if it is, then we've
59  // done something wrong, since we
60  // shouldn't have initialized an
61  // iterator for an empty line (what
62  // would it point to?)
63  Assert (ncols != 0, ExcInternalError());
64  colnum_cache.reset (new std::vector<size_type> (colnums, colnums+ncols));
65  value_cache.reset (new std::vector<PetscScalar> (values, values+ncols));
66 
67  // and finally restore the matrix
68  ierr = MatRestoreRow(*matrix, this->a_row, &ncols, &colnums, &values);
69  AssertThrow (ierr == 0, MatrixBase::ExcPETScError(ierr));
70  }
71  }
72 
73 
74 
76  :
77  last_action (VectorOperation::unknown)
78  {}
79 
80 
81 
83  {
84 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
85  const int ierr = MatDestroy (matrix);
86 #else
87  const int ierr = MatDestroy (&matrix);
88 #endif
89  AssertThrow (ierr == 0, ExcPETScError(ierr));
90  }
91 
92 
93 
94  void
96  {
97  // destroy the matrix...
98 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
99  int ierr = MatDestroy (matrix);
100 #else
101  int ierr = MatDestroy (&matrix);
102 #endif
103  AssertThrow (ierr == 0, ExcPETScError(ierr));
104  // ...and replace it by an empty
105  // sequential matrix
106  const int m=0, n=0, n_nonzero_per_row=0;
107  ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, m, n, n_nonzero_per_row,
108  0, &matrix);
109  AssertThrow (ierr == 0, ExcPETScError(ierr));
110  }
111 
112 
113 
114  MatrixBase &
116  {
117  (void)d;
118  Assert (d==value_type(), ExcScalarAssignmentOnlyForZeroValue());
119 
121 
122  const int ierr = MatZeroEntries (matrix);
123  AssertThrow (ierr == 0, ExcPETScError(ierr));
124 
125  return *this;
126  }
127 
128 
129 
130  void
132  const PetscScalar new_diag_value)
133  {
135 
136  // now set all the entries of this row to zero
137  const PetscInt petsc_row = row;
138 
139  IS index_set;
140 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
141  ISCreateGeneral (get_mpi_communicator(), 1, &petsc_row, &index_set);
142 #else
143  ISCreateGeneral (get_mpi_communicator(), 1, &petsc_row, PETSC_COPY_VALUES, &index_set);
144 #endif
145 
146 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
147  const int ierr
148  = MatZeroRowsIS(matrix, index_set, new_diag_value);
149 #else
150  const int ierr
151  = MatZeroRowsIS(matrix, index_set, new_diag_value, PETSC_NULL, PETSC_NULL);
152 #endif
153  AssertThrow (ierr == 0, ExcPETScError(ierr));
154 
155 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
156  ISDestroy (index_set);
157 #else
158  ISDestroy (&index_set);
159 #endif
160  }
161 
162 
163 
164  void
165  MatrixBase::clear_rows (const std::vector<size_type> &rows,
166  const PetscScalar new_diag_value)
167  {
169 
170  // now set all the entries of these rows
171  // to zero
172  const std::vector<PetscInt> petsc_rows (rows.begin(), rows.end());
173 
174  // call the functions. note that we have
175  // to call them even if #rows is empty,
176  // since this is a collective operation
177  IS index_set;
178 
179 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
180  ISCreateGeneral (get_mpi_communicator(), rows.size(),
181  &petsc_rows[0], &index_set);
182 #else
183  ISCreateGeneral (get_mpi_communicator(), rows.size(),
184  &petsc_rows[0], PETSC_COPY_VALUES, &index_set);
185 #endif
186 
187 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
188  const int ierr
189  = MatZeroRowsIS(matrix, index_set, new_diag_value);
190 #else
191  const int ierr
192  = MatZeroRowsIS(matrix, index_set, new_diag_value, PETSC_NULL, PETSC_NULL);
193 #endif
194  AssertThrow (ierr == 0, ExcPETScError(ierr));
195 
196 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
197  ISDestroy (index_set);
198 #else
199  ISDestroy (&index_set);
200 #endif
201  }
202 
203 
204 
205  PetscScalar
207  const size_type j) const
208  {
209  PetscInt petsc_i = i, petsc_j = j;
210 
211  PetscScalar value;
212 
213  const int ierr
214  = MatGetValues (matrix, 1, &petsc_i, 1, &petsc_j,
215  &value);
216  AssertThrow (ierr == 0, ExcPETScError(ierr));
217 
218  return value;
219  }
220 
221 
222 
223  PetscScalar
225  {
226  Assert (m() == n(), ExcNotQuadratic());
227 
228  // this doesn't seem to work any
229  // different than any other element
230  return el(i,i);
231  }
232 
233 
234 
235  void
236  MatrixBase::compress (const VectorOperation::values operation)
237  {
238 #ifdef DEBUG
239 #ifdef DEAL_II_WITH_MPI
240  // Check that all processors agree that last_action is the same (or none!)
241 
242  int my_int_last_action = last_action;
243  int all_int_last_action;
244 
245  MPI_Allreduce(&my_int_last_action, &all_int_last_action, 1, MPI_INT,
246  MPI_BOR, get_mpi_communicator());
247 
248  AssertThrow(all_int_last_action != (VectorOperation::add | VectorOperation::insert),
249  ExcMessage("Error: not all processors agree on the last VectorOperation before this compress() call."));
250 #endif
251 #endif
252 
253  AssertThrow(last_action == VectorOperation::unknown
254  || last_action == operation,
255  ExcMessage("Missing compress() or calling with wrong VectorOperation argument."));
256 
257  // flush buffers
258  int ierr;
259  ierr = MatAssemblyBegin (matrix,MAT_FINAL_ASSEMBLY);
260  AssertThrow (ierr == 0, ExcPETScError(ierr));
261 
262  ierr = MatAssemblyEnd (matrix,MAT_FINAL_ASSEMBLY);
263  AssertThrow (ierr == 0, ExcPETScError(ierr));
264 
265  last_action = VectorOperation::unknown;
266  }
267 
268 
269 
271  MatrixBase::m () const
272  {
273  PetscInt n_rows, n_cols;
274 
275  int ierr = MatGetSize (matrix, &n_rows, &n_cols);
276  AssertThrow (ierr == 0, ExcPETScError(ierr));
277 
278  return n_rows;
279  }
280 
281 
282 
284  MatrixBase::n () const
285  {
286  PetscInt n_rows, n_cols;
287 
288  int ierr = MatGetSize (matrix, &n_rows, &n_cols);
289  AssertThrow (ierr == 0, ExcPETScError(ierr));
290 
291  return n_cols;
292  }
293 
294 
295 
298  {
299  PetscInt n_rows, n_cols;
300 
301  int ierr = MatGetLocalSize (matrix, &n_rows, &n_cols);
302  AssertThrow (ierr == 0, ExcPETScError(ierr));
303 
304  return n_rows;
305  }
306 
307 
308 
309  std::pair<MatrixBase::size_type, MatrixBase::size_type>
311  {
312  PetscInt begin, end;
313 
314  const int ierr = MatGetOwnershipRange (static_cast<const Mat &>(matrix),
315  &begin, &end);
316  AssertThrow (ierr == 0, ExcPETScError(ierr));
317 
318  return std::make_pair (begin, end);
319  }
320 
321 
322 
325  {
326  MatInfo mat_info;
327  const int ierr
328  = MatGetInfo (matrix, MAT_GLOBAL_SUM, &mat_info);
329  AssertThrow (ierr == 0, ExcPETScError(ierr));
330 
331  return static_cast<size_type>(mat_info.nz_used);
332  }
333 
334 
335 
338  row_length (const size_type row) const
339  {
340 //TODO: this function will probably only work if compress() was called on the
341 //matrix previously. however, we can't do this here, since it would impose
342 //global communication and one would have to make sure that this function is
343 //called the same number of times from all processors, something that is
344 //unreasonable. there should simply be a way in PETSc to query the number of
345 //entries in a row bypassing the call to compress(), but I can't find one
346  Assert (row < m(), ExcInternalError());
347 
348  // get a representation of the present
349  // row
350  PetscInt ncols;
351  const PetscInt *colnums;
352  const PetscScalar *values;
353 
354 //TODO: this is probably horribly inefficient; we should lobby for a way to
355 //query this information from PETSc
356  int ierr;
357  ierr = MatGetRow(*this, row, &ncols, &colnums, &values);
358  AssertThrow (ierr == 0, MatrixBase::ExcPETScError(ierr));
359 
360  // then restore the matrix and return the number of columns in this row as
361  // queried previously. Starting with PETSc 3.4, MatRestoreRow actually
362  // resets the last three arguments to zero/NULL, to avoid abuse of pointers
363  // now dangling. as a consequence, we need to save the size of the array
364  // and return the saved value.
365  const PetscInt ncols_saved = ncols;
366  ierr = MatRestoreRow(*this, row, &ncols, &colnums, &values);
367  AssertThrow (ierr == 0, MatrixBase::ExcPETScError(ierr));
368 
369  return ncols_saved;
370  }
371 
372 
373  PetscReal
375  {
376  PetscReal result;
377 
378  const int ierr
379  = MatNorm (matrix, NORM_1, &result);
380  AssertThrow (ierr == 0, ExcPETScError(ierr));
381 
382  return result;
383  }
384 
385 
386 
387  PetscReal
389  {
390  PetscReal result;
391 
392  const int ierr
393  = MatNorm (matrix, NORM_INFINITY, &result);
394  AssertThrow (ierr == 0, ExcPETScError(ierr));
395 
396  return result;
397  }
398 
399 
400 
401  PetscReal
403  {
404  PetscReal result;
405 
406  const int ierr
407  = MatNorm (matrix, NORM_FROBENIUS, &result);
408  AssertThrow (ierr == 0, ExcPETScError(ierr));
409 
410  return result;
411  }
412 
413 
414  PetscScalar
416  {
417  Vector tmp(v.size());
418  vmult (tmp, v);
419  return tmp*v;
420  }
421 
422 
423  PetscScalar
425  const VectorBase &v) const
426  {
427  Vector tmp(v.size());
428  vmult (tmp, v);
429  return u*tmp;
430  }
431 
432 
433 #if DEAL_II_PETSC_VERSION_GTE(3,1,0)
434  PetscScalar
435  MatrixBase::trace () const
436  {
437  PetscScalar result;
438 
439  const int ierr
440  = MatGetTrace (matrix, &result);
441  AssertThrow (ierr == 0, ExcPETScError(ierr));
442 
443  return result;
444  }
445 #endif
446 
447 
448 
449  MatrixBase &
450  MatrixBase::operator *= (const PetscScalar a)
451  {
452  const int ierr = MatScale (matrix, a);
453  AssertThrow (ierr == 0, ExcPETScError(ierr));
454 
455  return *this;
456  }
457 
458 
459 
460  MatrixBase &
461  MatrixBase::operator /= (const PetscScalar a)
462  {
463  const PetscScalar factor = 1./a;
464  const int ierr = MatScale (matrix, factor);
465 
466  AssertThrow (ierr == 0, ExcPETScError(ierr));
467 
468  return *this;
469  }
470 
471 
472  MatrixBase &
474  const PetscScalar factor)
475  {
476  const int ierr = MatAXPY (matrix, factor,
477  other, DIFFERENT_NONZERO_PATTERN);
478  (void)ierr;
479 
480  Assert (ierr == 0, ExcPETScError(ierr));
481 
482  return *this;
483  }
484 
485 
486  void
488  const VectorBase &src) const
489  {
490  Assert (&src != &dst, ExcSourceEqualsDestination());
491 
492  const int ierr = MatMult (matrix, src, dst);
493  (void)ierr;
494  AssertThrow (ierr == 0, ExcPETScError(ierr));
495  }
496 
497 
498 
499  void
501  const VectorBase &src) const
502  {
503  Assert (&src != &dst, ExcSourceEqualsDestination());
504 
505  const int ierr = MatMultTranspose (matrix, src, dst);
506  (void)ierr;
507  AssertThrow (ierr == 0, ExcPETScError(ierr));
508  }
509 
510 
511 
512  void
514  const VectorBase &src) const
515  {
516  Assert (&src != &dst, ExcSourceEqualsDestination());
517 
518  const int ierr = MatMultAdd (matrix, src, dst, dst);
519  (void)ierr;
520  AssertThrow (ierr == 0, ExcPETScError(ierr));
521  }
522 
523 
524 
525  void
527  const VectorBase &src) const
528  {
529  Assert (&src != &dst, ExcSourceEqualsDestination());
530 
531  const int ierr = MatMultTransposeAdd (matrix, src, dst, dst);
532  (void)ierr;
533  AssertThrow (ierr == 0, ExcPETScError(ierr));
534  }
535 
536 
537  PetscScalar
539  const VectorBase &x,
540  const VectorBase &b) const
541  {
542  // avoid the use of a temporary, and
543  // rather do one negation pass more than
544  // necessary
545  vmult (dst, x);
546  dst -= b;
547  dst *= -1;
548 
549  return dst.l2_norm();
550  }
551 
552 
553 
554  MatrixBase::operator Mat () const
555  {
556  return matrix;
557  }
558 
559  void
561  {
562  int ierr = MatTranspose(matrix, MAT_REUSE_MATRIX, &matrix);
563  (void)ierr;
564  AssertThrow (ierr == 0, ExcPETScError(ierr));
565  }
566 
567 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
568  PetscTruth
569 #else
570  PetscBool
571 #endif
572  MatrixBase::is_symmetric (const double tolerance)
573  {
574 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
575  PetscTruth
576 #else
577  PetscBool
578 #endif
579  truth;
581  MatIsSymmetric (matrix, tolerance, &truth);
582  return truth;
583  }
584 
585 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
586  PetscTruth
587 #else
588  PetscBool
589 #endif
590  MatrixBase::is_hermitian (const double tolerance)
591  {
592 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
593  PetscTruth
594 #else
595  PetscBool
596 #endif
597  truth;
598 
600  MatIsHermitian (matrix, tolerance, &truth);
601 
602  return truth;
603  }
604 
605  void
606  MatrixBase::write_ascii (const PetscViewerFormat format)
607  {
609 
610  // Set options
611  PetscViewerSetFormat (PETSC_VIEWER_STDOUT_WORLD,
612  format);
613 
614  // Write to screen
615  MatView (matrix, PETSC_VIEWER_STDOUT_WORLD);
616  }
617 
618  void
619  MatrixBase::print (std::ostream &out,
620  const bool /*alternative_output*/) const
621  {
622  std::pair<MatrixBase::size_type, MatrixBase::size_type>
623  loc_range = local_range();
624 
625  PetscInt ncols;
626  const PetscInt *colnums;
627  const PetscScalar *values;
628 
630  for (row = loc_range.first; row < loc_range.second; ++row)
631  {
632  int ierr = MatGetRow(*this, row, &ncols, &colnums, &values);
633  (void)ierr;
634  AssertThrow (ierr == 0, MatrixBase::ExcPETScError(ierr));
635 
636  for (PetscInt col = 0; col < ncols; ++col)
637  {
638  out << "(" << row << "," << colnums[col] << ") " << values[col] << std::endl;
639  }
640 
641  ierr = MatRestoreRow(*this, row, &ncols, &colnums, &values);
642  AssertThrow (ierr == 0, MatrixBase::ExcPETScError(ierr));
643  }
644 
645  AssertThrow (out, ExcIO());
646  }
647 
648 
649 
650  std::size_t
652  {
653  MatInfo info;
654  MatGetInfo(matrix, MAT_LOCAL, &info);
655 
656  return sizeof(*this) + static_cast<size_type>(info.memory);
657  }
658 
659 }
660 
661 DEAL_II_NAMESPACE_CLOSE
662 
663 #endif // DEAL_II_WITH_PETSC
PetscScalar matrix_scalar_product(const VectorBase &u, const VectorBase &v) const
void add(const size_type i, const size_type j, const PetscScalar value)
void vmult(VectorBase &dst, const VectorBase &src) const
const_iterator begin() const
PetscReal frobenius_norm() const
::ExceptionBase & ExcMessage(std::string arg1)
size_type row_length(const size_type row) const
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
PetscBool is_symmetric(const double tolerance=1.e-12)
PetscBool is_hermitian(const double tolerance=1.e-12)
void clear_rows(const std::vector< size_type > &rows, const PetscScalar new_diag_value=0)
void compress(const VectorOperation::values operation)
size_type n_nonzero_elements() const
#define Assert(cond, exc)
Definition: exceptions.h:294
MatrixBase & operator=(const value_type d)
void vmult_add(VectorBase &dst, const VectorBase &src) const
PetscScalar residual(VectorBase &dst, const VectorBase &x, const VectorBase &b) const
types::global_dof_index size_type
PetscScalar diag_element(const size_type i) const
void Tvmult(VectorBase &dst, const VectorBase &src) const
PetscScalar el(const size_type i, const size_type j) const
void Tvmult_add(VectorBase &dst, const VectorBase &src) const
const_iterator end() const
void print(std::ostream &out, const bool alternative_output=false) const
void write_ascii(const PetscViewerFormat format=PETSC_VIEWER_DEFAULT)
VectorOperation::values last_action
MatrixBase & operator*=(const PetscScalar factor)
std::pair< size_type, size_type > local_range() const
PetscScalar matrix_norm_square(const VectorBase &v) const
MatrixBase & operator/=(const PetscScalar factor)
virtual const MPI_Comm & get_mpi_communicator() const =0
void clear_row(const size_type row, const PetscScalar new_diag_value=0)
std::size_t memory_consumption() const