Loading...
Searching...
No Matches
system_solvers_split.cpp
1/*---------------------------------------------------------------------------*\
2 *
3 * bitpit
4 *
5 * Copyright (C) 2015-2021 OPTIMAD engineering Srl
6 *
7 * -------------------------------------------------------------------------
8 * License
9 * This file is part of bitpit.
10 *
11 * bitpit is free software: you can redistribute it and/or modify it
12 * under the terms of the GNU Lesser General Public License v3 (LGPL)
13 * as published by the Free Software Foundation.
14 *
15 * bitpit is distributed in the hope that it will be useful, but WITHOUT
16 * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
17 * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
18 * License for more details.
19 *
20 * You should have received a copy of the GNU Lesser General Public License
21 * along with bitpit. If not, see <http://www.gnu.org/licenses/>.
22 *
23\*---------------------------------------------------------------------------*/
24
25#include "matrix_utilities.hpp"
26#include "system_solvers_split.hpp"
27
28#include "bitpit_common.hpp"
29#include "bitpit_containers.hpp"
30#include "bitpit_IO.hpp"
31#include "bitpit_operators.hpp"
32
33#include "petscksp.h"
34#include "petscmat.h"
35#include "petscvec.h"
36
37#include <fstream>
38#include <stdexcept>
39#include <string>
40#include <numeric>
41#include <unordered_set>
42
43#ifndef PETSC_NULLPTR
44#define PETSC_NULLPTR PETSC_NULL
45#endif
46
47namespace bitpit {
48
63SplitSystemMatrixAssembler::SplitSystemMatrixAssembler(SystemSplitStrategy splitStrategy, const std::vector<int> &splitSizes)
65 m_splitStrategy(splitStrategy), m_splitSizes(splitSizes)
66{
67}
68
75{
76 return m_splitStrategy;
77}
78
85{
86 return m_splitSizes.size();
87}
88
108const std::vector<int> & SplitSystemMatrixAssembler::getSplitSizes() const
109{
110 return m_splitSizes;
111}
112
129 const std::vector<int> &splitSizes)
131 SplitSystemMatrixAssembler(splitStrategy, splitSizes)
132{
133}
134
212#if BITPIT_ENABLE_MPI==1
220#else
225#endif
250#if BITPIT_ENABLE_MPI==1
258#else
263#endif
286#if BITPIT_ENABLE_MPI==1
294#else
299#endif
300
307 : SplitSystemSolver("", false, false, debug)
308{
309}
310
317SplitSystemSolver::SplitSystemSolver(bool transpose, bool debug)
318 : SplitSystemSolver("", false, transpose, debug)
319{
320}
321
334SplitSystemSolver::SplitSystemSolver(bool flatten, bool transpose, bool debug)
335 : SplitSystemSolver("", flatten, transpose, debug)
336{
337}
338
345SplitSystemSolver::SplitSystemSolver(const std::string &prefix, bool debug)
346 : SplitSystemSolver(prefix, false, debug)
347{
348}
349
356SplitSystemSolver::SplitSystemSolver(const std::string &prefix, bool transpose, bool debug)
357 : SplitSystemSolver(prefix, false, transpose, debug)
358{
359}
360
374SplitSystemSolver::SplitSystemSolver(const std::string &prefix, bool flatten, bool transpose, bool debug)
375 : SystemSolver(prefix, flatten, transpose, debug),
376 m_rhsSplitPermutation(PETSC_NULLPTR), m_solutionSplitPermutation(PETSC_NULLPTR)
377{
378}
379
386{
387 if (m_A == PETSC_NULLPTR) {
388 return 0;
389 }
390
391 std::vector<int> splitBlockSizes = getSplitSizes();
392 int blockSize = std::accumulate(splitBlockSizes.begin(), splitBlockSizes.end(), 0);
393
394 return blockSize;
395}
396
402SystemSplitStrategy SplitSystemSolver::getSplitStrategy() const
403{
404 return m_splitStrategy;
405}
406
413{
414 if (m_A == PETSC_NULLPTR) {
415 return -1;
416 }
417
418 PetscInt nNestedRows;
419 MatNestGetSize(m_A, &nNestedRows, PETSC_NULLPTR);
420
421 return nNestedRows;
422}
423
443std::vector<int> SplitSystemSolver::getSplitSizes() const
444{
445 if (m_A == PETSC_NULLPTR) {
446 return std::vector<int>();
447 }
448
449 int nSplits = getSplitCount();
450 std::vector<int> splitSizes(nSplits, 0);
451 for (int k = 0; k < nSplits; ++k) {
452 int kk = getBlockSplitLinearIndex(k, k, nSplits);
453
454 PetscInt splitBlockSize;
455 MatGetBlockSize(m_splitAs[kk], &splitBlockSize);
456 splitSizes[k] = splitBlockSize;
457 }
458
459 return splitSizes;
460}
461
468{
469 std::vector<int> splitBlockSizes = getSplitSizes();
470 std::vector<int> splitBlockOffsets(splitBlockSizes.size());
471 splitBlockOffsets[0] = 0;
472 std::partial_sum(splitBlockSizes.begin(), std::prev(splitBlockSizes.end()), std::next(splitBlockOffsets.begin()));
473
474 return splitBlockOffsets;
475}
476
481{
482 // Base actions
484
485 // Fill statuses of split KSP
487}
488
498{
499 if (!isAssembled()) {
500 throw std::runtime_error("The options associated with the Krylov solver can only be accessed after assembling the system.");
501 }
502
503 return m_splitKSPOptions[split];
504}
505
516{
517 if (!isAssembled()) {
518 throw std::runtime_error("The options associated with the Krylov solver can only be accessed after assembling the system.");
519 }
520
521 return m_splitKSPOptions[split];
522}
523
532
537{
538 int nSplits = getSplitCount();
539 m_splitKSPOptions.resize(nSplits);
540 for (KSPOptions &splitOptions : m_splitKSPOptions) {
541 resetKSPOptions(&splitOptions);
542 }
543}
544
553
558{
559 m_splitKSPOptions.clear();
560}
561
571{
572 if (!isAssembled()) {
573 throw std::runtime_error("The status of the Krylov solver can only be accessed after assembling the system.");
574 }
575
576 return m_splitKSPStatuses[split];
577}
586
591{
592 int nSplits = getSplitCount();
593 for (int k = static_cast<int>(m_splitKSPStatuses.size()); k < nSplits; ++k) {
594 m_splitKSPStatuses.emplace_back();
595 KSPStatus *splitStatus = &(m_splitKSPStatuses.back());
596 resetKSPStatus(splitStatus);
597 }
598 m_splitKSPStatuses.resize(nSplits);
599}
600
609
614{
615 PC pc;
616 KSPGetPC(m_KSP, &pc);
617
618 PetscInt nSplits;
619 KSP *splitKSPs;
620 PCFieldSplitGetSubKSP(pc, &nSplits, &splitKSPs);
621 for (PetscInt k = 0; k < nSplits; ++k) {
622 KSPStatus *splitKSPStatus = &(m_splitKSPStatuses[k]);
623 fillKSPStatus(splitKSPs[k], splitKSPStatus);
624 }
625}
626
635
640{
641 for (KSPStatus &splitStatus : m_splitKSPStatuses) {
642 resetKSPStatus(&splitStatus);
643 }
644}
645
654
659{
660 m_splitKSPStatuses.clear();
661}
662
663
669void SplitSystemSolver::dumpInfo(std::ostream &systemStream) const
670{
671 SystemSolver::dumpInfo(systemStream);
672
673 if (systemStream.good()) {
674 utils::binary::write(systemStream, m_splitStrategy);
675 }
676}
677
683void SplitSystemSolver::restoreInfo(std::istream &systemStream)
684{
685 SystemSolver::restoreInfo(systemStream);
686
687 utils::binary::read(systemStream, m_splitStrategy);
688}
689
699void SplitSystemSolver::assembly(const SparseMatrix &matrix, SystemSplitStrategy splitStrategy,
700 const std::vector<int> &splitSizes)
701{
702 assembly(matrix, splitStrategy, splitSizes, NaturalSystemMatrixOrdering());
703}
704
713void SplitSystemSolver::assembly(const SparseMatrix &matrix, SystemSplitStrategy splitStrategy,
714 const std::vector<int> &splitSizes,
715 const SystemMatrixOrdering &reordering)
716{
717 // Check if the matrix is assembled
718 if (!matrix.isAssembled()) {
719 throw std::runtime_error("Unable to assembly the system. The matrix is not yet assembled.");
720 }
721
722 // Update matrix
723 SplitSystemSparseMatrixAssembler assembler(&matrix, splitStrategy, splitSizes);
724 assembly<SplitSystemSolver>(static_cast<const Assembler &>(assembler), reordering);
725}
726
735{
737}
738
745void SplitSystemSolver::assembly(const Assembler &assembler, const SystemMatrixOrdering &reordering)
746{
747 assembly<SplitSystemSolver>(assembler, reordering);
748}
749
759{
760 update(getRowCount(), nullptr, elements);
761}
762
773void SplitSystemSolver::update(long nRows, const long *rows, const SparseMatrix &elements)
774{
775 // Check if the element storage is assembled
776 if (!elements.isAssembled()) {
777 throw std::runtime_error("Unable to update the system. The element storage is not yet assembled.");
778 }
779
780 // Update matrix
782 update<SplitSystemSolver>(nRows, rows, static_cast<const Assembler &>(assembler));
783}
784
794{
795 update(getRowCount(), nullptr, assembler);
796}
797
808void SplitSystemSolver::update(long nRows, const long *rows, const Assembler &assembler)
809{
810 update<SplitSystemSolver>(nRows, rows, assembler);
811}
812
819{
820 const PetscInt *rowReordering = PETSC_NULLPTR;
821 if (m_rowReordering) {
822 ISGetIndices(m_rowReordering, &rowReordering);
823 }
824
825 // Matrix information
826 long nRows = assembler.getRowCount();
827
828 // Split information
829 m_splitStrategy = assembler.getSplitStrategy();
830
831 const std::vector<int> &splitBlockSizes = assembler.getSplitSizes();
832 int nSplits = splitBlockSizes.size();
833
834 // Create split matrices
835 m_splitAs.assign(nSplits * nSplits, PETSC_NULLPTR);
836 for (int splitRow = 0; splitRow < nSplits; ++splitRow) {
837 for (int splitCol = 0; splitCol < nSplits; ++splitCol) {
838 // Create matrix
839 //
840 // In LOWER mode, only the lower triangular portion of the splits need to
841 // be created, whereas in DIAGONAL mode only the diagonal splits need to
842 // be created.
843 if (m_splitStrategy == SystemSplitStrategy::SYSTEM_SPLIT_STRATEGY_LOWER) {
844 if (splitCol > splitRow) {
845 continue;
846 }
847 } else if (m_splitStrategy == SystemSplitStrategy::SYSTEM_SPLIT_STRATEGY_DIAGONAL) {
848 if (splitCol != splitRow) {
849 continue;
850 }
851 }
852
853 int splitIndex = getBlockSplitLinearIndex(splitRow, splitCol, nSplits);
854 Mat *splitMatrix = m_splitAs.data() + splitIndex;
855 createMatrix(splitBlockSizes[splitRow], splitBlockSizes[splitCol], splitMatrix);
856
857 MatType splitMatrixType;
858 MatGetType(*splitMatrix, &splitMatrixType);
859
860 // Set matrix sizes
861 long nSplitRowsElements = assembler.getRowCount() * splitBlockSizes[splitRow];
862 long nSplitColsElements = assembler.getColCount() * splitBlockSizes[splitCol];
863
864 long nGlobalSplitRowsElements;
865 long nGlobalSplitColsElements;
866#if BITPIT_ENABLE_MPI == 1
867 nGlobalSplitRowsElements = assembler.getRowGlobalCount() * splitBlockSizes[splitRow];
868 nGlobalSplitColsElements = assembler.getColGlobalCount() * splitBlockSizes[splitCol];
869#else
870 nGlobalSplitRowsElements = nSplitRowsElements;
871 nGlobalSplitColsElements = nSplitColsElements;
872#endif
873
874 MatSetSizes(*splitMatrix, nSplitRowsElements, nSplitColsElements, nGlobalSplitRowsElements, nGlobalSplitColsElements);
875
876 // Allocate matrix storage
877 //
878 // When the internal storage of the system matrix was created without taking into account
879 // block information, preallocation information should be provided for each row of each
880 // block.
881 int rowAllocationExpansion;
882 int colAllocationExpansion;
883 if (strcmp(splitMatrixType, MATSEQAIJ) == 0) {
884 rowAllocationExpansion = splitBlockSizes[splitRow];
885 colAllocationExpansion = splitBlockSizes[splitCol];
886#if BITPIT_ENABLE_MPI == 1
887 } else if (strcmp(splitMatrixType, MATMPIAIJ) == 0) {
888 rowAllocationExpansion = splitBlockSizes[splitRow];
889 colAllocationExpansion = splitBlockSizes[splitCol];
890#endif
891 } else {
892 rowAllocationExpansion = 1;
893 colAllocationExpansion = 1;
894 }
895
896 long nAllocatedRowElements = rowAllocationExpansion * nRows;
897
898 std::vector<int> d_nnz(nAllocatedRowElements, 0);
899 for (long n = 0; n < nRows; ++n) {
900 long matrixRow = n;
901 if (rowReordering) {
902 matrixRow = rowReordering[matrixRow];
903 }
904
905 int nAssemblerRowNZ = assembler.getRowNZCount(n);
906
907 long matrixRowOffset = matrixRow * rowAllocationExpansion;
908 for (int n = 0; n < rowAllocationExpansion; ++n) {
909 d_nnz[matrixRowOffset + n] = colAllocationExpansion * nAssemblerRowNZ;
910 }
911 }
912
913#if BITPIT_ENABLE_MPI == 1
914 std::vector<int> o_nnz(nAllocatedRowElements, 0);
915 if (isPartitioned()) {
916 long nAssemblerCols = assembler.getColCount();
917
918 long assemblerDiagonalBegin = assembler.getColGlobalOffset();
919 long assemblerDiagonalEnd = assemblerDiagonalBegin + nAssemblerCols;
920
921 ConstProxyVector<long> assemblerRowPattern(static_cast<std::size_t>(0), assembler.getMaxRowNZCount());
922 for (long n = 0; n < nRows; ++n) {
923 long matrixRow = n;
924 if (rowReordering) {
925 matrixRow = rowReordering[matrixRow];
926 }
927
928 assembler.getRowPattern(n, &assemblerRowPattern);
929 int nAssemblerRowNZ = assemblerRowPattern.size();
930
931 long matrixRowOffset = matrixRow * rowAllocationExpansion;
932 for (int k = 0; k < nAssemblerRowNZ; ++k) {
933 long id = assemblerRowPattern[k];
934 if (id < assemblerDiagonalBegin || id >= assemblerDiagonalEnd) {
935 for (int n = 0; n < rowAllocationExpansion; ++n) {
936 o_nnz[matrixRowOffset + n] += colAllocationExpansion;
937 d_nnz[matrixRowOffset + n] -= colAllocationExpansion;
938 }
939 }
940 }
941 }
942 }
943#endif
944
945 if (strcmp(splitMatrixType, MATSEQAIJ) == 0) {
946 MatSeqAIJSetPreallocation(*splitMatrix, 0, d_nnz.data());
947 } else if (strcmp(splitMatrixType, MATSEQBAIJ) == 0) {
948 MatSeqBAIJSetPreallocation(*splitMatrix, splitBlockSizes[splitRow], 0, d_nnz.data());
949#if BITPIT_ENABLE_MPI == 1
950 } else if (strcmp(splitMatrixType, MATMPIAIJ) == 0) {
951 MatMPIAIJSetPreallocation(*splitMatrix, 0, d_nnz.data(), 0, o_nnz.data());
952 } else if (strcmp(splitMatrixType, MATMPIBAIJ) == 0) {
953 MatMPIBAIJSetPreallocation(*splitMatrix, splitBlockSizes[splitRow], 0, d_nnz.data(), 0, o_nnz.data());
954#endif
955 } else {
956 throw std::runtime_error("Matrix format not supported.");
957 }
958
959 // Each process will only set values for its own rows
960 MatSetOption(*splitMatrix, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE);
961
962#if PETSC_VERSION_GE(3, 12, 0)
963 // The first assembly will set a superset of the off-process entries
964 // required for all subsequent assemblies. This avoids a rendezvous
965 // step in the MatAssembly functions.
966 MatSetOption(*splitMatrix, MAT_SUBSET_OFF_PROC_ENTRIES, PETSC_TRUE);
967#endif
968 }
969 }
970
971 // Create nest matrix
972 createMatrix(1, 1, nSplits, nSplits, m_splitAs.data(), &m_A);
973
974 // Cleanup
975 if (m_rowReordering) {
976 ISRestoreIndices(m_rowReordering, &rowReordering);
977 }
978
979 // Fill matrix
980 matrixUpdate(assembler.getRowCount(), nullptr, assembler);
981
982 // No new allocations are now allowed
983 //
984 // When updating the matrix it will not be possible to alter the pattern,
985 // it will be possible to change only the values.
986 for (Mat splitMatrix : m_splitAs) {
987 if (!splitMatrix) {
988 continue;
989 }
990
991 MatSetOption(splitMatrix, MAT_NEW_NONZERO_LOCATIONS, PETSC_FALSE);
992 }
993 MatSetOption(m_A, MAT_NEW_NONZERO_LOCATIONS, PETSC_FALSE);
994}
995
1005void SplitSystemSolver::matrixFill(const std::string &filePath)
1006{
1007 // Check if the matrix exists
1008 if (!m_A) {
1009 throw std::runtime_error("Matrix should be created before filling it.");
1010 }
1011
1012 // Fill the matrix
1013 int nSplits = getSplitCount();
1014 for (PetscInt splitRow = 0; splitRow < nSplits; ++splitRow) {
1015 for (PetscInt splitCol = 0; splitCol < nSplits; ++splitCol) {
1016 int splitIndex = getBlockSplitLinearIndex(splitRow, splitCol, nSplits);
1017 std::string splitFilePath = generateSplitPath(filePath, splitRow, splitCol);
1018 fillMatrix(m_splitAs[splitIndex], filePath);
1019 }
1020 }
1021}
1022
1039void SplitSystemSolver::matrixUpdate(long nRows, const long *rows, const Assembler &assembler)
1040{
1041 // Updating the matrix invalidates the KSP
1042 m_KSPDirty = true;
1043
1044 // Get split information
1045 std::vector<int> splitBlockSizes = getSplitSizes();
1046 std::vector<int> splitBlockOffsets = getSplitOffsets();
1047 int nSplits = splitBlockSizes.size();
1048
1049 // Get matrix information
1050 int blockSize = getBlockSize();
1051
1052 std::vector<int> blockSplitRows(blockSize);
1053 for (int i = 0; i < nSplits; ++i) {
1054 for (int k = 0; k < splitBlockSizes[i]; ++k) {
1055 blockSplitRows[splitBlockOffsets[i] + k] = i;
1056 }
1057 }
1058
1059 // Get assembler information
1060 int assemblerBlockSize = assembler.getBlockSize();
1061 if (assemblerBlockSize != blockSize) {
1062 std::string message = "Unable to update the matrix.";
1063 message += " The block size of the assembler is not equal to the block size of the system matrix.";
1064 throw std::runtime_error(message);
1065 }
1066
1067 const long assemblerMaxRowNZ = assembler.getMaxRowNZCount();
1068
1069 PetscInt colGlobalBegin;
1070 PetscInt colGlobalEnd;
1071 MatGetOwnershipRangeColumn(m_A, &colGlobalBegin, &colGlobalEnd);
1072 colGlobalBegin /= assemblerBlockSize;
1073 colGlobalEnd /= assemblerBlockSize;
1074
1075 PetscInt rowGlobalOffset;
1076 MatGetOwnershipRange(m_A, &rowGlobalOffset, PETSC_NULLPTR);
1077 rowGlobalOffset /= assemblerBlockSize;
1078
1079 // Initialize reordering
1080 const PetscInt *rowReordering = PETSC_NULLPTR;
1081 if (m_rowReordering) {
1082 ISGetIndices(m_rowReordering, &rowReordering);
1083 }
1084
1085 const PetscInt *colReordering = PETSC_NULLPTR;
1086 if (m_colReordering) {
1087 ISGetIndices(m_colReordering, &colReordering);
1088 }
1089
1090 // Generate the split indexes
1091 std::vector<std::vector<std::size_t>> scatterIndexes(nSplits);
1092 for (int split = 0; split < nSplits; ++split) {
1093 generateSplitIndexes(split, assemblerMaxRowNZ, scatterIndexes.data() + split);
1094 }
1095
1096 // Get the options for assembling the matrix
1097 SystemMatrixAssembler::AssemblyOptions assemblyOptions = assembler.getOptions();
1098
1099#if PETSC_VERSION_GE(3, 12, 0)
1100 // Check if it is possible to speedup insertion of values
1101 //
1102 // The option MAT_SORTED_FULL means that each process provides exactly its
1103 // local rows; all column indices for a given row are passed in a single call
1104 // to MatSetValues(), preallocation is perfect, row oriented, INSERT_VALUES
1105 // is used. If this options is set to PETSC_TRUE, the function MatSetValues
1106 // will be faster.
1107 //
1108 // This options needs at least PETSc 3.12.
1109 PetscBool matrixSortedFull = (assemblyOptions.full && assemblyOptions.sorted) ? PETSC_TRUE : PETSC_FALSE;
1110 MatSetOption(m_A, MAT_SORTED_FULL, matrixSortedFull);
1111#endif
1112
1113 // Insert matrix values
1114 //
1115 // If the sizes of PETSc data types match the sizes of data types expected by bitpit
1116 // and no column reordering is needed, a direct update of the pattern can be performed,
1117 // otherwise an intermediate pattern storage is required.
1118 bool patternDirectUpdate = !colReordering && (sizeof(long) == sizeof(PetscInt));
1119
1120 ConstProxyVector<long> rowPattern;
1121 std::vector<PetscInt> petscRowPatternStorage;
1122 const PetscInt *petscRowPattern;
1123 if (!patternDirectUpdate) {
1124 rowPattern.set(ConstProxyVector<long>::INTERNAL_STORAGE, 0, assemblerMaxRowNZ);
1125 petscRowPatternStorage.resize(assemblerMaxRowNZ);
1126 petscRowPattern = petscRowPatternStorage.data();
1127 }
1128
1129 ConstProxyVector<double> rowValues;
1130 std::vector<std::vector<PetscScalar>> petscSplitRowValues(m_splitAs.size());
1131 for (PetscInt splitRow = 0; splitRow < nSplits; ++splitRow) {
1132 for (PetscInt splitCol = 0; splitCol < nSplits; ++splitCol) {
1133 int splitIndex = getBlockSplitLinearIndex(splitRow, splitCol, nSplits);
1134 if (!m_splitAs[splitIndex]) {
1135 continue;
1136 }
1137
1138 int splitMatriBlockSize = splitBlockSizes[splitRow] * splitBlockSizes[splitCol];
1139 petscSplitRowValues[splitIndex].resize(splitMatriBlockSize * assemblerMaxRowNZ);
1140 }
1141 }
1142
1143 for (long n = 0; n < nRows; ++n) {
1144 // Get row information
1145 long row;
1146 if (rows) {
1147 row = rows[n];
1148 } else {
1149 row = n;
1150 }
1151
1152 if (rowReordering) {
1153 row = rowReordering[row];
1154 }
1155
1156 const PetscInt globalRow = rowGlobalOffset + row;
1157
1158 // Get row data
1159 assembler.getRowData(n, &rowPattern, &rowValues);
1160 if (rowValues.size() == 0) {
1161 continue;
1162 }
1163
1164 // Get pattern in PETSc format
1165 const std::size_t rowPatternSize = rowPattern.size();
1166 if (patternDirectUpdate) {
1167 petscRowPattern = reinterpret_cast<const PetscInt *>(rowPattern.data());
1168 } else {
1169 for (std::size_t k = 0; k < rowPatternSize; ++k) {
1170 long globalCol = rowPattern[k];
1171 if (colReordering) {
1172 if (globalCol >= colGlobalBegin && globalCol < colGlobalEnd) {
1173 long col = globalCol - colGlobalBegin;
1174 col = colReordering[col];
1175 globalCol = colGlobalBegin + col;
1176 }
1177 }
1178
1179 petscRowPatternStorage[k] = globalCol;
1180 }
1181 }
1182
1183 // Get values in PETSc format
1184 //
1185 // The assembler returns the block row values of the matrix in a row-major order.
1186 // We need to scatter these values to obtain the values of the block rows of the
1187 // split matrices.
1188 for (int blockRow = 0; blockRow < blockSize; ++blockRow) {
1189 int splitRow = blockSplitRows[blockRow];
1190 int splitBlockRow = blockRow - splitBlockOffsets[splitRow];
1191
1192 std::size_t nRowElements = blockSize * rowPatternSize;
1193 std::size_t blockRowValuesOffset = nRowElements * blockRow;
1194 const double *blockRowValues = rowValues.data() + blockRowValuesOffset;
1195
1196 for (int splitCol = 0; splitCol < nSplits; ++splitCol) {
1197 int splitIndex = getBlockSplitLinearIndex(splitRow, splitCol, nSplits);
1198 if (!m_splitAs[splitIndex]) {
1199 continue;
1200 }
1201
1202 const std::vector<std::size_t> &splitScatterIndexes = scatterIndexes[splitCol];
1203 std::size_t nSplitRowElements = splitBlockSizes[splitCol] * rowPatternSize;
1204 std::size_t splitBlockRowValuesOffset = nSplitRowElements * splitBlockRow;
1205 double *splitBlockRowValues = petscSplitRowValues[splitIndex].data() + splitBlockRowValuesOffset;
1206 for (std::size_t k = 0; k < nSplitRowElements; ++k) {
1207 splitBlockRowValues[k] = blockRowValues[splitScatterIndexes[k]];
1208 }
1209 }
1210 }
1211
1212 // Insert values
1213 for (std::size_t k = 0; k < m_splitAs.size(); ++k) {
1214 if (!m_splitAs[k]) {
1215 continue;
1216 }
1217
1218 MatSetValuesBlocked(m_splitAs[k], 1, &globalRow, rowPatternSize, petscRowPattern,
1219 petscSplitRowValues[k].data(), INSERT_VALUES);
1220 }
1221 }
1222
1223 // Let PETSc assembly the matrix after the update
1224 for (Mat splitMatrix : m_splitAs) {
1225 if (!splitMatrix) {
1226 continue;
1227 }
1228
1229 MatAssemblyBegin(splitMatrix, MAT_FINAL_ASSEMBLY);
1230 }
1231
1232 for (Mat splitMatrix : m_splitAs) {
1233 if (!splitMatrix) {
1234 continue;
1235 }
1236
1237 MatAssemblyEnd(splitMatrix, MAT_FINAL_ASSEMBLY);
1238 }
1239
1240 MatAssemblyBegin(m_A, MAT_FINAL_ASSEMBLY);
1241 MatAssemblyEnd(m_A, MAT_FINAL_ASSEMBLY);
1242
1243 // Cleanup
1244 if (rowReordering) {
1245 ISRestoreIndices(m_rowReordering, &rowReordering);
1246 }
1247
1248 if (colReordering) {
1249 ISRestoreIndices(m_colReordering, &colReordering);
1250 }
1251}
1252
1260void SplitSystemSolver::matrixDump(std::ostream &systemStream, const std::string &directory,
1261 const std::string &prefix) const
1262{
1263 // Dump split matrices
1264 int nSplits = getSplitCount();
1265 if (systemStream.good()) {
1266 utils::binary::write(systemStream, nSplits);
1267 }
1268
1269 for (PetscInt splitRow = 0; splitRow < nSplits; ++splitRow) {
1270 for (PetscInt splitCol = 0; splitCol < nSplits; ++splitCol) {
1271 int splitIndex = getBlockSplitLinearIndex(splitRow, splitCol, nSplits);
1272 std::string splitMatrixName = generateSplitPath(prefix + "A", splitRow, splitCol);
1273 dumpMatrix(m_splitAs[splitIndex], directory, splitMatrixName);
1274 }
1275 }
1276
1277 // Dump information for creating the main matrix
1278 if (systemStream.good()) {
1279 PetscInt rowBlockSize;
1280 PetscInt colBlockSize;
1281 MatGetBlockSizes(m_A, &rowBlockSize, &colBlockSize);
1282 utils::binary::write(systemStream, static_cast<int>(rowBlockSize));
1283 utils::binary::write(systemStream, static_cast<int>(colBlockSize));
1284 }
1285}
1286
1294#if BITPIT_ENABLE_MPI==1
1300void SplitSystemSolver::matrixRestore(std::istream &systemStream, const std::string &directory,
1301 const std::string &prefix, bool redistribute)
1302#else
1303void SplitSystemSolver::matrixRestore(std::istream &systemStream, const std::string &directory,
1304 const std::string &prefix)
1305#endif
1306{
1307 // Restore split matrices
1308 int nSplits;
1309 utils::binary::read(systemStream, nSplits);
1310
1311 m_splitAs.assign(nSplits * nSplits, PETSC_NULLPTR);
1312 for (PetscInt splitRow = 0; splitRow < nSplits; ++splitRow) {
1313 for (PetscInt splitCol = 0; splitCol < nSplits; ++splitCol) {
1314 int splitIndex = getBlockSplitLinearIndex(splitRow, splitCol, nSplits);
1315 std::string splitMatrixName = generateSplitPath(prefix + "A", splitRow, splitCol);
1316#if BITPIT_ENABLE_MPI==1
1317 restoreMatrix(directory, splitMatrixName, redistribute, m_splitAs.data() + splitIndex);
1318#else
1319 restoreMatrix(directory, splitMatrixName, m_splitAs.data() + splitIndex);
1320#endif
1321 }
1322 }
1323
1324 // Restore main matrix
1325 int rowBlockSize;
1326 utils::binary::read(systemStream, rowBlockSize);
1327 int colBlockSize;
1328 utils::binary::read(systemStream, colBlockSize);
1329
1330 createMatrix(rowBlockSize, colBlockSize, nSplits, nSplits, m_splitAs.data(), &m_A);
1331}
1332
1337{
1338 for (Mat &splitMatrix : m_splitAs) {
1339 destroyMatrix(&splitMatrix);
1340 }
1341
1342 destroyMatrix(&m_A);
1343}
1344
1354{
1355 // Create vectors
1357
1358 // Create the split permutations
1360}
1361
1368{
1370
1371 reorderVector(m_rhs, m_solutionSplitPermutation, !invert);
1372 reorderVector(m_solution, m_rhsSplitPermutation, !invert);
1373}
1374
1382void SplitSystemSolver::vectorsRestore(std::istream &systemStream, const std::string &directory,
1383 const std::string &prefix)
1384{
1385 // Create vectors
1386 SystemSolver::vectorsRestore(systemStream, directory, prefix);
1387
1388 // Create the split permutations
1390}
1391
1399{
1400 int blockSize = getBlockSize();
1401
1402#if BITPIT_ENABLE_MPI == 1
1403 PetscInt rhsGlobalBegin;
1404 PetscInt rhsGlobalEnd;
1405 VecGetOwnershipRange(m_rhs, &rhsGlobalBegin, &rhsGlobalEnd);
1406 rhsGlobalBegin /= blockSize;
1407 rhsGlobalEnd /= blockSize;
1408 PetscInt rhsSize = rhsGlobalEnd - rhsGlobalBegin;
1409
1410 PetscInt solutionGlobalBegin;
1411 PetscInt solutionGlobalEnd;
1412 VecGetOwnershipRange(m_solution, &solutionGlobalBegin, &solutionGlobalEnd);
1413 solutionGlobalBegin /= blockSize;
1414 solutionGlobalEnd /= blockSize;
1415 PetscInt solutionSize = solutionGlobalEnd - solutionGlobalBegin;
1416
1417 generateSplitPermutation(rhsSize, static_cast<long>(rhsGlobalBegin), &m_rhsSplitPermutation);
1418 generateSplitPermutation(solutionSize, static_cast<long>(solutionGlobalBegin), &m_solutionSplitPermutation);
1419#else
1420 PetscInt rhsSize;
1421 VecGetSize(m_rhs, &rhsSize);
1422 rhsSize /= blockSize;
1423
1424 PetscInt solutionSize;
1425 VecGetSize(m_solution, &solutionSize);
1426 solutionSize /= blockSize;
1427
1428 generateSplitPermutation(rhsSize, &m_rhsSplitPermutation);
1429 generateSplitPermutation(solutionSize, &m_solutionSplitPermutation);
1430#endif
1431}
1432
1437{
1438 // Destroy split permutations
1439 if (m_rhsSplitPermutation) {
1440 ISDestroy(&m_rhsSplitPermutation);
1441 }
1442
1443 if (m_solutionSplitPermutation) {
1444 ISDestroy(&m_solutionSplitPermutation);
1445 }
1446
1447 // Destroy RHS and solution vectors
1449}
1450
1458void SplitSystemSolver::exportMatrix(const std::string &filePath, FileFormat fileFormat) const
1459{
1460 // Export main matrix
1461 exportMatrix(m_A, filePath, fileFormat);
1462
1463 // Export split matrices
1464 int nSplits = getSplitCount();
1465 for (int splitRow = 0; splitRow < nSplits; ++splitRow) {
1466 for (int splitCol = 0; splitCol < nSplits; ++splitCol) {
1467 int splitIndex = getBlockSplitLinearIndex(splitRow, splitCol, nSplits);
1468 std::string splitFilePath = generateSplitPath(filePath, splitRow, splitCol);
1469 exportMatrix(m_splitAs[splitIndex], splitFilePath, fileFormat);
1470 }
1471 }
1472}
1473
1478{
1479 PC pc;
1480 KSPGetPC(m_KSP, &pc);
1481
1482 // Setup main preconditioner
1483 PCSetType(pc, PCFIELDSPLIT);
1484 if (m_splitStrategy == SystemSplitStrategy::SYSTEM_SPLIT_STRATEGY_LOWER) {
1485 PCFieldSplitSetType(pc, PC_COMPOSITE_MULTIPLICATIVE);
1486 } else if (m_splitStrategy == SystemSplitStrategy::SYSTEM_SPLIT_STRATEGY_DIAGONAL) {
1487 PCFieldSplitSetType(pc, PC_COMPOSITE_ADDITIVE);
1488 } else {
1489 PCFieldSplitSetType(pc, PC_COMPOSITE_ADDITIVE);
1490 }
1491
1492 int nSplits = getSplitCount();
1493 std::vector<IS> splitIndexSets(nSplits);
1494 MatNestGetISs(m_A, splitIndexSets.data(), PETSC_NULLPTR);
1495 for (int k = 0; k < nSplits; ++k) {
1496 PCFieldSplitSetIS(pc, PETSC_NULLPTR, splitIndexSets[k]);
1497 }
1498
1499 PCSetUp(pc);
1500
1501 // Set preconditioners of the split
1503}
1504
1509{
1510 PC pc;
1511 KSPGetPC(m_KSP, &pc);
1512
1513 PetscInt nSplits;
1514 KSP *splitKSPs;
1515 PCFieldSplitGetSubKSP(pc, &nSplits, &splitKSPs);
1516
1517 for (PetscInt k = 0; k < nSplits; ++k) {
1518 PC splitPc;
1519 KSPGetPC(splitKSPs[k], &splitPc);
1520 const KSPOptions &splitOptions = getSplitKSPOptions(k);
1521
1522 setupPreconditioner(splitPc, splitOptions);
1523 }
1524}
1525
1530{
1531 // Setup main Krylov subspace method
1532 //
1533 // In LOWER and DIAGONAL modes, the solution of the system is performed by the Krylov subspace
1534 // methods of the single splits. Thus, in these modes, the main Krylov subspace method only
1535 // needs to apply the preconditioner.
1536 if (m_splitStrategy == SystemSplitStrategy::SYSTEM_SPLIT_STRATEGY_LOWER) {
1537 KSPSetType(m_KSP, KSPPREONLY);
1538 } else if (m_splitStrategy == SystemSplitStrategy::SYSTEM_SPLIT_STRATEGY_DIAGONAL) {
1539 KSPSetType(m_KSP, KSPPREONLY);
1540 } else {
1541 setupKrylov(m_KSP, getKSPOptions());
1542 }
1543
1544 // Setup Krylov subspace methods of the splits
1546}
1547
1552{
1553 PC pc;
1554 KSPGetPC(m_KSP, &pc);
1555
1556 PetscInt nSplits;
1557 KSP *splitKSPs;
1558 PCFieldSplitGetSubKSP(pc, &nSplits, &splitKSPs);
1559
1560 for (PetscInt k = 0; k < nSplits; ++k) {
1561 // In FULL mode, the solution of the system is performed by the Krylov subspace
1562 // methods of the single splits. Thus, in these modes, the Krylov subspaces of
1563 // the single blocks only need to apply the preconditioner.
1564 const KSPOptions &splitOptions = getSplitKSPOptions(k);
1565 if (m_splitStrategy == SystemSplitStrategy::SYSTEM_SPLIT_STRATEGY_LOWER) {
1566 setupKrylov(splitKSPs[k], splitOptions);
1567 } else if (m_splitStrategy == SystemSplitStrategy::SYSTEM_SPLIT_STRATEGY_DIAGONAL) {
1568 setupKrylov(splitKSPs[k], splitOptions);
1569 } else {
1570 KSPSetType(splitKSPs[k], KSPPREONLY);
1571 }
1572 }
1573}
1574
1579{
1580 // Execute base actions
1582
1583 // Enable convergence monitor
1584 if (m_convergenceMonitorEnabled) {
1585 int nSplits = getSplitCount();
1586 for (int split = 0; split < nSplits; ++split) {
1587 std::string prefix = "fieldsplit_" + std::to_string(split) + "_";
1588
1589#if (PETSC_VERSION_MAJOR >= 3 && PETSC_VERSION_MINOR >= 7)
1590 PetscOptionsSetValue(PETSC_NULLPTR, ("-" + prefix + "ksp_monitor_true_residual").c_str(), "");
1591 PetscOptionsSetValue(PETSC_NULLPTR, ("-" + prefix + "ksp_converged_reason").c_str(), "");
1592 PetscOptionsSetValue(PETSC_NULLPTR, ("-" + prefix + "ksp_monitor_singular_value").c_str(), "");
1593#else
1594 PetscOptionsSetValue(("-" + prefix + "ksp_monitor_true_residual").c_str(), "");
1595 PetscOptionsSetValue(("-" + prefix + "ksp_converged_reason").c_str(), "");
1596 PetscOptionsSetValue(("-" + prefix + "ksp_monitor_singular_value").c_str(), "");
1597#endif
1598 }
1599 }
1600}
1601
1610#if BITPIT_ENABLE_MPI == 1
1614#endif
1619#if BITPIT_ENABLE_MPI == 1
1620void SplitSystemSolver::generateSplitPermutation(long nItems, long itemsGlobalOffset, IS *permutation) const
1621#else
1622void SplitSystemSolver::generateSplitPermutation(long nItems, IS *permutation) const
1623#endif
1624{
1625 std::vector<std::size_t> workspace;
1626
1627 // Get split information
1628 std::vector<int> splitBlockSizes = getSplitSizes();
1629 std::vector<int> splitBlockOffsets = getSplitOffsets();
1630 int nSplits = splitBlockSizes.size();
1631
1632 // Matrix information
1633 int blockSize = getBlockSize();
1634 std::size_t nElements = blockSize * nItems;
1635
1636 // Create permutations
1637 PetscInt elementsGlobalOffset;
1638#if BITPIT_ENABLE_MPI == 1
1639 elementsGlobalOffset = static_cast<PetscInt>(blockSize * itemsGlobalOffset);
1640#else
1641 elementsGlobalOffset = 0;
1642#endif
1643
1644 PetscInt *permutationStorage;
1645 PetscMalloc(nElements * sizeof(PetscInt), &permutationStorage);
1646
1647 PetscInt *splitPermutationStorage = permutationStorage;
1648 for (int split = 0; split < nSplits; ++split) {
1649 int splitBlockSize = splitBlockSizes[split];
1650 std::size_t nSplitElements = splitBlockSize * nItems;
1651
1652 workspace.resize(nSplitElements);
1653 generateSplitIndexes(split, nItems, &workspace);
1654 for (std::size_t k = 0; k < nSplitElements; ++k) {
1655 splitPermutationStorage[k] = elementsGlobalOffset + static_cast<PetscInt>(workspace[k]);
1656
1657 }
1658 splitPermutationStorage += nSplitElements;
1659 }
1660
1661#if BITPIT_ENABLE_MPI == 1
1662 ISCreateGeneral(getCommunicator(), nItems, permutationStorage, PETSC_OWN_POINTER, permutation);
1663#else
1664 ISCreateGeneral(PETSC_COMM_SELF, nItems, permutationStorage, PETSC_OWN_POINTER, permutation);
1665#endif
1666 ISSetPermutation(*permutation);
1667}
1668
1680void SplitSystemSolver::generateSplitIndexes(int split, long nItems, std::vector<std::size_t> *indexes) const
1681{
1682 // Get split information
1683 int splitBlockSize = getSplitSizes()[split];
1684 int splitBlockOffset = getSplitOffsets()[split];
1685
1686 // Get matrix information
1687 int blockSize = getBlockSize();
1688
1689 // Initialize index storage
1690 indexes->resize(nItems * splitBlockSize);
1691
1692 // Evaluate the indexes
1693 for (long k = 0; k < nItems; ++k) {
1694 for (int n = 0; n < splitBlockSize; ++n) {
1695 (*indexes)[splitBlockSize * k + n] = blockSize * k + splitBlockOffset + n;
1696 }
1697 }
1698}
1699
1710std::string SplitSystemSolver::generateSplitPath(const std::string &path, int i) const
1711{
1712 return generateSplitPath(path, std::to_string(i));
1713}
1714
1726std::string SplitSystemSolver::generateSplitPath(const std::string &path, int i, int j) const
1727{
1728 return generateSplitPath(path, std::to_string(i) + std::to_string(j));
1729}
1730
1741std::string SplitSystemSolver::generateSplitPath(const std::string &path, const std::string &index) const
1742{
1743 std::string splitLabel = "_split_" + index;
1744
1745 std::size_t dotPosition = path.find_last_of(".");
1746
1747 std::string splitPath;
1748 std::string extension;
1749 std::string basePath;
1750 if (dotPosition != std::string::npos) {
1751 extension = path.substr(dotPosition + 1);
1752 basePath = path.substr(0, dotPosition);
1753 splitPath = basePath + splitLabel + "." + extension;
1754 } else {
1755 splitPath = path + splitLabel;
1756 }
1757
1758 return splitPath;
1759}
1760
1770int SplitSystemSolver::getBlockSplitLinearIndex(int i, int j) const
1771{
1772 int nSplits = getSplitCount();
1773
1774 return getBlockSplitLinearIndex(i, j, nSplits);
1775}
1776
1785int SplitSystemSolver::getBlockSplitLinearIndex(int i, int j, int nSplits) const
1786{
1787 return linearalgebra::linearIndexRowMajor(i, j, nSplits, nSplits);
1788}
1789
1790}
The NaturalSystemMatrixOrdering class defines allows to use a matrix natural ordering.
Metafunction for generating a list of elements that can be either stored in an external vectror or,...
std::size_t size() const
__PXV_POINTER__ data() noexcept
void set(__PXV_POINTER__ data, std::size_t size)
The SplitSystemMatrixAssembler class is the base class for defining assemblers for split system solve...
SystemSplitStrategy getSplitStrategy() const
SplitSystemMatrixAssembler(SystemSplitStrategy splitStrategy, const std::vector< int > &splitSizes)
const std::vector< int > & getSplitSizes() const
The SplitSystemSolver class allows to solve a split linear system.
std::vector< int > getSplitSizes() const
SystemSplitStrategy getSplitStrategy() const
void matrixUpdate(long nRows, const long *rows, const Assembler &assembler)
std::vector< int > getSplitOffsets() const
void generateSplitIndexes(int split, long nItems, std::vector< std::size_t > *indexes) const
void vectorsReorder(bool invert) override
const KSPStatus & getSplitKSPStatus(int split) const
void matrixAssembly(const Assembler &assembler)
virtual void resetKSPOptions(KSPOptions *options) const
void matrixRestore(std::istream &systemStream, const std::string &directory, const std::string &prefix, bool redistribute) override
KSPOptions & getSplitKSPOptions(int split)
void vectorsRestore(std::istream &systemStream, const std::string &directory, const std::string &prefix) override
void assembly(const SparseMatrix &matrix, SystemSplitStrategy splitStrategy, const std::vector< int > &splitSizes)
void restoreInfo(std::istream &systemStream) override
void matrixFill(const std::string &filePath) override
std::string generateSplitPath(const std::string &path, int i) const
void generateSplitPermutation(long nItems, long itemGlobalOffset, IS *splitReordering) const
void update(const SparseMatrix &elements)
void matrixDump(std::ostream &systemStream, const std::string &directory, const std::string &prefix) const override
void exportMatrix(const std::string &filePath, FileFormat exportFormat=FILE_BINARY) const override
void dumpInfo(std::ostream &systemStream) const override
The SplitSystemSparseMatrixAssembler class allows to assembly a split system solver from a sparse mat...
SplitSystemSparseMatrixAssembler(const SparseMatrix *matrix, SystemSplitStrategy splitStrategy, const std::vector< int > &splitSizes)
The SystemMatrixAssembler class provides an interface for defining system matrix assemblers.
The SystemMatrixOrdering class provides an interface for defining classes that allows to reorder the ...
The SystemSolver class provides methods for building and solving large linear systems.
void dumpMatrix(Mat matrix, const std::string &directory, const std::string &name) const
void fillMatrix(Mat matrix, const std::string &filePath) const
void reorderVector(Vec vector, IS permutations, bool invert) const
virtual void restoreInfo(std::istream &stream)
void destroyMatrix(Mat *matrix) const
virtual void dumpInfo(std::ostream &stream) const
void createMatrix(int rowBlockSize, int colBlockSize, Mat *matrix) const
virtual void vectorsReorder(bool invert)
virtual void vectorsRestore(std::istream &systemStream, const std::string &directory, const std::string &prefix)
The SystemSparseMatrixAssembler class defines an assembler for building the system matrix form a spar...
void write(std::ostream &stream, const std::vector< bool > &container)
void read(std::istream &stream, std::vector< bool > &container)
int linearIndexRowMajor(int row, int col, int nRows, int nCols)
bool full
Controls if the assembler is providing all the non-zero values of a row.
--- layout: doxygen_footer ---