Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
DistributedCheckDirect.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#if defined(ARCCORE_ALINA_HAVE_EIGEN)
20// To remove warnings about deprecated Eigen usage.
21#pragma GCC diagnostic ignored "-Wdeprecated-copy"
22#pragma GCC diagnostic ignored "-Wint-in-bool-context"
23#include "arccore/alina/DistributedEigenSparseLUDirectSolver.h"
24#endif
25
26#include "arccore/trace/ITraceMng.h"
27
28#include "arccore/alina/DistributedDirectSolverRuntime.h"
29#include "arccore/alina/Profiler.h"
30
31#include "AlinaSamplesCommon.h"
32
33#include <boost/program_options.hpp>
34
35using namespace Arcane;
36
37int main2(const Alina::SampleMainContext& ctx, int argc, char* argv[])
38{
39 ITraceMng* tm = ctx.traceMng();
40
41 auto& prof = Alina::Profiler::globalProfiler();
42
43 Alina::mpi_communicator comm(MPI_COMM_WORLD);
44
45 tm->info() << "World size: " << comm.size;
46
47 namespace po = boost::program_options;
48 po::options_description desc("Options");
49
50 desc.add_options()("help,h", "show help")("size,n",
51 po::value<int>()->default_value(128),
52 "domain size");
53 desc.add_options()("prm-file,P",
54 po::value<std::string>(),
55 "Parameter file in json format. ");
56 desc.add_options()("prm,p",
57 po::value<std::vector<std::string>>()->multitoken(),
58 "Parameters specified as name=value pairs. "
59 "May be provided multiple times. Examples:\n"
60 " -p solver.tol=1e-3\n"
61 " -p precond.coarse_enough=300");
62
63 po::variables_map vm;
64 po::store(po::parse_command_line(argc, argv, desc), vm);
65 po::notify(vm);
66
67 if (vm.count("help")) {
68 if (comm.rank == 0)
69 std::cout << desc << std::endl;
70 return 0;
71 }
72
74 if (vm.count("prm-file")) {
75 prm.read_json(vm["prm-file"].as<std::string>());
76 }
77
78 if (vm.count("prm")) {
79 for (const std::string& v : vm["prm"].as<std::vector<std::string>>()) {
80 prm.putKeyValue(v);
81 }
82 }
83
84 const int n = vm["size"].as<int>();
85 const int n2 = n * n;
86
87 int chunk = (n2 + comm.size - 1) / comm.size;
88 int chunk_start = comm.rank * chunk;
89 int chunk_end = std::min(chunk_start + chunk, n2);
90
91 chunk = chunk_end - chunk_start;
92
93 std::vector<int> domain(comm.size + 1);
94 MPI_Allgather(&chunk, 1, MPI_INT, &domain[1], 1, MPI_INT, comm);
95 std::partial_sum(domain.begin(), domain.end(), domain.begin());
96
97 prof.tic("assemble");
99 A.set_size(chunk, domain.back(), true);
100 A.set_nonzeros(chunk * 5);
101 std::vector<double> rhs(chunk, 1);
102
103 const double h2i = (n - 1) * (n - 1);
104 for (int idx = chunk_start, row = 0, head = 0; idx < chunk_end; ++idx, ++row) {
105 int j = idx / n;
106 int i = idx % n;
107
108 if (j > 0) {
109 A.col[head] = idx - n;
110 A.val[head] = -h2i;
111 ++head;
112 }
113
114 if (i > 0) {
115 A.col[head] = idx - 1;
116 A.val[head] = -h2i;
117 ++head;
118 }
119
120 A.col[head] = idx;
121 A.val[head] = 4 * h2i;
122 ++head;
123
124 if (i + 1 < n) {
125 A.col[head] = idx + 1;
126 A.val[head] = -h2i;
127 ++head;
128 }
129
130 if (j + 1 < n) {
131 A.col[head] = idx + n;
132 A.val[head] = -h2i;
133 ++head;
134 }
135
136 A.ptr[row + 1] = head;
137 }
138 A.setNbNonZero(A.ptr[chunk]);
139 prof.toc("assemble");
140
141 std::vector<double> x(chunk);
142
143 {
144 auto t = prof.scoped_tic("skyline");
145
146 prof.tic("setup");
148 tm->info() << "Solver is = " << solve.type();
149 prof.toc("setup");
150
151 prof.tic("solve1");
152 solve(rhs, x);
153 prof.toc("solve1");
154 }
155
156#if defined(ARCCORE_ALINA_HAVE_EIGEN)
157 {
158 auto t = prof.scoped_tic("eigen");
159
160 prof.tic("setup");
162 prof.toc("setup");
163
164 prof.tic("solve");
165 std::vector<double> x2(chunk);
166 solve(rhs, x2);
167 prof.toc("solve");
168 // Compare values between skyline and eigen
169 Int32 nb_error = 0;
170 for (Int32 i = 0; i < chunk; ++i) {
171 Real x_ref = x[i];
172 Real x_eigen = x2[i];
173 Real abs_sum = std::abs(x_ref);
174
175 Real diff = std::abs(x_eigen - x_ref);
176 if (abs_sum != 0.0)
177 diff /= abs_sum;
178
179 if (diff > 1.0e-12) {
180 tm->info() << "I=" << i << " compare skyline=" << x[i] << " eigen=" << x2[i] << " diff=" << diff;
181 ++nb_error;
182 }
183 }
184 if (nb_error != 0)
185 ARCCORE_FATAL("Error when comparing Eigen and SkylineLU nb_error={0}", nb_error);
186 }
187#endif
188
189 tm->info() << prof;
190 return 0;
191}
192
193int main(int argc, char* argv[])
194{
195 return Arcane::Alina::SampleMainContext::execMain(main2, argc, argv);
196}
#define ARCCORE_FATAL(...)
Macro envoyant une exception FatalErrorException.
Runtime wrapper for distributed direct solvers.
Interface du gestionnaire de traces.
virtual TraceMessage info()=0
Flot pour un message d'information.
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
Sparse matrix stored in CSR (Compressed Sparse Row) format.
Definition CSRMatrix.h:98
Convenience wrapper around MPI_Comm.