Reference documentation for deal.II version 8.5.1
petsc_parallel_sparse_matrix.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2017 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_parallel_sparse_matrix.h>
17 
18 #ifdef DEAL_II_WITH_PETSC
19 
20 # include <deal.II/base/mpi.h>
21 # include <deal.II/lac/exceptions.h>
22 # include <deal.II/lac/petsc_compatibility.h>
23 # include <deal.II/lac/petsc_vector.h>
24 # include <deal.II/lac/sparsity_pattern.h>
25 # include <deal.II/lac/dynamic_sparsity_pattern.h>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 namespace PETScWrappers
30 {
31  namespace MPI
32  {
34  : communicator(MPI_COMM_SELF)
35  {
36  // just like for vectors: since we
37  // create an empty matrix, we can as
38  // well make it sequential
39  const int m=0, n=0, n_nonzero_per_row=0;
40  const PetscErrorCode ierr
41  = MatCreateSeqAIJ(PETSC_COMM_SELF, m, n, n_nonzero_per_row,
42  0, &matrix);
43  AssertThrow (ierr == 0, ExcPETScError(ierr));
44  }
45 
46 
48  {
50  }
51 
52  SparseMatrix::SparseMatrix (const MPI_Comm &communicator,
53  const size_type m,
54  const size_type n,
55  const size_type local_rows,
56  const size_type local_columns,
57  const size_type n_nonzero_per_row,
58  const bool is_symmetric,
59  const size_type n_offdiag_nonzero_per_row)
60  :
61  communicator (communicator)
62  {
63  do_reinit (m, n, local_rows, local_columns,
64  n_nonzero_per_row, is_symmetric,
65  n_offdiag_nonzero_per_row);
66  }
67 
68 
69 
70  SparseMatrix::SparseMatrix (const MPI_Comm &communicator,
71  const size_type m,
72  const size_type n,
73  const size_type local_rows,
74  const size_type local_columns,
75  const std::vector<size_type> &row_lengths,
76  const bool is_symmetric,
77  const std::vector<size_type> &offdiag_row_lengths)
78  :
79  communicator (communicator)
80  {
81  do_reinit (m, n, local_rows, local_columns,
82  row_lengths, is_symmetric, offdiag_row_lengths);
83  }
84 
85 
86 
87  template <typename SparsityPatternType>
89  SparseMatrix (const MPI_Comm &communicator,
90  const SparsityPatternType &sparsity_pattern,
91  const std::vector<size_type> &local_rows_per_process,
92  const std::vector<size_type> &local_columns_per_process,
93  const unsigned int this_process,
94  const bool preset_nonzero_locations)
95  :
96  communicator (communicator)
97  {
98  do_reinit (sparsity_pattern, local_rows_per_process,
99  local_columns_per_process, this_process,
100  preset_nonzero_locations);
101  }
102 
103 
104  void
106  reinit (const SparseMatrix &other)
107  {
108  if (&other == this)
109  return;
110 
111  this->communicator = other.communicator;
112 
113  PetscErrorCode ierr = destroy_matrix (matrix);
114  AssertThrow (ierr == 0, ExcPETScError(ierr));
115 
116  ierr = MatDuplicate (other.matrix, MAT_DO_NOT_COPY_VALUES, &matrix);
117  AssertThrow (ierr == 0, ExcPETScError(ierr));
118  }
119 
120 
121  SparseMatrix &
123  {
125  return *this;
126  }
127 
128  void
130  {
131  if (&other == this)
132  return;
133 
134  this->communicator = other.communicator;
135 
136  const PetscErrorCode ierr = MatCopy(other.matrix, matrix,
137  SAME_NONZERO_PATTERN);
138  AssertThrow (ierr == 0, ExcPETScError(ierr));
139  }
140 
141  void
142  SparseMatrix::reinit (const MPI_Comm &communicator,
143  const size_type m,
144  const size_type n,
145  const size_type local_rows,
146  const size_type local_columns,
147  const size_type n_nonzero_per_row,
148  const bool is_symmetric,
149  const size_type n_offdiag_nonzero_per_row)
150  {
151  this->communicator = communicator;
152 
153  // get rid of old matrix and generate a new one
154  const PetscErrorCode ierr = destroy_matrix (matrix);
155  AssertThrow (ierr == 0, ExcPETScError (ierr));
156 
157  do_reinit (m, n, local_rows, local_columns,
158  n_nonzero_per_row, is_symmetric,
159  n_offdiag_nonzero_per_row);
160  }
161 
162 
163 
164  void
165  SparseMatrix::reinit (const MPI_Comm &communicator,
166  const size_type m,
167  const size_type n,
168  const size_type local_rows,
169  const size_type local_columns,
170  const std::vector<size_type> &row_lengths,
171  const bool is_symmetric,
172  const std::vector<size_type> &offdiag_row_lengths)
173  {
174  this->communicator = communicator;
175 
176  // get rid of old matrix and generate a
177  // new one
178  const PetscErrorCode ierr = destroy_matrix (matrix);
179  AssertThrow (ierr == 0, ExcPETScError (ierr));
180 
181  do_reinit (m, n, local_rows, local_columns,
182  row_lengths, is_symmetric, offdiag_row_lengths);
183  }
184 
185 
186 
187  template <typename SparsityPatternType>
188  void
190  reinit (const MPI_Comm &communicator,
191  const SparsityPatternType &sparsity_pattern,
192  const std::vector<size_type> &local_rows_per_process,
193  const std::vector<size_type> &local_columns_per_process,
194  const unsigned int this_process,
195  const bool preset_nonzero_locations)
196  {
197  this->communicator = communicator;
198 
199  // get rid of old matrix and generate a new one
200  const PetscErrorCode ierr = destroy_matrix (matrix);
201  AssertThrow (ierr == 0, ExcPETScError (ierr));
202 
203 
204  do_reinit (sparsity_pattern, local_rows_per_process,
205  local_columns_per_process, this_process,
206  preset_nonzero_locations);
207  }
208 
209  template <typename SparsityPatternType>
210  void
212  reinit (const IndexSet &local_rows,
213  const IndexSet &local_columns,
214  const SparsityPatternType &sparsity_pattern,
215  const MPI_Comm &communicator)
216  {
217  this->communicator = communicator;
218 
219  // get rid of old matrix and generate a new one
220  const PetscErrorCode ierr = destroy_matrix (matrix);
221  AssertThrow(ierr == 0, ExcPETScError (ierr));
222 
223  do_reinit (local_rows, local_columns, sparsity_pattern);
224  }
225 
226  void
228  const size_type n,
229  const size_type local_rows,
230  const size_type local_columns,
231  const size_type n_nonzero_per_row,
232  const bool is_symmetric,
233  const size_type n_offdiag_nonzero_per_row)
234  {
235  Assert (local_rows <= m, ExcLocalRowsTooLarge (local_rows, m));
236 
237  // use the call sequence indicating only
238  // a maximal number of elements per row
239  // for all rows globally
240 #if DEAL_II_PETSC_VERSION_LT(3,3,0)
241  const PetscErrorCode ierr = MatCreateMPIAIJ
242  (communicator,
243  local_rows, local_columns,
244  m, n,
245  n_nonzero_per_row, 0,
246  n_offdiag_nonzero_per_row, 0,
247  &matrix);
248 #else
249  const PetscErrorCode ierr = MatCreateAIJ
250  (communicator,
251  local_rows, local_columns,
252  m, n,
253  n_nonzero_per_row, 0,
254  n_offdiag_nonzero_per_row, 0,
255  &matrix);
256  set_matrix_option (matrix, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
257 #endif
258  AssertThrow (ierr == 0, ExcPETScError(ierr));
259 
260  // set symmetric flag, if so requested
261  if (is_symmetric == true)
262  {
263  set_matrix_option (matrix, MAT_SYMMETRIC, PETSC_TRUE);
264  }
265  }
266 
267 
268 
269  void
271  const size_type n,
272  const size_type local_rows,
273  const size_type local_columns,
274  const std::vector<size_type> &row_lengths,
275  const bool is_symmetric,
276  const std::vector<size_type> &offdiag_row_lengths)
277  {
278  Assert (local_rows <= m, ExcLocalRowsTooLarge (local_rows, m));
279 
280  Assert (row_lengths.size() == m,
281  ExcDimensionMismatch (row_lengths.size(), m));
282 
283  // For the case that
284  // local_columns is smaller
285  // than one of the row lengths
286  // MatCreateMPIAIJ throws an
287  // error. In this case use a
288  // PETScWrappers::SparseMatrix
289  for (size_type i=0; i<row_lengths.size(); ++i)
290  Assert(row_lengths[i]<=local_columns,
291  ExcIndexRange(row_lengths[i], 1, local_columns+1));
292 
293  // use the call sequence indicating a
294  // maximal number of elements for each
295  // row individually. annoyingly, we
296  // always use unsigned ints for cases
297  // like this, while PETSc wants to see
298  // signed integers. so we have to
299  // convert, unless we want to play dirty
300  // tricks with conversions of pointers
301  const std::vector<PetscInt> int_row_lengths (row_lengths.begin(),
302  row_lengths.end());
303  const std::vector<PetscInt> int_offdiag_row_lengths (offdiag_row_lengths.begin(),
304  offdiag_row_lengths.end());
305 
306 //TODO: There must be a significantly better way to provide information about the off-diagonal blocks of the matrix. this way, petsc keeps allocating tiny chunks of memory, and gets completely hung up over this
307 #if DEAL_II_PETSC_VERSION_LT(3,3,0)
308  const PetscErrorCode ierr = MatCreateMPIAIJ
309  (communicator,
310  local_rows, local_columns,
311  m, n,
312  0, &int_row_lengths[0],
313  0, offdiag_row_lengths.size() ? &int_offdiag_row_lengths[0] : 0,
314  &matrix);
315 #else
316  const PetscErrorCode ierr = MatCreateAIJ
317  (communicator,
318  local_rows, local_columns,
319  m, n,
320  0, &int_row_lengths[0],
321  0,
322  offdiag_row_lengths.size() ? &int_offdiag_row_lengths[0] : 0,
323  &matrix);
324 
325 //TODO: Sometimes the actual number of nonzero entries allocated is greater than the number of nonzero entries, which petsc will complain about unless explicitly disabled with MatSetOption. There is probably a way to prevent a different number nonzero elements being allocated in the first place. (See also previous TODO).
326  set_matrix_option (matrix, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
327 #endif
328  AssertThrow (ierr == 0, ExcPETScError(ierr));
329 
330  // set symmetric flag, if so requested
331  if (is_symmetric == true)
332  {
333  set_matrix_option (matrix, MAT_SYMMETRIC, PETSC_TRUE);
334  }
335  }
336 
337 
338  template <typename SparsityPatternType>
339  void
341  do_reinit (const IndexSet &local_rows,
342  const IndexSet &local_columns,
343  const SparsityPatternType &sparsity_pattern)
344  {
345  Assert(sparsity_pattern.n_rows()==local_rows.size(),
346  ExcMessage("SparsityPattern and IndexSet have different number of rows"));
347  Assert(sparsity_pattern.n_cols()==local_columns.size(),
348  ExcMessage("SparsityPattern and IndexSet have different number of columns"));
349  Assert(local_rows.is_contiguous() && local_columns.is_contiguous(),
350  ExcMessage("PETSc only supports contiguous row/column ranges"));
352 
353 #ifdef DEBUG
354  {
355  // check indexsets
358  Assert(row_owners == sparsity_pattern.n_rows(),
359  ExcMessage(std::string("Each row has to be owned by exactly one owner (n_rows()=")
360  + Utilities::to_string(sparsity_pattern.n_rows())
361  + " but sum(local_rows.n_elements())="
362  + Utilities::to_string(row_owners)
363  + ")"));
364  Assert(col_owners == sparsity_pattern.n_cols(),
365  ExcMessage(std::string("Each column has to be owned by exactly one owner (n_cols()=")
366  + Utilities::to_string(sparsity_pattern.n_cols())
367  + " but sum(local_columns.n_elements())="
368  + Utilities::to_string(col_owners)
369  + ")"));
370  }
371 #endif
372 
373 
374  // create the matrix. We do not set row length but set the
375  // correct SparsityPattern later.
376  PetscErrorCode ierr = MatCreate(communicator,&matrix);
377  AssertThrow (ierr == 0, ExcPETScError(ierr));
378 
379  ierr = MatSetSizes(matrix,
380  local_rows.n_elements(),
381  local_columns.n_elements(),
382  sparsity_pattern.n_rows(),
383  sparsity_pattern.n_cols());
384  AssertThrow (ierr == 0, ExcPETScError(ierr));
385 
386  ierr = MatSetType(matrix,MATMPIAIJ);
387  AssertThrow (ierr == 0, ExcPETScError(ierr));
388 
389 
390  // next preset the exact given matrix
391  // entries with zeros. this doesn't avoid any
392  // memory allocations, but it at least
393  // avoids some searches later on. the
394  // key here is that we can use the
395  // matrix set routines that set an
396  // entire row at once, not a single
397  // entry at a time
398  //
399  // for the usefulness of this option
400  // read the documentation of this
401  // class.
402  //if (preset_nonzero_locations == true)
403  if (local_rows.n_elements()>0)
404  {
405  Assert(local_columns.n_elements()>0, ExcInternalError());
406  // MatMPIAIJSetPreallocationCSR
407  // can be used to allocate the sparsity
408  // pattern of a matrix
409 
410  const PetscInt local_row_start = local_rows.nth_index_in_set(0);
411  const PetscInt
412  local_row_end = local_row_start + local_rows.n_elements();
413 
414 
415  // first set up the column number
416  // array for the rows to be stored
417  // on the local processor. have one
418  // dummy entry at the end to make
419  // sure petsc doesn't read past the
420  // end
421  std::vector<PetscInt>
422 
423  rowstart_in_window (local_row_end - local_row_start + 1, 0),
424  colnums_in_window;
425  {
426  unsigned int n_cols = 0;
427  for (PetscInt i=local_row_start; i<local_row_end; ++i)
428  {
429  const PetscInt row_length = sparsity_pattern.row_length(i);
430  rowstart_in_window[i+1-local_row_start]
431  = rowstart_in_window[i-local_row_start] + row_length;
432  n_cols += row_length;
433  }
434  colnums_in_window.resize (n_cols+1, -1);
435  }
436 
437  // now copy over the information
438  // from the sparsity pattern.
439  {
440  PetscInt *ptr = & colnums_in_window[0];
441  for (PetscInt i=local_row_start; i<local_row_end; ++i)
442  for (typename SparsityPatternType::iterator p=sparsity_pattern.begin(i);
443  p != sparsity_pattern.end(i); ++p, ++ptr)
444  *ptr = p->column();
445  }
446 
447 
448  // then call the petsc function
449  // that summarily allocates these
450  // entries:
451  ierr = MatMPIAIJSetPreallocationCSR (matrix,
452  &rowstart_in_window[0],
453  &colnums_in_window[0],
454  0);
455  AssertThrow (ierr == 0, ExcPETScError(ierr));
456  }
457  else
458  {
459  PetscInt i=0;
460  ierr = MatMPIAIJSetPreallocationCSR (matrix,
461  &i,
462  &i,
463  0);
464  AssertThrow (ierr == 0, ExcPETScError(ierr));
465  }
467 
468  {
471  }
472  }
473 
474 
475  template <typename SparsityPatternType>
476  void
478  do_reinit (const SparsityPatternType &sparsity_pattern,
479  const std::vector<size_type> &local_rows_per_process,
480  const std::vector<size_type> &local_columns_per_process,
481  const unsigned int this_process,
482  const bool preset_nonzero_locations)
483  {
484  Assert (local_rows_per_process.size() == local_columns_per_process.size(),
485  ExcDimensionMismatch (local_rows_per_process.size(),
486  local_columns_per_process.size()));
487  Assert (this_process < local_rows_per_process.size(),
488  ExcInternalError());
490 
491  // for each row that we own locally, we
492  // have to count how many of the
493  // entries in the sparsity pattern lie
494  // in the column area we have locally,
495  // and how many arent. for this, we
496  // first have to know which areas are
497  // ours
498  size_type local_row_start = 0;
499  size_type local_col_start = 0;
500  for (unsigned int p=0; p<this_process; ++p)
501  {
502  local_row_start += local_rows_per_process[p];
503  local_col_start += local_columns_per_process[p];
504  }
505  const size_type
506  local_row_end = local_row_start + local_rows_per_process[this_process];
507 
508 #if DEAL_II_PETSC_VERSION_LT(2,3,3)
509  //old version to create the matrix, we
510  //can skip calculating the row length
511  //at least starting from 2.3.3 (tested,
512  //see below)
513 
514  const size_type
515  local_col_end = local_col_start + local_columns_per_process[this_process];
516 
517  // then count the elements in- and
518  // out-of-window for the rows we own
519  std::vector<PetscInt>
520 
521  row_lengths_in_window (local_row_end - local_row_start),
522  row_lengths_out_of_window (local_row_end - local_row_start);
523  for (size_type row = local_row_start; row<local_row_end; ++row)
524  for (size_type c=0; c<sparsity_pattern.row_length(row); ++c)
525  {
526  const size_type column = sparsity_pattern.column_number(row,c);
527 
528  if ((column >= local_col_start) &&
529  (column < local_col_end))
530  ++row_lengths_in_window[row-local_row_start];
531  else
532  ++row_lengths_out_of_window[row-local_row_start];
533  }
534 
535 
536  // create the matrix. completely
537  // confusingly, PETSc wants us to pass
538  // arrays for the local number of
539  // elements that starts with zero for
540  // the first _local_ row, i.e. it
541  // doesn't index into an array for
542  // _all_ rows.
543  const PetscErrorCode ierr = MatCreateMPIAIJ(communicator,
544  local_rows_per_process[this_process],
545  local_columns_per_process[this_process],
546  sparsity_pattern.n_rows(),
547  sparsity_pattern.n_cols(),
548  0, &row_lengths_in_window[0],
549  0, &row_lengths_out_of_window[0],
550  &matrix);
551  AssertThrow (ierr == 0, ExcPETScError(ierr));
552 
553 #else //PETSC_VERSION>=2.3.3
554  // create the matrix. We
555  // do not set row length but set the
556  // correct SparsityPattern later.
557  PetscErrorCode ierr = MatCreate(communicator,&matrix);
558  AssertThrow (ierr == 0, ExcPETScError(ierr));
559 
560  ierr = MatSetSizes(matrix,
561  local_rows_per_process[this_process],
562  local_columns_per_process[this_process],
563  sparsity_pattern.n_rows(),
564  sparsity_pattern.n_cols());
565  AssertThrow (ierr == 0, ExcPETScError(ierr));
566 
567  ierr = MatSetType(matrix,MATMPIAIJ);
568  AssertThrow (ierr == 0, ExcPETScError(ierr));
569 #endif
570 
571 
572  // next preset the exact given matrix
573  // entries with zeros, if the user
574  // requested so. this doesn't avoid any
575  // memory allocations, but it at least
576  // avoids some searches later on. the
577  // key here is that we can use the
578  // matrix set routines that set an
579  // entire row at once, not a single
580  // entry at a time
581  //
582  // for the usefulness of this option
583  // read the documentation of this
584  // class.
585  if (preset_nonzero_locations == true)
586  {
587  // MatMPIAIJSetPreallocationCSR
588  // can be used to allocate the sparsity
589  // pattern of a matrix if it is already
590  // available:
591 
592  // first set up the column number
593  // array for the rows to be stored
594  // on the local processor. have one
595  // dummy entry at the end to make
596  // sure petsc doesn't read past the
597  // end
598  std::vector<PetscInt>
599 
600  rowstart_in_window (local_row_end - local_row_start + 1, 0),
601  colnums_in_window;
602  {
603  size_type n_cols = 0;
604  for (size_type i=local_row_start; i<local_row_end; ++i)
605  {
606  const size_type row_length = sparsity_pattern.row_length(i);
607  rowstart_in_window[i+1-local_row_start]
608  = rowstart_in_window[i-local_row_start] + row_length;
609  n_cols += row_length;
610  }
611  colnums_in_window.resize (n_cols+1, -1);
612  }
613 
614  // now copy over the information
615  // from the sparsity pattern.
616  {
617  PetscInt *ptr = & colnums_in_window[0];
618  for (size_type i=local_row_start; i<local_row_end; ++i)
619  for (typename SparsityPatternType::iterator p=sparsity_pattern.begin(i);
620  p != sparsity_pattern.end(i); ++p, ++ptr)
621  *ptr = p->column();
622  }
623 
624 
625  // then call the petsc function
626  // that summarily allocates these
627  // entries:
628  ierr = MatMPIAIJSetPreallocationCSR (matrix,
629  &rowstart_in_window[0],
630  &colnums_in_window[0],
631  0);
632  AssertThrow (ierr == 0, ExcPETScError(ierr));
633 
634 #if DEAL_II_PETSC_VERSION_LT(2,3,3)
635  // this is only needed for old
636  // PETSc versions:
637 
638  // for some reason, it does not
639  // seem to be possible to force
640  // actual allocation of actual
641  // entries by using the last
642  // arguments to the call above. if
643  // we don't initialize the entries
644  // like in the following loop, then
645  // the program is unbearably slow
646  // because elements are allocated
647  // and accessed in random order,
648  // which is not what PETSc likes
649  //
650  // note that we actually have to
651  // set the entries to something
652  // non-zero! do the allocation one
653  // row at a time
654  {
655  const std::vector<PetscScalar>
656  values (sparsity_pattern.max_entries_per_row(),
657  1.);
658 
659  for (size_type i=local_row_start; i<local_row_end; ++i)
660  {
661  PetscInt petsc_i = i;
662  ierr = MatSetValues (matrix, 1, &petsc_i,
663  sparsity_pattern.row_length(i),
664  &colnums_in_window[rowstart_in_window[i-local_row_start]],
665  &values[0], INSERT_VALUES);
666  AssertThrow (ierr == 0, ExcPETScError(ierr));
667  }
668  }
669 
671 
672  // set the dummy entries set above
673  // back to zero
674  *this = 0;
675 #endif // version <=2.3.3
676 
679  }
680  }
681 
682  // explicit instantiations
683  //
684  template
685  SparseMatrix::SparseMatrix (const MPI_Comm &,
686  const SparsityPattern &,
687  const std::vector<size_type> &,
688  const std::vector<size_type> &,
689  const unsigned int,
690  const bool);
691  template
692  SparseMatrix::SparseMatrix (const MPI_Comm &,
693  const DynamicSparsityPattern &,
694  const std::vector<size_type> &,
695  const std::vector<size_type> &,
696  const unsigned int,
697  const bool);
698 
699  template void
700  SparseMatrix::reinit (const MPI_Comm &,
701  const SparsityPattern &,
702  const std::vector<size_type> &,
703  const std::vector<size_type> &,
704  const unsigned int,
705  const bool);
706  template void
707  SparseMatrix::reinit (const MPI_Comm &,
708  const DynamicSparsityPattern &,
709  const std::vector<size_type> &,
710  const std::vector<size_type> &,
711  const unsigned int,
712  const bool);
713 
714  template void
716  reinit (const IndexSet &,
717  const IndexSet &,
718  const SparsityPattern &,
719  const MPI_Comm &);
720 
721  template void
723  reinit (const IndexSet &,
724  const IndexSet &,
725  const DynamicSparsityPattern &,
726  const MPI_Comm &);
727 
728  template void
730  const std::vector<size_type> &,
731  const std::vector<size_type> &,
732  const unsigned int ,
733  const bool);
734  template void
736  const std::vector<size_type> &,
737  const std::vector<size_type> &,
738  const unsigned int ,
739  const bool);
740 
741  template void
743  do_reinit (const IndexSet &,
744  const IndexSet &,
745  const SparsityPattern &);
746 
747  template void
749  do_reinit (const IndexSet &,
750  const IndexSet &,
751  const DynamicSparsityPattern &);
752 
753 
754  PetscScalar
756  {
757  Vector tmp (v);
758  vmult (tmp, v);
759  // note, that v*tmp returns sum_i conjugate(v)_i * tmp_i
760  return v*tmp;
761  }
762 
763  PetscScalar
765  const Vector &v) const
766  {
767  Vector tmp (v);
768  vmult (tmp, v);
769  // note, that v*tmp returns sum_i conjugate(v)_i * tmp_i
770  return u*tmp;
771  }
772 
773  IndexSet
775  {
776 #if DEAL_II_PETSC_VERSION_LT(3,3,0)
777  Assert(false,ExcNotImplemented());
778  return IndexSet();
779 #else
780  PetscInt n_rows, n_cols, n_loc_rows, n_loc_cols, min, max;
781  PetscErrorCode ierr;
782 
783  ierr = MatGetSize (matrix, &n_rows, &n_cols);
784  AssertThrow (ierr == 0, ExcPETScError(ierr));
785 
786  ierr = MatGetLocalSize(matrix, &n_loc_rows, &n_loc_cols);
787  AssertThrow (ierr == 0, ExcPETScError(ierr));
788 
789  ierr = MatGetOwnershipRangeColumn(matrix, &min, &max);
790  AssertThrow (ierr == 0, ExcPETScError(ierr));
791 
792  Assert(n_loc_cols==max-min, ExcMessage("PETSc is requiring non contiguous memory allocation."));
793 
794  IndexSet indices(n_cols);
795  indices.add_range(min, max);
796  indices.compress();
797 
798  return indices;
799 #endif
800  }
801 
802  IndexSet
804  {
805 #if DEAL_II_PETSC_VERSION_LT(3,3,0)
806  Assert(false,ExcNotImplemented());
807  return IndexSet();
808 #else
809  PetscInt n_rows, n_cols, n_loc_rows, n_loc_cols, min, max;
810  PetscErrorCode ierr;
811 
812  ierr = MatGetSize (matrix, &n_rows, &n_cols);
813  AssertThrow (ierr == 0, ExcPETScError(ierr));
814 
815  ierr = MatGetLocalSize(matrix, &n_loc_rows, &n_loc_cols);
816  AssertThrow (ierr == 0, ExcPETScError(ierr));
817 
818  ierr = MatGetOwnershipRange(matrix, &min, &max);
819  AssertThrow (ierr == 0, ExcPETScError(ierr));
820 
821  Assert(n_loc_rows==max-min, ExcMessage("PETSc is requiring non contiguous memory allocation."));
822 
823  IndexSet indices(n_rows);
824  indices.add_range(min, max);
825  indices.compress();
826 
827  return indices;
828 #endif
829  }
830 
831  }
832 }
833 
834 
835 DEAL_II_NAMESPACE_CLOSE
836 
837 #endif // DEAL_II_WITH_PETSC
static ::ExceptionBase & ExcLocalRowsTooLarge(int arg1, int arg2)
PetscErrorCode destroy_matrix(Mat &matrix)
void copy_from(const SparseMatrix &other)
size_type nth_index_in_set(const unsigned int local_index) const
Definition: index_set.h:1612
void vmult(VectorBase &dst, const VectorBase &src) const
void do_reinit(const size_type m, const size_type n, const size_type local_rows, const size_type local_columns, const size_type n_nonzero_per_row, const bool is_symmetric=false, const size_type n_offdiag_nonzero_per_row=0)
size_type row_length(const size_type row) const
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:92
bool is_ascending_and_one_to_one(const MPI_Comm &communicator) const
Definition: index_set.cc:604
size_type size() const
Definition: index_set.h:1419
void compress(const VectorOperation::values operation)
static ::ExceptionBase & ExcMessage(std::string arg1)
SparseMatrix & operator=(const value_type d)
PetscScalar matrix_scalar_product(const Vector &u, const Vector &v) const
T sum(const T &t, const MPI_Comm &mpi_communicator)
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
MatrixBase & operator=(const value_type d)
types::global_dof_index size_type
void set_keep_zero_rows(Mat &matrix)
void reinit(const MPI_Comm &communicator, const size_type m, const size_type n, const size_type local_rows, const size_type local_columns, const size_type n_nonzero_per_row, const bool is_symmetric=false, const size_type n_offdiag_nonzero_per_row=0)
void add_range(const size_type begin, const size_type end)
Definition: index_set.cc:88
PetscBooleanType is_symmetric(const double tolerance=1.e-12)
void compress() const
Definition: index_set.h:1428
bool is_contiguous() const
Definition: index_set.h:1542
PetscScalar matrix_norm_square(const Vector &v) const
void set_matrix_option(Mat &matrix, const MatOption option_name, const PetscBooleanType option_value=PETSC_FALSE)
static ::ExceptionBase & ExcNotImplemented()
size_type n_elements() const
Definition: index_set.h:1560
static ::ExceptionBase & ExcInternalError()
void close_matrix(Mat &matrix)