18 #ifdef DEAL_II_WITH_PETSC
35 : communicator(MPI_COMM_SELF)
40 const int m = 0,
n = 0, n_nonzero_per_row = 0;
41 const PetscErrorCode ierr = MatCreateSeqAIJ(
42 PETSC_COMM_SELF,
m,
n, n_nonzero_per_row,
nullptr, &
matrix);
58 const bool is_symmetric,
59 const size_type n_offdiag_nonzero_per_row)
60 : communicator(communicator)
68 n_offdiag_nonzero_per_row);
79 const std::vector<size_type> &row_lengths,
80 const bool is_symmetric,
81 const std::vector<size_type> &offdiag_row_lengths)
82 : communicator(communicator)
95 template <
typename SparsityPatternType>
98 const SparsityPatternType & sparsity_pattern,
99 const std::vector<size_type> &local_rows_per_process,
100 const std::vector<size_type> &local_columns_per_process,
101 const unsigned int this_process,
102 const bool preset_nonzero_locations)
103 : communicator(communicator)
106 local_rows_per_process,
107 local_columns_per_process,
109 preset_nonzero_locations);
124 ierr = MatDuplicate(other.
matrix, MAT_DO_NOT_COPY_VALUES, &
matrix);
144 const PetscErrorCode ierr =
156 const bool is_symmetric,
157 const size_type n_offdiag_nonzero_per_row)
171 n_offdiag_nonzero_per_row);
182 const std::vector<size_type> &row_lengths,
183 const bool is_symmetric,
184 const std::vector<size_type> &offdiag_row_lengths)
199 offdiag_row_lengths);
204 template <
typename SparsityPatternType>
208 const SparsityPatternType & sparsity_pattern,
209 const std::vector<size_type> &local_rows_per_process,
210 const std::vector<size_type> &local_columns_per_process,
211 const unsigned int this_process,
212 const bool preset_nonzero_locations)
222 local_rows_per_process,
223 local_columns_per_process,
225 preset_nonzero_locations);
228 template <
typename SparsityPatternType>
232 const SparsityPatternType &sparsity_pattern,
241 do_reinit(local_rows, local_columns, sparsity_pattern);
250 const bool is_symmetric,
251 const size_type n_offdiag_nonzero_per_row)
265 n_offdiag_nonzero_per_row,
285 const std::vector<size_type> &row_lengths,
286 const bool is_symmetric,
287 const std::vector<size_type> &offdiag_row_lengths)
291 Assert(row_lengths.size() ==
m,
312 const std::vector<PetscInt> int_row_lengths(row_lengths.begin(),
314 const std::vector<PetscInt> int_offdiag_row_lengths(
315 offdiag_row_lengths.begin(), offdiag_row_lengths.end());
320 const PetscErrorCode ierr =
327 int_row_lengths.data(),
329 offdiag_row_lengths.size() ?
330 int_offdiag_row_lengths.data() :
350 template <
typename SparsityPatternType>
354 const SparsityPatternType &sparsity_pattern)
356 Assert(sparsity_pattern.n_rows() == local_rows.
size(),
358 "SparsityPattern and IndexSet have different number of rows"));
360 sparsity_pattern.n_cols() == local_columns.
size(),
362 "SparsityPattern and IndexSet have different number of columns"));
364 ExcMessage(
"PETSc only supports contiguous row/column ranges"));
375 Assert(row_owners == sparsity_pattern.n_rows(),
378 "Each row has to be owned by exactly one owner (n_rows()=") +
380 " but sum(local_rows.n_elements())=" +
383 col_owners == sparsity_pattern.n_cols(),
386 "Each column has to be owned by exactly one owner (n_cols()=") +
388 " but sum(local_columns.n_elements())=" +
399 ierr = MatSetSizes(
matrix,
402 sparsity_pattern.n_rows(),
403 sparsity_pattern.n_cols());
406 ierr = MatSetType(
matrix, MATMPIAIJ);
430 const PetscInt local_row_end =
440 std::vector<PetscInt>
442 rowstart_in_window(local_row_end - local_row_start + 1, 0),
445 unsigned int n_cols = 0;
446 for (PetscInt i = local_row_start; i < local_row_end; ++i)
448 const PetscInt
row_length = sparsity_pattern.row_length(i);
449 rowstart_in_window[i + 1 - local_row_start] =
450 rowstart_in_window[i - local_row_start] +
row_length;
453 colnums_in_window.resize(n_cols + 1, -1);
459 PetscInt *ptr = colnums_in_window.data();
460 for (PetscInt i = local_row_start; i < local_row_end; ++i)
461 for (
typename SparsityPatternType::iterator p =
462 sparsity_pattern.begin(i);
463 p != sparsity_pattern.end(i);
472 ierr = MatMPIAIJSetPreallocationCSR(
matrix,
473 rowstart_in_window.data(),
474 colnums_in_window.data(),
481 ierr = MatMPIAIJSetPreallocationCSR(
matrix, &i, &i,
nullptr);
493 template <
typename SparsityPatternType>
496 const SparsityPatternType & sparsity_pattern,
497 const std::vector<size_type> &local_rows_per_process,
498 const std::vector<size_type> &local_columns_per_process,
499 const unsigned int this_process,
500 const bool preset_nonzero_locations)
502 Assert(local_rows_per_process.size() == local_columns_per_process.size(),
504 local_columns_per_process.size()));
517 for (
unsigned int p = 0; p < this_process; ++p)
519 local_row_start += local_rows_per_process[p];
520 local_col_start += local_columns_per_process[p];
523 local_row_start + local_rows_per_process[this_process];
531 ierr = MatSetSizes(
matrix,
532 local_rows_per_process[this_process],
533 local_columns_per_process[this_process],
534 sparsity_pattern.n_rows(),
535 sparsity_pattern.n_cols());
538 ierr = MatSetType(
matrix, MATMPIAIJ);
554 if (preset_nonzero_locations ==
true)
567 std::vector<PetscInt>
569 rowstart_in_window(local_row_end - local_row_start + 1, 0),
573 for (
size_type i = local_row_start; i < local_row_end; ++i)
576 rowstart_in_window[i + 1 - local_row_start] =
577 rowstart_in_window[i - local_row_start] +
row_length;
580 colnums_in_window.resize(n_cols + 1, -1);
586 PetscInt *ptr = colnums_in_window.data();
587 for (
size_type i = local_row_start; i < local_row_end; ++i)
588 for (
typename SparsityPatternType::iterator p =
589 sparsity_pattern.begin(i);
590 p != sparsity_pattern.end(i);
599 ierr = MatMPIAIJSetPreallocationCSR(
matrix,
600 rowstart_in_window.data(),
601 colnums_in_window.data(),
615 const std::vector<size_type> &,
616 const std::vector<size_type> &,
621 const std::vector<size_type> &,
622 const std::vector<size_type> &,
629 const std::vector<size_type> &,
630 const std::vector<size_type> &,
636 const std::vector<size_type> &,
637 const std::vector<size_type> &,
655 const std::vector<size_type> &,
656 const std::vector<size_type> &,
661 const std::vector<size_type> &,
662 const std::vector<size_type> &,
699 PetscInt n_rows, n_cols, n_loc_rows, n_loc_cols,
min,
max;
702 ierr = MatGetSize(
matrix, &n_rows, &n_cols);
705 ierr = MatGetLocalSize(
matrix, &n_loc_rows, &n_loc_cols);
713 "PETSc is requiring non contiguous memory allocation."));
725 PetscInt n_rows, n_cols, n_loc_rows, n_loc_cols,
min,
max;
728 ierr = MatGetSize(
matrix, &n_rows, &n_cols);
731 ierr = MatGetLocalSize(
matrix, &n_loc_rows, &n_loc_cols);
739 "PETSc is requiring non contiguous memory allocation."));
776 #endif // DEAL_II_WITH_PETSC