Loading...
Searching...
No Matches
LA_example_00001.cpp
1/*---------------------------------------------------------------------------*\
2 *
3 * bitpit
4 *
5 * Copyright (C) 2015-2023 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 <array>
26#if BITPIT_ENABLE_MPI==1
27#include <mpi.h>
28#endif
29
30#include "bitpit_common.hpp"
31#include "bitpit_LA.hpp"
32
33using namespace bitpit;
34
74
75public:
83 AMGSystemSolver(bool transpose, bool multigrid, bool debug)
84 : SystemSolver(multigrid, transpose, debug),
85 m_multigrid(multigrid)
86 {
87 }
88
89protected:
90 const static PetscInt GAMG_MAX_LEVELS;
91 const static PetscReal GAMG_THRESHOLD;
92 const static PetscReal GAMG_THRESHOLD_SCALE;
93
94 bool m_multigrid;
95
102 void setupPreconditioner(PC pc, const KSPOptions &options) const override
103 {
104 // Early return for non-multigrid
105 if (!m_multigrid) {
107 return;
108 }
109
110 // Set multigrid options
111 KSPGetPC(this->m_KSP, &pc);
112 PCSetType(pc, PCGAMG);
113
114 std::vector<PetscReal> thresholds(GAMG_MAX_LEVELS);
115 for (int i = 0; i < GAMG_MAX_LEVELS; ++i) {
116 thresholds[i] = GAMG_THRESHOLD;
117 }
118
119 PCGAMGSetNSmooths(pc, 0);
120#if (PETSC_VERSION_MAJOR >= 3 && PETSC_VERSION_MINOR <= 17)
121#if BITPIT_ENABLE_MPI == 1
122 if (this->isPartitioned()) {
123 PCGAMGSetSymGraph(pc, PETSC_TRUE);
124 }
125#endif
126#endif
127
128 PCGAMGSetNlevels(pc, GAMG_MAX_LEVELS);
129#if (PETSC_VERSION_MAJOR >= 3 && PETSC_VERSION_MINOR >= 16)
130 PCGAMGSetThreshold(pc, thresholds.data(), GAMG_MAX_LEVELS);
131 PCGAMGSetThresholdScale(pc, GAMG_THRESHOLD_SCALE);
132#else
133 PCGAMGSetThreshold(pc, GAMG_THRESHOLD);
134#endif
135
136 // Setup
137 PCSetUp(pc);
138
139 // Set options for the multigrid coarse level
140 KSP coarseKSP;
141 PCMGGetCoarseSolve(pc, &coarseKSP);
142
143 PC coarsePreconditioner;
144 KSPGetPC(coarseKSP, &coarsePreconditioner);
145 PCSetType(coarsePreconditioner, PCSVD);
146 }
147};
148
149const PetscInt AMGSystemSolver::GAMG_MAX_LEVELS = 10;
150const PetscReal AMGSystemSolver::GAMG_THRESHOLD = 0.02;
151const PetscReal AMGSystemSolver::GAMG_THRESHOLD_SCALE = 1.0;
152
159int run_dump(int rank, int nProcs)
160{
161 log::cout() << std::endl;
162 log::cout() << "Running dump example..." << std::endl;
163
164 int nRows;
165 int nCols;
166 int nNZ;
167 if (nProcs == 1) {
168 nRows = 10;
169 nCols = 10;
170 nNZ = 10;
171 } else if (nProcs == 2) {
172 if (rank == 0) {
173 nRows = 6;
174 nCols = 6;
175 nNZ = 6;
176 } if (rank == 1) {
177 nRows = 4;
178 nCols = 4;
179 nNZ = 4;
180 }
181 } else {
182 if (rank == 0) {
183 nRows = 2;
184 nCols = 2;
185 nNZ = 2;
186 } if (rank == 1) {
187 nRows = 5;
188 nCols = 5;
189 nNZ = 5;
190 } if (rank == 2) {
191 nRows = 3;
192 nCols = 3;
193 nNZ = 3;
194 } else {
195 nRows = 0;
196 nCols = 0;
197 nNZ = 0;
198 }
199 }
200
201 int blockSize = 5;
202
203 // Build matrix
204 log::cout() << "Building matrix..." << std::endl;
205
206 std::vector<long> rowPattern(1);
207 std::vector<double> rowValues(blockSize * blockSize);
208
209#if BITPIT_ENABLE_MPI==1
210 SparseMatrix matrix(MPI_COMM_WORLD, true, blockSize, nRows, nCols, nNZ);
211
212 int rowOffset = matrix.getRowGlobalOffset();
213 int colOffset = matrix.getColGlobalOffset();
214#else
215 SparseMatrix matrix(blockSize, nRows, nCols, nNZ);
216
217 int rowOffset = 0;
218 int colOffset = 0;
219#endif
220
221 for (int i = 0; i < nRows; ++i) {
222 rowPattern[0] = colOffset + i;
223 for (int k = 0; k < blockSize; ++k) {
224 int kk = linearalgebra::linearIndexRowMajor(k, k, blockSize, blockSize);
225
226 rowValues[kk] = 1. / (double) (blockSize * (rowOffset + i) + k + 1);
227 }
228
229 matrix.addRow(rowPattern, rowValues);
230 }
231
232 matrix.assembly();
233
234 // Build solver
235 log::cout() << "Building solver..." << std::endl;
236
237 bool multigrid = true;
238 bool debug = false;
239 bool transpose = false;
240 AMGSystemSolver solver(transpose, multigrid, debug);
241
242 solver.assembly(matrix);
243
244 double *rhs = solver.getRHSRawPtr();
245 for (int i = 0; i < blockSize * nCols; ++i) {
246 rhs[i] = 1.;
247 }
248 solver.restoreRHSRawPtr(rhs);
249
250 double *initialSolution = solver.getSolutionRawPtr();
251 for (int i = 0; i < nRows; ++i) {
252 initialSolution[i] = 0;
253 }
254 solver.restoreSolutionRawPtr(initialSolution);
255
256 // Export system
257 solver.exportMatrix("LA_example_0001_matrix.txt", SystemSolver::FILE_ASCII);
258 solver.exportRHS("LA_example_0001_rhs.txt", SystemSolver::FILE_ASCII);
259 solver.exportSolution("LA_example_0001_solution.txt", SystemSolver::FILE_ASCII);
260
261 // Dump system
262 log::cout() << "Dumping initialized linear system..." << std::endl;
263
264 solver.dumpSystem("bitpit linear system", ".", "LA_example_0001_linear_system_");
265
266 // Set KSP
267 KSPOptions &options = solver.getKSPOptions();
268 options.atol = 1e-50;
269 options.rtol = 1e-7;
270 options.maxits = 100;
271 options.initial_non_zero = PETSC_FALSE;
272
273 // Solve System
274 log::cout() << "Solve linear system..." << std::endl;
275 solver.solve();
276
277 const KSPStatus &status = solver.getKSPStatus();
278 log::cout() << "Reason: " << status.convergence << " in its: " << status.its << std::endl;
279 log::cout() << "Linear system solved." << std::endl;
280
281 // Export solution
282 log::cout() << "Export solution..." << std::endl;
283 solver.exportSolution("LA_example_0001_linear_system_solution.dat", bitpit::SystemSolver::FILE_BINARY);
284
285 // Compare computed solution with expected one
286 //
287 // Set comparison tolerance to the linear system relative tolerance
288 log::cout() << "Comparing solutions..." << std::endl;
289
290 const double TOLERANCE = options.rtol;
291 log::cout() << " Tolerance = " << options.rtol << std::endl;
292
293 std::size_t nRowElements = solver.getRowElementCount();
294
295 bool isSolutionCorrect = true;
296 const double *solution = solver.getSolutionRawReadPtr();
297 for (std::size_t i = 0; i < nRowElements; ++i) {
298 std::size_t globalDOF = rowOffset * blockSize + i;
299
300 // Check for absolute tolerance only
301 double expectedSolution = static_cast<double>(globalDOF + 1);
302 if (!utils::DoubleFloatingEqual()(expectedSolution - solution[i], 0.0, TOLERANCE, TOLERANCE)) {
303 std::stringstream message;
304 message << "Solutions do not match for global DOF #" << globalDOF << ".";
305 message << " Expected solution is " << expectedSolution;
306 message << ", current solution is " << solution[i] << ".";
307 log::cout() << message.str() << std::endl;
308 isSolutionCorrect = false;
309 }
310 }
311 solver.restoreSolutionRawReadPtr(solution);
312 if (isSolutionCorrect) {
313 log::cout() << "Solutions match." << std::endl;
314 }
315
316 return 0;
317}
318
325int run_restore(int rank, int nProcs)
326{
327 BITPIT_UNUSED(rank);
328 BITPIT_UNUSED(nProcs);
329
330 log::cout() << std::endl;
331 log::cout() << "Running restore example..." << std::endl;
332
333 // Restore solver with initial guess
334 log::cout() << "Restoring linear system ..." << std::endl;
335 bool multigrid = true;
336 bool debug = false;
337 bool transpose = false;
338 AMGSystemSolver solver(transpose, multigrid, debug);
339
340#if BITPIT_ENABLE_MPI==1
341 solver.restoreSystem(MPI_COMM_WORLD, false, ".", "LA_example_0001_linear_system_");
342#else
343 solver.restoreSystem(".", "LA_example_0001_linear_system_");
344#endif
345 log::cout() << "Linear system restored." << std::endl;
346
347 // Export system
348 solver.exportMatrix("LA_example_0001_restored_matrix.txt", SystemSolver::FILE_ASCII);
349 solver.exportRHS("LA_example_0001_restored_rhs.txt", SystemSolver::FILE_ASCII);
350 solver.exportSolution("LA_example_0001_restored_solution.txt", SystemSolver::FILE_ASCII);
351
352 // Set KSP options
353 KSPOptions &options = solver.getKSPOptions();
354 options.atol = 1e-50;
355 options.rtol = 1e-7;
356 options.maxits = 100;
357 options.initial_non_zero = PETSC_FALSE;
358
359 // Solve linear system
360 log::cout() << "Solving linear system..." << std::endl;
361 solver.solve();
362
363 const KSPStatus &status = solver.getKSPStatus();
364 log::cout() << "Reason: " << status.convergence << " in its: " << status.its << std::endl;
365 log::cout() << "Linear system solved." << std::endl;
366
367 // Turn solution comparison on
368 bool compareSolutions = true;
369 if (compareSolutions) {
370 // Store expected solution for comparison
371 log::cout() << "Storing computed solution for comparison..." << std::endl;
372 std::size_t nRowElements = solver.getRowElementCount();
373 std::vector<double> expectedSolution(nRowElements, 0.0);
374 double *computedSolution = solver.getSolutionRawPtr();
375 for (std::size_t i = 0; i < nRowElements; ++i) {
376 expectedSolution[i] = computedSolution[i];
377 }
378 solver.restoreSolutionRawPtr(computedSolution);
379
380 // Import solution with solution from outside
381 log::cout() << "Restoring solution..." << std::endl;
382 solver.importSolution("LA_example_0001_linear_system_solution.dat");
383 log::cout() << "Solution restored." << std::endl;
384
385 // Compare computed solution with restored one
386 //
387 // Set comparison tolerance to the linear system relative tolerance
388 log::cout() << "Comparing solutions..." << std::endl;
389
390 const double TOLERANCE = options.rtol;
391 log::cout() << " Tolerance = " << options.rtol << std::endl;
392
393 bool isSolutionCorrect = true;
394 const double *solution = solver.getSolutionRawReadPtr();
395 for (std::size_t i = 0; i < nRowElements; ++i) {
396 if (!utils::DoubleFloatingEqual()(expectedSolution[i] - solution[i], 0.0, TOLERANCE, TOLERANCE)) {
397 std::stringstream message;
398 message << "Solutions do not match for local DOF #" << i << ".";
399 message << " Expected solution is " << expectedSolution[i];
400 message << ", current solution is " << solution[i] << ".";
401 log::cout() << message.str() << std::endl;
402 isSolutionCorrect = false;
403 }
404 }
405 solver.restoreSolutionRawReadPtr(solution);
406 if (isSolutionCorrect) {
407 log::cout() << "Solutions match." << std::endl;
408 }
409 }
410
411 return 0;
412}
413
417int main(int argc, char *argv[])
418{
419#if BITPIT_ENABLE_MPI==1
420 MPI_Init(&argc, &argv);
421#else
422 BITPIT_UNUSED(argc);
423 BITPIT_UNUSED(argv);
424#endif
425
426 // Initialize the logger
427 int nProcs;
428 int rank;
429#if BITPIT_ENABLE_MPI==1
430 MPI_Comm_size(MPI_COMM_WORLD, &nProcs);
431 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
432#else
433 rank = 0;
434 nProcs = 1;
435#endif
436 log::manager().initialize(log::COMBINED, true, nProcs, rank);
437 log::cout().setDefaultVisibility(log::GLOBAL);
438
439 // Run the subtests
440 log::cout() << "Build, dump, solve and restore linear system" << std::endl;
441
442 int status;
443 try {
444 status = run_dump(rank, nProcs);
445 if (status != 0) {
446 return status;
447 }
448
449 status = run_restore(rank, nProcs);
450 if (status != 0) {
451 return status;
452 }
453 } catch (const std::exception &exception) {
454 log::cout() << exception.what();
455 exit(1);
456 }
457
458#if BITPIT_ENABLE_MPI==1
459 MPI_Finalize();
460#endif
461}
The class AMGSystemSolver is derived form SystemSolver and it allows to solve the linear system using...
AMGSystemSolver(bool transpose, bool multigrid, bool debug)
void setupPreconditioner(PC pc, const KSPOptions &options) const override
void initialize(log::Mode mode, bool reset, int nProcesses, int rank)
Definition logger.cpp:1268
void setDefaultVisibility(log::Visibility visibility)
Definition logger.cpp:579
The SystemSolver class provides methods for building and solving large linear systems.
#define BITPIT_UNUSED(variable)
Definition compiler.hpp:63
int linearIndexRowMajor(int row, int col, int nRows, int nCols)
void transpose(std::vector< std::vector< T > > &, std::vector< std::vector< T > > &)
Logger & cout(log::Level defaultSeverity, log::Visibility defaultVisibility)
Definition logger.cpp:1705
Logger & debug(log::Visibility defaultVisibility)
Definition logger.cpp:1882
LoggerManager & manager()
Definition logger.cpp:1685
PetscInt maxits
Number of iterations at which the GMRES method restarts.
PetscScalar atol
Relative convergence tolerance, relative decrease in the preconditioned residual norm.
PetscScalar rtol
Maximum number of iterations.
PetscBool initial_non_zero
Number of levels of fill used by the preconditioner.
--- layout: doxygen_footer ---