Arcane  4.1.12.0
Developer documentation
Loading...
Searching...
No Matches
AlephTestModule.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/* AlephUnitTest.cc (C) 2000-2024 */
9/* */
10/* Aleph service testing. */
11/*---------------------------------------------------------------------------*/
12
13#include "arcane/aleph/tests/AlephTest.h"
14
15/*---------------------------------------------------------------------------*/
16/*---------------------------------------------------------------------------*/
17
18namespace ArcaneTest
19{
20
21/*---------------------------------------------------------------------------*/
22/*---------------------------------------------------------------------------*/
23
24using namespace Arcane;
25
26AlephTestModule::
27AlephTestModule(const ModuleBuildInfo& mbi)
29, m_total_nb_cell(0)
30, m_local_nb_cell(0)
31, m_rows_nb_element(0)
32, m_vector_indexs(0)
33, m_vector_zeroes(0)
34, m_vector_rhs(0)
35, m_vector_lhs(0)
36, m_aleph_kernel(NULL)
37, m_aleph_factory(new AlephFactory(subDomain()->application(), traceMng()))
38, m_aleph_params(NULL)
39, m_aleph_mat(0)
40, m_aleph_rhs(0)
41, m_aleph_sol(0)
42, m_get_solution_idx(0)
43, m_fake_nb_iteration(0)
44{
45 debug() << "\33[37m[AlephTestModule::AlephTestModule] Loading AlephTestModule\33[0m";
46}
47
48AlephTestModule::
49~AlephTestModule()
50{
51 delete m_aleph_kernel;
52 delete m_aleph_params;
53 delete m_aleph_factory;
54}
55
56/***************************************************************************
57 * AlephTestModule::init *
58 ***************************************************************************/
59void AlephTestModule::
60init()
61{
62 ItacFunction(AlephTestModule);
63 m_aleph_mat.resize(options()->alephNumberOfSolvers);
64 m_aleph_rhs.resize(options()->alephNumberOfSolvers);
65 m_aleph_sol.resize(options()->alephNumberOfSolvers);
66
67 // We retrieve the number of cells
68 m_local_nb_cell = MESH_OWN_ACTIVE_CELLS(mesh()).size();
69 m_total_nb_cell = subDomain()->parallelMng()->reduce(Parallel::ReduceSum, m_local_nb_cell);
70 debug() << "\33[37m[AlephTestModule::init] m_total_nb_cell=" << m_total_nb_cell
71 << ", local=" << m_local_nb_cell << "\33[0m";
72
73 // We set our time step
74 m_global_deltat = options()->deltaT;
75
76 // initialization of temperature on all cells and faces
77 ENUMERATE_CELL (icell, allCells())
78 m_cell_temperature[icell] = options()->initTemperature();
79 ENUMERATE_FACE (iFace, allFaces())
80 m_face_temperature[iFace] = options()->initTemperature();
81
82 // Update of UID and SDID variables
83 ENUMERATE_CELL (icell, allCells())
84 m_sub_domain_id[icell] = subDomain()->subDomainId();
85 ENUMERATE_CELL (icell, allCells())
86 m_unique_id[icell] = icell->uniqueId().asInteger();
87
88 // initialization of temperature at boundaries
89 options()->schema()->boundaries(options());
90
91 // We refine the cell ratio from the dataset
92 if ((options()->initAmr() < 0) || (options()->initAmr() > 1.0))
93 throw FatalErrorException("AlephTestModule::init()", "AMR ratio out of range");
94
95 initAmrRefineMesh((Integer)(mesh()->cellFamily()->allItems().own().size() * options()->initAmr()));
96 mesh()->checkValidMeshFull();
97
98 initAlgebra();
99}
100
101/***************************************************************************
102 *
103 ***************************************************************************/
104
105void AlephTestModule::
106initAlgebra()
107{
108 const Integer nb_row_size = m_total_nb_cell;
109 const Integer nb_row_rank = m_local_nb_cell;
110 ItacFunction(AlephTestModule);
111
112 Int32 underlying_solver = options()->alephUnderlyingSolver;
113
114 info() << "[AlephTestModule::initAlgebra] ALEPH init, " << nb_row_size
115 << " lines in total, " << nb_row_rank << " for me";
116 info() << "[AlephTestModule::initAlgebra] alephUnderlyingSolver="
117 << underlying_solver;
118 info() << "[AlephTestModule::initAlgebra] alephNumberOfCores="
119 << options()->alephNumberOfCores;
120
121 delete m_aleph_factory;
122 m_aleph_factory = new AlephFactory(subDomain()->application(), traceMng());
123
124 // We create our scheduling policy in place relative to the described topology
125 if (m_aleph_kernel)
126 delete m_aleph_kernel;
127
128 m_aleph_kernel = new AlephKernel(traceMng(),
129 subDomain(),
130 m_aleph_factory,
131 underlying_solver,
132 options()->alephNumberOfCores,
133 options()->alephCellOrdering);
134 // To test the specific PETSc initialization, add arguments
135 // if the number of PE is even
136
137 const bool use_is_even = (subDomain()->parallelMng()->commSize() % 2) == 0;
138 const bool use_petsc_args = (underlying_solver == AlephKernel::SOLVER_PETSC && use_is_even);
139 if (use_petsc_args) {
140 StringList slist1;
141 slist1.add("-trmalloc");
142 slist1.add("-log_trace");
143 slist1.add("-help");
144 CommandLineArguments args(slist1);
145 m_aleph_kernel->solverInitializeArgs().setCommandLineArguments(args);
146 }
147 m_aleph_kernel->initialize(nb_row_size, nb_row_rank);
148
149 delete m_aleph_params;
150 m_aleph_params = new AlephParams(traceMng(),
151 1.0e-10, // m_param_epsilon convergence epsilon
152 2000, // m_param_max_iteration max iterations
153 TypesSolver::DIAGONAL, // m_param_preconditioner_method: DIAGONAL, AMG, IC
154 TypesSolver::PCG, // m_param_solver_method solution method
155 -1, // m_param_gamma
156 -1.0, // m_param_alpha
157 false, // m_param_xo_user by default Xo is not equal to 0
158 false, // m_param_check_real_residue
159 false, // m_param_print_real_residue
160 false, // m_param_debug_info
161 1.e-20, // m_param_min_rhs_norm
162 false, // m_param_convergence_analyse
163 true, // m_param_stop_error_strategy
164 false, // m_param_write_matrix_to_file_error_strategy
165 "SolveErrorAlephMatrix.dbg", // m_param_write_matrix_name_error_strategy
166 false, // m_param_listing_output
167 0., // m_param_threshold
168 false, // m_param_print_cpu_time_resolution
169 0, // m_param_amg_coarsening_method: by default Sloop,
170 0, // m_param_output_level
171 1, // m_param_amg_cycle: 1=V, 2= W, 3=Full Multigrid V
172 1, // m_param_amg_solver_iterations
173 1, // m_param_amg_smoother_iterations
174 TypesSolver::SymHybGSJ_smoother, // m_param_amg_smootherOption
175 TypesSolver::ParallelRugeStuben, // m_param_amg_coarseningOption
176 TypesSolver::CG_coarse_solver, // m_param_amg_coarseSolverOption
177 false, // m_param_keep_solver_structure
178 false, // m_param_sequential_solver
179 TypesSolver::RB); // m_param_criteria_stop
180
181 // Calculate the global indices of the cell x cell matrix
182 // And we prepare the vector indices once and for all
183 m_cell_matrix_idx.fill(-1);
184 {
185 AlephInt idx = 0;
186 // Very important, because we return to it during reconfiguration!!
187 m_vector_indexs.resize(0);
188 m_vector_zeroes.resize(0);
189 m_vector_lhs.resize(0);
190 const AlephInt row_offset = m_aleph_kernel->topology()->part()[m_aleph_kernel->rank()];
191 ENUMERATE_CELL (iCell, MESH_OWN_ACTIVE_CELLS(mesh())) {
192 m_cell_matrix_idx[iCell] = row_offset + idx;
193 m_vector_indexs.add(m_cell_matrix_idx[iCell] = row_offset + idx);
194 m_vector_zeroes.add(0.0);
195 m_vector_lhs.add(0.0);
196 idx += 1;
197 }
198 }
199 m_cell_matrix_idx.synchronize();
200
201 // Now, we will count the number of elements per row
202 // We perform a scan as we will have to do in the loop
203 //debug() << "\33[37m[AlephTestModule::initAlgebra] Counting the number of elements of local rows"<<"\33[0m";
204 m_rows_nb_element.resize(nb_row_rank);
205 // Can now be forgotten
206 options()->schema()->preFetchNumElementsForEachRow(m_rows_nb_element,
207 m_aleph_kernel->topology()->part()[m_aleph_kernel->rank()]);
208 debug() << "\33[37m[AlephTestModule::initAlgebra] done\33[0m";
209}
210
211/***************************************************************************
212 * Entry points leads us here
213 ***************************************************************************/
214void AlephTestModule::
215compute()
216{
217 Integer delta_amr = 0;
218 Integer nb_iteration;
219 Real residual_norm[4];
220 debug() << "\33[37m[AlephTestModule::compute]\33[0m";
221 const Integer rank_offset = m_aleph_kernel->topology()->part()[m_aleph_kernel->rank()];
222 ItacFunction(AlephTestModule);
223
224 for (int i = 0; i < options()->alephNumberOfSolvers; i += 1)
225 postSolver(i);
226
227 // And should be able to get the solutions after the last resolution was fired
228 for (int i = options()->alephNumberOfSolvers - 1; i >= 0; i -= 1) {
229 debug() << "\33[37m[AlephTestModule::compute] Getting solution #" << i << "\33[0m";
230 AlephVector* solution = m_aleph_kernel->syncSolver(i, nb_iteration, &residual_norm[0]);
231 //if (i!=m_get_solution_idx) continue;
232 info() << "\33[37mSolved in \33[7m" << nb_iteration << "\33[m iterations,"
233 << "residuals=[\33[1m" << residual_norm[0] << "\33[m," << residual_norm[3] << "]";
234 debug() << "\33[37m[AlephTestModule::compute] Applying solution #" << m_get_solution_idx << "\33[0m";
235 solution->getLocalComponents(m_vector_indexs.size(), m_vector_indexs.view(), m_vector_lhs.view());
236 m_get_solution_idx += 1;
237 }
238 m_get_solution_idx %= options()->alephNumberOfSolvers;
239
241 // AMR COARSEN & REFINE //
243 /*if (options()->schema()->amrRefine(m_vector_lhs, options()->trigRefine)==true){
244 debug()<<"\33[37m[AlephTestModule::compute] AMR REFINE"<<"\33[0m";
245 initAlgebra();
246 delta_amr+=1;
247 return;
248 }
249 if (options()->schema()->amrCoarsen(m_vector_lhs, options()->trigCoarse)==true){
250 debug()<<"\33[37m[AlephTestModule::compute] AMR COARSEN"<<"\33[0m";
251 initAlgebra();
252 delta_amr+=1;
253 return;
254 }*/
255
257 // DELETE or NOT the KERNEL between each iteration //
259 if (options()->alephDeleteKernel() == true) {
260 debug() << "\33[37m[AlephTestModule::compute] DELETE the KERNEL between each iteration"
261 << "\33[0m";
262 initAlgebra();
263 }
264
265 // Otherwise, we copy the results
266 debug() << "\33[37m[AlephTestModule::compute] Now get our results"
267 << "\33[0m";
268 ENUMERATE_CELL (iCell, MESH_OWN_ACTIVE_CELLS(mesh())) {
269 m_cell_temperature[iCell] = m_vector_lhs[m_cell_matrix_idx[iCell] - rank_offset];
270 }
271
272 // If we have reached our maximum number of iterations, we exit
273 if (subDomain()->commonVariables().globalIteration() >= options()->iterations + delta_amr)
274 subDomain()->timeLoopMng()->stopComputeLoop(true);
275
276 debug() << "\33[37m[AlephTestModule::compute] done"
277 << "\33[0m";
278}
279
280/***************************************************************************
281 * job
282 ***************************************************************************/
283void AlephTestModule::
284postSolver(const Integer i)
285{
286 ItacFunction(AlephTestModule);
287 debug() << "\33[37m[AlephTestModule::postSolver] #" << i << "\33[0m";
288 const Integer rank_offset = m_aleph_kernel->topology()->part()[m_aleph_kernel->rank()];
289
290 // Filling the right-hand side: boundary conditions + right-hand side
291 debug() << "\33[37m[AlephTestModule::postSolver] Filling the right-hand side"
292 << "\33[0m";
293 m_vector_rhs.resize(0);
294 ENUMERATE_CELL (iCell, MESH_OWN_ACTIVE_CELLS(mesh())) {
295 m_vector_rhs.add(m_cell_temperature[iCell]);
296 }
297
298 debug() << "\33[37m[AlephTestModule::postSolver] ENUMERATE_FACE"
299 << "\33[0m";
300 ENUMERATE_FACE (iFace, OUTER_ACTIVE_FACE_GROUP(allCells())) {
301 if (!iFace->cell(0).isOwn())
302 continue;
303 m_vector_rhs[m_cell_matrix_idx[iFace->cell(0)] - rank_offset] +=
304 options()->deltaT * (m_face_temperature[iFace]) / geoFaceSurface(*iFace, nodesCoordinates());
305 }
306
307 // Creation of the MatVec matrix and Aleph requirements
308 m_aleph_mat.setAt(i, m_aleph_kernel->createSolverMatrix());
309 m_aleph_rhs.setAt(i, m_aleph_kernel->createSolverVector()); // First vector returned IS the rhs
310 m_aleph_sol.setAt(i, m_aleph_kernel->createSolverVector()); // Next one IS the solution
311
312 //m_aleph_mat.at(i)->create(m_rows_nb_element);
313 m_aleph_mat.at(i)->create();
314 m_aleph_rhs.at(i)->create();
315 m_aleph_sol.at(i)->create();
316
317 // Filling the matrix and assembly
318 options()->schema()->setValues(options()->deltaT, m_aleph_mat.at(i));
319 m_aleph_mat.at(i)->assemble();
320
321 debug() << "\33[37m[AlephTestModule::postSolver] setLocalComponents"
322 << "\33[0m";
323 m_aleph_rhs.at(i)->setLocalComponents(m_vector_indexs.size(), m_vector_indexs.view(), m_vector_rhs.view());
324 m_aleph_rhs.at(i)->assemble();
325
326 m_aleph_sol.at(i)->setLocalComponents(m_vector_indexs.size(), m_vector_indexs.view(), m_vector_zeroes.view());
327 m_aleph_sol.at(i)->assemble();
328
329 // Now solve with Aleph
330 debug() << "\33[37m[AlephTestModule::postSolver] Now solve with Aleph"
331 << "\33[0m";
332 m_aleph_mat.at(i)->solve(m_aleph_sol.at(i),
333 m_aleph_rhs.at(i),
334 m_fake_nb_iteration,
335 &m_fake_residual_norm[0],
336 m_aleph_params,
337 true); // We wish to post asynchronously
338}
339
340/***************************************************************************
341 * geoFaceSurface
342 ***************************************************************************/
343Real AlephTestModule::
344geoFaceSurface(Face face, VariableItemReal3& nodes_coords)
345{
346 if (face.nbNode() == 4) {
347 const Real3 xyz0 = nodes_coords[face.node(0)];
348 const Real3 xyz1 = nodes_coords[face.node(1)];
349 const Real3 xyz3 = nodes_coords[face.node(3)];
350 return (xyz0 - xyz1).normL2() * (xyz0 - xyz3).normL2();
351 }
352 if (face.nbNode() == 2) {
353 const Real3 xyz0 = nodes_coords[face.node(0)];
354 const Real3 xyz1 = nodes_coords[face.node(1)];
355 return (xyz0 - xyz1).squareNormL2();
356 }
357 throw FatalErrorException("geoFace", "Nb nodes != 4 !=2");
358}
359
360/***************************************************************************
361 * amrRefineMesh *
362 ***************************************************************************/
363void AlephTestModule::
364initAmrRefineMesh(Integer nb_to_refine)
365{
366 ItacFunction(AlephTestModule);
367 if (nb_to_refine == 0)
368 return;
369
370 info() << "Refining nb_to_refine=" << nb_to_refine;
371
372 Int32UniqueArray cells_local_id;
373 ENUMERATE_CELL (iCell, allCells()) {
374 Cell cell = *iCell;
375 if ((cell.type() == IT_Hexaedron8) || (cell.type() == IT_Quad4)) {
376 if (nb_to_refine-- <= 0)
377 break;
378 info(5) << "\t[amrRefineMesh] refine cell.uid=" << cell.uniqueId();
379 cells_local_id.add(cell.localId());
380 //iItem->setFlags(iItem->flags() | ItemInternal::II_Refine);
381 }
382 }
383
384 info() << "now refine";
385 mesh()->modifier()->flagCellToRefine(cells_local_id);
386 mesh()->modifier()->adapt();
387
388 //MESH_MODIFIER_REFINE_ITEMS(mesh());
389
390 // Now call back the values
391 CellInfoListView cells(mesh()->cellFamily());
392 //ItemInternalList faces = mesh()->faceFamily()->itemsInternal();
393 for (Integer i = 0, is = cells_local_id.size(); i < is; ++i) {
394 Int32 lid = cells_local_id[i];
395 Cell cell = cells[lid];
396 //debug()<<"[amrRefineMesh] focus on cell #"<<lid<<", nbHChildren="<<cell.nbHChildren();
397 for (Integer j = 0, js = CELL_NB_H_CHILDREN(cell); j < js; ++j) {
398 //debug()<<"\t\t[amrRefineMesh] child cell #"<<cell.hChild(j).localId();
399 m_cell_temperature[cells[CELL_H_CHILD(cell, j).localId()]] = m_cell_temperature[cells[lid]];
400 auto faces = allCells().view()[CELL_H_CHILD(cell, j).localId()].toCell().faces();
401 Integer index = 0;
402 for (Face face : faces) {
403 if (face.isSubDomainBoundary()) {
404 //debug() << "\t\t\t[amrRefineMesh] outer face #"<< (*iFace).localId()<<", index="<<iFace.index()<<", T="<<m_face_temperature[cell.face(iFace.index())];
405 m_face_temperature[face] = m_face_temperature[cell.face(index)];
406 }
407 else {
408 //debug() << "\t\t\t[amrRefineMesh] inner face #"<< (*iFace).localId();//<<", T="<<m_face_temperature[face.toFace()];
409 m_face_temperature[face] = 0;
410 }
411 ++index;
412 }
413 }
414 }
415}
416
417/*---------------------------------------------------------------------------*/
418/*---------------------------------------------------------------------------*/
419
420ARCANE_DEFINE_STANDARD_MODULE(AlephTestModule, AlephTest);
421
422/*---------------------------------------------------------------------------*/
423/*---------------------------------------------------------------------------*/
424
425} // namespace ArcaneTest
426
427/*---------------------------------------------------------------------------*/
428/*---------------------------------------------------------------------------*/
#define ENUMERATE_FACE(name, group)
Generic enumerator for a face group.
#define ENUMERATE_CELL(name, group)
Generic enumerator for a cell group.
Test module for Aleph.
Generation de la classe de base du Module.
CaseOptionsAlephTestModule * options() const
Options du jeu de données du module.
Arcane::VariableCellReal m_cell_temperature
Variables du module.
Integer size() const
Number of elements in the vector.
ISubDomain * subDomain() const override
Sub-domain associated with the module.
Parameters of a linear system.
Definition AlephParams.h:32
Vector of a linear system.
Definition AlephVector.h:33
void add(ConstReferenceType val)
Adds element val to the end of the array.
View of cell information.
Cell of a mesh.
Definition Item.h:1300
Face face(Int32 i) const
i-th face of the cell
Definition Item.h:1400
VariableScalarReal m_global_deltat
Global Delta T.
Face of a cell.
Definition Item.h:1032
bool isSubDomainBoundary() const
Indicates if the face is on the subdomain boundary (i.e nbCell()==1).
Definition Item.h:1148
Node node(Int32 i) const
i-th node of the entity
Definition Item.h:840
Int32 nbNode() const
Number of nodes of the entity.
Definition Item.h:837
constexpr Int32 localId() const
Local identifier of the entity in the processor subdomain.
Definition Item.h:233
ItemUniqueId uniqueId() const
Unique identifier across all domains.
Definition Item.h:239
Int16 type() const
Entity type.
Definition Item.h:255
CellGroup allCells() const
Returns the group containing all cells.
FaceGroup allFaces() const
Returns the group containing all faces.
Information for building a module.
Class managing a 3-dimensional real vector.
Definition Real3.h:132
TraceMessageDbg debug(Trace::eDebugLevel=Trace::Medium) const
Flow for a debug message.
ItemVariableScalarRefT< Real3 > VariableItemReal3
3D coordinate type quantity
@ ReduceSum
Sum of values.
Int32 Integer
Type representing an integer.
List< String > StringList
Unicode string list.
Definition UtilsTypes.h:509
UniqueArray< Int32 > Int32UniqueArray
Dynamic 1D array of 32-bit integers.
Definition UtilsTypes.h:341
int AlephInt
Default type for indexing rows and columns of matrices and vectors.
Definition AlephGlobal.h:50
std::int32_t Int32
Signed integer type of 32 bits.