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