Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
coupcons3d_spc.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 <iostream>
33#include <string>
34
35#include "arccore/alina/BuiltinBackend.h"
36#include "arccore/alina/Adapters.h"
37#include "arccore/alina/StaticMatrix.h"
38#include "arccore/alina/SchurPressureCorrectionPreconditioner.h"
39#include "arccore/alina/PreconditionedSolver.h"
40#include "arccore/alina/make_block_solver.h"
41#include "arccore/alina/AMG.h"
42#include "arccore/alina/BiCGStabSolver.h"
43#include "arccore/alina/PreconditionerOnlySolver.h"
44#include "arccore/alina/Coarsening.h"
45#include "arccore/alina/Relaxation.h"
46
47#include "arccore/alina/IO.h"
48#include "arccore/alina/Profiler.h"
49
50using namespace Arcane;
51
52//---------------------------------------------------------------------------
53int main(int argc, char* argv[])
54{
55 // The command line should contain the matrix file name:
56 if (argc < 3) {
57 std::cerr << "Usage: " << argv[0] << " <matrix.mtx> <nu>" << std::endl;
58 return 1;
59 }
60
61 // The profiler:
62 auto& prof = Alina::Profiler::globalProfiler();
63
64 // Read the system matrix:
65 ptrdiff_t rows, cols;
66 std::vector<ptrdiff_t> ptr, col;
67 std::vector<double> val;
68
69 prof.tic("read");
70 std::tie(rows, cols) = Alina::IO::mm_reader(argv[1])(ptr, col, val);
71 std::cout << "Matrix " << argv[1] << ": " << rows << "x" << cols << std::endl;
72 prof.toc("read");
73
74 // The RHS is filled with ones:
75 std::vector<double> f(rows, 1.0);
76
77 // The number of unknowns in the U subsystem
78 ptrdiff_t nu = std::stoi(argv[2]);
79
80 // We use the tuple of CRS arrays to represent the system matrix.
81 // Note that std::tie creates a tuple of references, so no data is actually
82 // copied here:
83 auto A = std::tie(rows, ptr, col, val);
84
85 // Compose the solver type
86 typedef Alina::BuiltinBackend<double> SBackend; // the outer iterative solver backend
87 typedef Alina::BuiltinBackend<float> PBackend; // the PSolver backend
88 typedef Alina::BuiltinBackend<Alina::StaticMatrix<float, 4, 4>> UBackend; // the USolver backend
89
94 UBackend,
100 PBackend,
104 Solver;
105
106 // Solver parameters
107 Solver::params prm;
108 prm.precond.pmask.resize(rows);
109 for (ptrdiff_t i = 0; i < rows; ++i)
110 prm.precond.pmask[i] = (i >= nu);
111
112 // Initialize the solver with the system matrix.
113 prof.tic("setup");
114 Solver solve(A, prm);
115 prof.toc("setup");
116
117 // Show the mini-report on the constructed solver:
118 std::cout << solve << std::endl;
119
120 // Solve the system with the zero initial approximation:
121 std::vector<double> x(rows, 0.0);
122 prof.tic("solve");
123 Alina::SolverResult r = solve(A, f, x);
124 prof.toc("solve");
125
126 // Output the number of iterations, the relative error,
127 // and the profiling data:
128 std::cout << "Iters: " << r.nbIteration() << std::endl
129 << "Error: " << r.residual() << std::endl
130 << prof << std::endl;
131}
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.
Solver which only apply preconditioner once.
Allows to use an AMG smoother as standalone preconditioner.
Definition Relaxation.h:144
Result of a solving.
Definition AlinaUtils.h:52
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
Non-smoothed aggregation.
Definition Coarsening.h:639
Sparse approximate interface smoother.