Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
AdapterCSRBuilder.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/*---------------------------------------------------------------------------*/
9/*
10 * This file is based on the work on AMGCL library (version march 2026)
11 * which can be found at https://github.com/ddemidov/amgcl.
12 *
13 * Copyright (c) 2012-2022 Denis Demidov <dennis.demidov@gmail.com>
14 * SPDX-License-Identifier: MIT
15 */
16/*---------------------------------------------------------------------------*/
17/*---------------------------------------------------------------------------*/
18
19#include <iostream>
20#include <vector>
21#include <algorithm>
22
23#include "arccore/alina/AMG.h"
24#include "arccore/alina/PreconditionedSolver.h"
25#include "arccore/alina/BuiltinBackend.h"
26#include "arccore/alina/Adapters.h"
27#include "arccore/alina/Coarsening.h"
28#include "arccore/alina/Relaxation.h"
29#include "arccore/alina/ConjugateGradientSolver.h"
30#include "arccore/alina/Profiler.h"
31
32#include "SampleProblemCommon.h"
33
34using namespace Arcane;
35using namespace Arcane::Alina;
36
37//---------------------------------------------------------------------------
38struct poisson_2d
39{
40 typedef double val_type;
41 typedef ptrdiff_t col_type;
42
43 size_t n;
44 double h2i;
45
46 poisson_2d(size_t n)
47 : n(n)
48 , h2i((n - 1) * (n - 1))
49 {}
50
51 size_t rows() const { return n * n; }
52 size_t nonzeros() const { return 5 * rows(); }
53
54 void operator()(size_t row,
55 std::vector<col_type>& col,
56 std::vector<val_type>& val) const
57 {
58 size_t i = row % n;
59 size_t j = row / n;
60
61 if (j > 0) {
62 col.push_back(row - n);
63 val.push_back(-h2i);
64 }
65
66 if (i > 0) {
67 col.push_back(row - 1);
68 val.push_back(-h2i);
69 }
70
71 col.push_back(row);
72 val.push_back(4 * h2i);
73
74 if (i + 1 < n) {
75 col.push_back(row + 1);
76 val.push_back(-h2i);
77 }
78
79 if (j + 1 < n) {
80 col.push_back(row + n);
81 val.push_back(-h2i);
82 }
83 }
84};
85
86//---------------------------------------------------------------------------
87
88template <class Vec>
89double norm(const Vec& v)
90{
91 return sqrt(Arcane::Alina::backend::inner_product(v, v));
92}
93
94//---------------------------------------------------------------------------
95int main(int argc, char* argv[])
96{
97 auto& prof = Alina::Profiler::globalProfiler();
98 int m = argc > 1 ? atoi(argv[1]) : 1024;
99 int n = m * m;
100
101 // Create iterative solver preconditioned by AMG.
102 // The use of make_matrix() from crs_builder.hpp allows to construct the
103 // system matrix on demand row by row.
104 prof.tic("build");
110
111 Solver solve(Alina::adapter::make_matrix(poisson_2d(m)));
112 prof.toc("build");
113
114 std::cout << solve.precond() << std::endl;
115
116 std::vector<double> f(n, 1);
117 std::vector<double> x(n, 0);
118
119 prof.tic("solve");
120 SolverResult r = solve(f, x);
121 prof.toc("solve");
122
123 std::cout << "Solver:" << std::endl
124 << " Iterations: " << r.nbIteration() << std::endl
125 << " Error: " << r.residual() << std::endl
126 << std::endl;
127
128 // Use the constructed solver as a preconditioner for another iterative
129 // solver.
130 //
131 // Iterative methods use estimated residual for exit condition. For some
132 // problems the value of estimated residual can get too far from true
133 // residual due to round-off errors.
134 //
135 // Nesting iterative solvers in this way allows to shave last bits off the
136 // error.
138 std::fill(x.begin(), x.end(), 0);
139
140 prof.tic("nested solver");
141 r = S(solve.system_matrix(), solve, f, x);
142 prof.toc("nested solver");
143
144 std::cout << "Nested solver:" << std::endl
145 << " Iterations: " << r.nbIteration() << std::endl
146 << " Error: " << r.residual() << std::endl
147 << std::endl;
148
149 std::cout << prof << std::endl;
150}
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 -*-
Gauss-Seidel relaxation.
Definition Relaxation.h:510
Smoothed aggregation coarsening.