Loading...
Searching...
No Matches
system_solvers_large.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 "system_solvers_large.hpp"
26
27#include "bitpit_common.hpp"
28#include "bitpit_containers.hpp"
29#include "bitpit_IO.hpp"
30#include "bitpit_operators.hpp"
31
32#include "petscksp.h"
33#include "petscmat.h"
34#include "petscsystypes.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
79{
80 return row;
81}
82
94{
95 return col;
96}
97
123
124#if BITPIT_ENABLE_MPI==1
131{
132 return m_matrix->isPartitioned();
133}
134
141{
142 return m_matrix->getCommunicator();
143}
144#endif
145
152{
153 AssemblyOptions options;
154 options.full = true;
155 options.sorted = false;
156
157 return options;
158}
159
166{
167 return m_transpose;
168}
169
181void SystemSolver::setTranspose(bool transpose)
182{
183 // Early return if the transpose flag is already correctly set.
184 if (m_transpose == transpose) {
185 return;
186 }
187
188 // Clear the workspace
190
191 // Re-create the vectors
192 if (isAssembled()) {
193 VecDestroy(&m_rhs);
194 VecDestroy(&m_solution);
195
197 }
198
199 // Set the flag
200 m_transpose = transpose;
201}
202
209{
210 return m_matrix->getBlockSize();
211}
212
223{
224 return m_matrix->getRowCount();
225}
226
237{
238 return m_matrix->getColCount();
239}
240
250{
251 long nRowElements = getBlockSize() * getRowCount();
252
253 return nRowElements;
254}
255
265{
266 long nColElements = getBlockSize() * getColCount();
267
268 return nColElements;
269}
270
271#if BITPIT_ENABLE_MPI==1
282{
283 return m_matrix->getRowGlobalCount();
284}
285
297 return m_matrix->getColGlobalCount();
298}
299
309{
310 long nElements = getBlockSize() * getRowGlobalCount();
311
312 return nElements;
313}
314
324{
325 long nElements = getBlockSize() * getColGlobalCount();
326
327 return nElements;
328}
329
339
346{
347 return m_matrix->getColGlobalOffset();
348}
349
359{
360 long offset = getBlockSize() * getRowGlobalOffset();
361
362 return offset;
363}
364
374{
375 long offset = getBlockSize() * getColGlobalOffset();
376
377 return offset;
378}
379#endif
380
388{
389 return m_matrix->getRowNZCount(rowIndex);
390}
391
398{
399 return m_matrix->getMaxRowNZCount();
400}
401
413{
414 m_matrix->getRowPattern(rowIndex, pattern);
415}
416
432{
433 m_matrix->getRowValues(rowIndex, values);
434}
435
456 m_matrix->getRowPattern(rowIndex, pattern);
457 m_matrix->getRowValues(rowIndex, values);
458}
459
472PetscErrorCode PetscManager::displayLogView()
473{
474 return PetscLogView(PETSC_VIEWER_STDOUT_WORLD);
475}
476
482void PetscManager::addInitOption(const std::string &option)
483{
484 if (!areOptionsEditable()) {
485 throw std::runtime_error("Initialization options can be set only before initializing the solver.");
486 }
487
488 m_options.push_back(option);
489}
490
505void PetscManager::addInitOptions(int argc, char **argv)
506{
507 if (!areOptionsEditable()) {
508 throw std::runtime_error("Initialization options can be set only before initializing the solver.");
509 }
510
511 for (int i = 1; i < argc; ++i) {
512 m_options.push_back(argv[i]);
513 }
514}
515
521void PetscManager::addInitOptions(const std::vector<std::string> &options)
522{
523 if (!areOptionsEditable()) {
524 throw std::runtime_error("Initialization options can be set only before initializing the solver.");
525 }
526
527 for (const std::string &option : options) {
528 m_options.push_back(option);
529 }
530}
531
536{
537 m_options.clear();
538}
539
545void PetscManager::enableLogView()
546{
547 if (m_logViewEnabled) {
548 return;
549 }
550
551#if (PETSC_VERSION_MAJOR >= 3 && PETSC_VERSION_MINOR >= 7)
552 PetscLogDefaultBegin();
553#else
554 PetscLogBegin();
555#endif
556 PetscRegisterFinalize(&(PetscManager::displayLogView));
557 m_logViewEnabled = true;
558}
559
564 : m_externalMPIInitialization(true), m_externalPETScInitialization(true),
565 m_options(std::vector<std::string>(1, "bitpit")), m_logViewEnabled(false)
566{
567}
568
576
583{
584 PetscBool isPetscInitialized;
585 PetscInitialized(&isPetscInitialized);
586
587 return (!isPetscInitialized);
588}
589
595#if BITPIT_ENABLE_MPI==1
600#endif
610{
611 // Early return if PETSc is already initialized
612 PetscBool isPetscInitialized;
613 PetscInitialized(&isPetscInitialized);
614 if (isPetscInitialized) {
615 return false;
616 }
617
618#if BITPIT_ENABLE_MPI==1
619 // Early return if MPI is already finalized
620 int isMPIFinalized;
621 MPI_Finalized(&isMPIFinalized);
622 if (isMPIFinalized) {
623 throw std::runtime_error("PETSc finalization cannot be called after MPI finaliation.");
624 }
625#endif
626
627 // Generate command line arguments
628 //
629 // The first argument is the executable name and it is set to a
630 // dummy value.
631 std::string help = "None";
632 std::string programName = "bitpit_petsc_manager";
633
634 int argc = 1 + m_options.size();
635 char **argv = new char*[argc + 1];
636 argv[0] = strdup(programName.data());
637 for (std::size_t i = 0; i < m_options.size(); i++) {
638 argv[1 + i] = strdup(m_options[i].data());
639 }
640 argv[argc] = nullptr;
641
642#if BITPIT_ENABLE_MPI==1
643 // Initialize MPI
644 int isMPIInitialized;
645 MPI_Initialized(&isMPIInitialized);
646 if (!isMPIInitialized) {
647 MPI_Init(&argc, &argv);
648 m_externalMPIInitialization = false;
649 }
650#endif
651
652 // Initialize PETSc
653 PetscInitialize(&argc, &argv, 0, help.data());
654 m_externalPETScInitialization = false;
655
656 // Enable log view
657 if (debug) {
658 enableLogView();
659 }
660
661 // Clean-up command line arguments
662 for (int i = 0; i < argc; ++i) {
663 free(argv[i]);
664 }
665
666 delete[] argv;
667
668 // Initialization completed
669 return true;
670}
671
677#if BITPIT_ENABLE_MPI==1
683#else
688#endif
694{
695#if BITPIT_ENABLE_MPI==1
696 // Early return if MPI is already finalized
697 int isMPIFinalized;
698 MPI_Finalized(&isMPIFinalized);
699 if (isMPIFinalized) {
700 return false;
701 }
702#endif
703
704 // Early return if PETSc was initialized externally
705 if (m_externalPETScInitialization) {
706 return false;
707 }
708
709 // Finalize PETSc
710 PetscBool isPetscFinalized;
711 PetscFinalized(&isPetscFinalized);
712 if (!isPetscFinalized) {
713 PetscFinalize();
714 }
715
716#if BITPIT_ENABLE_MPI==1
717 // Finalize MPI
718 if (!m_externalMPIInitialization) {
719 MPI_Finalize();
720 }
721#endif
722
723 return (!isPetscFinalized);
724}
725
745#if BITPIT_ENABLE_MPI==1
754#else
759#endif
760
761PetscManager SystemSolver::m_petscManager = PetscManager();
762
763int SystemSolver::m_nInstances = 0;
764
771 : SystemSolver("", false, false, debug)
772{
773}
774
781SystemSolver::SystemSolver(bool transpose, bool debug)
782 : SystemSolver("", false, transpose, debug)
783{
784}
785
786
799SystemSolver::SystemSolver(bool flatten, bool transpose, bool debug)
800 : SystemSolver("", flatten, transpose, debug)
801{
802}
803
810SystemSolver::SystemSolver(const std::string &prefix, bool debug)
811 : SystemSolver(prefix, false, debug)
812{
813}
814
821SystemSolver::SystemSolver(const std::string &prefix, bool transpose, bool debug)
822 : SystemSolver(prefix, false, transpose, debug)
823{
824}
825
839SystemSolver::SystemSolver(const std::string &prefix, bool flatten, bool transpose, bool debug)
840 : m_flatten(flatten), m_transpose(transpose),
841 m_A(PETSC_NULLPTR), m_rhs(PETSC_NULLPTR), m_solution(PETSC_NULLPTR),
842 m_rowReordering(PETSC_NULLPTR), m_colReordering(PETSC_NULLPTR),
843 m_convergenceMonitorEnabled(debug),
844 m_KSP(PETSC_NULLPTR), m_KSPDirty(true),
845 m_prefix(prefix), m_assembled(false),
846#if BITPIT_ENABLE_MPI==1
847 m_communicator(MPI_COMM_SELF), m_partitioned(false),
848#endif
849 m_forceConsistency(false)
850{
851 // Initialize PETSc
852 if (m_nInstances == 0) {
853 m_petscManager.initialize(debug);
854 }
855
856 // Increase the number of instances
857 ++m_nInstances;
858}
859
864{
865 // Clear the solver
866 clear();
867
868 // Decrease the number of instances
869 --m_nInstances;
870}
871
876{
878
881
884
886
887#if BITPIT_ENABLE_MPI==1
888 freeCommunicator();
889#endif
890
891 m_assembled = false;
892}
893
900{
901 // Early return if the KSP has not been created
902 if (!m_KSP) {
903 return;
904 }
905
906 // Destroy the KSP
907 destroyKSP();
908}
909
918{
920}
921
930void SystemSolver::assembly(const SparseMatrix &matrix, const SystemMatrixOrdering &reordering)
931{
932 // Check if the matrix is assembled
933 if (!matrix.isAssembled()) {
934 throw std::runtime_error("Unable to assembly the system. The matrix is not yet assembled.");
935 }
936
937 // Assembly the system matrix
938 SystemSparseMatrixAssembler assembler(&matrix);
939 assembly<SystemSolver>(static_cast<const Assembler &>(assembler), reordering);
940}
941
949void SystemSolver::assembly(const Assembler &assembler)
950{
952}
953
962void SystemSolver::assembly(const Assembler &assembler, const SystemMatrixOrdering &reordering)
963{
964 assembly<SystemSolver>(assembler, reordering);
965}
966
977void SystemSolver::update(long nRows, const long *rows, const SparseMatrix &elements)
978{
979 // Check if the element storage is assembled
980 if (!elements.isAssembled()) {
981 throw std::runtime_error("Unable to update the system. The element storage is not yet assembled.");
982 }
983
984 // Update matrix
985 SystemSparseMatrixAssembler assembler(&elements);
986 update<SystemSolver>(nRows, rows, static_cast<const Assembler &>(assembler));
987}
988
997void SystemSolver::update(const Assembler &assembler)
998{
999 update(getRowCount(), nullptr, assembler);
1000}
1001
1012void SystemSolver::update(long nRows, const long *rows, const Assembler &assembler)
1013{
1014 update<SystemSolver>(nRows, rows, assembler);
1015}
1016
1023{
1024 if (m_A == PETSC_NULLPTR) {
1025 return 0;
1026 }
1027
1028 PetscInt blockSize;
1029 MatGetBlockSize(m_A, &blockSize);
1030
1031 return blockSize;
1032}
1033
1044{
1045 if (m_A == PETSC_NULLPTR) {
1046 return 0;
1047 }
1048
1049 PetscInt nRows;
1050 MatGetLocalSize(m_A, &nRows, NULL);
1051 nRows /= getBlockSize();
1052
1053 return nRows;
1054}
1055
1066{
1067 if (m_A == PETSC_NULLPTR) {
1068 return 0;
1069 }
1070
1071 PetscInt nCols;
1072 MatGetLocalSize(m_A, NULL, &nCols);
1073 nCols /= getBlockSize();
1074
1075 return nCols;
1076}
1077
1088{
1089 if (m_A == PETSC_NULLPTR) {
1090 return 0;
1091 }
1092
1093 PetscInt nRows;
1094 MatGetLocalSize(m_A, &nRows, NULL);
1095
1096 return nRows;
1097}
1098
1109{
1110 if (m_A == PETSC_NULLPTR) {
1111 return 0;
1112 }
1113
1114 PetscInt nCols;
1115 MatGetLocalSize(m_A, NULL, &nCols);
1116
1117 return nCols;
1118}
1119
1120#if BITPIT_ENABLE_MPI==1
1131{
1132 if (m_A == PETSC_NULLPTR) {
1133 return 0;
1134 }
1135
1136 PetscInt nRows;
1137 MatGetSize(m_A, &nRows, NULL);
1138 nRows /= getBlockSize();
1139
1140 return nRows;
1141}
1142
1153{
1154 if (m_A == PETSC_NULLPTR) {
1155 return 0;
1156 }
1157
1158 PetscInt nCols;
1159 MatGetSize(m_A, NULL, &nCols);
1160 nCols /= getBlockSize();
1161
1162 return nCols;
1163}
1164
1174{
1175 if (m_A == PETSC_NULLPTR) {
1176 return 0;
1177 }
1178
1179 PetscInt nRows;
1180 MatGetSize(m_A, &nRows, NULL);
1181
1182 return nRows;
1183}
1184
1194{
1195 if (m_A == PETSC_NULLPTR) {
1196 return 0;
1197 }
1198
1199 PetscInt nCols;
1200 MatGetSize(m_A, NULL, &nCols);
1201
1202 return nCols;
1203}
1204
1211{
1212 return m_partitioned;
1213}
1214#endif
1215
1222{
1223 return m_assembled;
1224}
1225
1234{
1235 return true;
1236}
1237
1242{
1243 // Check if the system is assembled
1244 if (!isAssembled()) {
1245 throw std::runtime_error("Unable to solve the system. The system is not yet assembled.");
1246 }
1247
1248 // Perform actions before KSP solution
1250
1251 // Solve KSP
1252 solveKSP();
1253
1254 // Perform actions after KSP solution
1256}
1257
1265void SystemSolver::solve(const std::vector<double> &rhs, std::vector<double> *solution)
1266{
1267 // Fills the vectors
1268 vectorsFill(rhs, *solution);
1269
1270 // Solve the system
1271 solve();
1272
1273 // Export the solution
1274 exportVector(m_solution, solution);
1275}
1276
1281{
1282 PetscErrorCode solverError;
1283 if (!m_transpose) {
1284 solverError = KSPSolve(m_KSP, m_rhs, m_solution);
1285 } else {
1286 solverError = KSPSolveTranspose(m_KSP, m_rhs, m_solution);
1287 }
1288
1289 if (solverError) {
1290 const char *petscMessage = nullptr;
1291 PetscErrorMessage(solverError, &petscMessage, PETSC_NULLPTR);
1292 std::string message = "Unable to solver the system. " + std::string(petscMessage);
1293 throw std::runtime_error(message);
1294 }
1295}
1296
1301{
1302 // Prepare KSP
1303 prepareKSP();
1304
1305 // Reorder vectors
1306 vectorsReorder(true);
1307
1308 // Force consistency
1309 if (m_forceConsistency) {
1310 removeNullSpaceFromRHS();
1311 }
1312}
1313
1318{
1319 // Fill status of KSP
1320 fillKSPStatus();
1321
1322 // Reorder vectors
1323 vectorsReorder(false);
1324
1325 // Finalize KSP
1326 finalizeKSP();
1327}
1328
1335{
1336 const int DUMP_VERSION = 1;
1337
1338 return DUMP_VERSION;
1339}
1340
1347{
1348 const PetscInt *rowReordering = PETSC_NULLPTR;
1349 if (m_rowReordering) {
1350 ISGetIndices(m_rowReordering, &rowReordering);
1351 }
1352
1353 // Create the matrix
1354 int blockSize = assembler.getBlockSize();
1355 createMatrix(blockSize, blockSize, &m_A);
1356
1357 MatType matrixType;
1358 MatGetType(m_A, &matrixType);
1359
1360 // Get sizes
1361 long nAssemblerRows = assembler.getRowCount();
1362
1363 long nRowsElements = assembler.getRowElementCount();
1364 long nColsElements = assembler.getColElementCount();
1365
1366 long nGlobalRowsElements;
1367 long nGlobalColsElements;
1368#if BITPIT_ENABLE_MPI == 1
1369 nGlobalRowsElements = assembler.getRowGlobalElementCount();
1370 nGlobalColsElements = assembler.getColGlobalElementCount();
1371#else
1372 nGlobalRowsElements = nRowsElements;
1373 nGlobalColsElements = nColsElements;
1374#endif
1375
1376 MatSetSizes(m_A, nRowsElements, nColsElements, nGlobalRowsElements, nGlobalColsElements);
1377
1378 // Allocate storage
1379 //
1380 // When the internal storage of the system matrix was created without taking into account
1381 // block information, preallocation information should be provided for each row of each
1382 // block.
1383 int allocationExpansion;
1384 if (strcmp(matrixType, MATSEQAIJ) == 0) {
1385 allocationExpansion = blockSize;
1386#if BITPIT_ENABLE_MPI == 1
1387 } else if (strcmp(matrixType, MATMPIAIJ) == 0) {
1388 allocationExpansion = blockSize;
1389#endif
1390 } else {
1391 allocationExpansion = 1;
1392 }
1393
1394 long nAllocatedRows = allocationExpansion * nAssemblerRows;
1395
1396 std::vector<int> d_nnz(nAllocatedRows, 0);
1397 for (long n = 0; n < nAssemblerRows; ++n) {
1398 long matrixRow = n;
1399 if (rowReordering) {
1400 matrixRow = rowReordering[matrixRow];
1401 }
1402
1403 int nAssemblerRowNZ = assembler.getRowNZCount(n);
1404
1405 long matrixRowOffset = matrixRow * allocationExpansion;
1406 for (int i = 0; i < allocationExpansion; ++i) {
1407 d_nnz[matrixRowOffset + i] = allocationExpansion * nAssemblerRowNZ;
1408 }
1409 }
1410
1411#if BITPIT_ENABLE_MPI == 1
1412 std::vector<int> o_nnz(nAllocatedRows, 0);
1413 if (m_partitioned) {
1414 long nAssemblerCols = assembler.getColCount();
1415
1416 long assemblerDiagonalBegin = assembler.getColGlobalOffset();
1417 long assemblerDiagonalEnd = assemblerDiagonalBegin + nAssemblerCols;
1418
1419 ConstProxyVector<long> assemblerRowPattern(static_cast<std::size_t>(0), assembler.getMaxRowNZCount());
1420 for (long n = 0; n < nAssemblerRows; ++n) {
1421 long matrixRow = n;
1422 if (rowReordering) {
1423 matrixRow = rowReordering[matrixRow];
1424 }
1425
1426 assembler.getRowPattern(n, &assemblerRowPattern);
1427 int nAssemblerRowNZ = assemblerRowPattern.size();
1428
1429 long matrixRowOffset = matrixRow * allocationExpansion;
1430 for (int k = 0; k < nAssemblerRowNZ; ++k) {
1431 long id = assemblerRowPattern[k];
1432 if (id < assemblerDiagonalBegin || id >= assemblerDiagonalEnd) {
1433 for (int i = 0; i < allocationExpansion; ++i) {
1434 o_nnz[matrixRowOffset + i] += allocationExpansion;
1435 d_nnz[matrixRowOffset + i] -= allocationExpansion;
1436 }
1437 }
1438 }
1439 }
1440 }
1441#endif
1442
1443 if (strcmp(matrixType, MATSEQAIJ) == 0) {
1444 MatSeqAIJSetPreallocation(m_A, 0, d_nnz.data());
1445 } else if (strcmp(matrixType, MATSEQBAIJ) == 0) {
1446 MatSeqBAIJSetPreallocation(m_A, blockSize, 0, d_nnz.data());
1447#if BITPIT_ENABLE_MPI == 1
1448 } else if (strcmp(matrixType, MATMPIAIJ) == 0) {
1449 MatMPIAIJSetPreallocation(m_A, 0, d_nnz.data(), 0, o_nnz.data());
1450 } else if (strcmp(matrixType, MATMPIBAIJ) == 0) {
1451 MatMPIBAIJSetPreallocation(m_A, blockSize, 0, d_nnz.data(), 0, o_nnz.data());
1452#endif
1453 } else {
1454 throw std::runtime_error("Matrix format not supported.");
1455 }
1456
1457 // Each process will only set values for its own rows
1458 MatSetOption(m_A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE);
1459
1460#if PETSC_VERSION_GE(3, 12, 0)
1461 // The first assembly will set a superset of the off-process entries
1462 // required for all subsequent assemblies. This avoids a rendezvous
1463 // step in the MatAssembly functions.
1464 MatSetOption(m_A, MAT_SUBSET_OFF_PROC_ENTRIES, PETSC_TRUE);
1465#endif
1466
1467 // Cleanup
1468 if (m_rowReordering) {
1469 ISRestoreIndices(m_rowReordering, &rowReordering);
1470 }
1471
1472 // Fill matrix
1473 matrixUpdate(assembler.getRowCount(), nullptr, assembler);
1474
1475 // No new allocations are now allowed
1476 //
1477 // When updating the matrix it will not be possible to alter the pattern,
1478 // it will be possible to change only the values.
1479 MatSetOption(m_A, MAT_NEW_NONZERO_LOCATIONS, PETSC_FALSE);
1480}
1481
1491void SystemSolver::matrixFill(const std::string &filePath)
1492{
1493 // Check if the matrix exists
1494 if (!m_A) {
1495 throw std::runtime_error("Matrix should be created before filling it.");
1496 }
1497
1498 // Fill the matrix
1499 fillMatrix(m_A, filePath);
1500}
1501
1518void SystemSolver::matrixUpdate(long nRows, const long *rows, const Assembler &assembler)
1519{
1520 // Updating the matrix invalidates the KSP
1521 m_KSPDirty = true;
1522
1523 // Get block size
1524 const int blockSize = getBlockSize();
1525 if (assembler.getBlockSize() != blockSize) {
1526 std::string message = "Unable to update the matrix.";
1527 message += " The block size of the assembler is not equal to the block size of the system matrix.";
1528 throw std::runtime_error(message);
1529 }
1530
1531 // Initialize reordering
1532 const PetscInt *rowReordering = PETSC_NULLPTR;
1533 if (m_rowReordering) {
1534 ISGetIndices(m_rowReordering, &rowReordering);
1535 }
1536
1537 const PetscInt *colReordering = PETSC_NULLPTR;
1538 if (m_colReordering) {
1539 ISGetIndices(m_colReordering, &colReordering);
1540 }
1541
1542 // Global information
1543 PetscInt colGlobalBegin;
1544 PetscInt colGlobalEnd;
1545 MatGetOwnershipRangeColumn(m_A, &colGlobalBegin, &colGlobalEnd);
1546 colGlobalBegin /= blockSize;
1547 colGlobalEnd /= blockSize;
1548
1549 PetscInt rowGlobalOffset;
1550 MatGetOwnershipRange(m_A, &rowGlobalOffset, PETSC_NULLPTR);
1551 rowGlobalOffset /= blockSize;
1552
1553 // Get the options for assembling the matrix
1554 SystemMatrixAssembler::AssemblyOptions assemblyOptions = assembler.getOptions();
1555
1556#if PETSC_VERSION_GE(3, 12, 0)
1557 // Check if it is possible to speedup insertion of values
1558 //
1559 // The option MAT_SORTED_FULL means that each process provides exactly its
1560 // local rows; all column indices for a given row are passed in a single call
1561 // to MatSetValues(), preallocation is perfect, row oriented, INSERT_VALUES
1562 // is used. If this options is set to PETSC_TRUE, the function MatSetValues
1563 // will be faster.
1564 //
1565 // This options needs at least PETSc 3.12.
1566 PetscBool matrixSortedFull = (assemblyOptions.full && assemblyOptions.sorted) ? PETSC_TRUE : PETSC_FALSE;
1567 MatSetOption(m_A, MAT_SORTED_FULL, matrixSortedFull);
1568#endif
1569
1570 // Check if it possible to perform a fast update
1571 //
1572 // A fast update allows to set all the values of a row at once (without
1573 // the need to get the row pattern), it can be performed if:
1574 // - the system matrix has already been assembled;
1575 // - the system matrix has a unitary block size;
1576 // - the assembler is providing all the values of the row;
1577 // - values provided by the assembler are sorted by ascending column.
1578 //
1579 // If fast update is used, row values will be set using a special PETSc
1580 // function (MatSetValuesRow) that allows to set all the values of a
1581 // row at once, without requiring the pattern of the row.
1582 //
1583 // Fast update is not related to the option MAT_SORTED_FULL, that option
1584 // is used to speedup the standard function MatSetValues (which still
1585 // requires the pattern of the row).
1586 bool fastUpdate = isAssembled() && (blockSize == 1) && assemblyOptions.full && assemblyOptions.sorted;
1587
1588 // Update element values
1589 //
1590 // If the sizes of PETSc data types match the sizes of data types expected by
1591 // bitpit a direct update can be performed, otherwise the matrix is updated
1592 // using intermediate data storages.
1593 const long assemblerMaxRowNZ = std::max(assembler.getMaxRowNZCount(), 0L);
1594
1595 bool patternDirectUpdate = !colReordering && (sizeof(long) == sizeof(PetscInt));
1596 bool valuesDirectUpdate = (sizeof(double) == sizeof(PetscScalar));
1597
1598 ConstProxyVector<long> rowPattern;
1599 std::vector<PetscInt> petscRowPatternStorage;
1600 const PetscInt *petscRowPattern;
1601 if (!patternDirectUpdate) {
1602 rowPattern.set(ConstProxyVector<long>::INTERNAL_STORAGE, 0, assemblerMaxRowNZ);
1603 petscRowPatternStorage.resize(assemblerMaxRowNZ);
1604 petscRowPattern = petscRowPatternStorage.data();
1605 }
1606
1607 ConstProxyVector<double> rowValues;
1608 std::vector<PetscScalar> petscRowValuesStorage;
1609 const PetscScalar *petscRowValues;
1610 if (!valuesDirectUpdate) {
1611 long assemblerMaxRowNZElements = blockSize * blockSize * assemblerMaxRowNZ;
1612 rowValues.set(ConstProxyVector<double>::INTERNAL_STORAGE, 0, assemblerMaxRowNZElements);
1613 petscRowValuesStorage.resize(assemblerMaxRowNZElements);
1614 petscRowValues = petscRowValuesStorage.data();
1615 }
1616
1617 for (long n = 0; n < nRows; ++n) {
1618 // Get row information
1619 long row;
1620 if (rows) {
1621 row = rows[n];
1622 } else {
1623 row = n;
1624 }
1625
1626 if (rowReordering) {
1627 row = rowReordering[row];
1628 }
1629
1630 const PetscInt globalRow = rowGlobalOffset + row;
1631
1632 // Get row data
1633 if (fastUpdate) {
1634 assembler.getRowValues(n, &rowValues);
1635 } else {
1636 assembler.getRowData(n, &rowPattern, &rowValues);
1637 }
1638
1639 if (rowValues.size() == 0) {
1640 continue;
1641 }
1642
1643 // Get values in PETSc format
1644 const std::size_t rowPatternSize = rowPattern.size();
1645
1646 if (valuesDirectUpdate) {
1647 petscRowValues = reinterpret_cast<const PetscScalar *>(rowValues.data());
1648 } else {
1649 std::copy(rowValues.cbegin(), rowValues.cend(), petscRowValuesStorage.begin());
1650 }
1651
1652 if (fastUpdate) {
1653 MatSetValuesRow(m_A, globalRow, petscRowValues);
1654 } else {
1655 // Get pattern in PETSc format
1656 if (patternDirectUpdate) {
1657 petscRowPattern = reinterpret_cast<const PetscInt *>(rowPattern.data());
1658 } else {
1659 for (std::size_t k = 0; k < rowPatternSize; ++k) {
1660 long globalCol = rowPattern[k];
1661 if (colReordering) {
1662 if (globalCol >= colGlobalBegin && globalCol < colGlobalEnd) {
1663 long col = globalCol - colGlobalBegin;
1664 col = colReordering[col];
1665 globalCol = colGlobalBegin + col;
1666 }
1667 }
1668
1669 petscRowPatternStorage[k] = globalCol;
1670 }
1671 }
1672
1673 // Set data
1674 if (blockSize > 1) {
1675 MatSetValuesBlocked(m_A, 1, &globalRow, rowPatternSize, petscRowPattern, petscRowValues, INSERT_VALUES);
1676 } else {
1677 MatSetValues(m_A, 1, &globalRow, rowPatternSize, petscRowPattern, petscRowValues, INSERT_VALUES);
1678 }
1679 }
1680 }
1681
1682 // Let petsc assembly the matrix after the update
1683 MatAssemblyBegin(m_A, MAT_FINAL_ASSEMBLY);
1684 MatAssemblyEnd(m_A, MAT_FINAL_ASSEMBLY);
1685
1686 // Cleanup
1687 if (rowReordering) {
1688 ISRestoreIndices(m_rowReordering, &rowReordering);
1689 }
1690
1691 if (colReordering) {
1692 ISRestoreIndices(m_colReordering, &colReordering);
1693 }
1694}
1695
1703void SystemSolver::matrixDump(std::ostream &systemStream, const std::string &directory,
1704 const std::string &prefix) const
1705{
1706 BITPIT_UNUSED(systemStream);
1707
1708 dumpMatrix(m_A, directory, prefix + "A");
1709}
1710
1718#if BITPIT_ENABLE_MPI==1
1724void SystemSolver::matrixRestore(std::istream &systemStream, const std::string &directory,
1725 const std::string &prefix, bool redistribute)
1726#else
1727void SystemSolver::matrixRestore(std::istream &systemStream, const std::string &directory,
1728 const std::string &prefix)
1729#endif
1730{
1731 BITPIT_UNUSED(systemStream);
1732
1733#if BITPIT_ENABLE_MPI==1
1734 restoreMatrix(directory, prefix + "A", redistribute, &m_A);
1735#else
1736 restoreMatrix(directory, prefix + "A", &m_A);
1737#endif
1738}
1739
1744{
1745 destroyMatrix(&m_A);
1746}
1747
1754{
1755 if (!m_transpose) {
1756 MatCreateVecs(m_A, &m_solution, &m_rhs);
1757 } else {
1758 MatCreateVecs(m_A, &m_rhs, &m_solution);
1759 }
1760}
1761
1770void SystemSolver::vectorsFill(const std::vector<double> &rhs, const std::vector<double> &solution)
1771{
1772 if (!rhs.empty()) {
1773 fillVector(m_rhs, rhs);
1774 }
1775
1776 if (!solution.empty()) {
1777 fillVector(m_solution, solution);
1778 }
1779}
1780
1789void SystemSolver::vectorsFill(const std::string &rhsFilePath, const std::string &solutionFilePath)
1790{
1791 if (!rhsFilePath.empty()) {
1792 fillVector(m_rhs, rhsFilePath);
1793 }
1794
1795 if (!solutionFilePath.empty()) {
1796 fillVector(m_solution, solutionFilePath);
1797 }
1798}
1799
1806{
1807 if (!m_transpose) {
1808 reorderVector(m_rhs, m_colReordering, invert);
1809 reorderVector(m_solution, m_rowReordering, invert);
1810 } else {
1811 reorderVector(m_rhs, m_rowReordering, invert);
1812 reorderVector(m_solution, m_colReordering, invert);
1813 }
1814}
1815
1823void SystemSolver::vectorsDump(std::ostream &systemStream, const std::string &directory,
1824 const std::string &prefix) const
1825{
1826 BITPIT_UNUSED(systemStream);
1827
1828 dumpVector(m_rhs, directory, prefix + "rhs");
1829 dumpVector(m_solution, directory, prefix + "solution");
1830}
1831
1839void SystemSolver::vectorsRestore(std::istream &systemStream, const std::string &directory,
1840 const std::string &prefix)
1841{
1842 BITPIT_UNUSED(systemStream);
1843
1844 restoreVector(directory, prefix + "rhs", m_A, VECTOR_SIDE_RIGHT, &m_rhs);
1845 restoreVector(directory, prefix + "solution", m_A, VECTOR_SIDE_LEFT, &m_solution);
1846}
1847
1852{
1853 destroyVector(&m_rhs);
1854 destroyVector(&m_solution);
1855}
1856
1863{
1864 PetscScalar *raw_rhs;
1865 VecGetArray(m_rhs, &raw_rhs);
1866
1867 return raw_rhs;
1868}
1869
1875const double * SystemSolver::getRHSRawPtr() const
1876{
1877 return getRHSRawReadPtr();
1878}
1879
1886{
1887 const PetscScalar *raw_rhs;
1888 VecGetArrayRead(m_rhs, &raw_rhs);
1889
1890 return raw_rhs;
1891}
1892
1900{
1901 VecRestoreArray(m_rhs, &raw_rhs);
1902}
1903
1910void SystemSolver::restoreRHSRawReadPtr(const double *raw_rhs) const
1911{
1912 VecRestoreArrayRead(m_rhs, &raw_rhs);
1913}
1914
1921{
1922 PetscScalar *raw_solution;
1923 VecGetArray(m_solution, &raw_solution);
1924
1925 return raw_solution;
1926}
1927
1934{
1935 return getSolutionRawReadPtr();
1936}
1937
1944{
1945 const PetscScalar *raw_solution;
1946 VecGetArrayRead(m_solution, &raw_solution);
1947
1948 return raw_solution;
1949}
1950
1957void SystemSolver::restoreSolutionRawPtr(double *raw_solution)
1958{
1959 VecRestoreArray(m_solution, &raw_solution);
1960}
1961
1968void SystemSolver::restoreSolutionRawReadPtr(const double *raw_solution) const
1969{
1970 VecRestoreArrayRead(m_solution, &raw_solution);
1971}
1972
1980void SystemSolver::exportMatrix(const std::string &filePath, FileFormat fileFormat) const
1981{
1982 exportMatrix(m_A, filePath, fileFormat);
1983}
1984
1994void SystemSolver::importMatrix(const std::string &filePath)
1995{
1996 // Check if the matrix exists
1997 if (!m_A) {
1998 throw std::runtime_error("Matrix should be created before filling it.");
1999 }
2000
2001 // Clear workspace
2003
2004 // Clear reordering
2006
2007 // Fill the matrix
2008 matrixFill(filePath);
2009
2010 // Re-create vectors
2012 vectorsCreate();
2013}
2014
2022void SystemSolver::exportRHS(const std::string &filePath, FileFormat fileFormat) const
2023{
2024 exportVector(m_rhs, filePath, fileFormat);
2025}
2026
2038void SystemSolver::importRHS(const std::string &filePath)
2039{
2040 // Check if the system is assembled
2041 if (!m_assembled) {
2042 throw std::runtime_error("The RHS vector can be loaded only after assembling the system.");
2043 }
2044
2045 // Fill the RHS vector
2046 vectorsFill(filePath, "");
2047
2048 // Check if the imported RHS is compatible with the matrix
2049 PetscInt size;
2050 VecGetLocalSize(m_rhs, &size);
2051
2052 PetscInt expectedSize;
2053 if (!m_transpose) {
2054 expectedSize = getColCount();
2055 } else {
2056 expectedSize = getRowCount();
2057 }
2058 expectedSize *= getBlockSize();
2059
2060 if (size != expectedSize) {
2061 log::cout() << "The imported RHS vector is not compatible with the matrix" << std::endl;
2062 log::cout() << "The size of the imported RHS vector is " << size << std::endl;
2063 log::cout() << "The expected size of RHS vector is " << expectedSize << std::endl;
2064 throw std::runtime_error("The imported RHS vector is not compatible with the matrix");
2065 }
2066}
2067
2075void SystemSolver::exportSolution(const std::string &filePath, FileFormat fileFormat) const
2076{
2077 exportVector(m_solution, filePath, fileFormat);
2078}
2079
2091void SystemSolver::importSolution(const std::string &filePath)
2092{
2093 // Check if the system is assembled
2094 if (!m_assembled) {
2095 throw std::runtime_error("The solution vector can be loaded only after assembling the system.");
2096 }
2097
2098 // Fill the solution vector
2099 vectorsFill("", filePath);
2100
2101 // Check if the imported solution is compatible with the matrix
2102 PetscInt size;
2103 VecGetLocalSize(m_solution, &size);
2104
2105 PetscInt expectedSize;
2106 if (!m_transpose) {
2107 expectedSize = getRowCount();
2108 } else {
2109 expectedSize = getColCount();
2110 }
2111 expectedSize *= getBlockSize();
2112
2113 if (size != expectedSize) {
2114 log::cout() << "The imported solution vector is not compatible with the matrix" << std::endl;
2115 log::cout() << "The size of the imported solution vector is " << size << std::endl;
2116 log::cout() << "The expected size of solution vector is " << expectedSize << std::endl;
2117 throw std::runtime_error("The imported solution vector is not compatible with the matrix");
2118 }
2119}
2120
2132void SystemSolver::dumpSystem(const std::string &header, const std::string &directory,
2133 const std::string &prefix) const
2134{
2135 int partitioningBlock = getBinaryArchiveBlock();
2136
2137 // Open stream that will contain system information
2138 OBinaryArchive systemArchive;
2139 if (partitioningBlock <= 0) {
2140 openBinaryArchive(header, directory, prefix + "info", -1, &systemArchive);
2141 }
2142 std::ostream &systemStream = systemArchive.getStream();
2143
2144 // Dump system information
2145 dumpInfo(systemStream);
2146
2147 // Dump matrix
2148 matrixDump(systemStream, directory, prefix);
2149
2150 // Dump vectors
2151 vectorsDump(systemStream, directory, prefix);
2152
2153 // Open stream with system information
2154 closeBinaryArchive(&systemArchive);
2155}
2156
2167#if BITPIT_ENABLE_MPI==1
2174#endif
2179#if BITPIT_ENABLE_MPI==1
2180void SystemSolver::restoreSystem(MPI_Comm communicator, bool redistribute,
2181 const std::string &directory, const std::string &prefix)
2182#else
2183void SystemSolver::restoreSystem(const std::string &directory, const std::string &prefix)
2184#endif
2185{
2186 // Clear the system
2187 clear();
2188
2189#if BITPIT_ENABLE_MPI == 1
2190 // Set the communicator
2191 setCommunicator(communicator);
2192#endif
2193
2194 // Open stream with system information
2195 IBinaryArchive systemArchive;
2196 openBinaryArchive(directory, prefix + "info", -1, &systemArchive);
2197 std::istream &systemStream = systemArchive.getStream();
2198
2199 // Dump system information
2200 restoreInfo(systemStream);
2201
2202 // Restore the matrix
2203#if BITPIT_ENABLE_MPI==1
2204 matrixRestore(systemStream, directory, prefix, redistribute);
2205#else
2206 matrixRestore(systemStream, directory, prefix);
2207#endif
2208
2209 // Restore RHS and solution vectors
2210 vectorsRestore(systemStream, directory, prefix);
2211
2212 // Close stream with system information
2213 closeBinaryArchive(&systemArchive);
2214
2215 // The system is now assembled
2216 m_assembled = true;
2217
2218 // Initialize KSP options
2219 initializeKSPOptions();
2220
2221 // Initialize KSP statuses
2222 initializeKSPStatus();
2223}
2224
2230void SystemSolver::dumpInfo(std::ostream &systemStream) const
2231{
2232 if (systemStream.good()) {
2233#if BITPIT_ENABLE_MPI==1
2234 utils::binary::write(systemStream, m_partitioned);
2235#endif
2236 utils::binary::write(systemStream, m_transpose);
2237 }
2238}
2239
2245void SystemSolver::restoreInfo(std::istream &systemStream)
2246{
2247#if BITPIT_ENABLE_MPI == 1
2248 // Detect if the system is partitioned
2249 utils::binary::read(systemStream, m_partitioned);
2250#endif
2251
2252 // Set transpose flag
2253 bool transpose;
2254 utils::binary::read(systemStream, transpose);
2255 setTranspose(transpose);
2256
2257}
2258
2266void SystemSolver::createMatrix(int rowBlockSize, int colBlockSize, Mat *matrix) const
2267{
2268 // Create the matrix
2269#if BITPIT_ENABLE_MPI == 1
2270 MatCreate(m_communicator, matrix);
2271#else
2272 MatCreate(PETSC_COMM_SELF, matrix);
2273#endif
2274
2275 // Set matrix type
2276 bool creatBlockMatrix = false;
2277 if (m_flatten) {
2278 creatBlockMatrix = false;
2279 } else {
2280 creatBlockMatrix = (rowBlockSize == colBlockSize) && (rowBlockSize != 1);
2281 }
2282
2283#if BITPIT_ENABLE_MPI == 1
2284 if (m_partitioned) {
2285 if (creatBlockMatrix) {
2286 MatSetType(*matrix, MATMPIBAIJ);
2287 } else {
2288 MatSetType(*matrix, MATMPIAIJ);
2289 }
2290 } else
2291#endif
2292 {
2293 if (creatBlockMatrix) {
2294 MatSetType(*matrix, MATSEQBAIJ);
2295 } else {
2296 MatSetType(*matrix, MATSEQAIJ);
2297 }
2298 }
2299
2300 // Set block size
2301 if (rowBlockSize == colBlockSize) {
2302 MatSetBlockSize(*matrix, rowBlockSize);
2303 } else if (rowBlockSize != 1 || colBlockSize != 1) {
2304 MatSetBlockSizes(*matrix, rowBlockSize, colBlockSize);
2305 }
2306}
2307
2321void SystemSolver::createMatrix(int rowBlockSize, int colBlockSize, int nNestRows, int nNestCols,
2322 Mat *subMatrices, Mat *matrix) const
2323{
2324 // Create the matrix
2325#if BITPIT_ENABLE_MPI == 1
2326 MatCreateNest(getCommunicator(), nNestRows, PETSC_NULLPTR, nNestCols, PETSC_NULLPTR, subMatrices, matrix);
2327#else
2328 MatCreateNest(PETSC_COMM_SELF, nNestRows, PETSC_NULLPTR, nNestCols, PETSC_NULLPTR, subMatrices, matrix);
2329#endif
2330
2331 // Set block size
2332 if (rowBlockSize == colBlockSize) {
2333 MatSetBlockSize(*matrix, rowBlockSize);
2334 } else if (rowBlockSize != 1 || colBlockSize != 1) {
2335 MatSetBlockSizes(*matrix, rowBlockSize, colBlockSize);
2336 }
2337}
2338
2348void SystemSolver::fillMatrix(Mat matrix, const std::string &filePath) const
2349{
2350 // Check if the file exists
2351 std::ifstream fileStream(filePath.c_str());
2352 if (!fileStream.good()) {
2353 throw std::runtime_error("The PETSc matrix file \"" + filePath + "\" doesn't exists.");
2354 }
2355 fileStream.close();
2356
2357 // Fill matrix content
2358 PetscViewer viewer;
2359#if BITPIT_ENABLE_MPI==1
2360 PetscViewerCreate(m_communicator, &viewer);
2361#else
2362 PetscViewerCreate(PETSC_COMM_SELF, &viewer);
2363#endif
2364 PetscViewerSetType(viewer, PETSCVIEWERBINARY);
2365 PetscViewerFileSetMode(viewer, FILE_MODE_READ);
2366 PetscViewerFileSetName(viewer, filePath.c_str());
2367 MatLoad(matrix, viewer);
2368 PetscViewerDestroy(&viewer);
2369}
2370
2378#if BITPIT_ENABLE_MPI==1
2383#endif
2390void SystemSolver::dumpMatrix(Mat matrix, const std::string &directory, const std::string &name) const
2391{
2392 int partitioningBlock = getBinaryArchiveBlock();
2393
2394 // Store information needed to create the matrix
2395 OBinaryArchive infoArchive;
2396 if (partitioningBlock <= 0) {
2397 openBinaryArchive("", directory, name + ".info", -1, &infoArchive);
2398 std::ostream &infoStream = infoArchive.getStream();
2399
2400 bool matrixExists = matrix;
2401 utils::binary::write(infoStream, matrixExists);
2402 if (!matrixExists) {
2403 closeBinaryArchive(&infoArchive);
2404 return;
2405 }
2406
2407 PetscInt rowBlockSize;
2408 PetscInt colBlockSize;
2409 MatGetBlockSizes(matrix, &rowBlockSize, &colBlockSize);
2410 utils::binary::write(infoStream, static_cast<int>(rowBlockSize));
2411 utils::binary::write(infoStream, static_cast<int>(colBlockSize));
2412
2413 closeBinaryArchive(&infoArchive);
2414 }
2415
2416 if (!matrix) {
2417 return;
2418 }
2419
2420#if BITPIT_ENABLE_MPI==1
2421 // Store partitioning information
2422 if (isPartitioned()) {
2423 OBinaryArchive partitioningArchive;
2424 openBinaryArchive("", directory, name + ".partitioning", partitioningBlock, &partitioningArchive);
2425 std::ostream &partitioningStream = partitioningArchive.getStream();
2426
2427 PetscInt nLocalRows;
2428 PetscInt nLocalCols;
2429 MatGetLocalSize(matrix, &nLocalRows, &nLocalCols);
2430 utils::binary::write(partitioningStream, static_cast<std::size_t>(nLocalRows));
2431 utils::binary::write(partitioningStream, static_cast<std::size_t>(nLocalCols));
2432
2433 PetscInt nGlobalRows;
2434 PetscInt nGlobalCols;
2435 MatGetSize(matrix, &nGlobalRows, &nGlobalCols);
2436 utils::binary::write(partitioningStream, static_cast<std::size_t>(nGlobalRows));
2437 utils::binary::write(partitioningStream, static_cast<std::size_t>(nGlobalCols));
2438
2439 closeBinaryArchive(&partitioningArchive);
2440 }
2441#endif
2442
2443 // Store matrix data
2444 std::string filePath = getDataFilePath(directory, name);
2445 exportMatrix(matrix, filePath, FILE_BINARY);
2446}
2447
2454#if BITPIT_ENABLE_MPI==1
2460#endif
2464#if BITPIT_ENABLE_MPI==1
2465void SystemSolver::restoreMatrix(const std::string &directory, const std::string &name,
2466 bool redistribute, Mat *matrix) const
2467{
2468#else
2469void SystemSolver::restoreMatrix(const std::string &directory, const std::string &name,
2470 Mat *matrix) const
2471{
2472#endif
2473
2474 // Create matrix
2475 IBinaryArchive infoArchive;
2476 openBinaryArchive(directory, name + ".info", -1, &infoArchive);
2477 std::istream &infoStream = infoArchive.getStream();
2478
2479 bool matrixExists;
2480 utils::binary::read(infoStream, matrixExists);
2481 if (!matrixExists) {
2482 *matrix = PETSC_NULLPTR;
2483 closeBinaryArchive(&infoArchive);
2484 return;
2485 }
2486
2487 int rowBlockSize;
2488 int colBlockSize;
2489 utils::binary::read(infoStream, rowBlockSize);
2490 utils::binary::read(infoStream, colBlockSize);
2491 createMatrix(rowBlockSize, colBlockSize, matrix);
2492
2493 closeBinaryArchive(&infoArchive);
2494
2495#if BITPIT_ENABLE_MPI==1
2496 // Set partitioning information
2497 if (isPartitioned() && !redistribute) {
2498 int partitioningBlock = getBinaryArchiveBlock();
2499
2500 IBinaryArchive partitioningArchive;
2501 openBinaryArchive(directory, name + ".partitioning", partitioningBlock, &partitioningArchive);
2502 std::istream &partitioningStream = partitioningArchive.getStream();
2503
2504 std::size_t nLocalRows;
2505 std::size_t nLocalCols;
2506 utils::binary::read(partitioningStream, nLocalRows);
2507 utils::binary::read(partitioningStream, nLocalCols);
2508 std::size_t nGlobalRows;
2509 std::size_t nGlobalCols;
2510 utils::binary::read(partitioningStream, nGlobalRows);
2511 utils::binary::read(partitioningStream, nGlobalCols);
2512 MatSetSizes(*matrix, nLocalRows, nLocalCols, nGlobalRows, nGlobalCols);
2513
2514 closeBinaryArchive(&partitioningArchive);
2515 }
2516#endif
2517
2518 // Fill matrix
2519 std::string filePath = getDataFilePath(directory, name);
2520 fillMatrix(*matrix, filePath);
2521}
2522
2531void SystemSolver::exportMatrix(Mat matrix, const std::string &filePath, FileFormat fileFormat) const
2532{
2533 // Early return if the matrix doesn't exist
2534 bool matrixExists = matrix;
2535 if (!matrixExists) {
2536 std::ofstream dataFile(filePath);
2537 dataFile.close();
2538
2539 std::ofstream infoFile(filePath + ".info");
2540 infoFile.close();
2541
2542 return;
2543 }
2544
2545 // Create the viewer
2546 PetscViewerType viewerType;
2547 PetscViewerFormat viewerFormat;
2548 if (fileFormat == FILE_BINARY) {
2549 viewerType = PETSCVIEWERBINARY;
2550 viewerFormat = PETSC_VIEWER_DEFAULT;
2551 } else {
2552 viewerType = PETSCVIEWERASCII;
2553 viewerFormat = PETSC_VIEWER_ASCII_MATLAB;
2554 }
2555
2556 PetscViewer matViewer;
2557#if BITPIT_ENABLE_MPI==1
2558 PetscViewerCreate(m_communicator, &matViewer);
2559#else
2560 PetscViewerCreate(PETSC_COMM_SELF, &matViewer);
2561#endif
2562 PetscViewerSetType(matViewer, viewerType);
2563 PetscViewerFileSetMode(matViewer, FILE_MODE_WRITE);
2564 PetscViewerPushFormat(matViewer, viewerFormat);
2565 PetscViewerFileSetName(matViewer, filePath.c_str());
2566
2567 // Export the matrix
2568 MatView(matrix, matViewer);
2569
2570 // Destroy the viewer
2571 PetscViewerDestroy(&matViewer);
2572}
2573
2579void SystemSolver::destroyMatrix(Mat *matrix) const
2580{
2581 if (matrix) {
2582 MatDestroy(matrix);
2583 *matrix = PETSC_NULLPTR;
2584 }
2585}
2586
2593void SystemSolver::createVector(int blockSize, Vec *vector) const
2594{
2595 // Create the vector
2596#if BITPIT_ENABLE_MPI == 1
2597 VecCreate(m_communicator, vector);
2598#else
2599 VecCreate(PETSC_COMM_SELF, vector);
2600#endif
2601
2602 // Set vector type
2603#if BITPIT_ENABLE_MPI == 1
2604 if (m_partitioned) {
2605 VecSetType(*vector, VECMPI);
2606 } else
2607#endif
2608 {
2609 VecSetType(*vector, VECSEQ);
2610 }
2611
2612 // Set block size
2613 if (blockSize != 1) {
2614 VecSetBlockSize(*vector, blockSize);
2615 }
2616}
2617
2626void SystemSolver::createVector(int blockSize, int nestSize, Vec *subVectors, Vec *vector) const
2627{
2628 // Create the vector
2629#if BITPIT_ENABLE_MPI == 1
2630 VecCreateNest(getCommunicator(), nestSize, PETSC_NULLPTR, subVectors, vector);
2631#else
2632 VecCreateNest(PETSC_COMM_SELF, nestSize, PETSC_NULLPTR, subVectors, vector);
2633#endif
2634
2635 // Set block size
2636 if (blockSize != 1) {
2637 VecSetBlockSize(*vector, blockSize);
2638 }
2639}
2640
2650void SystemSolver::fillVector(Vec vector, const std::string &filePath) const
2651{
2652 // Check if the file exists
2653 std::ifstream fileStream(filePath);
2654 if (!fileStream.good()) {
2655 throw std::runtime_error("The file \"" + filePath + "\" cannot be read.");
2656 }
2657 fileStream.close();
2658
2659 // Fill the vector
2660 PetscViewer viewer;
2661#if BITPIT_ENABLE_MPI==1
2662 PetscViewerCreate(m_communicator, &viewer);
2663#else
2664 PetscViewerCreate(PETSC_COMM_SELF, &viewer);
2665#endif
2666 PetscViewerSetType(viewer, PETSCVIEWERBINARY);
2667 PetscViewerFileSetMode(viewer, FILE_MODE_READ);
2668 PetscViewerFileSetName(viewer, filePath.c_str());
2669 VecLoad(vector, viewer);
2670 PetscViewerDestroy(&viewer);
2671}
2672
2679void SystemSolver::fillVector(Vec vector, const std::vector<double> &data) const
2680{
2681 int size;
2682 VecGetLocalSize(vector, &size);
2683
2684 PetscScalar *petscData;
2685 VecGetArray(vector, &petscData);
2686 for (int i = 0; i < size; ++i) {
2687 petscData[i] = data[i];
2688 }
2689 VecRestoreArray(vector, &petscData);
2690}
2691
2699#if BITPIT_ENABLE_MPI==1
2704#endif
2711void SystemSolver::dumpVector(Vec vector, const std::string &directory, const std::string &name) const
2712{
2713 int partitioningBlock = getBinaryArchiveBlock();
2714
2715 // Store information needed to create the matrix
2716 if (partitioningBlock <= 0) {
2717 OBinaryArchive infoArchive;
2718 openBinaryArchive("", directory, name + ".info", -1, &infoArchive);
2719 std::ostream &infoStream = infoArchive.getStream();
2720
2721 bool vectorExists = vector;
2722 utils::binary::write(infoStream, vectorExists);
2723 if (!vectorExists) {
2724 closeBinaryArchive(&infoArchive);
2725 return;
2726 }
2727
2728 PetscInt blockSize;
2729 VecGetBlockSize(vector, &blockSize);
2730 utils::binary::write(infoStream, static_cast<int>(blockSize));
2731
2732 closeBinaryArchive(&infoArchive);
2733 }
2734
2735 if (!vector) {
2736 return;
2737 }
2738
2739#if BITPIT_ENABLE_MPI==1
2740 // Store partitioning information
2741 if (isPartitioned()) {
2742 int partitioningBlock = getBinaryArchiveBlock();
2743
2744 OBinaryArchive partitioningArchive;
2745 openBinaryArchive("", directory, name + ".partitioning", partitioningBlock, &partitioningArchive);
2746 std::ostream &partitioningStream = partitioningArchive.getStream();
2747
2748 PetscInt localSize;
2749 VecGetLocalSize(vector, &localSize);
2750 utils::binary::write(partitioningStream, static_cast<std::size_t>(localSize));
2751
2752 PetscInt globalSize;
2753 VecGetSize(vector, &globalSize);
2754 utils::binary::write(partitioningStream, static_cast<std::size_t>(globalSize));
2755
2756 closeBinaryArchive(&partitioningArchive);
2757 }
2758#endif
2759
2760 // Store vector data
2761 std::string filePath = getDataFilePath(directory, name);
2762 exportVector(vector, filePath, FILE_BINARY);
2763}
2764
2775void SystemSolver::restoreVector(const std::string &directory, const std::string &name,
2776 Mat matrix, VectorSide side, Vec *vector) const
2777{
2778 // Get size information
2779 PetscInt nLocalRows;
2780 PetscInt nLocalCols;
2781 MatGetLocalSize(matrix, &nLocalRows, &nLocalCols);
2782
2783 PetscInt nGlobalRows;
2784 PetscInt nGlobalCols;
2785 MatGetSize(matrix, &nGlobalRows, &nGlobalCols);
2786
2787 std::size_t localSize = std::numeric_limits<std::size_t>::max();
2788 std::size_t globalSize = std::numeric_limits<std::size_t>::max();
2789 if (!m_transpose) {
2790 if (side == VECTOR_SIDE_RIGHT) {
2791 localSize = nLocalCols;
2792 globalSize = nGlobalCols;
2793 } else if (side == VECTOR_SIDE_LEFT) {
2794 localSize = nLocalRows;
2795 globalSize = nGlobalRows;
2796 }
2797 } else {
2798 if (side == VECTOR_SIDE_RIGHT) {
2799 localSize = nLocalRows;
2800 globalSize = nGlobalRows;
2801 } else if (side == VECTOR_SIDE_LEFT) {
2802 localSize = nLocalCols;
2803 globalSize = nGlobalCols;
2804 }
2805 }
2806
2807 // Restore the vector
2808 restoreVector(directory, name, localSize, globalSize, vector);
2809}
2810
2817#if BITPIT_ENABLE_MPI==1
2823#endif
2827#if BITPIT_ENABLE_MPI==1
2828void SystemSolver::restoreVector(const std::string &directory, const std::string &name,
2829 bool redistribute, Vec *vector) const
2830#else
2831void SystemSolver::restoreVector(const std::string &directory, const std::string &name,
2832 Vec *vector) const
2833#endif
2834{
2835
2836 // Get size information
2837 std::size_t localSize = std::numeric_limits<std::size_t>::max();
2838 std::size_t globalSize = std::numeric_limits<std::size_t>::max();
2839#if BITPIT_ENABLE_MPI==1
2840 if (isPartitioned()) {
2841 if (!redistribute) {
2842 int partitioningBlock = getBinaryArchiveBlock();
2843
2844 IBinaryArchive partitioningArchive;
2845 openBinaryArchive(directory, name + ".partitioning", partitioningBlock, &partitioningArchive);
2846 std::istream &partitioningStream = partitioningArchive.getStream();
2847
2848 utils::binary::read(partitioningStream, localSize);
2849 utils::binary::read(partitioningStream, globalSize);
2850
2851 closeBinaryArchive(&partitioningArchive);
2852 }
2853 }
2854#endif
2855
2856 restoreVector(directory, name, localSize, globalSize, vector);
2857}
2858
2868void SystemSolver::restoreVector(const std::string &directory, const std::string &name,
2869 std::size_t localSize, std::size_t globalSize, Vec *vector) const
2870{
2871 // Create vector
2872 IBinaryArchive infoArchive;
2873 openBinaryArchive(directory, name + ".info", -1, &infoArchive);
2874 std::istream &infoStream = infoArchive.getStream();
2875
2876 bool vectorExists;;
2877 utils::binary::read(infoStream, vectorExists);
2878 if (!vectorExists) {
2879 *vector = PETSC_NULLPTR;
2880 closeBinaryArchive(&infoArchive);
2881 return;
2882 }
2883
2884 int blockSize;
2885 utils::binary::read(infoStream, blockSize);
2886 createVector(blockSize, vector);
2887
2888 closeBinaryArchive(&infoArchive);
2889
2890 // Set size information
2891 PetscInt petscLocalSize;
2892 if (localSize != std::numeric_limits<std::size_t>::max()) {
2893 petscLocalSize = static_cast<PetscInt>(localSize);
2894 } else {
2895 petscLocalSize = PETSC_DECIDE;
2896 }
2897
2898 PetscInt petscGlobalSize;
2899 if (globalSize != std::numeric_limits<std::size_t>::max()) {
2900 petscGlobalSize = static_cast<PetscInt>(globalSize);
2901 } else {
2902 petscGlobalSize = PETSC_DETERMINE;
2903 }
2904
2905 if (petscLocalSize != PETSC_DECIDE || petscGlobalSize != PETSC_DETERMINE) {
2906 VecSetSizes(*vector, petscLocalSize, petscGlobalSize);
2907 }
2908
2909 // Fill vector data
2910 std::string filePath = getDataFilePath(directory, name);
2911 fillVector(*vector, filePath);
2912}
2913
2919void SystemSolver::reorderVector(Vec vector, IS permutations, bool invert) const
2920{
2921 if (!permutations) {
2922 return;
2923 }
2924
2925 PetscBool petscInvert;
2926 if (invert) {
2927 petscInvert = PETSC_TRUE;
2928 } else {
2929 petscInvert = PETSC_FALSE;
2930 }
2931
2932 VecPermute(vector, permutations, petscInvert);
2933}
2934
2943void SystemSolver::exportVector(Vec vector, const std::string &filePath, FileFormat fileFormat) const
2944{
2945 // Early return if the matrix doesn't exist
2946 bool vectorExists = vector;
2947 if (!vectorExists) {
2948 std::ofstream dataFile(filePath);
2949 dataFile.close();
2950
2951 std::ofstream infoFile(filePath + ".info");
2952 infoFile.close();
2953
2954 return;
2955 }
2956
2957 // Create the viewer
2958 PetscViewerType viewerType;
2959 PetscViewerFormat viewerFormat;
2960 if (fileFormat == FILE_BINARY) {
2961 viewerType = PETSCVIEWERBINARY;
2962 viewerFormat = PETSC_VIEWER_DEFAULT;
2963 } else {
2964 viewerType = PETSCVIEWERASCII;
2965 viewerFormat = PETSC_VIEWER_ASCII_MATLAB;
2966 }
2967
2968 PetscViewer viewer;
2969#if BITPIT_ENABLE_MPI==1
2970 PetscViewerCreate(m_communicator, &viewer);
2971#else
2972 PetscViewerCreate(PETSC_COMM_SELF, &viewer);
2973#endif
2974 PetscViewerSetType(viewer, viewerType);
2975 PetscViewerFileSetMode(viewer, FILE_MODE_WRITE);
2976 PetscViewerPushFormat(viewer, viewerFormat);
2977 PetscViewerFileSetName(viewer, filePath.c_str());
2978
2979 // Export the vector
2980 VecView(vector, viewer);
2981
2982 // Destroy the viewer
2983 PetscViewerDestroy(&viewer);
2984}
2985
2992void SystemSolver::exportVector(Vec vector, std::vector<double> *data) const
2993{
2994 int size;
2995 VecGetLocalSize(m_solution, &size);
2996
2997 const PetscScalar *vectorData;
2998 VecGetArrayRead(m_solution, &vectorData);
2999 for (int i = 0; i < size; ++i) {
3000 (*data)[i] = vectorData[i];
3001 }
3002 VecRestoreArrayRead(vector, &vectorData);
3003}
3004
3010void SystemSolver::destroyVector(Vec *vector) const
3011{
3012 if (vector) {
3013 VecDestroy(vector);
3014 *vector = PETSC_NULLPTR;
3015 }
3016}
3017
3025std::string SystemSolver::getInfoFilePath(const std::string &directory, const std::string &name) const
3026{
3027 return getFilePath(directory, name + ".info");
3028}
3029
3037std::string SystemSolver::getDataFilePath(const std::string &directory, const std::string &name) const
3038{
3039 return getFilePath(directory, name + ".dat");
3040}
3041
3049std::string SystemSolver::getFilePath(const std::string &directory, const std::string &name) const
3050{
3051 std::string path = directory + "/" + name;
3052
3053 return path;
3054}
3055
3060{
3061 MatNullSpace nullspace;
3062#if BITPIT_ENABLE_MPI==1
3063 MatNullSpaceCreate(m_communicator, PETSC_TRUE, 0, NULL, &nullspace);
3064#else
3065 MatNullSpaceCreate(PETSC_COMM_SELF, PETSC_TRUE, 0, NULL, &nullspace);
3066#endif
3067 MatSetNullSpace(m_A, nullspace);
3068 MatNullSpaceDestroy(&nullspace);
3069}
3070
3075{
3076 MatSetNullSpace(m_A, NULL);
3077}
3078
3096void SystemSolver::setReordering(long nRows, long nCols, const SystemMatrixOrdering &reordering)
3097{
3098 // Clear existing reordering
3100
3101 // Early return if natural ordering is used
3102 try {
3103 dynamic_cast<const NaturalSystemMatrixOrdering &>(reordering);
3104 return;
3105 } catch(const std::bad_cast &exception) {
3106 BITPIT_UNUSED(exception);
3107
3108 // A reordering other than the natural one has been passed
3109 }
3110
3111 // Set row reordering
3112 PetscInt *rowReorderingStorage;
3113 PetscMalloc(nRows * sizeof(PetscInt), &rowReorderingStorage);
3114 for (long i = 0; i < nRows; ++i) {
3115 rowReorderingStorage[reordering.getRowPermutationRank(i)] = i;
3116 }
3117
3118#if BITPIT_ENABLE_MPI == 1
3119 ISCreateGeneral(m_communicator, nRows, rowReorderingStorage, PETSC_OWN_POINTER, &m_rowReordering);
3120#else
3121 ISCreateGeneral(PETSC_COMM_SELF, nRows, rowReorderingStorage, PETSC_OWN_POINTER, &m_rowReordering);
3122#endif
3123 ISSetPermutation(m_rowReordering);
3124
3125 // Create new permutations
3126 PetscInt *colReorderingStorage;
3127 PetscMalloc(nCols * sizeof(PetscInt), &colReorderingStorage);
3128 for (long j = 0; j < nCols; ++j) {
3129 colReorderingStorage[reordering.getColPermutationRank(j)] = j;
3130 }
3131
3132#if BITPIT_ENABLE_MPI == 1
3133 ISCreateGeneral(m_communicator, nCols, colReorderingStorage, PETSC_OWN_POINTER, &m_colReordering);
3134#else
3135 ISCreateGeneral(PETSC_COMM_SELF, nCols, colReorderingStorage, PETSC_OWN_POINTER, &m_colReordering);
3136#endif
3137 ISSetPermutation(m_colReordering);
3138}
3139
3147{
3148 if (m_rowReordering) {
3149 ISDestroy(&m_rowReordering);
3150 }
3151
3152 if (m_colReordering) {
3153 ISDestroy(&m_colReordering);
3154 }
3155}
3156
3163{
3164 // Nothing to do
3165}
3166
3171{
3172 // Check if the system is assembled
3173 if (!isAssembled()) {
3174 throw std::runtime_error("Unable to solve the system. The system is not yet assembled.");
3175 }
3176
3177 // Early return if the KSP can be reused
3178 if (!m_KSPDirty) {
3179 return;
3180 }
3181
3182 // Create KSP
3183 bool setupNeeded = false;
3184 if (!m_KSP) {
3185 createKSP();
3186 setupNeeded = true;
3187 }
3188
3189 // Set the matrix associated with the linear system
3190 KSPSetOperators(m_KSP, m_A, m_A);
3191
3192 // Set up
3193 if (setupNeeded) {
3194 // Perform actions before preconditioner set up
3196
3197 // Initialize preconditioner from options
3198 PC pc;
3199 KSPGetPC(m_KSP, &pc);
3200 PCSetFromOptions(pc);
3201
3202 // Set up preconditioner
3204
3205 // Perform actions after preconditioner set up
3207
3208 // Perform actions before Krylov subspace method set up set up
3210
3211 // Initialize Krylov subspace from options
3212 KSPSetFromOptions(m_KSP);
3213
3214 // Set up the Krylov subspace method
3215 setupKrylov();
3216
3217 // Perform actions after Krylov subspace method set up
3219 }
3220
3221 // KSP is now ready
3222 m_KSPDirty = false;
3223
3224 // Reset KSP status
3226}
3227
3232{
3233 // Nothing to do
3234}
3235
3240{
3241 // Create KSP object
3242#if BITPIT_ENABLE_MPI==1
3243 KSPCreate(m_communicator, &m_KSP);
3244#else
3245 KSPCreate(PETSC_COMM_SELF, &m_KSP);
3246#endif
3247
3248 // Set options prefix
3249 if (!m_prefix.empty()) {
3250 KSPSetOptionsPrefix(m_KSP, m_prefix.c_str());
3251 }
3252}
3253
3258{
3259 m_KSPDirty = true;
3260 KSPDestroy(&m_KSP);
3261 m_KSP = PETSC_NULLPTR;
3262}
3263
3268{
3269 PC pc;
3270 KSPGetPC(m_KSP, &pc);
3272}
3273
3280void SystemSolver::setupPreconditioner(PC pc, const KSPOptions &options) const
3281{
3282 // Set preconditioner type
3283 PCType pcType;
3284#if BITPIT_ENABLE_MPI == 1
3285 if (isPartitioned()) {
3286 pcType = PCASM;
3287 } else {
3288 pcType = PCILU;
3289 }
3290#else
3291 pcType = PCILU;
3292#endif
3293
3294 PCSetType(pc, pcType);
3295
3296 // Configure preconditioner
3297 if (strcmp(pcType, PCASM) == 0) {
3298 if (options.overlap != PETSC_DEFAULT) {
3299 PCASMSetOverlap(pc, options.overlap);
3300 }
3301 } else if (strcmp(pcType, PCILU) == 0) {
3302 if (options.levels != PETSC_DEFAULT) {
3303 PCFactorSetLevels(pc, options.levels);
3304 }
3305 }
3306
3307 PCSetUp(pc);
3308
3309 if (strcmp(pcType, PCASM) == 0) {
3310 KSP *subKSPs;
3311 PetscInt nSubKSPs;
3312 PCASMGetSubKSP(pc, &nSubKSPs, PETSC_NULLPTR, &subKSPs);
3313
3314 for (PetscInt i = 0; i < nSubKSPs; ++i) {
3315 KSP subKSP = subKSPs[i];
3316 KSPSetType(subKSP, KSPPREONLY);
3317
3318 PC subPC;
3319 KSPGetPC(subKSP, &subPC);
3320 PCSetType(subPC, PCILU);
3321 if (options.sublevels != PETSC_DEFAULT) {
3322 log::warning() << " Setting ASM ILU levels using the member \"sublevels\" is deprecated."
3323 << " ASM ILU levels should be set using the member \"levels\"." << std::endl;
3324
3325 PCFactorSetLevels(subPC, options.sublevels);
3326 } else if (options.levels != PETSC_DEFAULT) {
3327 PCFactorSetLevels(subPC, options.levels);
3328 }
3329
3330 if (options.subrtol != PETSC_DEFAULT) {
3331 log::warning() << " The member \"subrtol\" is deprecated. Since ASM is only use a a preconditioner"
3332 << " setting the tolerance of the Krylov subspace doens't have any effect on the"
3333 << " solution of the system."
3334 << std::endl;
3335 }
3336 }
3337 }
3338}
3339
3344{
3345 // Nothing to do
3346}
3347
3352{
3353 // Nothing to do
3354}
3355
3360{
3361 setupKrylov(m_KSP, getKSPOptions());
3362}
3363
3374void SystemSolver::setupKrylov(KSP ksp, const KSPOptions &options) const
3375{
3376 KSPSetType(ksp, KSPFGMRES);
3377 if (options.restart != PETSC_DEFAULT) {
3378 KSPGMRESSetRestart(ksp, options.restart);
3379 }
3380 if (options.rtol != PETSC_DEFAULT || options.atol != PETSC_DEFAULT || options.maxits != PETSC_DEFAULT) {
3381 KSPSetTolerances(ksp, options.rtol, options.atol, PETSC_DEFAULT, options.maxits);
3382 }
3383 KSPSetInitialGuessNonzero(ksp, options.initial_non_zero);
3384}
3385
3390{
3391 // Enable convergence monitor
3392 if (m_convergenceMonitorEnabled) {
3393#if (PETSC_VERSION_MAJOR >= 3 && PETSC_VERSION_MINOR >= 7)
3394 PetscOptionsSetValue(PETSC_NULLPTR, ("-" + m_prefix + "ksp_monitor_true_residual").c_str(), "");
3395 PetscOptionsSetValue(PETSC_NULLPTR, ("-" + m_prefix + "ksp_monitor_singular_value").c_str(), "");
3396 PetscOptionsSetValue(PETSC_NULLPTR, ("-" + m_prefix + "ksp_converged_reason").c_str(), "");
3397#else
3398 PetscOptionsSetValue(("-" + m_prefix + "ksp_monitor_true_residual").c_str(), "");
3399 PetscOptionsSetValue(("-" + m_prefix + "ksp_monitor_singular_value").c_str(), "");
3400 PetscOptionsSetValue(("-" + m_prefix + "ksp_converged_reason").c_str(), "");
3401#endif
3402 }
3403}
3404
3411
3420{
3421 if (!isAssembled()) {
3422 throw std::runtime_error("The options associated with the Krylov solver can only be accessed after assembling the system.");
3423 }
3424
3425 return m_KSPOptions;
3426}
3427
3436{
3437 if (!isAssembled()) {
3438 throw std::runtime_error("The options associated with the Krylov solver can only be accessed after assembling the system.");
3439 }
3440
3441 return m_KSPOptions;
3442}
3443
3448{
3449 resetKSPOptions(&m_KSPOptions);
3450}
3451
3458{
3459 *options = KSPOptions();
3460}
3461
3466{
3467 resetKSPOptions(&m_KSPOptions);
3468}
3469
3478{
3479 if (!isAssembled()) {
3480 throw std::runtime_error("The status of the Krylov solver can only be accessed after assembling the system.");
3481 }
3482
3483 return m_KSPStatus;
3484}
3485
3493
3498{
3499 fillKSPStatus(m_KSP, &m_KSPStatus);
3500}
3501
3508void SystemSolver::fillKSPStatus(KSP ksp, KSPStatus *status) const
3509{
3510 status->error = 0;
3511 KSPGetIterationNumber(ksp, &(status->its));
3512 KSPGetConvergedReason(ksp, &(status->convergence));
3513}
3514
3519{
3520 resetKSPStatus(&m_KSPStatus);
3521}
3522
3529{
3530 status->error = 0;
3531 status->its = -1;
3532 status->convergence = KSP_CONVERGED_ITERATING;
3533}
3534
3542
3543#if BITPIT_ENABLE_MPI==1
3549const MPI_Comm & SystemSolver::getCommunicator() const
3550{
3551 return m_communicator;
3552}
3553
3560void SystemSolver::setCommunicator(MPI_Comm communicator)
3561{
3562 if ((communicator != MPI_COMM_NULL) && (communicator != MPI_COMM_SELF)) {
3563 MPI_Comm_dup(communicator, &m_communicator);
3564 } else {
3565 m_communicator = MPI_COMM_SELF;
3566 }
3567}
3568
3572void SystemSolver::freeCommunicator()
3573{
3574 if (m_communicator != MPI_COMM_SELF) {
3575 int finalizedCalled;
3576 MPI_Finalized(&finalizedCalled);
3577 if (!finalizedCalled) {
3578 MPI_Comm_free(&m_communicator);
3579 }
3580 }
3581}
3582#endif
3583
3590void SystemSolver::removeNullSpaceFromRHS()
3591{
3592 // Get the null space from the matrix
3593 MatNullSpace nullspace;
3594 MatGetNullSpace(m_A, &nullspace);
3595
3596 // Remove null space components from right hand side
3597 MatNullSpaceRemove(nullspace, m_rhs);
3598}
3599
3607{
3608 return m_forceConsistency;
3609}
3610
3618{
3619 m_forceConsistency = enable;
3620}
3621
3634void SystemSolver::openBinaryArchive(const std::string &directory, const std::string &prefix,
3635 int block, IBinaryArchive *archive) const
3636{
3637 int expectedVersion = getDumpVersion();
3638
3639 std::string path = getFilePath(directory, prefix);
3640 archive->open(path, "dat", block);
3641 bool versionsMatch = archive->checkVersion(expectedVersion);
3642 if (!versionsMatch) {
3643 std::string message = "Version stored in binary archive " + path + " does not match";
3644 message += " the expected version " + std::to_string(expectedVersion) + ".";
3645
3646 throw std::runtime_error(message);
3647 }
3648}
3649
3663void SystemSolver::openBinaryArchive(const std::string &header, const std::string &directory,
3664 const std::string &prefix, int block, OBinaryArchive *archive) const
3665{
3666 int version = getDumpVersion();
3667 std::string path = getFilePath(directory, prefix);
3668 archive->open(path, "dat", version, header, block);
3669}
3670
3676void SystemSolver::closeBinaryArchive(BinaryArchive *archive) const
3677{
3678 archive->close();
3679}
3680
3686int SystemSolver::getBinaryArchiveBlock() const
3687{
3688#if BITPIT_ENABLE_MPI == 1
3689 if (isPartitioned()) {
3690 int nProcesses;
3691 MPI_Comm_size(getCommunicator(), &nProcesses);
3692 if (nProcesses <= 1) {
3693 return -1;
3694 }
3695
3696 int archiveBlock;
3697 MPI_Comm_rank(getCommunicator(), &archiveBlock);
3698
3699 return archiveBlock;
3700 } else {
3701 return -1;
3702 }
3703#else
3704 return -1;
3705#endif
3706}
3707
3708}
Input binary archive.
bool checkVersion(int version)
void open(const std::string &name, int block=-1)
std::istream & getStream()
The NaturalSystemMatrixOrdering class defines allows to use a matrix natural ordering.
long getColPermutationRank(long col) const override
long getRowPermutationRank(long row) const override
Output binary archive.
std::ostream & getStream()
The PetscManager class handles the interaction with PETSc library.
void addInitOptions(int argc, char **args)
void addInitOption(const std::string &option)
Metafunction for generating a list of elements that can be either stored in an external vectror or,...
__PXV_CONST_ITERATOR__ cend()
std::size_t size() const
__PXV_POINTER__ data() noexcept
__PXV_CONST_ITERATOR__ cbegin()
void set(__PXV_POINTER__ data, std::size_t size)
const MPI_Comm & getCommunicator() const
long getColGlobalCount() const
ConstProxyVector< double > getRowValues(long row) const
long getColGlobalOffset() const
long getMaxRowNZCount() const
long getRowGlobalCount() const
ConstProxyVector< long > getRowPattern(long row) const
long getRowNZCount(long row) const
long getRowGlobalOffset() const
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 ...
virtual long getColPermutationRank(long col) const =0
virtual long getRowPermutationRank(long row) const =0
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 matrixUpdate(long nRows, const long *rows, const Assembler &assembler)
virtual int getDumpVersion() const
virtual void exportMatrix(const std::string &filePath, FileFormat exportFormat=FILE_BINARY) const
virtual void matrixFill(const std::string &filePath)
virtual void prePreconditionerSetupActions()
void restoreSolutionRawPtr(double *raw_solution)
virtual void matrixRestore(std::istream &systemStream, const std::string &directory, const std::string &prefix, bool redistribute)
void fillMatrix(Mat matrix, const std::string &filePath) const
void restoreRHSRawReadPtr(const double *raw_rhs) const
void dumpVector(Vec vector, const std::string &directory, const std::string &name) const
virtual void postPreconditionerSetupActions()
virtual void importMatrix(const std::string &filePath)
void reorderVector(Vec vector, IS permutations, bool invert) const
std::string getInfoFilePath(const std::string &directory, const std::string &name) const
virtual void vectorsDump(std::ostream &systemStream, const std::string &directory, const std::string &prefix) const
virtual void resetKSPOptions(KSPOptions *options) const
void exportVector(Vec vector, const std::string &filePath, FileFormat fileFormat) const
void dumpSystem(const std::string &header, const std::string &directory, const std::string &prefix="") const
void setReordering(long nRows, long nCols, const SystemMatrixOrdering &reordering)
virtual void restoreInfo(std::istream &stream)
void restoreSystem(MPI_Comm communicator, bool redistribute, const std::string &directory, const std::string &prefix="")
const double * getSolutionRawReadPtr() const
virtual void exportRHS(const std::string &filePath, FileFormat exportFormat=FILE_BINARY) const
virtual int getBlockSize() const
void restoreRHSRawPtr(double *raw_rhs)
void setTranspose(bool transpose)
const KSPStatus & getKSPStatus() const
void restoreMatrix(const std::string &directory, const std::string &name, bool redistribute, Mat *matrix) const
void matrixAssembly(const Assembler &assembler)
void createVector(int blockSize, Vec *vector) const
void destroyMatrix(Mat *matrix) const
std::string getFilePath(const std::string &directory, const std::string &name) const
virtual void matrixDump(std::ostream &systemStream, const std::string &directory, const std::string &prefix) const
virtual void dumpInfo(std::ostream &stream) const
void assembly(const SparseMatrix &matrix)
const MPI_Comm & getCommunicator() const
std::string getDataFilePath(const std::string &directory, const std::string &name) const
const double * getRHSRawReadPtr() const
virtual void importSolution(const std::string &filePath)
virtual void importRHS(const std::string &filePath)
void destroyVector(Vec *vector) const
void createMatrix(int rowBlockSize, int colBlockSize, Mat *matrix) const
void fillVector(Vec vector, const std::string &filePath) const
virtual void vectorsFill(const std::vector< double > &rhs, const std::vector< double > &solution)
void enableForceConsistency(bool enable)
void restoreSolutionRawReadPtr(const double *raw_solution) const
virtual void vectorsReorder(bool invert)
virtual void vectorsRestore(std::istream &systemStream, const std::string &directory, const std::string &prefix)
void restoreVector(const std::string &directory, const std::string &name, Mat matrix, VectorSide side, Vec *vector) const
virtual void exportSolution(const std::string &filePath, FileFormat exportFormat=FILE_BINARY) const
The SystemSparseMatrixAssembler class defines an assembler for building the system matrix form a spar...
void getRowValues(long rowIndex, ConstProxyVector< double > *values) const override
const MPI_Comm & getCommunicator() const override
void getRowData(long rowIndex, ConstProxyVector< long > *pattern, ConstProxyVector< double > *values) const override
long getRowNZCount(long rowIndex) const override
AssemblyOptions getOptions() const override
void getRowPattern(long rowIndex, ConstProxyVector< long > *pattern) const override
SystemSparseMatrixAssembler(const SparseMatrix *matrix)
void write(std::ostream &stream, const std::vector< bool > &container)
void read(std::istream &stream, std::vector< bool > &container)
#define BITPIT_UNUSED(variable)
Definition compiler.hpp:63
Logger & cout(log::Level defaultSeverity, log::Visibility defaultVisibility)
Definition logger.cpp:1705
Logger & warning(log::Visibility defaultVisibility)
Definition logger.cpp:1812
PetscScalar subrtol
Deprecated, ASM ILU levels should be set using the "levels" member.
PetscInt levels
Overlap between a pair of subdomains for the partitioned preconditioner.
PetscInt restart
Tells the iterative solver that the initial guess is nonzero.
PetscInt maxits
Number of iterations at which the GMRES method restarts.
PetscScalar atol
Relative convergence tolerance, relative decrease in the preconditioned residual norm.
PetscInt sublevels
Absolute convergence tolerance, absolute size of the preconditioned residual norm.
PetscScalar rtol
Maximum number of iterations.
PetscBool initial_non_zero
Number of levels of fill used by the preconditioner.
bool full
Controls if the assembler is providing all the non-zero values of a row.
--- layout: doxygen_footer ---