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
sparse_direct.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2001 - 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
17
21#include <deal.II/lac/vector.h>
22
23#include <complex>
24#include <vector>
25
26#ifdef DEAL_II_WITH_UMFPACK
27# include <umfpack.h>
28#endif
29
30#ifdef DEAL_II_WITH_MUMPS
31# include <dmumps_c.h>
32#endif
33
35
36namespace TrilinosWrappers
37{
38 class SparseMatrix;
39 namespace MPI
40 {
41 class SparseMatrix;
42 class Vector;
43 } // namespace MPI
44} // namespace TrilinosWrappers
45namespace PETScWrappers
46{
47 namespace MPI
48 {
49 class SparseMatrix;
50 class Vector;
51 } // namespace MPI
52} // namespace PETScWrappers
53
54namespace
55{
63 template <typename SparseMatrixType>
64 unsigned int
65 parallel_grainsize(const SparseMatrixType &matrix)
66 {
67 const unsigned int avg_entries_per_row =
68 matrix.n_nonzero_elements() / matrix.m();
69 return std::max(1000 / avg_entries_per_row, 1u);
70 }
71} // namespace
72
73
78
79
80void
83
84
85#ifdef DEAL_II_WITH_UMFPACK
86
88 : n_rows(0)
89 , n_cols(0)
90 , symbolic_decomposition(nullptr)
91 , numeric_decomposition(nullptr)
92 , control(UMFPACK_CONTROL)
93{
94 umfpack_dl_defaults(control.data());
95}
96
97
98
99void
101{
102 // delete objects that haven't been deleted yet
103 if (symbolic_decomposition != nullptr)
104 {
105 umfpack_dl_free_symbolic(&symbolic_decomposition);
106 symbolic_decomposition = nullptr;
107 }
108
109 if (numeric_decomposition != nullptr)
110 {
111 umfpack_dl_free_numeric(&numeric_decomposition);
112 numeric_decomposition = nullptr;
113 }
114
115 {
116 std::vector<types::suitesparse_index> tmp;
117 tmp.swap(Ap);
118 }
119
120 {
121 std::vector<types::suitesparse_index> tmp;
122 tmp.swap(Ai);
123 }
124
125 {
126 std::vector<double> tmp;
127 tmp.swap(Ax);
128 }
129
130 {
131 std::vector<double> tmp;
132 tmp.swap(Az);
133 }
134
135 umfpack_dl_defaults(control.data());
136}
137
138
139
140template <typename number>
141void
143{
144 // do the copying around of entries so that the diagonal entry is in the
145 // right place. note that this is easy to detect: since all entries apart
146 // from the diagonal entry are sorted, we know that the diagonal entry is
147 // in the wrong place if and only if its column index is larger than the
148 // column index of the second entry in a row
149 //
150 // ignore rows with only one or no entry
152 0,
153 matrix.m(),
154 [this](const size_type row_begin, const size_type row_end) {
155 for (size_type row = row_begin; row < row_end; ++row)
156 {
157 // we may have to move some elements that are left of the diagonal
158 // but presently after the diagonal entry to the left, whereas the
159 // diagonal entry has to move to the right. we could first figure out
160 // where to move everything to, but for simplicity we just make a
161 // series of swaps instead (this is kind of a single run of
162 // bubble-sort, which gives us the desired result since the array is
163 // already "almost" sorted)
164 //
165 // in the first loop, the condition in the while-header also checks
166 // that the row has at least two entries and that the diagonal entry
167 // is really in the wrong place
168 long int cursor = Ap[row];
169 while ((cursor < Ap[row + 1] - 1) && (Ai[cursor] > Ai[cursor + 1]))
170 {
171 std::swap(Ai[cursor], Ai[cursor + 1]);
172
173 std::swap(Ax[cursor], Ax[cursor + 1]);
174 if (numbers::NumberTraits<number>::is_complex == true)
175 std::swap(Az[cursor], Az[cursor + 1]);
176
177 ++cursor;
178 }
179 }
180 },
181 parallel_grainsize(matrix));
182}
183
184
185
186template <typename number>
187void
189{
190 // same thing for SparseMatrixEZ
192 0,
193 matrix.m(),
194 [this](const size_type row_begin, const size_type row_end) {
195 for (size_type row = row_begin; row < row_end; ++row)
196 {
197 long int cursor = Ap[row];
198 while ((cursor < Ap[row + 1] - 1) && (Ai[cursor] > Ai[cursor + 1]))
199 {
200 std::swap(Ai[cursor], Ai[cursor + 1]);
201
202 std::swap(Ax[cursor], Ax[cursor + 1]);
203 if (numbers::NumberTraits<number>::is_complex == true)
204 std::swap(Az[cursor], Az[cursor + 1]);
205
206 ++cursor;
207 }
208 }
209 },
210 parallel_grainsize(matrix));
211}
212
213
214
215template <typename number>
216void
218{
219 // the case for block matrices is a bit more difficult, since all we know
220 // is that *within each block*, the diagonal of that block may come
221 // first. however, that means that there may be as many entries per row
222 // in the wrong place as there are block columns. we can do the same
223 // thing as above, but we have to do it multiple times
225 0,
226 matrix.m(),
227 [this, &matrix](const size_type row_begin, const size_type row_end) {
228 for (size_type row = row_begin; row < row_end; ++row)
229 {
230 long int cursor = Ap[row];
231 for (size_type block = 0; block < matrix.n_block_cols(); ++block)
232 {
233 // find the next out-of-order element
234 while ((cursor < Ap[row + 1] - 1) &&
235 (Ai[cursor] < Ai[cursor + 1]))
236 ++cursor;
237
238 // if there is none, then just go on
239 if (cursor == Ap[row + 1] - 1)
240 break;
241
242 // otherwise swap this entry with successive ones as long as
243 // necessary
244 long int element = cursor;
245 while ((element < Ap[row + 1] - 1) &&
246 (Ai[element] > Ai[element + 1]))
247 {
248 std::swap(Ai[element], Ai[element + 1]);
249
250 std::swap(Ax[element], Ax[element + 1]);
251 if (numbers::NumberTraits<number>::is_complex == true)
252 std::swap(Az[element], Az[element + 1]);
253
254 ++element;
255 }
256 }
257 }
258 },
259 parallel_grainsize(matrix));
260}
261
262
263
264template <class Matrix>
265void
267{
268 Assert(matrix.m() == matrix.n(), ExcNotQuadratic());
269
270 clear();
271
272 using number = typename Matrix::value_type;
273
274 n_rows = matrix.m();
275 n_cols = matrix.n();
276
277 const size_type N = matrix.m();
278
279 // copy over the data from the matrix to the data structures UMFPACK
280 // wants. note two things: first, UMFPACK wants compressed column storage
281 // whereas we always do compressed row storage; we work around this by,
282 // rather than shuffling things around, copy over the data we have, but
283 // then call the umfpack_dl_solve function with the UMFPACK_At argument,
284 // meaning that we want to solve for the transpose system
285 //
286 // second: the data we have in the sparse matrices is "almost" right
287 // already; UMFPACK wants the entries in each row (i.e. really: column)
288 // to be sorted in ascending order. we almost have that, except that we
289 // usually store the diagonal first in each row to allow for some
290 // optimizations. thus, we have to resort things a little bit, but only
291 // within each row
292 //
293 // final note: if the matrix has entries in the sparsity pattern that are
294 // actually occupied by entries that have a zero numerical value, then we
295 // keep them anyway. people are supposed to provide accurate sparsity
296 // patterns.
297 Ap.resize(N + 1);
298 Ai.resize(matrix.n_nonzero_elements());
299 Ax.resize(matrix.n_nonzero_elements());
301 Az.resize(matrix.n_nonzero_elements());
302
303 // first fill row lengths array
304 Ap[0] = 0;
305 for (size_type row = 1; row <= N; ++row)
306 Ap[row] = Ap[row - 1] + matrix.get_row_length(row - 1);
307 Assert(static_cast<size_type>(Ap.back()) == Ai.size(), ExcInternalError());
308
309 // then copy over matrix elements. note that for sparse matrices,
310 // iterators are sorted so that they traverse each row from start to end
311 // before moving on to the next row.
313 0,
314 matrix.m(),
315 [this, &matrix](const size_type row_begin, const size_type row_end) {
316 for (size_type row = row_begin; row < row_end; ++row)
317 {
318 long int index = Ap[row];
319 for (typename Matrix::const_iterator p = matrix.begin(row);
320 p != matrix.end(row);
321 ++p)
322 {
323 // write entry into the first free one for this row
324 Ai[index] = p->column();
325 Ax[index] = std::real(p->value());
326 if (numbers::NumberTraits<number>::is_complex == true)
327 Az[index] = std::imag(p->value());
328
329 // then move pointer ahead
330 ++index;
331 }
332 Assert(index == Ap[row + 1], ExcInternalError());
333 }
334 },
335 parallel_grainsize(matrix));
336
337 // make sure that the elements in each row are sorted. we have to be more
338 // careful for block sparse matrices, so ship this task out to a
339 // different function
340 sort_arrays(matrix);
341
342 int status;
344 status = umfpack_dl_symbolic(N,
345 N,
346 Ap.data(),
347 Ai.data(),
348 Ax.data(),
349 &symbolic_decomposition,
350 control.data(),
351 nullptr);
352 else
353 status = umfpack_zl_symbolic(N,
354 N,
355 Ap.data(),
356 Ai.data(),
357 Ax.data(),
358 Az.data(),
359 &symbolic_decomposition,
360 control.data(),
361 nullptr);
362 AssertThrow(status == UMFPACK_OK,
363 ExcUMFPACKError("umfpack_dl_symbolic", status));
364
366 status = umfpack_dl_numeric(Ap.data(),
367 Ai.data(),
368 Ax.data(),
369 symbolic_decomposition,
370 &numeric_decomposition,
371 control.data(),
372 nullptr);
373 else
374 status = umfpack_zl_numeric(Ap.data(),
375 Ai.data(),
376 Ax.data(),
377 Az.data(),
378 symbolic_decomposition,
379 &numeric_decomposition,
380 control.data(),
381 nullptr);
382
383 // Clean up before we deal with the error code from the calls above:
384 umfpack_dl_free_symbolic(&symbolic_decomposition);
385 if (status == UMFPACK_WARNING_singular_matrix)
386 {
387 // UMFPACK sometimes warns that the matrix is singular, but that a
388 // factorization was successful nonetheless. Report this by
389 // throwing an exception that can be caught at a higher level:
390 AssertThrow(false,
392 "UMFPACK reports that the matrix is singular, "
393 "but that the factorization was successful anyway. "
394 "You can try and see whether you can still "
395 "solve a linear system with such a factorization "
396 "by catching and ignoring this exception, "
397 "though in practice this will typically not "
398 "work."));
399 }
400 else
401 AssertThrow(status == UMFPACK_OK,
402 ExcUMFPACKError("umfpack_dl_numeric", status));
403}
404
405
406
407void
409 const bool transpose /*=false*/) const
410{
411 // make sure that some kind of factorize() call has happened before
412 Assert(Ap.size() != 0, ExcNotInitialized());
413 Assert(Ai.size() != 0, ExcNotInitialized());
414 Assert(Ai.size() == Ax.size(), ExcNotInitialized());
415
416 Assert(Az.empty(),
417 ExcMessage("You have previously factored a matrix using this class "
418 "that had complex-valued entries. This then requires "
419 "applying the factored matrix to a complex-valued "
420 "vector, but you are only providing a real-valued vector "
421 "here."));
422
423 Vector<double> rhs(rhs_and_solution.size());
424 rhs = rhs_and_solution;
425
426 // solve the system. note that since UMFPACK wants compressed column
427 // storage instead of the compressed row storage format we use in
428 // deal.II's SparsityPattern classes, we solve for UMFPACK's A^T instead
429
430 // Conversely, if we solve for the transpose, we have to use UMFPACK_A
431 // instead.
432 const int status = umfpack_dl_solve(transpose ? UMFPACK_A : UMFPACK_At,
433 Ap.data(),
434 Ai.data(),
435 Ax.data(),
436 rhs_and_solution.begin(),
437 rhs.begin(),
439 control.data(),
440 nullptr);
441 AssertThrow(status == UMFPACK_OK,
442 ExcUMFPACKError("umfpack_dl_solve", status));
443}
444
445
446
447void
448SparseDirectUMFPACK::solve(Vector<std::complex<double>> &rhs_and_solution,
449 const bool transpose /*=false*/) const
450{
451# ifdef DEAL_II_WITH_COMPLEX_VALUES
452 // make sure that some kind of factorize() call has happened before
453 Assert(Ap.size() != 0, ExcNotInitialized());
454 Assert(Ai.size() != 0, ExcNotInitialized());
455 Assert(Ai.size() == Ax.size(), ExcNotInitialized());
456
457 // First see whether the matrix that was factorized was complex-valued.
458 // If so, just apply the complex factorization to the vector.
459 if (Az.size() != 0)
460 {
461 Assert(Ax.size() == Az.size(), ExcInternalError());
462
463 // It would be nice if we could just present a pointer to the
464 // first element of the complex-valued solution vector and let
465 // UMFPACK fill both the real and imaginary parts of the solution
466 // vector at that address. UMFPACK calls this the 'packed' format,
467 // and in those cases it only takes one pointer to the entire
468 // vector, rather than a pointer to the real and one pointer to
469 // the imaginary parts of the vector. The problem is that if we
470 // want to pack, then we also need to pack the matrix, and the
471 // functions above have already decided that we don't want to pack
472 // the matrix but instead deal with split format for the matrix,
473 // and then UMFPACK decides that it can't deal with a split matrix
474 // and a packed vector. We have to choose one or the other, not
475 // mix.
476 //
477 // So create four vectors, one each for the real and imaginary parts
478 // of the right hand side and solution.
479 Vector<double> rhs_re(rhs_and_solution.size());
480 Vector<double> rhs_im(rhs_and_solution.size());
481 for (unsigned int i = 0; i < rhs_and_solution.size(); ++i)
482 {
483 rhs_re(i) = std::real(rhs_and_solution(i));
484 rhs_im(i) = std::imag(rhs_and_solution(i));
485 }
486
487 Vector<double> solution_re(rhs_and_solution.size());
488 Vector<double> solution_im(rhs_and_solution.size());
489
490 // Solve the system. note that since UMFPACK wants compressed column
491 // storage instead of the compressed row storage format we use in
492 // deal.II's SparsityPattern classes, we solve for UMFPACK's A^T instead
493 //
494 // Conversely, if we solve for the transpose, we have to use UMFPACK_A
495 // instead.
496 //
497 // Note that for the complex case here, the transpose is selected using
498 // UMFPACK_Aat, not UMFPACK_At.
499 const int status = umfpack_zl_solve(transpose ? UMFPACK_A : UMFPACK_Aat,
500 Ap.data(),
501 Ai.data(),
502 Ax.data(),
503 Az.data(),
504 solution_re.data(),
505 solution_im.data(),
506 rhs_re.data(),
507 rhs_im.data(),
509 control.data(),
510 nullptr);
511 AssertThrow(status == UMFPACK_OK,
512 ExcUMFPACKError("umfpack_dl_solve", status));
513
514 // Now put things back together into the output vector
515 for (unsigned int i = 0; i < rhs_and_solution.size(); ++i)
516 rhs_and_solution(i) = {solution_re(i), solution_im(i)};
517 }
518 else
519 {
520 // We have factorized a real-valued matrix, but the rhs and solution
521 // vectors are complex-valued. UMFPACK does not natively support this
522 // case, but we can just apply the factorization to real and imaginary
523 // parts of the right hand side separately
524 const Vector<std::complex<double>> rhs = rhs_and_solution;
525
526 // Get the real part of the right hand side, solve with it, and copy the
527 // results into the result vector by just copying the real output
528 // into the complex-valued result vector (which implies setting the
529 // imaginary parts to zero):
530 Vector<double> rhs_real_or_imag(rhs_and_solution.size());
531 for (unsigned int i = 0; i < rhs.size(); ++i)
532 rhs_real_or_imag(i) = std::real(rhs(i));
533
534 solve(rhs_real_or_imag, transpose);
535
536 rhs_and_solution = rhs_real_or_imag;
537
538 // Then repeat the whole thing with the imaginary part. The copying step
539 // is more complicated now because we can only touch the imaginary
540 // component of the output vector:
541 for (unsigned int i = 0; i < rhs.size(); ++i)
542 rhs_real_or_imag(i) = std::imag(rhs(i));
543
544 solve(rhs_real_or_imag, transpose);
545
546 for (unsigned int i = 0; i < rhs.size(); ++i)
547 rhs_and_solution(i).imag(rhs_real_or_imag(i));
548 }
549
550# else
551
552 (void)rhs_and_solution;
553 (void)transpose;
554 Assert(false,
556 "This function can't be called if deal.II has been configured "
557 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
558# endif
559}
560
561
562void
564 const bool transpose /*=false*/) const
565{
566 // the UMFPACK functions want a contiguous array of elements, so
567 // there is no way around copying data around. thus, just copy the
568 // data into a regular vector and back
569 Vector<double> tmp(rhs_and_solution.size());
570 tmp = rhs_and_solution;
571 solve(tmp, transpose);
572 rhs_and_solution = tmp;
573}
574
575
576
577void
578SparseDirectUMFPACK::solve(BlockVector<std::complex<double>> &rhs_and_solution,
579 const bool transpose /*=false*/) const
580{
581# ifdef DEAL_II_WITH_COMPLEX_VALUES
582 // the UMFPACK functions want a contiguous array of elements, so
583 // there is no way around copying data around. thus, just copy the
584 // data into a regular vector and back
585 Vector<std::complex<double>> tmp(rhs_and_solution.size());
586 tmp = rhs_and_solution;
587 solve(tmp, transpose);
588 rhs_and_solution = tmp;
589
590# else
591 (void)rhs_and_solution;
592 (void)transpose;
593 Assert(false,
595 "This function can't be called if deal.II has been configured "
596 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
597# endif
598}
599
600
601
602template <class Matrix>
603void
604SparseDirectUMFPACK::solve(const Matrix &matrix,
605 Vector<double> &rhs_and_solution,
606 const bool transpose /*=false*/)
607{
608 factorize(matrix);
609 solve(rhs_and_solution, transpose);
610}
611
612
613
614template <class Matrix>
615void
616SparseDirectUMFPACK::solve(const Matrix &matrix,
617 Vector<std::complex<double>> &rhs_and_solution,
618 const bool transpose /*=false*/)
619{
620# ifdef DEAL_II_WITH_COMPLEX_VALUES
621 factorize(matrix);
622 solve(rhs_and_solution, transpose);
623
624# else
625
626 (void)matrix;
627 (void)rhs_and_solution;
628 (void)transpose;
629 Assert(false,
631 "This function can't be called if deal.II has been configured "
632 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
633# endif
634}
635
636
637
638template <class Matrix>
639void
640SparseDirectUMFPACK::solve(const Matrix &matrix,
641 BlockVector<double> &rhs_and_solution,
642 const bool transpose /*=false*/)
643{
644 factorize(matrix);
645 solve(rhs_and_solution, transpose);
646}
647
648
649
650template <class Matrix>
651void
652SparseDirectUMFPACK::solve(const Matrix &matrix,
653 BlockVector<std::complex<double>> &rhs_and_solution,
654 const bool transpose /*=false*/)
655{
656# ifdef DEAL_II_WITH_COMPLEX_VALUES
657 factorize(matrix);
658 solve(rhs_and_solution, transpose);
659
660# else
661
662 (void)matrix;
663 (void)rhs_and_solution;
664 (void)transpose;
665 Assert(false,
667 "This function can't be called if deal.II has been configured "
668 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
669# endif
670}
671
672
673#else
674
675
677 : n_rows(0)
678 , n_cols(0)
679 , symbolic_decomposition(nullptr)
680 , numeric_decomposition(nullptr)
681 , control(0)
682{}
683
684
685void
687{}
688
689
690template <class Matrix>
691void
692SparseDirectUMFPACK::factorize(const Matrix &)
693{
695 false,
697 "To call this function you need UMFPACK, but you configured deal.II "
698 "without passing the necessary switch to 'cmake'. Please consult the "
699 "installation instructions at https://dealii.org/current/readme.html"));
700}
701
702
703void
705{
707 false,
709 "To call this function you need UMFPACK, but you configured deal.II "
710 "without passing the necessary switch to 'cmake'. Please consult the "
711 "installation instructions at https://dealii.org/current/readme.html"));
712}
713
714
715
716void
717SparseDirectUMFPACK::solve(Vector<std::complex<double>> &, const bool) const
718{
720 false,
722 "To call this function you need UMFPACK, but you configured deal.II "
723 "without passing the necessary switch to 'cmake'. Please consult the "
724 "installation instructions at https://dealii.org/current/readme.html"));
725}
726
727
728
729void
731{
733 false,
735 "To call this function you need UMFPACK, but you configured deal.II "
736 "without passing the necessary switch to 'cmake'. Please consult the "
737 "installation instructions at https://dealii.org/current/readme.html"));
738}
739
740
741
742void
743SparseDirectUMFPACK::solve(BlockVector<std::complex<double>> &,
744 const bool) const
745{
747 false,
749 "To call this function you need UMFPACK, but you configured deal.II "
750 "without passing the necessary switch to 'cmake'. Please consult the "
751 "installation instructions at https://dealii.org/current/readme.html"));
752}
753
754
755
756template <class Matrix>
757void
758SparseDirectUMFPACK::solve(const Matrix &, Vector<double> &, const bool)
759{
761 false,
763 "To call this function you need UMFPACK, but you configured deal.II "
764 "without passing the necessary switch to 'cmake'. Please consult the "
765 "installation instructions at https://dealii.org/current/readme.html"));
766}
767
768
769
770template <class Matrix>
771void
772SparseDirectUMFPACK::solve(const Matrix &,
773 Vector<std::complex<double>> &,
774 const bool)
775{
777 false,
779 "To call this function you need UMFPACK, but you configured deal.II "
780 "without passing the necessary switch to 'cmake'. Please consult the "
781 "installation instructions at https://dealii.org/current/readme.html"));
782}
783
784
785
786template <class Matrix>
787void
788SparseDirectUMFPACK::solve(const Matrix &, BlockVector<double> &, const bool)
789{
791 false,
793 "To call this function you need UMFPACK, but you configured deal.II "
794 "without passing the necessary switch to 'cmake'. Please consult the "
795 "installation instructions at https://dealii.org/current/readme.html"));
796}
797
798
799
800template <class Matrix>
801void
802SparseDirectUMFPACK::solve(const Matrix &,
803 BlockVector<std::complex<double>> &,
804 const bool)
805{
807 false,
809 "To call this function you need UMFPACK, but you configured deal.II "
810 "without passing the necessary switch to 'cmake'. Please consult the "
811 "installation instructions at https://dealii.org/current/readme.html"));
812}
813
814#endif
815
816
817template <class Matrix>
818void
820{
821 this->factorize(M);
822}
823
824
825void
827{
828 dst = src;
829 this->solve(dst);
830}
831
832
833
834void
836 const BlockVector<double> &src) const
837{
838 dst = src;
839 this->solve(dst);
840}
841
842
843void
845 const Vector<double> &src) const
846{
847 dst = src;
848 this->solve(dst, /*transpose=*/true);
849}
850
851
852
853void
855 const BlockVector<double> &src) const
856{
857 dst = src;
858 this->solve(dst, /*transpose=*/true);
859}
860
863{
865 return n_rows;
866}
867
870{
872 return n_cols;
873}
874
875
876
877#ifdef DEAL_II_WITH_MUMPS
878
880 const MPI_Comm &communicator)
882 , mpi_communicator(communicator)
883{
884 // Initialize MUMPS instance:
885 id.job = -1;
886 id.par = 1;
887
888 Assert(!(additional_data.symmetric == false &&
889 additional_data.posdef == true),
891 "You can't have a positive definite matrix that is not symmetric."));
892
893 if (additional_data.symmetric == true)
894 {
895 if (additional_data.posdef == true)
896 id.sym = 1;
897 else
898 id.sym = 2;
899 }
900 else
901 id.sym = 0;
902
903 id.comm_fortran = (MUMPS_INT)MPI_Comm_c2f(mpi_communicator);
904 dmumps_c(&id);
905
906 if (additional_data.output_details == false)
907 {
908 // No outputs
909 id.icntl[0] = -1;
910 id.icntl[1] = -1;
911 id.icntl[2] = -1;
912 id.icntl[3] = 0;
913 }
914
915 if (additional_data.error_statistics == true)
916 id.icntl[10] = 2;
917
918 if (additional_data.blr_factorization == true)
919 {
920 id.icntl[34] = 2;
921 if (additional_data.blr.blr_ucfs == true)
922 id.icntl[35] = 1;
923 Assert(additional_data.blr.lowrank_threshold > 0,
924 ExcMessage("Lowrank threshold must be positive."));
925 id.cntl[6] = additional_data.blr.lowrank_threshold;
926 }
927}
928
929
930
932{
933 // MUMPS destructor
934 id.job = -2;
935 dmumps_c(&id);
936}
937
938
939
940template <class Matrix>
941void
943{
944 Assert(matrix.n() == matrix.m(), ExcMessage("Matrix needs to be square."));
945
946 n = matrix.n();
947 id.n = n;
948
949 if constexpr (std::is_same_v<Matrix, SparseMatrix<double>>)
950 {
951 // Serial matrix: hand over matrix to MUMPS as centralized assembled
952 // matrix
954 {
955 // number of nonzero elements in matrix
956 nnz = matrix.n_actually_nonzero_elements();
957
958 // representation of the matrix
959 a = std::make_unique<double[]>(nnz);
960
961 // matrix indices pointing to the row and column dimensions
962 // respectively of the matrix representation above (a): ie. a[k] is
963 // the matrix element (irn[k], jcn[k])
964 irn = std::make_unique<MUMPS_INT[]>(nnz);
965 jcn = std::make_unique<MUMPS_INT[]>(nnz);
966
967 size_type n_non_zero_elements = 0;
968
969 // loop over the elements of the matrix row by row, as suggested in
970 // the documentation of the sparse matrix iterator class
971 if (additional_data.symmetric == true)
972 {
973 for (size_type row = 0; row < matrix.m(); ++row)
974 {
975 for (typename Matrix::const_iterator ptr = matrix.begin(row);
976 ptr != matrix.end(row);
977 ++ptr)
978 if (std::abs(ptr->value()) > 0.0 && ptr->column() >= row)
979 {
980 a[n_non_zero_elements] = ptr->value();
981 irn[n_non_zero_elements] = row + 1;
982 jcn[n_non_zero_elements] = ptr->column() + 1;
983
984 ++n_non_zero_elements;
985 }
986 }
987 }
988 else
989 {
990 for (size_type row = 0; row < matrix.m(); ++row)
991 {
992 for (typename Matrix::const_iterator ptr = matrix.begin(row);
993 ptr != matrix.end(row);
994 ++ptr)
995 if (std::abs(ptr->value()) > 0.0)
996 {
997 a[n_non_zero_elements] = ptr->value();
998 irn[n_non_zero_elements] = row + 1;
999 jcn[n_non_zero_elements] = ptr->column() + 1;
1000 ++n_non_zero_elements;
1001 }
1002 }
1003 }
1004 id.n = n;
1005 id.nnz = n_non_zero_elements;
1006 id.irn = irn.get();
1007 id.jcn = jcn.get();
1008 id.a = a.get();
1009 }
1010 }
1011 else if constexpr (std::is_same_v<Matrix, TrilinosWrappers::SparseMatrix> ||
1012 std::is_same_v<Matrix, PETScWrappers::MPI::SparseMatrix>)
1013 {
1014 int result;
1015 MPI_Comm_compare(mpi_communicator,
1016 matrix.get_mpi_communicator(),
1017 &result);
1018 AssertThrow(result == MPI_IDENT,
1019 ExcMessage("The matrix communicator must match the MUMPS "
1020 "communicator."));
1021
1022 // Distributed matrix case
1023 id.icntl[17] = 3; // distributed matrix assembly
1024 id.nnz = matrix.n_nonzero_elements();
1025 nnz = id.nnz;
1026 size_type n_non_zero_local = 0;
1027
1028 // Get the range of rows owned by this process
1029 locally_owned_rows = matrix.locally_owned_range_indices();
1030 size_type local_non_zeros = 0;
1031
1032 if constexpr (std::is_same_v<Matrix, TrilinosWrappers::SparseMatrix>)
1033 {
1034 const auto &trilinos_matrix = matrix.trilinos_matrix();
1035 local_non_zeros = trilinos_matrix.NumMyNonzeros();
1036 }
1037 else if constexpr (std::is_same_v<Matrix,
1039 {
1040# ifdef DEAL_II_WITH_PETSC
1041 Mat &petsc_matrix =
1042 const_cast<PETScWrappers::MPI::SparseMatrix &>(matrix)
1043 .petsc_matrix();
1044 MatInfo info;
1045 MatGetInfo(petsc_matrix, MAT_LOCAL, &info);
1046 local_non_zeros = (size_type)info.nz_used;
1047# endif
1048 }
1049
1050
1051 // We allocate enough entries for the general, nonsymmetric, case. If case
1052 // of a symmetric matrix, we will end up with fewer entries.
1053 irn = std::make_unique<MUMPS_INT[]>(local_non_zeros);
1054 jcn = std::make_unique<MUMPS_INT[]>(local_non_zeros);
1055 a = std::make_unique<double[]>(local_non_zeros);
1056 irhs_loc.resize(locally_owned_rows.n_elements());
1057
1058 if (additional_data.symmetric == true)
1059 {
1060 if constexpr (std::is_same_v<Matrix,
1062 {
1063# ifdef DEAL_II_WITH_PETSC
1064 Mat &petsc_matrix =
1065 const_cast<PETScWrappers::MPI::SparseMatrix &>(matrix)
1066 .petsc_matrix();
1067
1068 PetscInt rstart, rend;
1069 MatGetOwnershipRange(petsc_matrix, &rstart, &rend);
1070 for (PetscInt i = rstart; i < rend; i++)
1071 {
1072 PetscInt n_cols;
1073 const PetscInt *cols;
1074 const PetscScalar *values;
1075 MatGetRow(petsc_matrix, i, &n_cols, &cols, &values);
1076
1077 for (PetscInt j = 0; j < n_cols; j++)
1078 {
1079 if (cols[j] >= i)
1080 {
1081 irn[n_non_zero_local] = i + 1;
1082 jcn[n_non_zero_local] = cols[j] + 1;
1083 a[n_non_zero_local] = values[j];
1084
1085 // Count local non-zeros
1086 n_non_zero_local++;
1087 }
1088 }
1089 MatRestoreRow(petsc_matrix, i, &n_cols, &cols, &values);
1090
1091 // Store the row index for the rhs vector
1092 const types::global_cell_index local_index =
1093 locally_owned_rows.index_within_set(i);
1094 irhs_loc[local_index] = i + 1;
1095 }
1096
1097 id.a_loc = a.get();
1098# endif
1099 }
1100 else if constexpr (std::is_same_v<Matrix,
1102 {
1103 const auto &trilinos_matrix = matrix.trilinos_matrix();
1104 for (int local_row = 0; local_row < trilinos_matrix.NumMyRows();
1105 ++local_row)
1106 {
1107 int num_entries;
1108 double *values;
1109 int *local_cols;
1110 int ierr = trilinos_matrix.ExtractMyRowView(local_row,
1111 num_entries,
1112 values,
1113 local_cols);
1114 (void)ierr;
1115 Assert(
1116 ierr == 0,
1117 ExcMessage(
1118 "Error extracting global row view from Trilinos matrix. Error code " +
1119 std::to_string(ierr) + "."));
1120
1121 int global_row = trilinos_matrix.GRID(local_row);
1122
1123 for (int j = 0; j < num_entries; ++j)
1124 {
1125 if (trilinos_matrix.GCID(local_cols[j]) >= global_row)
1126 {
1127 irn[n_non_zero_local] = global_row + 1;
1128 jcn[n_non_zero_local] =
1129 trilinos_matrix.GCID(local_cols[j]) + 1;
1130 a[n_non_zero_local] = values[j];
1131
1132 // Count local non-zeros
1133 n_non_zero_local++;
1134 }
1135 }
1136
1137 // Store the row index for the rhs vector
1138 irhs_loc[local_row] = global_row + 1;
1139 }
1140 id.a_loc = a.get();
1141 }
1142 else
1143 {
1145 }
1146 }
1147 else
1148 {
1149 // Unsymmetric case
1150 if constexpr (std::is_same_v<Matrix,
1152 {
1153# ifdef DEAL_II_WITH_PETSC
1154 Mat &petsc_matrix =
1155 const_cast<PETScWrappers::MPI::SparseMatrix &>(matrix)
1156 .petsc_matrix();
1157
1158 PetscInt rstart, rend;
1159 MatGetOwnershipRange(petsc_matrix, &rstart, &rend);
1160 for (PetscInt i = rstart; i < rend; i++)
1161 {
1162 PetscInt n_cols;
1163 const PetscInt *cols;
1164 const PetscScalar *values;
1165 MatGetRow(petsc_matrix, i, &n_cols, &cols, &values);
1166
1167 for (PetscInt j = 0; j < n_cols; j++)
1168 {
1169 irn[n_non_zero_local] = i + 1;
1170 jcn[n_non_zero_local] = cols[j] + 1;
1171 a[n_non_zero_local] = values[j];
1172
1173 // Count local non-zeros
1174 n_non_zero_local++;
1175 }
1176 MatRestoreRow(petsc_matrix, i, &n_cols, &cols, &values);
1177
1178 // Store the row index for the rhs vector
1179 const types::global_cell_index local_index =
1180 locally_owned_rows.index_within_set(i);
1181 irhs_loc[local_index] = i + 1;
1182 }
1183
1184 id.a_loc = a.get();
1185# endif
1186 }
1187 else if constexpr (std::is_same_v<Matrix,
1189 {
1190 const auto &trilinos_matrix = matrix.trilinos_matrix();
1191 for (int local_row = 0; local_row < trilinos_matrix.NumMyRows();
1192 ++local_row)
1193 {
1194 int num_entries;
1195 double *values;
1196 int *local_cols;
1197 int ierr = trilinos_matrix.ExtractMyRowView(local_row,
1198 num_entries,
1199 values,
1200 local_cols);
1201 (void)ierr;
1202 Assert(
1203 ierr == 0,
1204 ExcMessage(
1205 "Error extracting global row view from Trilinos matrix. Error code " +
1206 std::to_string(ierr) + "."));
1207
1208 int global_row = trilinos_matrix.GRID(local_row);
1209
1210 for (int j = 0; j < num_entries; ++j)
1211 {
1212 irn[n_non_zero_local] = global_row + 1;
1213 jcn[n_non_zero_local] =
1214 trilinos_matrix.GCID(local_cols[j]) + 1;
1215 a[n_non_zero_local] = values[j];
1216
1217 // Count local non-zeros
1218 n_non_zero_local++;
1219 }
1220
1221 // Store the row index for the rhs vector
1222 irhs_loc[local_row] = global_row + 1;
1223 }
1224 id.a_loc = a.get();
1225 }
1226 else
1227 {
1229 }
1230 }
1231
1232 // Hand over local arrays to MUMPS
1233 id.nnz_loc = n_non_zero_local;
1234 id.irn_loc = irn.get();
1235 id.jcn_loc = jcn.get();
1236 id.a_loc = a.get();
1237 id.irhs_loc = irhs_loc.data();
1238
1239 // rhs parameters
1240 id.icntl[19] = 10; // distributed rhs
1241 id.icntl[20] = 0; // centralized solution, stored on rank 0 by MUMPS
1242 id.nrhs = 1;
1243 id.lrhs_loc = n;
1244 id.nloc_rhs = locally_owned_rows.n_elements();
1245 }
1246 else
1247 {
1249 }
1250}
1251
1252
1253
1254void
1256{
1257 Assert(n == new_rhs.size(),
1258 ExcMessage("Matrix size and rhs length must be equal."));
1259
1261 {
1262 rhs.resize(n);
1263 for (size_type i = 0; i < n; ++i)
1264 rhs[i] = new_rhs(i);
1265
1266 id.rhs = &rhs[0];
1267 }
1268}
1269
1270
1271
1272void
1274{
1275 Assert(n == vector.size(),
1276 ExcMessage("Matrix size and solution vector length must be equal."));
1277 Assert(n == rhs.size(),
1278 ExcMessage("Class not initialized with a rhs vector."));
1279
1280 // Copy solution into the given vector
1282 {
1283 for (size_type i = 0; i < n; ++i)
1284 vector(i) = rhs[i];
1285
1286 rhs.resize(0); // remove rhs again
1287 }
1288}
1289
1290
1291
1292template <class Matrix>
1293void
1295{
1296 // Initialize MUMPS instance:
1297 initialize_matrix(matrix);
1298
1299 // Start analysis + factorization
1300 id.job = 4;
1301 dmumps_c(&id);
1302}
1303
1304
1305
1306template <typename VectorType>
1307void
1308SparseDirectMUMPS::vmult(VectorType &dst, const VectorType &src) const
1309{
1310 // Check that the matrix has at least one nonzero element:
1311 Assert(nnz != 0, ExcNotInitialized());
1312 Assert(n == dst.size(), ExcMessage("Destination vector has the wrong size."));
1313 Assert(n == src.size(), ExcMessage("Source vector has the wrong size."));
1314
1315
1316 if constexpr (std::is_same_v<VectorType, Vector<double>>)
1317 {
1318 // Centralized assembly for serial vectors.
1319
1320 // Hand over right-hand side
1321 copy_rhs_to_mumps(src);
1322
1323 // Start solver
1324 id.job = 3;
1325 dmumps_c(&id);
1326 copy_solution(dst);
1327 }
1328 else if constexpr (std::is_same_v<VectorType,
1330 std::is_same_v<VectorType, PETScWrappers::MPI::Vector>)
1331 {
1332 if constexpr (std::is_same_v<VectorType, TrilinosWrappers::MPI::Vector>)
1333 id.rhs_loc = const_cast<double *>(src.begin());
1334 else if constexpr (std::is_same_v<VectorType, PETScWrappers::MPI::Vector>)
1335 {
1336# ifdef DEAL_II_WITH_PETSC
1337 PetscScalar *local_array;
1338 VecGetArray(
1339 const_cast<PETScWrappers::MPI::Vector &>(src).petsc_vector(),
1340 &local_array);
1341 id.rhs_loc = local_array;
1342 VecRestoreArray(
1343 const_cast<PETScWrappers::MPI::Vector &>(src).petsc_vector(),
1344 &local_array);
1345# endif
1346 }
1347
1348
1350 {
1351 rhs.resize(id.lrhs_loc);
1352 id.rhs = rhs.data();
1353 }
1354
1355 // Start solver
1356 id.job = 3;
1357 dmumps_c(&id);
1358
1359 // Copy solution into the given vector
1360 // For MUMPS with centralized solution (icntl[20]=0), the solution is only
1361 // on the root process (0) and needs to be distributed to all processes
1363 {
1364 // Set all the values in the dst vector
1365 for (size_type i = 0; i < n; ++i)
1366 dst[i] = rhs[i];
1367 }
1368
1369 // Ensure the distributed vector has the correct values everywhere
1370 dst.compress(VectorOperation::insert);
1371
1372 rhs.resize(0); // remove rhs again
1373 }
1374 else
1375 {
1377 }
1378}
1379
1380
1381template <typename VectorType>
1382void
1383SparseDirectMUMPS::Tvmult(VectorType &dst, const VectorType &src) const
1384{
1385 // The matrix has at least one nonzero element:
1386 Assert(nnz != 0, ExcNotInitialized());
1387 Assert(n == dst.size(), ExcMessage("Destination vector has the wrong size."));
1388 Assert(n == src.size(), ExcMessage("Source vector has the wrong size."));
1389
1390 id.icntl[8] = 2; // transpose
1391 vmult(dst, src);
1392 id.icntl[8] = 1; // reset to default
1393}
1394
1395
1396
1397int *
1399{
1400 return id.icntl;
1401}
1402
1403
1404
1405#else
1406
1407
1408SparseDirectMUMPS::SparseDirectMUMPS(const AdditionalData &, const MPI_Comm &)
1409 : mpi_communicator(MPI_COMM_SELF)
1410{
1412 false,
1413 ExcMessage(
1414 "To call this function you need MUMPS, but you configured deal.II "
1415 "without passing the necessary switch to 'cmake'. Please consult the "
1416 "installation instructions at https://dealii.org/current/readme.html"));
1417}
1418
1419
1420
1422{}
1423
1424
1425#endif // DEAL_II_WITH_MUMPS
1426
1427
1428
1429// explicit instantiations for SparseMatrixUMFPACK
1430#define InstantiateUMFPACK(MatrixType) \
1431 template void SparseDirectUMFPACK::factorize(const MatrixType &); \
1432 template void SparseDirectUMFPACK::solve(const MatrixType &, \
1433 Vector<double> &, \
1434 const bool); \
1435 template void SparseDirectUMFPACK::solve(const MatrixType &, \
1436 Vector<std::complex<double>> &, \
1437 const bool); \
1438 template void SparseDirectUMFPACK::solve(const MatrixType &, \
1439 BlockVector<double> &, \
1440 const bool); \
1441 template void SparseDirectUMFPACK::solve( \
1442 const MatrixType &, BlockVector<std::complex<double>> &, const bool); \
1443 template void SparseDirectUMFPACK::initialize(const MatrixType &, \
1444 const AdditionalData)
1445
1446// Instantiate everything for real-valued matrices
1453
1454// Now also for complex-valued matrices
1455#ifdef DEAL_II_WITH_COMPLEX_VALUES
1456InstantiateUMFPACK(SparseMatrix<std::complex<double>>);
1457InstantiateUMFPACK(SparseMatrix<std::complex<float>>);
1460#endif
1461
1462// explicit instantiations for SparseDirectMUMPS
1463#ifdef DEAL_II_WITH_MUMPS
1464
1465# define InstantiateMUMPSMatVec(VECTOR) \
1466 template void SparseDirectMUMPS::vmult(VECTOR &, const VECTOR &) const; \
1467 template void SparseDirectMUMPS::Tvmult(VECTOR &, const VECTOR &) const;
1468# ifdef DEAL_II_WITH_TRILINOS
1470# endif
1471# ifdef DEAL_II_WITH_PETSC
1473# endif
1475
1476# define InstantiateMUMPS(MATRIX) \
1477 template void SparseDirectMUMPS::initialize(const MATRIX &);
1478
1481# ifdef DEAL_II_WITH_TRILINOS
1483# endif
1484# ifdef DEAL_II_WITH_PETSC
1487# endif
1488 // InstantiateMUMPS(SparseMatrixEZ<double>)
1489 // InstantiateMUMPS(SparseMatrixEZ<float>)
1492#endif
1493
virtual size_type size() const override
void copy_solution(Vector< double > &vector) const
void initialize_matrix(const Matrix &matrix)
void initialize(const Matrix &matrix)
std::vector< double > rhs
std::unique_ptr< double[]> a
void copy_rhs_to_mumps(const Vector< double > &rhs) const
types::global_dof_index n
const MPI_Comm mpi_communicator
SparseDirectMUMPS(const AdditionalData &additional_data=AdditionalData(), const MPI_Comm &communicator=MPI_COMM_WORLD)
std::unique_ptr< types::mumps_index[]> irn
std::unique_ptr< types::mumps_index[]> jcn
std::vector< types::mumps_index > irhs_loc
types::global_dof_index size_type
void vmult(VectorType &dst, const VectorType &src) const
types::mumps_nnz nnz
IndexSet locally_owned_rows
AdditionalData additional_data
void Tvmult(VectorType &, const VectorType &src) const
std::vector< double > Az
~SparseDirectUMFPACK() override
void initialize(const SparsityPattern &sparsity_pattern)
void Tvmult(Vector< double > &dst, const Vector< double > &src) const
size_type m() const
types::global_dof_index size_type
void solve(Vector< double > &rhs_and_solution, const bool transpose=false) const
std::vector< double > Ax
void sort_arrays(const SparseMatrixEZ< number > &)
size_type n() const
void factorize(const Matrix &matrix)
std::vector< double > control
void vmult(Vector< double > &dst, const Vector< double > &src) const
std::vector< types::suitesparse_index > Ap
std::vector< types::suitesparse_index > Ai
pointer data()
virtual size_type size() const override
iterator begin()
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcUMFPACKError(std::string arg1, int arg2)
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcNotQuadratic()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
@ matrix
Contents is actually a matrix.
SparseMatrix< double > SparseMatrix
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
Definition mpi.cc:120
void apply_to_subranges(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, const Function &f, const unsigned int grainsize)
Definition parallel.h:540
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int global_cell_index
Definition types.h:138
#define InstantiateUMFPACK(MatrixType)
#define InstantiateMUMPS(MATRIX)
#define InstantiateMUMPSMatVec(VECTOR)
static constexpr bool is_complex
Definition numbers.h:404