deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
petsc_matrix_base.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2004 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
16
17#ifdef DEAL_II_WITH_PETSC
18
24
26
27namespace PETScWrappers
28{
29 namespace MatrixIterators
30 {
31# ifndef DOXYGEN
32 void
33 MatrixBase::const_iterator::Accessor::visit_present_row()
34 {
35 // if we are asked to visit the past-the-end line (or a line that is not
36 // stored on the current processor), then simply release all our caches
37 // and go on with life
38 if (matrix->in_local_range(this->a_row) == false)
39 {
40 colnum_cache.reset();
41 value_cache.reset();
42
43 return;
44 }
45
46 // get a representation of the present row
47 PetscInt ncols;
48 const PetscInt *colnums;
49 const PetscScalar *values;
50
51 PetscErrorCode ierr =
52 MatGetRow(*matrix, this->a_row, &ncols, &colnums, &values);
53 AssertThrow(ierr == 0, ExcPETScError(ierr));
54
55 // copy it into our caches if the line
56 // isn't empty. if it is, then we've
57 // done something wrong, since we
58 // shouldn't have initialized an
59 // iterator for an empty line (what
60 // would it point to?)
61 Assert(ncols != 0, ExcInternalError());
62 if constexpr (running_in_debug_mode())
63 {
64 for (PetscInt j = 0; j < ncols; ++j)
65 {
66 const auto column = static_cast<PetscInt>(colnums[j]);
67 AssertIntegerConversion(column, colnums[j]);
68 }
69 }
70 colnum_cache =
71 std::make_shared<std::vector<size_type>>(colnums, colnums + ncols);
72 value_cache =
73 std::make_shared<std::vector<PetscScalar>>(values, values + ncols);
74
75 // and finally restore the matrix
76 ierr = MatRestoreRow(*matrix, this->a_row, &ncols, &colnums, &values);
77 AssertThrow(ierr == 0, ExcPETScError(ierr));
78 }
79# endif
80 } // namespace MatrixIterators
81
82
83
85 : matrix(nullptr)
87 {}
88
89
91 : matrix(A)
93 {
94 const PetscErrorCode ierr =
95 PetscObjectReference(reinterpret_cast<PetscObject>(matrix));
96 AssertThrow(ierr == 0, ExcPETScError(ierr));
97 }
98
99 void
101 {
103 ExcMessage("Cannot assign a new Mat."));
104 PetscErrorCode ierr =
105 PetscObjectReference(reinterpret_cast<PetscObject>(A));
106 AssertThrow(ierr == 0, ExcPETScError(ierr));
107 ierr = MatDestroy(&matrix);
108 AssertThrow(ierr == 0, ExcPETScError(ierr));
109 matrix = A;
110 }
111
113 {
114 PetscErrorCode ierr = MatDestroy(&matrix);
115 AssertNothrow(ierr == 0, ExcPETScError(ierr));
116 }
117
118 void
120 {
121 // destroy the matrix...
122 {
123 const PetscErrorCode ierr = MatDestroy(&matrix);
124 AssertThrow(ierr == 0, ExcPETScError(ierr));
125 }
126
127 // ...and replace it by an empty
128 // sequential matrix
129 const int m = 0, n = 0, n_nonzero_per_row = 0;
130 const PetscErrorCode ierr = MatCreateSeqAIJ(
131 PETSC_COMM_SELF, m, n, n_nonzero_per_row, nullptr, &matrix);
132 AssertThrow(ierr == 0, ExcPETScError(ierr));
133 }
134
135
136
137 MatrixBase &
139 {
141
143
144 const PetscErrorCode ierr = MatZeroEntries(matrix);
145 AssertThrow(ierr == 0, ExcPETScError(ierr));
146
147 return *this;
148 }
149
150
151
152 void
153 MatrixBase::clear_row(const size_type row, const PetscScalar new_diag_value)
154 {
155 clear_rows(ArrayView<const size_type>(row), new_diag_value);
156 }
157
158
159
160 void
162 const PetscScalar new_diag_value)
163 {
165
166 // now set all the entries of these rows to zero
167 if constexpr (running_in_debug_mode())
168 {
169 for (const auto &row : rows)
170 AssertIntegerConversion(static_cast<PetscInt>(row), row);
171 }
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 PetscErrorCode ierr;
180 ierr = ISCreateGeneral(get_mpi_communicator(),
181 rows.size(),
182 petsc_rows.data(),
183 PETSC_COPY_VALUES,
184 &index_set);
185 AssertThrow(ierr == 0, ExcPETScError(ierr));
186
187 ierr = MatZeroRowsIS(matrix, index_set, new_diag_value, nullptr, nullptr);
188 AssertThrow(ierr == 0, ExcPETScError(ierr));
189 ierr = ISDestroy(&index_set);
190 AssertThrow(ierr == 0, ExcPETScError(ierr));
191 }
192
193 void
194 MatrixBase::clear_rows_columns(const std::vector<size_type> &rows,
195 const PetscScalar new_diag_value)
196 {
198
199 // now set all the entries of these rows to zero
200 if constexpr (running_in_debug_mode())
201 {
202 for (const auto &row : rows)
203 AssertIntegerConversion(static_cast<PetscInt>(row), row);
204 }
205 const std::vector<PetscInt> petsc_rows(rows.begin(), rows.end());
206
207 // call the functions. note that we have
208 // to call them even if #rows is empty,
209 // since this is a collective operation
210 IS index_set;
211
212 PetscErrorCode ierr;
213 ierr = ISCreateGeneral(get_mpi_communicator(),
214 rows.size(),
215 petsc_rows.data(),
216 PETSC_COPY_VALUES,
217 &index_set);
218 AssertThrow(ierr == 0, ExcPETScError(ierr));
219
220 ierr =
221 MatZeroRowsColumnsIS(matrix, index_set, new_diag_value, nullptr, nullptr);
222 AssertThrow(ierr == 0, ExcPETScError(ierr));
223 ierr = ISDestroy(&index_set);
224 AssertThrow(ierr == 0, ExcPETScError(ierr));
225 }
226
227
228
229 PetscScalar
230 MatrixBase::el(const size_type i, const size_type j) const
231 {
232 const auto petsc_i = static_cast<PetscInt>(i);
233 AssertIntegerConversion(petsc_i, i);
234 const auto petsc_j = static_cast<PetscInt>(j);
235 AssertIntegerConversion(petsc_j, j);
236
237 PetscScalar value;
238
239 const PetscErrorCode ierr =
240 MatGetValues(matrix, 1, &petsc_i, 1, &petsc_j, &value);
241 AssertThrow(ierr == 0, ExcPETScError(ierr));
242
243 return value;
244 }
245
246
247
248 PetscScalar
250 {
251 Assert(m() == n(), ExcNotQuadratic());
252
253 // this doesn't seem to work any
254 // different than any other element
255 return el(i, i);
256 }
257
258
259
260 void
262 {
263 {
264 if constexpr (running_in_debug_mode())
265 {
266 // Check that all processors agree that last_action is the same (or
267 // none!)
268
269 int my_int_last_action = last_action;
270 int all_int_last_action;
271
272 const int ierr = MPI_Allreduce(&my_int_last_action,
273 &all_int_last_action,
274 1,
275 MPI_INT,
276 MPI_BOR,
278 AssertThrowMPI(ierr);
279
280 AssertThrow(all_int_last_action !=
283 "Error: not all processors agree on the last "
284 "VectorOperation before this compress() call."));
285 }
286 }
287
291 "Missing compress() or calling with wrong VectorOperation argument."));
292
293 // flush buffers
294 PetscErrorCode ierr = MatAssemblyBegin(matrix, MAT_FINAL_ASSEMBLY);
295 AssertThrow(ierr == 0, ExcPETScError(ierr));
296
297 ierr = MatAssemblyEnd(matrix, MAT_FINAL_ASSEMBLY);
298 AssertThrow(ierr == 0, ExcPETScError(ierr));
299
301 }
302
303
304
307 {
308 PetscInt n_rows, n_cols;
309
310 const PetscErrorCode ierr = MatGetSize(matrix, &n_rows, &n_cols);
311 AssertThrow(ierr == 0, ExcPETScError(ierr));
312
313 return n_rows;
314 }
315
316
317
320 {
321 PetscInt n_rows, n_cols;
322
323 const PetscErrorCode ierr = MatGetSize(matrix, &n_rows, &n_cols);
324 AssertThrow(ierr == 0, ExcPETScError(ierr));
325
326 return n_cols;
327 }
328
329
330
333 {
334 PetscInt n_rows;
335
336 const PetscErrorCode ierr = MatGetLocalSize(matrix, &n_rows, nullptr);
337 AssertThrow(ierr == 0, ExcPETScError(ierr));
338
339 return n_rows;
340 }
341
342
343
344 std::pair<MatrixBase::size_type, MatrixBase::size_type>
346 {
347 PetscInt begin, end;
348
349 const PetscErrorCode ierr =
350 MatGetOwnershipRange(static_cast<const Mat &>(matrix), &begin, &end);
351 AssertThrow(ierr == 0, ExcPETScError(ierr));
352
353 return {begin, end};
354 }
355
356
357
360 {
361 PetscInt n_cols;
362
363 const PetscErrorCode ierr = MatGetLocalSize(matrix, nullptr, &n_cols);
364 AssertThrow(ierr == 0, ExcPETScError(ierr));
365
366 return n_cols;
367 }
368
369
370
371 std::pair<MatrixBase::size_type, MatrixBase::size_type>
373 {
374 PetscInt begin, end;
375
376 const PetscErrorCode ierr =
377 MatGetOwnershipRangeColumn(static_cast<const Mat &>(matrix),
378 &begin,
379 &end);
380 AssertThrow(ierr == 0, ExcPETScError(ierr));
381
382 return {begin, end};
383 }
384
385
386
387 std::uint64_t
389 {
390 MatInfo mat_info;
391 const PetscErrorCode ierr = MatGetInfo(matrix, MAT_GLOBAL_SUM, &mat_info);
392 AssertThrow(ierr == 0, ExcPETScError(ierr));
393
394 // MatInfo logs quantities as PetscLogDouble. So we need to cast it to match
395 // our interface.
396 return static_cast<std::uint64_t>(mat_info.nz_used);
397 }
398
399
400
403 {
404 // TODO: this function will probably only work if compress() was called on
405 // the matrix previously. however, we can't do this here, since it would
406 // impose global communication and one would have to make sure that this
407 // function is called the same number of times from all processors,
408 // something that is unreasonable. there should simply be a way in PETSc to
409 // query the number of entries in a row bypassing the call to compress(),
410 // but I can't find one
411 AssertIndexRange(row, m());
412
413 // get a representation of the present
414 // row
415 PetscInt ncols;
416 const PetscInt *colnums;
417 const PetscScalar *values;
418
419 // TODO: this is probably horribly inefficient; we should lobby for a way to
420 // query this information from PETSc
421 PetscErrorCode ierr = MatGetRow(*this, row, &ncols, &colnums, &values);
422 AssertThrow(ierr == 0, ExcPETScError(ierr));
423
424 // then restore the matrix and return the number of columns in this row as
425 // queried previously. Starting with PETSc 3.4, MatRestoreRow actually
426 // resets the last three arguments to nullptr, to avoid abuse of pointers
427 // now dangling. as a consequence, we need to save the size of the array
428 // and return the saved value.
429 const PetscInt ncols_saved = ncols;
430 ierr = MatRestoreRow(*this, row, &ncols, &colnums, &values);
431 AssertThrow(ierr == 0, ExcPETScError(ierr));
432
433 return ncols_saved;
434 }
435
436
437 PetscReal
439 {
440 PetscReal result;
441
442 const PetscErrorCode ierr = MatNorm(matrix, NORM_1, &result);
443 AssertThrow(ierr == 0, ExcPETScError(ierr));
444
445 return result;
446 }
447
448
449
450 PetscReal
452 {
453 PetscReal result;
454
455 const PetscErrorCode ierr = MatNorm(matrix, NORM_INFINITY, &result);
456 AssertThrow(ierr == 0, ExcPETScError(ierr));
457
458 return result;
459 }
460
461
462
463 PetscReal
465 {
466 PetscReal result;
467
468 const PetscErrorCode ierr = MatNorm(matrix, NORM_FROBENIUS, &result);
469 AssertThrow(ierr == 0, ExcPETScError(ierr));
470
471 return result;
472 }
473
474
475 PetscScalar
477 {
478 AssertDimension(m(), v.size());
479
480 VectorBase tmp(v);
481 vmult(tmp, v);
482 return tmp * v;
483 }
484
485
486 PetscScalar
488 const VectorBase &v) const
489 {
490 AssertDimension(m(), u.size());
491 AssertDimension(m(), v.size());
492
493 VectorBase tmp(u);
494 vmult(tmp, v);
495 return u * tmp;
496 }
497
498
499 PetscScalar
501 {
502 PetscScalar result;
503
504 const PetscErrorCode ierr = MatGetTrace(matrix, &result);
505 AssertThrow(ierr == 0, ExcPETScError(ierr));
506
507 return result;
508 }
509
510
511
512 MatrixBase &
513 MatrixBase::operator*=(const PetscScalar a)
514 {
515 const PetscErrorCode ierr = MatScale(matrix, a);
516 AssertThrow(ierr == 0, ExcPETScError(ierr));
517
518 return *this;
519 }
520
521
522
523 MatrixBase &
524 MatrixBase::operator/=(const PetscScalar a)
525 {
526 const PetscScalar factor = 1. / a;
527 const PetscErrorCode ierr = MatScale(matrix, factor);
528 AssertThrow(ierr == 0, ExcPETScError(ierr));
529
530 return *this;
531 }
532
533
534
535 MatrixBase &
536 MatrixBase::add(const PetscScalar factor, const MatrixBase &other)
537 {
538 const PetscErrorCode ierr =
539 MatAXPY(matrix, factor, other, DIFFERENT_NONZERO_PATTERN);
540 AssertThrow(ierr == 0, ExcPETScError(ierr));
541
542 return *this;
543 }
544
545
546 void
548 {
549 Assert(&src != &dst, ExcSourceEqualsDestination());
550
551 const PetscErrorCode ierr = MatMult(matrix, src, dst);
552 AssertThrow(ierr == 0, ExcPETScError(ierr));
553 }
554
555
556
557 void
559 {
560 Assert(&src != &dst, ExcSourceEqualsDestination());
561
562 const PetscErrorCode ierr = MatMultTranspose(matrix, src, dst);
563 AssertThrow(ierr == 0, ExcPETScError(ierr));
564 }
565
566
567
568 void
570 {
571 Assert(&src != &dst, ExcSourceEqualsDestination());
572
573 const PetscErrorCode ierr = MatMultAdd(matrix, src, dst, dst);
574 AssertThrow(ierr == 0, ExcPETScError(ierr));
575 }
576
577
578
579 void
581 {
582 Assert(&src != &dst, ExcSourceEqualsDestination());
583
584 const PetscErrorCode ierr = MatMultTransposeAdd(matrix, src, dst, dst);
585 AssertThrow(ierr == 0, ExcPETScError(ierr));
586 }
587
588
589 namespace internals
590 {
591 void
592 perform_mmult(const MatrixBase &inputleft,
593 const MatrixBase &inputright,
594 MatrixBase &result,
595 const VectorBase &V,
596 const bool transpose_left)
597 {
598 const bool use_vector = (V.size() == inputright.m() ? true : false);
599 if (transpose_left == false)
600 {
601 Assert(inputleft.n() == inputright.m(),
602 ExcDimensionMismatch(inputleft.n(), inputright.m()));
603 }
604 else
605 {
606 Assert(inputleft.m() == inputright.m(),
607 ExcDimensionMismatch(inputleft.m(), inputright.m()));
608 }
609
610 result.clear();
611
612 PetscErrorCode ierr;
613
614 if (use_vector == false)
615 {
616 if (transpose_left)
617 {
618 ierr = MatTransposeMatMult(inputleft,
619 inputright,
620 MAT_INITIAL_MATRIX,
621 PETSC_DEFAULT,
622 &result.petsc_matrix());
623 AssertThrow(ierr == 0, ExcPETScError(ierr));
624 }
625 else
626 {
627 ierr = MatMatMult(inputleft,
628 inputright,
629 MAT_INITIAL_MATRIX,
630 PETSC_DEFAULT,
631 &result.petsc_matrix());
632 AssertThrow(ierr == 0, ExcPETScError(ierr));
633 }
634 }
635 else
636 {
637 Mat tmp;
638 ierr = MatDuplicate(inputleft, MAT_COPY_VALUES, &tmp);
639 AssertThrow(ierr == 0, ExcPETScError(ierr));
640 if (transpose_left)
641 {
642# if DEAL_II_PETSC_VERSION_LT(3, 8, 0)
643 ierr = MatTranspose(tmp, MAT_REUSE_MATRIX, &tmp);
644# else
645 ierr = MatTranspose(tmp, MAT_INPLACE_MATRIX, &tmp);
646# endif
647 AssertThrow(ierr == 0, ExcPETScError(ierr));
648 }
649 ierr = MatDiagonalScale(tmp, nullptr, V);
650 AssertThrow(ierr == 0, ExcPETScError(ierr));
651 ierr = MatMatMult(tmp,
652 inputright,
653 MAT_INITIAL_MATRIX,
654 PETSC_DEFAULT,
655 &result.petsc_matrix());
656 AssertThrow(ierr == 0, ExcPETScError(ierr));
657 ierr = MatDestroy(&tmp);
658 AssertThrow(ierr == 0, ExcPETScError(ierr));
659 }
660 }
661 } // namespace internals
662
663 void
665 const MatrixBase &B,
666 const VectorBase &V) const
667 {
668 internals::perform_mmult(*this, B, C, V, false);
669 }
670
671 void
673 const MatrixBase &B,
674 const VectorBase &V) const
675 {
676 internals::perform_mmult(*this, B, C, V, true);
677 }
678
679 PetscScalar
681 const VectorBase &x,
682 const VectorBase &b) const
683 {
684 // avoid the use of a temporary, and
685 // rather do one negation pass more than
686 // necessary
687 vmult(dst, x);
688 dst -= b;
689 dst *= -1;
690
691 return dst.l2_norm();
692 }
693
694
695
696 MatrixBase::operator Mat() const
697 {
698 return matrix;
699 }
700
701 Mat &
703 {
704 return matrix;
705 }
706
707 void
709 {
710# if DEAL_II_PETSC_VERSION_LT(3, 8, 0)
711 const PetscErrorCode ierr = MatTranspose(matrix, MAT_REUSE_MATRIX, &matrix);
712# else
713 const PetscErrorCode ierr =
714 MatTranspose(matrix, MAT_INPLACE_MATRIX, &matrix);
715# endif
716 AssertThrow(ierr == 0, ExcPETScError(ierr));
717 }
718
719 PetscBool
720 MatrixBase::is_symmetric(const double tolerance)
721 {
722 PetscBool truth;
724 const PetscErrorCode ierr = MatIsSymmetric(matrix, tolerance, &truth);
725 AssertThrow(ierr == 0, ExcPETScError(ierr));
726 return truth;
727 }
728
729 PetscBool
730 MatrixBase::is_hermitian(const double tolerance)
731 {
732 PetscBool truth;
733
735 const PetscErrorCode ierr = MatIsHermitian(matrix, tolerance, &truth);
736 AssertThrow(ierr == 0, ExcPETScError(ierr));
737
738 return truth;
739 }
740
741 void
742 MatrixBase::write_ascii(const PetscViewerFormat format)
743 {
745 MPI_Comm comm = PetscObjectComm(reinterpret_cast<PetscObject>(matrix));
746
747 // Set options
748 PetscErrorCode ierr =
749 PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(comm), format);
750 AssertThrow(ierr == 0, ExcPETScError(ierr));
751
752 // Write to screen
753 ierr = MatView(matrix, PETSC_VIEWER_STDOUT_(comm));
754 AssertThrow(ierr == 0, ExcPETScError(ierr));
755 ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(comm));
756 AssertThrow(ierr == 0, ExcPETScError(ierr));
757 }
758
759 void
760 MatrixBase::print(std::ostream &out, const bool /*alternative_output*/) const
761 {
762 PetscBool has;
763
764 PetscErrorCode ierr = MatHasOperation(matrix, MATOP_GET_ROW, &has);
765 AssertThrow(ierr == 0, ExcPETScError(ierr));
766
767 Mat vmatrix = matrix;
768 if (!has)
769 {
770 ierr = MatConvert(matrix, MATAIJ, MAT_INITIAL_MATRIX, &vmatrix);
771 AssertThrow(ierr == 0, ExcPETScError(ierr));
772 }
773
774 std::pair<MatrixBase::size_type, MatrixBase::size_type> loc_range =
775 local_range();
776
777 PetscInt ncols;
778 const PetscInt *colnums;
779 const PetscScalar *values;
780
782 for (row = loc_range.first; row < loc_range.second; ++row)
783 {
784 ierr = MatGetRow(vmatrix, row, &ncols, &colnums, &values);
785 AssertThrow(ierr == 0, ExcPETScError(ierr));
786
787 for (PetscInt col = 0; col < ncols; ++col)
788 {
789 out << "(" << row << "," << colnums[col] << ") " << values[col]
790 << std::endl;
791 }
792
793 ierr = MatRestoreRow(vmatrix, row, &ncols, &colnums, &values);
794 AssertThrow(ierr == 0, ExcPETScError(ierr));
795 }
796 if (vmatrix != matrix)
797 {
798 ierr = MatDestroy(&vmatrix);
799 AssertThrow(ierr == 0, ExcPETScError(ierr));
800 }
801 AssertThrow(out.fail() == false, ExcIO());
802 }
803
804
805
806 std::size_t
808 {
809 MatInfo info;
810 const PetscErrorCode ierr = MatGetInfo(matrix, MAT_LOCAL, &info);
811 AssertThrow(ierr == 0, ExcPETScError(ierr));
812
813 return sizeof(*this) + static_cast<size_type>(info.memory);
814 }
815
816} // namespace PETScWrappers
817
819
820#endif // DEAL_II_WITH_PETSC
iterator begin() const
Definition array_view.h:707
iterator end() const
Definition array_view.h:716
std::size_t size() const
Definition array_view.h:689
void add(const size_type i, const size_type j, const PetscScalar value)
size_type row_length(const size_type row) const
std::size_t memory_consumption() const
VectorOperation::values last_action
void vmult(VectorBase &dst, const VectorBase &src) const
MPI_Comm get_mpi_communicator() const
PetscScalar diag_element(const size_type i) const
size_type local_domain_size() const
const_iterator begin() const
MatrixBase & operator/=(const PetscScalar factor)
void mmult(MatrixBase &C, const MatrixBase &B, const VectorBase &V) const
PetscBool is_symmetric(const double tolerance=1.e-12)
types::global_dof_index size_type
std::pair< size_type, size_type > local_domain() const
const_iterator end() const
void Tvmult_add(VectorBase &dst, const VectorBase &src) const
void print(std::ostream &out, const bool alternative_output=false) const
PetscScalar el(const size_type i, const size_type j) const
MatrixBase & operator=(const MatrixBase &)=delete
PetscBool is_hermitian(const double tolerance=1.e-12)
PetscScalar matrix_scalar_product(const VectorBase &u, const VectorBase &v) const
void Tvmult(VectorBase &dst, const VectorBase &src) const
void clear_rows_columns(const std::vector< size_type > &row_and_column_indices, const PetscScalar new_diag_value=0)
MatrixBase & operator*=(const PetscScalar factor)
PetscScalar residual(VectorBase &dst, const VectorBase &x, const VectorBase &b) const
void Tmmult(MatrixBase &C, const MatrixBase &B, const VectorBase &V) const
void write_ascii(const PetscViewerFormat format=PETSC_VIEWER_DEFAULT)
void vmult_add(VectorBase &dst, const VectorBase &src) const
std::pair< size_type, size_type > local_range() const
void compress(const VectorOperation::values operation)
PetscScalar matrix_norm_square(const VectorBase &v) const
void clear_rows(const ArrayView< const size_type > &rows, const PetscScalar new_diag_value=0)
std::uint64_t n_nonzero_elements() const
void clear_row(const size_type row, const PetscScalar new_diag_value=0)
size_type size() const override
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
constexpr bool running_in_debug_mode()
Definition config.h:78
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define AssertIntegerConversion(index1, index2)
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
#define AssertNothrow(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcSourceEqualsDestination()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcNotQuadratic()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
const bool IsBlockVector< VectorType >::value
@ matrix
Contents is actually a matrix.
void perform_mmult(const MatrixBase &inputleft, const MatrixBase &inputright, MatrixBase &result, const VectorBase &V, const bool transpose_left)
*braid_SplitCommworld & comm