Reference documentation for deal.II version 9.6.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 - 2024 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 colnum_cache =
63 std::make_shared<std::vector<size_type>>(colnums, colnums + ncols);
64 value_cache =
65 std::make_shared<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, ExcPETScError(ierr));
70 }
71# endif
72 } // namespace MatrixIterators
73
74
75
77 : matrix(nullptr)
78 , last_action(VectorOperation::unknown)
79 {}
80
81
83 : matrix(A)
84 , last_action(VectorOperation::unknown)
85 {
86 const PetscErrorCode ierr =
87 PetscObjectReference(reinterpret_cast<PetscObject>(matrix));
88 AssertThrow(ierr == 0, ExcPETScError(ierr));
89 }
90
91 void
93 {
95 ExcMessage("Cannot assign a new Mat."));
96 PetscErrorCode ierr =
97 PetscObjectReference(reinterpret_cast<PetscObject>(A));
98 AssertThrow(ierr == 0, ExcPETScError(ierr));
99 ierr = MatDestroy(&matrix);
100 AssertThrow(ierr == 0, ExcPETScError(ierr));
101 matrix = A;
102 }
103
105 {
106 PetscErrorCode ierr = MatDestroy(&matrix);
107 (void)ierr;
108 AssertNothrow(ierr == 0, ExcPETScError(ierr));
109 }
110
111 void
113 {
114 // destroy the matrix...
115 {
116 const PetscErrorCode ierr = MatDestroy(&matrix);
117 AssertThrow(ierr == 0, ExcPETScError(ierr));
118 }
119
120 // ...and replace it by an empty
121 // sequential matrix
122 const int m = 0, n = 0, n_nonzero_per_row = 0;
123 const PetscErrorCode ierr = MatCreateSeqAIJ(
124 PETSC_COMM_SELF, m, n, n_nonzero_per_row, nullptr, &matrix);
125 AssertThrow(ierr == 0, ExcPETScError(ierr));
126 }
127
128
129
130 MatrixBase &
132 {
133 (void)d;
135
137
138 const PetscErrorCode ierr = MatZeroEntries(matrix);
139 AssertThrow(ierr == 0, ExcPETScError(ierr));
140
141 return *this;
142 }
143
144
145
146 void
147 MatrixBase::clear_row(const size_type row, const PetscScalar new_diag_value)
148 {
149 clear_rows(ArrayView<const size_type>(row), new_diag_value);
150 }
151
152
153
154 void
156 const PetscScalar new_diag_value)
157 {
159
160 // now set all the entries of these rows
161 // to zero
162 const std::vector<PetscInt> petsc_rows(rows.begin(), rows.end());
163
164 // call the functions. note that we have
165 // to call them even if #rows is empty,
166 // since this is a collective operation
167 IS index_set;
168
169 PetscErrorCode ierr;
170 ierr = ISCreateGeneral(get_mpi_communicator(),
171 rows.size(),
172 petsc_rows.data(),
173 PETSC_COPY_VALUES,
174 &index_set);
175 AssertThrow(ierr == 0, ExcPETScError(ierr));
176
177 ierr = MatZeroRowsIS(matrix, index_set, new_diag_value, nullptr, nullptr);
178 AssertThrow(ierr == 0, ExcPETScError(ierr));
179 ierr = ISDestroy(&index_set);
180 AssertThrow(ierr == 0, ExcPETScError(ierr));
181 }
182
183 void
184 MatrixBase::clear_rows_columns(const std::vector<size_type> &rows,
185 const PetscScalar new_diag_value)
186 {
188
189 // now set all the entries of these rows
190 // to zero
191 const std::vector<PetscInt> petsc_rows(rows.begin(), rows.end());
192
193 // call the functions. note that we have
194 // to call them even if #rows is empty,
195 // since this is a collective operation
196 IS index_set;
197
198 PetscErrorCode ierr;
199 ierr = ISCreateGeneral(get_mpi_communicator(),
200 rows.size(),
201 petsc_rows.data(),
202 PETSC_COPY_VALUES,
203 &index_set);
204 AssertThrow(ierr == 0, ExcPETScError(ierr));
205
206 ierr =
207 MatZeroRowsColumnsIS(matrix, index_set, new_diag_value, nullptr, nullptr);
208 AssertThrow(ierr == 0, ExcPETScError(ierr));
209 ierr = ISDestroy(&index_set);
210 AssertThrow(ierr == 0, ExcPETScError(ierr));
211 }
212
213
214
215 PetscScalar
216 MatrixBase::el(const size_type i, const size_type j) const
217 {
218 PetscInt petsc_i = i, petsc_j = j;
219
220 PetscScalar value;
221
222 const PetscErrorCode ierr =
223 MatGetValues(matrix, 1, &petsc_i, 1, &petsc_j, &value);
224 AssertThrow(ierr == 0, ExcPETScError(ierr));
225
226 return value;
227 }
228
229
230
231 PetscScalar
233 {
234 Assert(m() == n(), ExcNotQuadratic());
235
236 // this doesn't seem to work any
237 // different than any other element
238 return el(i, i);
239 }
240
241
242
243 void
245 {
246 {
247# ifdef DEBUG
248 // Check that all processors agree that last_action is the same (or none!)
249
250 int my_int_last_action = last_action;
251 int all_int_last_action;
252
253 const int ierr = MPI_Allreduce(&my_int_last_action,
254 &all_int_last_action,
255 1,
256 MPI_INT,
257 MPI_BOR,
259 AssertThrowMPI(ierr);
260
261 AssertThrow(all_int_last_action !=
263 ExcMessage("Error: not all processors agree on the last "
264 "VectorOperation before this compress() call."));
265# endif
266 }
267
271 "Missing compress() or calling with wrong VectorOperation argument."));
272
273 // flush buffers
274 PetscErrorCode ierr = MatAssemblyBegin(matrix, MAT_FINAL_ASSEMBLY);
275 AssertThrow(ierr == 0, ExcPETScError(ierr));
276
277 ierr = MatAssemblyEnd(matrix, MAT_FINAL_ASSEMBLY);
278 AssertThrow(ierr == 0, ExcPETScError(ierr));
279
281 }
282
283
284
287 {
288 PetscInt n_rows, n_cols;
289
290 const PetscErrorCode ierr = MatGetSize(matrix, &n_rows, &n_cols);
291 AssertThrow(ierr == 0, ExcPETScError(ierr));
292
293 return n_rows;
294 }
295
296
297
300 {
301 PetscInt n_rows, n_cols;
302
303 const PetscErrorCode ierr = MatGetSize(matrix, &n_rows, &n_cols);
304 AssertThrow(ierr == 0, ExcPETScError(ierr));
305
306 return n_cols;
307 }
308
309
310
313 {
314 PetscInt n_rows;
315
316 const PetscErrorCode ierr = MatGetLocalSize(matrix, &n_rows, nullptr);
317 AssertThrow(ierr == 0, ExcPETScError(ierr));
318
319 return n_rows;
320 }
321
322
323
324 std::pair<MatrixBase::size_type, MatrixBase::size_type>
326 {
327 PetscInt begin, end;
328
329 const PetscErrorCode ierr =
330 MatGetOwnershipRange(static_cast<const Mat &>(matrix), &begin, &end);
331 AssertThrow(ierr == 0, ExcPETScError(ierr));
332
333 return {begin, end};
334 }
335
336
337
340 {
341 PetscInt n_cols;
342
343 const PetscErrorCode ierr = MatGetLocalSize(matrix, nullptr, &n_cols);
344 AssertThrow(ierr == 0, ExcPETScError(ierr));
345
346 return n_cols;
347 }
348
349
350
351 std::pair<MatrixBase::size_type, MatrixBase::size_type>
353 {
354 PetscInt begin, end;
355
356 const PetscErrorCode ierr =
357 MatGetOwnershipRangeColumn(static_cast<const Mat &>(matrix),
358 &begin,
359 &end);
360 AssertThrow(ierr == 0, ExcPETScError(ierr));
361
362 return {begin, end};
363 }
364
365
366
367 std::uint64_t
369 {
370 MatInfo mat_info;
371 const PetscErrorCode ierr = MatGetInfo(matrix, MAT_GLOBAL_SUM, &mat_info);
372 AssertThrow(ierr == 0, ExcPETScError(ierr));
373
374 // MatInfo logs quantities as PetscLogDouble. So we need to cast it to match
375 // our interface.
376 return static_cast<std::uint64_t>(mat_info.nz_used);
377 }
378
379
380
383 {
384 // TODO: this function will probably only work if compress() was called on
385 // the matrix previously. however, we can't do this here, since it would
386 // impose global communication and one would have to make sure that this
387 // function is called the same number of times from all processors,
388 // something that is unreasonable. there should simply be a way in PETSc to
389 // query the number of entries in a row bypassing the call to compress(),
390 // but I can't find one
391 Assert(row < m(), ExcInternalError());
392
393 // get a representation of the present
394 // row
395 PetscInt ncols;
396 const PetscInt *colnums;
397 const PetscScalar *values;
398
399 // TODO: this is probably horribly inefficient; we should lobby for a way to
400 // query this information from PETSc
401 PetscErrorCode ierr = MatGetRow(*this, row, &ncols, &colnums, &values);
402 AssertThrow(ierr == 0, ExcPETScError(ierr));
403
404 // then restore the matrix and return the number of columns in this row as
405 // queried previously. Starting with PETSc 3.4, MatRestoreRow actually
406 // resets the last three arguments to nullptr, to avoid abuse of pointers
407 // now dangling. as a consequence, we need to save the size of the array
408 // and return the saved value.
409 const PetscInt ncols_saved = ncols;
410 ierr = MatRestoreRow(*this, row, &ncols, &colnums, &values);
411 AssertThrow(ierr == 0, ExcPETScError(ierr));
412
413 return ncols_saved;
414 }
415
416
417 PetscReal
419 {
420 PetscReal result;
421
422 const PetscErrorCode ierr = MatNorm(matrix, NORM_1, &result);
423 AssertThrow(ierr == 0, ExcPETScError(ierr));
424
425 return result;
426 }
427
428
429
430 PetscReal
432 {
433 PetscReal result;
434
435 const PetscErrorCode ierr = MatNorm(matrix, NORM_INFINITY, &result);
436 AssertThrow(ierr == 0, ExcPETScError(ierr));
437
438 return result;
439 }
440
441
442
443 PetscReal
445 {
446 PetscReal result;
447
448 const PetscErrorCode ierr = MatNorm(matrix, NORM_FROBENIUS, &result);
449 AssertThrow(ierr == 0, ExcPETScError(ierr));
450
451 return result;
452 }
453
454
455 PetscScalar
457 {
458 AssertDimension(m(), v.size());
459
460 VectorBase tmp(v);
461 vmult(tmp, v);
462 return tmp * v;
463 }
464
465
466 PetscScalar
468 const VectorBase &v) const
469 {
470 AssertDimension(m(), u.size());
471 AssertDimension(m(), v.size());
472
473 VectorBase tmp(u);
474 vmult(tmp, v);
475 return u * tmp;
476 }
477
478
479 PetscScalar
481 {
482 PetscScalar result;
483
484 const PetscErrorCode ierr = MatGetTrace(matrix, &result);
485 AssertThrow(ierr == 0, ExcPETScError(ierr));
486
487 return result;
488 }
489
490
491
492 MatrixBase &
493 MatrixBase::operator*=(const PetscScalar a)
494 {
495 const PetscErrorCode ierr = MatScale(matrix, a);
496 AssertThrow(ierr == 0, ExcPETScError(ierr));
497
498 return *this;
499 }
500
501
502
503 MatrixBase &
504 MatrixBase::operator/=(const PetscScalar a)
505 {
506 const PetscScalar factor = 1. / a;
507 const PetscErrorCode ierr = MatScale(matrix, factor);
508 AssertThrow(ierr == 0, ExcPETScError(ierr));
509
510 return *this;
511 }
512
513
514
515 MatrixBase &
516 MatrixBase::add(const PetscScalar factor, const MatrixBase &other)
517 {
518 const PetscErrorCode ierr =
519 MatAXPY(matrix, factor, other, DIFFERENT_NONZERO_PATTERN);
520 AssertThrow(ierr == 0, ExcPETScError(ierr));
521
522 return *this;
523 }
524
525
526 void
528 {
529 Assert(&src != &dst, ExcSourceEqualsDestination());
530
531 const PetscErrorCode ierr = MatMult(matrix, src, dst);
532 AssertThrow(ierr == 0, ExcPETScError(ierr));
533 }
534
535
536
537 void
539 {
540 Assert(&src != &dst, ExcSourceEqualsDestination());
541
542 const PetscErrorCode ierr = MatMultTranspose(matrix, src, dst);
543 AssertThrow(ierr == 0, ExcPETScError(ierr));
544 }
545
546
547
548 void
550 {
551 Assert(&src != &dst, ExcSourceEqualsDestination());
552
553 const PetscErrorCode ierr = MatMultAdd(matrix, src, dst, dst);
554 AssertThrow(ierr == 0, ExcPETScError(ierr));
555 }
556
557
558
559 void
561 {
562 Assert(&src != &dst, ExcSourceEqualsDestination());
563
564 const PetscErrorCode ierr = MatMultTransposeAdd(matrix, src, dst, dst);
565 AssertThrow(ierr == 0, ExcPETScError(ierr));
566 }
567
568
569 namespace internals
570 {
571 void
572 perform_mmult(const MatrixBase &inputleft,
573 const MatrixBase &inputright,
574 MatrixBase &result,
575 const VectorBase &V,
576 const bool transpose_left)
577 {
578 const bool use_vector = (V.size() == inputright.m() ? true : false);
579 if (transpose_left == false)
580 {
581 Assert(inputleft.n() == inputright.m(),
582 ExcDimensionMismatch(inputleft.n(), inputright.m()));
583 }
584 else
585 {
586 Assert(inputleft.m() == inputright.m(),
587 ExcDimensionMismatch(inputleft.m(), inputright.m()));
588 }
589
590 result.clear();
591
592 PetscErrorCode ierr;
593
594 if (use_vector == false)
595 {
596 if (transpose_left)
597 {
598 ierr = MatTransposeMatMult(inputleft,
599 inputright,
600 MAT_INITIAL_MATRIX,
601 PETSC_DEFAULT,
602 &result.petsc_matrix());
603 AssertThrow(ierr == 0, ExcPETScError(ierr));
604 }
605 else
606 {
607 ierr = MatMatMult(inputleft,
608 inputright,
609 MAT_INITIAL_MATRIX,
610 PETSC_DEFAULT,
611 &result.petsc_matrix());
612 AssertThrow(ierr == 0, ExcPETScError(ierr));
613 }
614 }
615 else
616 {
617 Mat tmp;
618 ierr = MatDuplicate(inputleft, MAT_COPY_VALUES, &tmp);
619 AssertThrow(ierr == 0, ExcPETScError(ierr));
620 if (transpose_left)
621 {
622# if DEAL_II_PETSC_VERSION_LT(3, 8, 0)
623 ierr = MatTranspose(tmp, MAT_REUSE_MATRIX, &tmp);
624# else
625 ierr = MatTranspose(tmp, MAT_INPLACE_MATRIX, &tmp);
626# endif
627 AssertThrow(ierr == 0, ExcPETScError(ierr));
628 }
629 ierr = MatDiagonalScale(tmp, nullptr, V);
630 AssertThrow(ierr == 0, ExcPETScError(ierr));
631 ierr = MatMatMult(tmp,
632 inputright,
633 MAT_INITIAL_MATRIX,
634 PETSC_DEFAULT,
635 &result.petsc_matrix());
636 AssertThrow(ierr == 0, ExcPETScError(ierr));
637 ierr = MatDestroy(&tmp);
638 AssertThrow(ierr == 0, ExcPETScError(ierr));
639 }
640 }
641 } // namespace internals
642
643 void
645 const MatrixBase &B,
646 const VectorBase &V) const
647 {
648 internals::perform_mmult(*this, B, C, V, false);
649 }
650
651 void
653 const MatrixBase &B,
654 const VectorBase &V) const
655 {
656 internals::perform_mmult(*this, B, C, V, true);
657 }
658
659 PetscScalar
661 const VectorBase &x,
662 const VectorBase &b) const
663 {
664 // avoid the use of a temporary, and
665 // rather do one negation pass more than
666 // necessary
667 vmult(dst, x);
668 dst -= b;
669 dst *= -1;
670
671 return dst.l2_norm();
672 }
673
674
675
676 MatrixBase::operator Mat() const
677 {
678 return matrix;
679 }
680
681 Mat &
683 {
684 return matrix;
685 }
686
687 void
689 {
690# if DEAL_II_PETSC_VERSION_LT(3, 8, 0)
691 const PetscErrorCode ierr = MatTranspose(matrix, MAT_REUSE_MATRIX, &matrix);
692# else
693 const PetscErrorCode ierr =
694 MatTranspose(matrix, MAT_INPLACE_MATRIX, &matrix);
695# endif
696 AssertThrow(ierr == 0, ExcPETScError(ierr));
697 }
698
699 PetscBool
700 MatrixBase::is_symmetric(const double tolerance)
701 {
702 PetscBool truth;
704 const PetscErrorCode ierr = MatIsSymmetric(matrix, tolerance, &truth);
705 AssertThrow(ierr == 0, ExcPETScError(ierr));
706 return truth;
707 }
708
709 PetscBool
710 MatrixBase::is_hermitian(const double tolerance)
711 {
712 PetscBool truth;
713
715 const PetscErrorCode ierr = MatIsHermitian(matrix, tolerance, &truth);
716 AssertThrow(ierr == 0, ExcPETScError(ierr));
717
718 return truth;
719 }
720
721 void
722 MatrixBase::write_ascii(const PetscViewerFormat format)
723 {
725 MPI_Comm comm = PetscObjectComm(reinterpret_cast<PetscObject>(matrix));
726
727 // Set options
728 PetscErrorCode ierr =
729 PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(comm), format);
730 AssertThrow(ierr == 0, ExcPETScError(ierr));
731
732 // Write to screen
733 ierr = MatView(matrix, PETSC_VIEWER_STDOUT_(comm));
734 AssertThrow(ierr == 0, ExcPETScError(ierr));
735 ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(comm));
736 AssertThrow(ierr == 0, ExcPETScError(ierr));
737 }
738
739 void
740 MatrixBase::print(std::ostream &out, const bool /*alternative_output*/) const
741 {
742 PetscBool has;
743
744 PetscErrorCode ierr = MatHasOperation(matrix, MATOP_GET_ROW, &has);
745 AssertThrow(ierr == 0, ExcPETScError(ierr));
746
747 Mat vmatrix = matrix;
748 if (!has)
749 {
750 ierr = MatConvert(matrix, MATAIJ, MAT_INITIAL_MATRIX, &vmatrix);
751 AssertThrow(ierr == 0, ExcPETScError(ierr));
752 }
753
754 std::pair<MatrixBase::size_type, MatrixBase::size_type> loc_range =
755 local_range();
756
757 PetscInt ncols;
758 const PetscInt *colnums;
759 const PetscScalar *values;
760
762 for (row = loc_range.first; row < loc_range.second; ++row)
763 {
764 ierr = MatGetRow(vmatrix, row, &ncols, &colnums, &values);
765 AssertThrow(ierr == 0, ExcPETScError(ierr));
766
767 for (PetscInt col = 0; col < ncols; ++col)
768 {
769 out << "(" << row << "," << colnums[col] << ") " << values[col]
770 << std::endl;
771 }
772
773 ierr = MatRestoreRow(vmatrix, row, &ncols, &colnums, &values);
774 AssertThrow(ierr == 0, ExcPETScError(ierr));
775 }
776 if (vmatrix != matrix)
777 {
778 ierr = MatDestroy(&vmatrix);
779 AssertThrow(ierr == 0, ExcPETScError(ierr));
780 }
781 AssertThrow(out.fail() == false, ExcIO());
782 }
783
784
785
786 std::size_t
788 {
789 MatInfo info;
790 const PetscErrorCode ierr = MatGetInfo(matrix, MAT_LOCAL, &info);
791 AssertThrow(ierr == 0, ExcPETScError(ierr));
792
793 return sizeof(*this) + static_cast<size_type>(info.memory);
794 }
795
796} // namespace PETScWrappers
797
799
800#endif // DEAL_II_WITH_PETSC
iterator begin() const
Definition array_view.h:702
iterator end() const
Definition array_view.h:711
std::size_t size() const
Definition array_view.h:684
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)
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:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
#define AssertNothrow(cond, exc)
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)
@ 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