Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
poisson3Db.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2026 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7
8/*
9The MIT License
10
11Copyright (c) 2012-2022 Denis Demidov <dennis.demidov@gmail.com>
12
13Permission is hereby granted, free of charge, to any person obtaining a copy
14of this software and associated documentation files (the "Software"), to deal
15in the Software without restriction, including without limitation the rights
16to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
17copies of the Software, and to permit persons to whom the Software is
18furnished to do so, subject to the following conditions:
19
20The above copyright notice and this permission notice shall be included in
21all copies or substantial portions of the Software.
22
23THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
24IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
26AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
28OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
29THE SOFTWARE.
30*/
31
32#include <vector>
33#include <iostream>
34
35#include "arccore/alina/BuiltinBackend.h"
36#include "arccore/alina/Adapters.h"
37#include "arccore/alina/PreconditionedSolver.h"
38#include "arccore/alina/AMG.h"
39#include "arccore/alina/Coarsening.h"
40#include "arccore/alina/Relaxation.h"
41#include "arccore/alina/BiCGStabSolver.h"
42
43#include "arccore/alina/IO.h"
44#include "arccore/alina/Profiler.h"
45
46using namespace Arcane;
47
48int main(int argc, char* argv[])
49{
50 // The matrix and the RHS file names should be in the command line options:
51 if (argc < 3) {
52 std::cerr << "Usage: " << argv[0] << " <matrix.mtx> <rhs.mtx>" << std::endl;
53 return 1;
54 }
55
56 // The profiler:
57 auto& prof = Alina::Profiler::globalProfiler();
58
59 // Read the system matrix and the RHS:
60 ptrdiff_t rows, cols;
61 std::vector<ptrdiff_t> ptr, col;
62 std::vector<double> val, rhs;
63
64 prof.tic("read");
65 std::tie(rows, cols) = Alina::IO::mm_reader(argv[1])(ptr, col, val);
66 std::cout << "Matrix " << argv[1] << ": " << rows << "x" << cols << std::endl;
67
68 std::tie(rows, cols) = Alina::IO::mm_reader(argv[2])(rhs);
69 std::cout << "RHS " << argv[2] << ": " << rows << "x" << cols << std::endl;
70 prof.toc("read");
71
72 // We use the tuple of CRS arrays to represent the system matrix.
73 // Note that std::tie creates a tuple of references, so no data is actually
74 // copied here:
75 auto A = std::tie(rows, ptr, col, val);
76
77 // Compose the solver type
78 // the solver backend:
79 typedef Alina::BuiltinBackend<double> SBackend;
80 // the preconditioner backend:
81#ifdef MIXED_PRECISION
82 typedef Alina::BuiltinBackend<float> PBackend;
83#else
84 typedef Alina::BuiltinBackend<double> PBackend;
85#endif
86
88 PBackend,
92
93 // Initialize the solver with the system matrix:
94 prof.tic("setup");
95 Solver solve(A);
96 prof.toc("setup");
97
98 // Show the mini-report on the constructed solver:
99 std::cout << solve << std::endl;
100
101 // Solve the system with the zero initial approximation:
102 std::vector<double> x(rows, 0.0);
103
104 prof.tic("solve");
105 Alina::SolverResult r = solve(A, rhs, x);
106 prof.toc("solve");
107
108 // Output the number of iterations, the relative error,
109 // and the profiling data:
110 std::cout << "Iters: " << r.nbIteration() << std::endl
111 << "Error: " << r.residual() << std::endl
112 << prof << std::endl;
113}
Algebraic multigrid method.
Definition AMG.h:71
BiConjugate Gradient Stabilized (BiCGSTAB) method.
Matrix market reader.
Definition IO.h:52
Convenience class that bundles together a preconditioner and an iterative solver.
Result of a solving.
Definition AlinaUtils.h:52
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
Sparse approximate interface smoother.
Smoothed aggregation coarsening.