Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
serena.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/ConjugateGradientSolver.h"
42#include "arccore/alina/StaticMatrix.h"
43#include "arccore/alina/Adapters.h"
44
45#include "arccore/alina/IO.h"
46#include "arccore/alina/Profiler.h"
47
48using namespace Arcane;
49
50int main(int argc, char *argv[]) {
51 // The command line should contain the matrix file name:
52 if (argc < 2) {
53 std::cerr << "Usage: " << argv[0] << " <matrix.mtx>" << std::endl;
54 return 1;
55 }
56
57 auto& prof = Alina::Profiler::globalProfiler();
58
59 // Read the system matrix:
60 ptrdiff_t rows, cols;
61 std::vector<ptrdiff_t> ptr, col;
62 std::vector<double> val;
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 prof.toc("read");
68
69 // The RHS is filled with ones:
70 std::vector<double> f(rows, 1.0);
71
72 // Scale the matrix so that it has the unit diagonal.
73 // First, find the diagonal values:
74 std::vector<double> D(rows, 1.0);
75 for(ptrdiff_t i = 0; i < rows; ++i) {
76 for(ptrdiff_t j = ptr[i], e = ptr[i+1]; j < e; ++j) {
77 if (col[j] == i) {
78 D[i] = 1 / sqrt(val[j]);
79 break;
80 }
81 }
82 }
83
84 // Then, apply the scaling in-place:
85 for(ptrdiff_t i = 0; i < rows; ++i) {
86 for(ptrdiff_t j = ptr[i], e = ptr[i+1]; j < e; ++j) {
87 val[j] *= D[i] * D[col[j]];
88 }
89 f[i] *= D[i];
90 }
91
92 // We use the tuple of CRS arrays to represent the system matrix.
93 // Note that std::tie creates a tuple of references, so no data is actually
94 // copied here:
95 auto A = std::tie(rows, ptr, col, val);
96
97 // Compose the solver type
98 typedef Alina::StaticMatrix<double, 3, 3> dmat_type; // matrix value type in double precision
99 typedef Alina::StaticMatrix<double, 3, 1> dvec_type; // the corresponding vector value type
100 typedef Alina::StaticMatrix<float, 3, 3> smat_type; // matrix value type in single precision
101
102 typedef Alina::BuiltinBackend<dmat_type> SBackend; // the solver backend
103 typedef Alina::BuiltinBackend<smat_type> PBackend; // the preconditioner backend
104
107 PBackend,
110 >,
112 > Solver;
113
114 // Solver parameters
115 Solver::params prm;
116 prm.solver.maxiter = 500;
117
118 // Initialize the solver with the system matrix.
119 // Use the block_matrix adapter to convert the matrix into
120 // the block format on the fly:
121 prof.tic("setup");
122 auto Ab = Alina::adapter::block_matrix<dmat_type>(A);
123 Solver solve(Ab, prm);
124 prof.toc("setup");
125
126 // Show the mini-report on the constructed solver:
127 std::cout << solve << std::endl;
128
129 // Solve the system with the zero initial approximation:
130 std::vector<double> x(rows, 0.0);
131
132 // Reinterpret both the RHS and the solution vectors as block-valued:
133 auto f_ptr = reinterpret_cast<dvec_type*>(f.data());
134 auto x_ptr = reinterpret_cast<dvec_type*>(x.data());
135 auto F = SmallSpan<dvec_type>(f_ptr, rows / 3);
136 auto X = SmallSpan<dvec_type>(x_ptr, rows / 3);
137
138 prof.tic("solve");
139 Alina::SolverResult r = solve(Ab, F, X);
140 prof.toc("solve");
141
142 // Output the number of iterations, the relative error,
143 // and the profiling data:
144 std::cout << "Iters: " << r.nbIteration() << std::endl
145 << "Error: " << r.residual() << std::endl
146 << prof << std::endl;
147}
Algebraic multigrid method.
Definition AMG.h:71
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
Vue d'un tableau d'éléments de type T.
Definition Span.h:801
__host__ __device__ double sqrt(double v)
Racine carrée de v.
Definition Math.h:135
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
Sparse approximate interface smoother.
Smoothed aggregation coarsening.