Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
nullspace_hybrid.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/HybridBuiltinBackend.h"
36#include "arccore/alina/StaticMatrix.h"
37#include "arccore/alina/Adapters.h"
38#include "arccore/alina/PreconditionedSolver.h"
39#include "arccore/alina/AMG.h"
40#include "arccore/alina/Coarsening.h"
41#include "arccore/alina/Relaxation.h"
42#include "arccore/alina/ConjugateGradientSolver.h"
43
44#include "arccore/alina/IO.h"
45#include "arccore/alina/Profiler.h"
46
47using namespace Arcane;
48
49int main(int argc, char* argv[])
50{
51 // The command line should contain the matrix, the RHS, and the coordinate files:
52 if (argc < 4) {
53 std::cerr << "Usage: " << argv[0] << " <A.mtx> <b.mtx> <coo.mtx>" << std::endl;
54 return 1;
55 }
56
57 // The profiler:
58 auto& prof = Alina::Profiler::globalProfiler();
59
60 // Read the system matrix, the RHS, and the coordinates:
61 ptrdiff_t rows, cols, ndim, ncoo;
62 std::vector<ptrdiff_t> ptr, col;
63 std::vector<double> val, rhs, coo;
64
65 prof.tic("read");
66 std::tie(rows, rows) = Alina::IO::mm_reader(argv[1])(ptr, col, val);
67 std::tie(rows, cols) = Alina::IO::mm_reader(argv[2])(rhs);
68 std::tie(ncoo, ndim) = Alina::IO::mm_reader(argv[3])(coo);
69 prof.toc("read");
70
71 Alina::precondition(ncoo * ndim == rows && (ndim == 2 || ndim == 3),
72 "The coordinate file has wrong dimensions");
73
74 std::cout << "Matrix " << argv[1] << ": " << rows << "x" << rows << std::endl;
75 std::cout << "RHS " << argv[2] << ": " << rows << "x" << cols << std::endl;
76 std::cout << "Coords " << argv[3] << ": " << ncoo << "x" << ndim << std::endl;
77
78 // Declare the solver type
81 typedef Alina::backend::HybridBuiltinBackend<DBlock> SBackend; // the solver backend
82 typedef Alina::backend::HybridBuiltinBackend<FBlock> PBackend; // the preconditioner backend
83
86 Solver;
87
88 // Solver parameters:
89 Solver::params prm;
90 prm.precond.coarsening.aggr.eps_strong = 0;
91
92 // Convert the coordinates to the rigid body modes.
93 // The function returns the number of near null-space vectors
94 // (3 in 2D case, 6 in 3D case) and writes the vectors to the
95 // std::vector<double> specified as the last argument:
96 prm.precond.coarsening.nullspace.cols = Alina::rigid_body_modes(ndim, coo, prm.precond.coarsening.nullspace.B);
97
98 // We use the tuple of CRS arrays to represent the system matrix.
99 auto A = std::tie(rows, ptr, col, val);
100
101 // Initialize the solver with the system matrix.
102 prof.tic("setup");
103 Solver solve(A, prm);
104 prof.toc("setup");
105
106 // Show the mini-report on the constructed solver:
107 std::cout << solve << std::endl;
108
109 // Solve the system with the zero initial approximation:
110 std::vector<double> x(rows, 0.0);
111
112 prof.tic("solve");
113 Alina::SolverResult r = solve(A, rhs, x);
114 prof.toc("solve");
115
116 // Output the number of iterations, the relative error,
117 // and the profiling data:
118 std::cout << "Iters: " << r.nbIteration() << std::endl
119 << "Error: " << r.residual() << std::endl
120 << prof << std::endl;
121}
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 -*-
Hybrid backend uses scalar matrices to build the hierarchy, but stores the computed matrices in the b...