Arcane  v3.14.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
AlephKernel.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2022 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/* AlephKernel.cc (C) 2000-2022 */
9/* */
10/* */
11/*---------------------------------------------------------------------------*/
12/*---------------------------------------------------------------------------*/
13
14#include "arcane/aleph/AlephArcane.h"
15#include "arcane/utils/StringList.h"
16
17/*---------------------------------------------------------------------------*/
18/*---------------------------------------------------------------------------*/
19
20namespace Arcane
21{
22
23/*---------------------------------------------------------------------------*/
24/*---------------------------------------------------------------------------*/
25
26AlephKernelSolverInitializeArguments::
27AlephKernelSolverInitializeArguments()
28: m_arguments(StringList{})
29{}
30
31/*---------------------------------------------------------------------------*/
32/*---------------------------------------------------------------------------*/
36AlephKernel::
37AlephKernel(IParallelMng* wpm,
38 Integer size,
39 IAlephFactory* factory,
40 Integer alephUnderlyingSolver,
41 Integer alephNumberOfCores,
42 bool alephOrdering)
43: TraceAccessor(wpm->traceMng())
44, m_isParallel(wpm->isParallel())
45, m_rank(wpm->commRank())
46, m_size(size)
47, m_world_size(wpm->commSize())
48, m_there_are_idles(true)
49, m_i_am_an_other(true)
50, m_parallel(wpm)
51, m_world_parallel(wpm)
52, m_factory(factory)
53, m_aleph_vector_idx(0)
54// Pour l'instant, on met Sloop par défaut
55, m_underlying_solver((alephUnderlyingSolver == 0 ? 1 : alephUnderlyingSolver))
56, m_reorder(alephOrdering)
57, m_solver_index(0)
58, m_solver_size(alephNumberOfCores)
59, m_solved(false)
60, m_has_been_initialized(false)
61, m_solver_ranks(0)
62, m_matrix_queue(0)
63, m_arguments_queue(0)
64, m_results_queue(0)
65{
66 setup();
67}
68
69/*---------------------------------------------------------------------------*/
70/*---------------------------------------------------------------------------*/
79 IAlephFactory* factory,
80 Integer alephUnderlyingSolver,
81 Integer alephNumberOfCores,
82 bool alephOrdering)
84, m_sub_domain(sd)
85, m_isParallel(sd->parallelMng()->isParallel())
86, m_rank(sd->parallelMng()->commRank())
87, m_size(sd->parallelMng()->commSize())
88, m_world_size(sd->parallelMng()->worldParallelMng()->commSize())
89, m_there_are_idles(m_size != m_world_size)
90, m_i_am_an_other(sd->parallelMng()->worldParallelMng()->commRank() > m_size)
91, m_parallel(sd->parallelMng())
92, m_world_parallel(sd->parallelMng()->worldParallelMng())
93, m_factory(factory)
94, m_aleph_vector_idx(0)
95// Pour l'instant, on met Sloop par défaut
96, m_underlying_solver((alephUnderlyingSolver == 0 ? 1 : alephUnderlyingSolver))
97, m_reorder(alephOrdering)
98, m_solver_index(0)
99, m_solver_size(alephNumberOfCores)
100, m_solved(false)
101, m_has_been_initialized(false)
102{
103 setup();
104}
105
106/*---------------------------------------------------------------------------*/
107/*---------------------------------------------------------------------------*/
116 Integer alephUnderlyingSolver,
117 Integer alephNumberOfCores)
118: TraceAccessor(sd->parallelMng()->traceMng())
119, m_sub_domain(sd)
120, m_isParallel(sd->parallelMng()->isParallel())
121, m_rank(sd->parallelMng()->commRank())
122, m_size(sd->parallelMng()->commSize())
123, m_world_size(sd->parallelMng()->worldParallelMng()->commSize())
124, m_there_are_idles(m_size != m_world_size)
125, m_i_am_an_other(sd->parallelMng()->worldParallelMng()->commRank() > m_size)
126, m_parallel(sd->parallelMng())
127, m_world_parallel(sd->parallelMng()->worldParallelMng())
128, m_factory(new AlephFactory(sd->application(), sd->parallelMng()->traceMng()))
129, m_topology(new AlephTopology(this))
130, m_ordering(new AlephOrdering(this))
131, m_indexing(new AlephIndexing(this))
132, m_aleph_vector_idx(0)
133, m_underlying_solver(alephUnderlyingSolver == 0 ? 1 : alephUnderlyingSolver)
134, m_reorder(false)
135, m_solver_index(0)
136, m_solver_size(alephNumberOfCores)
137, m_solved(false)
138, m_has_been_initialized(false)
139{
140 debug() << "\33[1;31m[AlephKernel] New kernel with indexing+init options!\33[0m";
141 setup();
142}
143
144/*---------------------------------------------------------------------------*/
145/*---------------------------------------------------------------------------*/
150setup(void)
151{
152 ItacFunction(AlephKernel);
153 if (m_sub_domain) {
154 debug() << "\33[1;31m[AlephKernel] thisParallelMng's size=" << m_size << "\33[0m";
155 debug() << "\33[1;31m[AlephKernel] worldParallelMng's size=" << m_world_size << "\33[0m";
156 }
157 else
158 debug() << "\33[1;31m[AlephKernel] I am an additional site #" << m_rank << " among " << m_world_size << "\33[0m";
159 // Solveur utilisé par défaut
160 debug() << "\33[1;31m[AlephKernel] Aleph underlying solver has been set to "
161 << m_underlying_solver << "\33[0m";
162 // Par défaut, on utilise tous les coeurs alloués au calcul pour chaque résolution
163 if (m_solver_size == 0) {
164 m_solver_size = m_world_size;
165 debug() << "\33[1;31m[AlephKernel] Aleph Number of Cores"
166 << " now matches world's number of processors: "
167 << m_solver_size << "\33[0m";
168 }
169 if (m_solver_size > m_size) {
170 m_solver_size = m_size;
171 debug() << "\33[1;31m[AlephKernel] Aleph Number of Cores"
172 << " exceeds in size, reverting to " << m_size << "\33[0m";
173 }
174 if ((m_size % m_solver_size) != 0)
175 throw FatalErrorException("AlephKernel", "Aleph Number of Cores modulo size");
176 debug() << "\33[1;31m[AlephKernel] Each solver takes "
177 << m_solver_size << " site(s)"
178 << "\33[0m";
179 // S'il y a des 'autres, on les tient au courant de la configuration
180 if (m_there_are_idles && !m_i_am_an_other) {
182 cfg.add(m_underlying_solver);
183 cfg.add(m_solver_size);
184 cfg.add((m_reorder == true) ? 1 : 0);
185 cfg.add(m_size);
186 debug() << "\33[1;31m[AlephKernel] Sending to others configuration: " << cfg << "\33[0m";
187 m_world_parallel->broadcast(cfg.view(), 0);
188 }
189}
190
191/*---------------------------------------------------------------------------*/
192/*---------------------------------------------------------------------------*/
193
194AlephKernel::
195~AlephKernel()
196{
197 info(4) << "Destroying ~AlephKernel";
198
199 delete m_topology;
200 delete m_indexing;
201 delete m_ordering;
202
203 for ( AlephKernelResults* rq : m_results_queue )
204 delete rq;
205 m_results_queue.clear();
206
207 for ( AlephMatrix* mq : m_matrix_queue )
208 delete mq;
209 m_matrix_queue.clear();
210
211 for ( AlephKernelArguments* aq : m_arguments_queue ){
212 // TODO: regarder pourquoi cela n'est pas fait dans le destructeur de AlephKernelArguments.
213 delete aq->m_x_vector;
214 delete aq->m_b_vector;
215 delete aq->m_tmp_vector;
216 delete aq;
217 }
218 m_arguments_queue.clear();
219
220 // PETSc seems not to like this too much but this is needed
221 //for ( IParallelMng* pm : m_sub_parallel_mng_queue )
222 //delete pm;
223
224 info(4) << "Destroying ~AlephKernel] done";
225}
226
227/*---------------------------------------------------------------------------*/
228/*---------------------------------------------------------------------------*/
229// d80dee82
230void AlephKernel::
231initialize(Integer global_nb_row,
232 Integer local_nb_row)
233{
234 ItacFunction(AlephKernel);
235 //Timer::Action ta(subDomain(),"AlephKernel::initialize");
236 if (m_there_are_idles && !m_i_am_an_other) {
237 m_world_parallel->broadcast(UniqueArray<unsigned long>(1, 0xd80dee82l).view(), 0);
238 UniqueArray<Integer> args(0);
239 args.add(global_nb_row);
240 args.add(local_nb_row);
241 m_world_parallel->broadcast(args.view(), 0);
242 }
243 debug() << "\33[1;31m[initialize] Geometry set to " << global_nb_row
244 << " lines, I see " << local_nb_row << " of them"
245 << "\33[0m";
246 m_topology = new AlephTopology(traceMng(), this, global_nb_row, local_nb_row);
247 m_ordering = new AlephOrdering(this, global_nb_row, local_nb_row, m_reorder);
248 debug() << "\33[1;31m[initialize] done"
249 << "\33[0m";
250 m_has_been_initialized = true;
251}
252
253/*---------------------------------------------------------------------------*/
254/*---------------------------------------------------------------------------*/
255// 4b97b15d
256void AlephKernel::
257break_and_return(void)
258{
259 if (m_there_are_idles && !m_i_am_an_other)
260 m_world_parallel->broadcast(UniqueArray<unsigned long>(1, 0x4b97b15dl).view(), 0);
261}
262
263/*---------------------------------------------------------------------------*/
264/*---------------------------------------------------------------------------*/
273{
274 debug() << "\33[1;31m[mapranks] mapranks starting @ "
275 << m_solver_index * m_size
276 << ", m_size=" << m_size
277 << ", m_solver_size=" << m_solver_size
278 << ", m_world_size=" << m_world_size << "\33[0m";
279 traceMng()->flush();
280 for (int rnk = m_solver_index * m_size; rnk < (m_solver_index + 1) * m_size; rnk += 1) {
281 const int map = (rnk / (m_size / m_solver_size)) % m_world_size;
282 debug() << "\33[1;31m[mapranks] map=" << map << "\33[0m";
283 ranks[rnk % m_size] = map;
284 debug() << "\33[1;31m[mapranks] mapped solver #" << m_solver_index
285 << ", core " << rnk % m_size << " --> site " << map << "\33[0m";
286 }
287}
288
289/*---------------------------------------------------------------------------*/
290/*---------------------------------------------------------------------------*/
291
292bool AlephKernel::
293hitranks(Integer rank, ArrayView<Integer> ranks)
294{
295 for (int rnk = ranks.size() - 1; rnk >= 0; rnk -= 1)
296 if (ranks[rnk] == rank)
297 return true;
298 return false;
299}
300
301/*---------------------------------------------------------------------------*/
302/*---------------------------------------------------------------------------*/
303
304Ref<IParallelMng> AlephKernel::
305_createUnderlyingParallelMng(Integer nb_wanted_sites)
306{
307 info(4) << "[createUnderlyingParallelMng] nb_wanted_sites=" << nb_wanted_sites;
308 UniqueArray<Integer> kept_ranks;
309 for (Int32 rank = 0; rank < m_world_size; ++rank) {
310 if (hitranks(rank, m_solver_ranks[m_solver_index]))
311 kept_ranks.add(rank);
312 }
313 info(4) << "[createUnderlyingParallelMng] Now createSubParallelMng of size=" << kept_ranks.size();
314 Ref<IParallelMng> upm = m_world_parallel->createSubParallelMngRef(kept_ranks.constView());
315 info(4) << "[createUnderlyingParallelMng] done: upm=" << upm.get();
316 traceMng()->flush();
317 return upm;
318}
319
320/*---------------------------------------------------------------------------*/
321/*---------------------------------------------------------------------------*/
322// BaseForm[Hash["createSolverMatrix", "CRC32"], 16] = ef162166
323AlephMatrix* AlephKernel::
324createSolverMatrix(void)
325{
326 ItacFunction(AlephService);
327
328 if (isInitialized() == false) {
329 debug() << "\33[1;31m[createSolverMatrix] has_NOT_been_initialized!\33[0m"
330 << "\33[0m";
331 return new AlephMatrix(this);
332 }
333
334 if (m_there_are_idles && !m_i_am_an_other)
335 m_world_parallel->broadcast(UniqueArray<unsigned long>(1, 0xef162166l).view(), 0);
336
337 debug() << "\33[1;31m[createSolverMatrix]\33[0m"
338 << "\33[0m";
339
340 m_solved = false;
341
342 if (!m_configured) {
343 debug() << "\33[1;31m[createSolverMatrix] UN configured, building Underlying Parallel Managers index="
344 << index() << "\33[0m";
345 traceMng()->flush();
346 m_solver_ranks.add(SharedArray<Integer>(m_world_size));
347 m_solver_ranks[m_solver_index].fill(-1);
348 mapranks(m_solver_ranks[m_solver_index]);
349 traceMng()->flush();
350 Ref<IParallelMng> upm = _createUnderlyingParallelMng(m_solver_size);
351 if (upm.get()) {
352 debug() << "\33[1;31m[createSolverMatrix] upm->isParallel()=" << upm->isParallel() << "\33[0m";
353 debug() << "\33[1;31m[createSolverMatrix] upm->commSize()=" << upm->commSize() << "\33[0m";
354 debug() << "\33[1;31m[createSolverMatrix] upm->commRank()=" << upm->commRank() << "\33[0m";
355 }
356 else {
357 debug() << "\33[1;31m[createSolverMatrix] upm NULL"
358 << "\33[0m";
359 }
360 m_sub_parallel_mng_queue.add(upm);
361 debug() << "\33[1;31m[createSolverMatrix] Queuing new kernel arguments: X, B and Tmp with their topolgy"
362 << "\33[0m";
363 // On va chercher la topologie avant toute autres choses afin que la bibliothèque
364 // sous-jacente la prenne en compte pour les vecteurs et la matrice à venir
365 IAlephTopology* underlying_topology = factory()->GetTopology(this, index(), topology()->nb_row_size());
366 // Dans le cas d'une bibliothèque qui possède une IAlephTopology
367 // On trig le prefix, on fera le postfix apres les solves
368 if (underlying_topology != NULL)
369 underlying_topology->backupAndInitialize();
370 m_arguments_queue.add(new AlephKernelArguments(traceMng(),
371 new AlephVector(this), // Vecteur X
372 new AlephVector(this), // Vecteur B
373 new AlephVector(this), // Vecteur tmp (pour l'isAlreadySolved)
374 underlying_topology));
375 // On initialise le vecteur temporaire qui n'est pas vu de l'API
376 // Pas besoin de l'assembler, il sert en output puis comme buffer
377 debug() << "\33[1;31m[createSolverMatrix] Creating Tmp vector for this set of arguments"
378 << "\33[0m";
379 m_arguments_queue.at(m_solver_index)->m_tmp_vector->create();
380 // On initialise la matrice apres la topologies et les vecteurs
381 debug() << "\33[1;31m[createSolverMatrix] Now queuing the matrix\33[0m";
382 m_matrix_queue.add(new AlephMatrix(this));
383 debug() << "\33[1;31m[createSolverMatrix] Now queuing the space for the resolution results\33[0m";
384 m_results_queue.add(new AlephKernelResults());
385 }
386 else {
387 if (getTopologyImplementation(m_solver_index) != NULL)
388 getTopologyImplementation(m_solver_index)->backupAndInitialize();
389 }
390 debug() << "\33[1;31m[createSolverMatrix] done!\33[0m";
391 traceMng()->flush();
392 return m_matrix_queue[m_solver_index];
393}
394
395/*---------------------------------------------------------------------------*/
396/*---------------------------------------------------------------------------*/
405{
406 ItacFunction(AlephKernel);
407 if (m_has_been_initialized == false) {
408 debug() << "\33[1;31m[createSolverVector] has_NOT_been_initialized!\33[0m";
409 return new AlephVector(this);
410 }
411 if (m_there_are_idles && !m_i_am_an_other)
412 m_world_parallel->broadcast(UniqueArray<unsigned long>(1, 0xc4b28f2l).view(), 0);
413 m_aleph_vector_idx++;
414 if ((m_aleph_vector_idx % 2) == 0) {
415 debug() << "\33[1;31m[createSolverVector] Get " << m_solver_index << "th X vector\33[0m";
416 return m_arguments_queue.at(m_solver_index)->m_x_vector;
417 }
418 else {
419 debug() << "\33[1;31m[createSolverVector] Get " << m_solver_index << "th B vector\33[0m";
420 return m_arguments_queue.at(m_solver_index)->m_b_vector;
421 }
422 traceMng()->flush();
423}
424
425/*---------------------------------------------------------------------------*/
426/*---------------------------------------------------------------------------*/
427// ba9488be
428
429void AlephKernel::
430postSolver(AlephParams* params,
434{
435 ItacFunction(AlephKernel);
436
437 if (!isInitialized()) {
438 debug() << "\33[1;31m[postSolver] Trying to post a solver to an uninitialized kernel!\33[0m";
439 debug() << "\33[1;31m[postSolver] Now telling Indexer to do its job!\33[0m";
440 indexing()->nowYouCanBuildTheTopology(fromThisMatrix, fromeThisX, fromThisB);
441 }
442
443 if (m_there_are_idles && !m_i_am_an_other) {
444 m_world_parallel->broadcast(UniqueArray<unsigned long>(1, 0xba9488bel).view(), 0);
445 UniqueArray<Real> real_args(0);
446 real_args.add(params->epsilon());
447 real_args.add(params->alpha());
448 real_args.add(params->minRHSNorm());
449 real_args.add(params->DDMCParameterAmgDiagonalThreshold());
450
451 UniqueArray<int> bool_args(0);
452 bool_args.add(params->xoUser());
453 bool_args.add(params->checkRealResidue());
454 bool_args.add(params->printRealResidue());
455 bool_args.add(params->debugInfo());
456 bool_args.add(params->convergenceAnalyse());
457 bool_args.add(params->stopErrorStrategy());
458 bool_args.add(params->writeMatrixToFileErrorStrategy());
459 bool_args.add(params->DDMCParameterListingOutput());
460 bool_args.add(params->printCpuTimeResolution());
461 bool_args.add(params->getKeepSolverStructure());
462 bool_args.add(params->getSequentialSolver());
463
464 UniqueArray<Integer> int_args(0);
465 int_args.add(params->maxIter());
466 int_args.add(params->gamma());
467 int_args.add((Integer)params->precond());
468 int_args.add((Integer)params->method());
469 int_args.add((Integer)params->amgCoarseningMethod());
470 int_args.add(params->getOutputLevel());
471 int_args.add(params->getAmgCycle());
472 int_args.add(params->getAmgSolverIter());
473 int_args.add(params->getAmgSmootherIter());
474 int_args.add((Integer)params->getAmgSmootherOption());
475 int_args.add((Integer)params->getAmgCoarseningOption());
476 int_args.add((Integer)params->getAmgCoarseSolverOption());
477 int_args.add((Integer)params->getCriteriaStop());
478
479 // not broadcasted writeMatrixNameErrorStrategy
480 m_world_parallel->broadcast(real_args.view(), 0);
481 m_world_parallel->broadcast(bool_args.view(), 0);
482 m_world_parallel->broadcast(int_args.view(), 0);
483 }
484
485 debug() << "\33[1;31m[postSolver] Queuing solver " << m_solver_index << "\33[0m";
486 debug() << "\33[1;31m[postSolver] Backuping its params @" << params << "\33[0m";
487 m_arguments_queue.at(m_solver_index)->m_params = params;
488
489 // Mis à jour des indices solvers et sites
490 debug() << "\33[1;31m[postSolver] Mis à jour des indices solvers et sites"
491 << "\33[0m";
492 m_solver_index += 1;
493 debug() << "\33[1;31m[postSolver] m_solver_index=" << m_solver_index << "\33[0m";
494 traceMng()->flush();
495}
496
497/*---------------------------------------------------------------------------*/
498/*---------------------------------------------------------------------------*/
504syncSolver(Integer gid, Integer& nb_iteration, Real* residual_norm)
505{
506 ItacFunction(AlephKernel);
507
508 if (m_there_are_idles && !m_i_am_an_other) {
509 m_world_parallel->broadcast(UniqueArray<unsigned long>(1, 0xbf8d3adfl).view(), 0);
510 m_world_parallel->broadcast(UniqueArray<Integer>(1, gid).view(), 0);
511 }
512
513 if (!m_solved) {
514 debug() << "\33[1;31m[syncSolver] Not solved, launching the work"
515 << "\33[0m";
516 workSolver();
517 m_solved = true;
518 }
519
520 debug() << "\33[1;31m[syncSolver] Syncing " << gid << "\33[0m";
521 AlephVector* aleph_vector_x = m_arguments_queue.at(gid)->m_x_vector;
522 aleph_vector_x->reassemble_waitAndFill();
523 m_matrix_queue.at(gid)->reassemble_waitAndFill(m_results_queue.at(gid)->m_nb_iteration,
524 m_results_queue.at(gid)->m_residual_norm);
525
526 debug() << "\33[1;31m[syncSolver] Finishing " << gid << "\33[0m";
527 nb_iteration = m_results_queue.at(gid)->m_nb_iteration;
528 residual_norm[0] = m_results_queue.at(gid)->m_residual_norm[0];
529 residual_norm[1] = m_results_queue.at(gid)->m_residual_norm[1];
530 residual_norm[2] = m_results_queue.at(gid)->m_residual_norm[2];
531 residual_norm[3] = m_results_queue.at(gid)->m_residual_norm[3];
532
533 if (getTopologyImplementation(gid) != NULL)
534 getTopologyImplementation(gid)->restore();
535
536 debug() << "\33[1;31m[syncSolver] Done " << gid << "\33[0m";
537 traceMng()->flush();
538 return m_arguments_queue.at(gid)->m_x_vector;
539}
540
541/*---------------------------------------------------------------------------*/
542/*---------------------------------------------------------------------------*/
543
544void AlephKernel::
545workSolver(void)
546{
547 ItacFunction(AlephKernel);
548
549 debug() << "\33[1;31m[workSolver] Now working"
550 << "\33[0m";
551 for (int gid = 0; gid < m_solver_index; ++gid) {
553 debug() << "\33[1;31m[workSolver] Waiting for assembling " << gid << "\33[0m";
554 AlephVector* aleph_vector_x = m_arguments_queue.at(gid)->m_x_vector;
555 AlephVector* aleph_vector_b = m_arguments_queue.at(gid)->m_b_vector;
556 AlephMatrix* aleph_matrix_A = m_matrix_queue.at(gid);
557 aleph_matrix_A->assemble_waitAndFill();
558 aleph_vector_b->assemble_waitAndFill();
559 aleph_vector_x->assemble_waitAndFill();
560 traceMng()->flush();
561 }
562
563 if (!m_configured) {
564 debug() << "\33[1;31m[workSolver] NOW CONFIGURED!"
565 << "\33[0m";
566 traceMng()->flush();
567 m_configured = true;
568 }
569
570 for (int gid = 0; gid < m_solver_index; ++gid) {
571 ItacRegion(gidSolving, AlephKernel);
572 debug() << "\33[1;31m[workSolver] Solving " << gid << " ?"
573 << "\33[0m";
574 if (getTopologyImplementation(gid) != NULL)
575 getTopologyImplementation(gid)->backupAndInitialize();
576 m_matrix_queue.at(gid)->solveNow(m_arguments_queue.at(gid)->m_x_vector,
577 m_arguments_queue.at(gid)->m_b_vector,
578 m_arguments_queue.at(gid)->m_tmp_vector,
579 m_results_queue.at(gid)->m_nb_iteration,
580 m_results_queue.at(gid)->m_residual_norm,
581 m_arguments_queue.at(gid)->m_params);
582 // Après le solve, on 'restore' la session
583 if (getTopologyImplementation(gid) != NULL)
584 getTopologyImplementation(gid)->restore();
585 }
586
587 for (int gid = 0; gid < m_solver_index; ++gid) {
588 ItacRegion(gidReAssemble, AlephKernel);
589 debug() << "\33[1;31m[workSolver] Posting re-assembling " << gid << "\33[0m";
590 m_arguments_queue.at(gid)->m_x_vector->reassemble();
591 m_matrix_queue.at(gid)->reassemble(m_results_queue.at(gid)->m_nb_iteration,
592 m_results_queue.at(gid)->m_residual_norm);
593 }
594 traceMng()->flush();
595 // Flush number of pending solvers
596 m_solver_index = 0;
597}
598
599/*---------------------------------------------------------------------------*/
600/*---------------------------------------------------------------------------*/
601
602ISubDomain* AlephKernel::
603subDomain(void)
604{
605 if (!m_sub_domain && !m_i_am_an_other)
606 ARCANE_FATAL("[AlephKernel::subDomain]", "No sub-domain to work on!");
607 return m_sub_domain;
608}
609
610/*---------------------------------------------------------------------------*/
611/*---------------------------------------------------------------------------*/
612
613} // End namespace Arcane
614
615/*---------------------------------------------------------------------------*/
616/*---------------------------------------------------------------------------*/
#define ARCANE_FATAL(...)
Macro envoyant une exception FatalErrorException.
Gestionaire d'indexing.
AlephKernel(IParallelMng *, Integer, IAlephFactory *, Integer=0, Integer=0, bool=false)
AlephVector * createSolverVector(void)
AlephVector * syncSolver(Integer, Integer &, Real *)
void mapranks(Array< Integer > &)
Matrice d'un système linéaire.
Definition AlephMatrix.h:33
Gestionnaire de reordering.
Paramètres d'un système linéraire.
Definition AlephParams.h:34
Informations sur l'environnement parallèle.
Vecteur d'un système linéaire.
Definition AlephVector.h:33
Interface du gestionnaire de parallélisme pour un sous-domaine.
virtual Ref< IParallelMng > createSubParallelMngRef(Int32ConstArrayView kept_ranks)=0
Créé un nouveau gestionnaire de parallélisme pour un sous-ensemble des rangs.
Interface du gestionnaire d'un sous-domaine.
Definition ISubDomain.h:74
Lecteur des fichiers de maillage via la bibliothèque LIMA.
Definition Lima.cc:120
Exception lorsqu'une erreur fatale est survenue.
Interface du gestionnaire de traces.
virtual void flush()=0
Flush tous les flots.
ITraceMng * traceMng() const
Gestionnaire de trace.
TraceMessageDbg debug(Trace::eDebugLevel=Trace::Medium) const
Flot pour un message de debug.
TraceMessage info() const
Flot pour un message d'information.
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
List< String > StringList
Tableau de chaînes de caractères unicode.
Definition UtilsTypes.h:667