Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
serena_mpi.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/StaticMatrix.h"
37#include "arccore/alina/Adapters.h"
38
39#include "arccore/alina/DistributedMatrix.h"
40#include "arccore/alina/DistributedPreconditionedSolver.h"
41#include "arccore/alina/DistributedAMG.h"
42#include "arccore/alina/DistributedCoarsening.h"
43#include "arccore/alina/DistributedRelaxation.h"
44#include "arccore/alina/DistributedSolver.h"
45
46#include "arccore/alina/IO.h"
47#include "arccore/alina/Profiler.h"
48
49#if defined(ARCCORE_ALINA_HAVE_PARMETIS)
50#include "arccore/alina/ParmetisMatrixPartitioner.h"
51#endif
52
53// Block size
54const int B = 3;
55
56using namespace Arcane;
57using namespace Arcane::Alina;
58
59//---------------------------------------------------------------------------
60int main(int argc, char* argv[])
61{
62 // The command line should contain the matrix file name:
63 if (argc < 2) {
64 std::cerr << "Usage: " << argv[0] << " <matrix.bin>" << std::endl;
65 return 1;
66 }
67
68 Alina::mpi_init mpi(&argc, &argv);
69 Alina::mpi_communicator world(MPI_COMM_WORLD);
70
71 auto& prof = Alina::Profiler::globalProfiler();
72
73 prof.tic("read");
74 // Get the global size of the matrix:
75 ptrdiff_t rows = Alina::IO::crs_size<ptrdiff_t>(argv[1]);
76
77 // Split the matrix into approximately equal chunks of rows, and
78 // make sure each chunk size is divisible by the block size.
79 ptrdiff_t chunk = (rows + world.size - 1) / world.size;
80 if (chunk % B)
81 chunk += B - chunk % B;
82
83 ptrdiff_t row_beg = std::min(rows, chunk * world.rank);
84 ptrdiff_t row_end = std::min(rows, row_beg + chunk);
85 chunk = row_end - row_beg;
86
87 // Read our part of the system matrix.
88 std::vector<ptrdiff_t> ptr, col;
89 std::vector<double> val;
90 Alina::IO::read_crs(argv[1], rows, ptr, col, val, row_beg, row_end);
91 prof.toc("read");
92
93 if (world.rank == 0)
94 std::cout
95 << "World size: " << world.size << std::endl
96 << "Matrix " << argv[1] << ": " << rows << "x" << rows << std::endl;
97
98 // Declare the backend and the solver types
99 typedef Alina::StaticMatrix<double, B, B> dmat_type;
100 typedef Alina::StaticMatrix<double, B, 1> dvec_type;
101 typedef Alina::StaticMatrix<float, B, B> fmat_type;
102 typedef Alina::BuiltinBackend<dmat_type> DBackend;
103 typedef Alina::BuiltinBackend<fmat_type> FBackend;
104
107 FBackend,
111 Solver;
112
113 // Solver parameters
114 Solver::params prm;
115 prm.solver.maxiter = 200;
116
117 // We need to scale the matrix, so that it has the unit diagonal.
118 // Since we only have the local rows for the matrix, and we may need the
119 // remote diagonal values, it is more convenient to represent the scaling
120 // with the matrix-matrix product (As = D^-1/2 A D^-1/2).
121 prof.tic("scale");
122 // Find the local diagonal values,
123 // and form the CRS arrays for a diagonal matrix.
124 std::vector<double> dia(chunk, 1.0);
125 std::vector<ptrdiff_t> d_ptr(chunk + 1), d_col(chunk);
126 for (ptrdiff_t i = 0, I = row_beg; i < chunk; ++i, ++I) {
127 d_ptr[i] = i;
128 d_col[i] = I;
129 for (ptrdiff_t j = ptr[i], e = ptr[i + 1]; j < e; ++j) {
130 if (col[j] == I) {
131 dia[i] = 1 / sqrt(val[j]);
132 break;
133 }
134 }
135 }
136 d_ptr.back() = chunk;
137
138 // Create the distributed diagonal matrix:
140 Alina::adapter::block_matrix<dmat_type>(
141 std::tie(chunk, d_ptr, d_col, dia)));
142
143 // The scaled matrix is formed as product D * A * D,
144 // where A is the local chunk of the matrix
145 // converted to the block format on the fly.
146 auto A = product(D, *product(Alina::DistributedMatrix<DBackend>(world, Alina::adapter::block_matrix<dmat_type>(std::tie(chunk, ptr, col, val))), D));
147 prof.toc("scale");
148
149 // Since the RHS in this case is filled with ones,
150 // the scaled RHS is equal to dia.
151 // Reinterpret the pointer to dia data to get the RHS in the block format:
152 auto f_ptr = reinterpret_cast<dvec_type*>(dia.data());
153 std::vector<dvec_type> rhs(f_ptr, f_ptr + chunk / B);
154
155 // Partition the matrix and the RHS vector.
156 // If neither ParMETIS not PT-SCOTCH are not available,
157 // just keep the current naive partitioning.
158#if defined(ARCCORE_ALINA_HAVE_PARMETIS)
160
161 if (world.size > 1) {
162 prof.tic("partition");
163 Partition part;
164
165 // part(A) returns the distributed permutation matrix:
166 auto P = part(*A);
167 auto R = transpose(*P);
168
169 // Reorder the matrix:
170 A = product(*R, *product(*A, *P));
171
172 // and the RHS vector:
173 std::vector<dvec_type> new_rhs(R->loc_rows());
174 R->move_to_backend();
175 Alina::backend::spmv(1, *R, rhs, 0, new_rhs);
176 rhs.swap(new_rhs);
177
178 // Update the number of the local rows
179 // (it may have changed as a result of permutation).
180 // Note that A->loc_rows() returns the number of blocks,
181 // as the matrix uses block values.
182 chunk = A->loc_rows();
183 prof.toc("partition");
184 }
185#endif
186
187 // Initialize the solver:
188 prof.tic("setup");
189 Solver solve(world, A, prm);
190 prof.toc("setup");
191
192 // Show the mini-report on the constructed solver:
193 if (world.rank == 0)
194 std::cout << solve << std::endl;
195
196 // Solve the system with the zero initial approximation:
197 std::vector<dvec_type> x(chunk, Alina::math::zero<dvec_type>());
198
199 prof.tic("solve");
200 Alina::SolverResult r = solve(*A, rhs, x);
201 prof.toc("solve");
202
203 // Output the number of iterations, the relative error,
204 // and the profiling data:
205 if (world.rank == 0)
206 std::cout << "Iters: " << r.nbIteration() << std::endl
207 << "Error: " << r.residual() << std::endl
208 << prof << std::endl;
209}
Distributed Matrix using message passing.
Iterative solver wrapper for distributed linear systems.
Result of a solving.
Definition AlinaUtils.h:52
__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 -*-
Distributed smoothed aggregation coarsening scheme.
Convenience wrapper around MPI_Comm.
Convenience wrapper around MPI_Init/MPI_Finalize.