14#include "arcane/utils/StringBuilder.h"
15#include "arcane/utils/Iostream.h"
16#include "arcane/utils/PlatformUtils.h"
17#include "arcane/utils/Convert.h"
18#include "arcane/utils/NotImplementedException.h"
19#include "arcane/utils/ArgumentException.h"
20#include "arcane/utils/FloatingPointExceptionSentry.h"
22#include "arcane/ISubDomain.h"
23#include "arcane/IParallelMng.h"
24#include "arcane/ItemEnumerator.h"
25#include "arcane/IMesh.h"
26#include "arcane/IMeshSubMeshTransition.h"
27#include "arcane/ItemGroup.h"
28#include "arcane/Service.h"
29#include "arcane/Timer.h"
30#include "arcane/FactoryService.h"
31#include "arcane/ItemPrinter.h"
32#include "arcane/IItemFamily.h"
33#include "arcane/MeshVariable.h"
34#include "arcane/VariableBuildInfo.h"
36#include "arcane/std/MeshPartitionerBase.h"
37#include "arcane/std/MetisMeshPartitioner_axl.h"
40#define MPICH_SKIP_MPICXX
41#define OMPI_SKIP_MPICXX
45#if (PARMETIS_MAJOR_VERSION < 4)
46#error "Your version of Parmetis is too old .Version 4.0+ is required"
49typedef real_t realtype;
52#include "arcane/std/PartitionConverter.h"
53#include "arcane/std/GraphDistributor.h"
54#include "arcane/std/internal/MetisWrapper.h"
64using MetisCallStrategy = TypesMetisMeshPartitioner::MetisCallStrategy;
65using MetisEmptyPartitionStrategy = TypesMetisMeshPartitioner::MetisEmptyPartitionStrategy;
72class MetisMeshPartitioner
94 bool m_disable_floatingexception =
false;
105 const String& Name)
const;
111MetisMeshPartitioner::
112MetisMeshPartitioner(
const ServiceBuildInfo& sbi)
113: ArcaneMetisMeshPartitionerObject(sbi)
115 m_parallel_mng = mesh()->parallelMng();
117 if (s==
"1" || s==
"TRUE")
118 m_disable_floatingexception =
true;
127 Int32 nb_part = m_parallel_mng->commSize();
140 if (m_disable_floatingexception){
159 bool force_partition =
false;
166 info() <<
"INFO: Unable to test load balancing on a single sub-domain";
174 bool force_full_repartition =
false;
176 force_full_repartition |= (m_nb_refine < 0);
178 if (nb_part != nb_rank) {
179 force_full_repartition =
true;
180 force_partition =
true;
183 info() <<
"WARNING: compensating the potential lack of 'Metis' options in case of manual instanciation";
184 Integer max_diffusion_count = 10;
185 Real imbalance_relative_tolerance = 4;
186 float tolerance_target = 1.05f;
187 bool dump_graph =
false;
188 MetisCallStrategy call_strategy = MetisCallStrategy::one_processor_per_node;
189 bool in_out_digest =
false;
194 if (call_strategy_env ==
"all-processors"){
195 call_strategy = MetisCallStrategy::all_processors;
196 }
else if (call_strategy_env ==
"one-processor-per-node") {
197 call_strategy = MetisCallStrategy::one_processor_per_node;
198 }
else if (call_strategy_env ==
"two-processors-two-nodes") {
199 call_strategy = MetisCallStrategy::two_processors_two_nodes;
200 }
else if (call_strategy_env ==
"two-gathered-processors") {
201 call_strategy = MetisCallStrategy::two_gathered_processors;
202 }
else if (call_strategy_env ==
"two-scattered-processors") {
203 call_strategy = MetisCallStrategy::two_scattered_processors;
204 }
else if (!call_strategy_env.
null()) {
205 ARCANE_FATAL(
"Invalid value '{0}' for ARCANE_METIS_CALL_STRATEGY environment variable",call_strategy_env);
209 in_out_digest = (v.value()!=0);
211 dump_graph = (v.value()!=0);
214 max_diffusion_count =
options()->maxDiffusiveCount();
215 imbalance_relative_tolerance =
options()->imbalanceRelativeTolerance();
216 tolerance_target = (float)(
options()->toleranceTarget());
217 dump_graph =
options()->dumpGraph();
218 call_strategy =
options()->metisCallStrategy();
219 in_out_digest =
options()->inputOutputDigest();
222 if (max_diffusion_count > 0)
223 force_full_repartition |= (m_nb_refine>=max_diffusion_count);
225 force_full_repartition |= (
imbalance()>1.0);
228 initConstraints(
false);
234 if (is_shared_memory)
235 call_strategy = MetisCallStrategy::one_processor_per_node;
237 idxtype nb_weight = nbCellWeight();
239 info() <<
"Load balancing with Metis nb_weight=" << nb_weight <<
" initial=" << initial_partition
240 <<
" call_strategy=" << (int)call_strategy
241 <<
" is_shared_memory?=" << is_shared_memory
242 <<
" disabling_fpe?=" << m_disable_floatingexception
243 <<
" sizeof(idxtype)==" <<
sizeof(idxtype);
246 initial_partition =
true;
257 Integer nb_own_cell = nbOwnCellsWithConstraints();
259 Int32 nb_empty_part = 0;
262 metis_vtkdist[0] = 0;
263 for(
Integer i=0; i<nb_rank; ++i ){
265 Int32 n = global_nb_own_cell[i];
269 metis_vtkdist[i+1] =
static_cast<idxtype
>(total_nb_cell);
274 info() <<
"Total nb_cell=" << total_nb_cell <<
" nb_empty_partition=" << nb_empty_part;
275 if (total_nb_cell==0){
276 info() <<
"INFO: mesh '" <<
mesh->name() <<
" has no cell. No partitioning is needed";
294 Integer nb_max_face_neighbour_cell = 0;
301 if (cellUsedWithConstraints(item)){
302 cell_metis_uid[item] = mid;
305 if (item.
hasFlags(ItemFlags::II_HasEdgeFor1DItems))
306 nb_max_face_neighbour_cell += item.
nbEdge();
308 nb_max_face_neighbour_cell += item.
nbFace();
310 cell_metis_uid.synchronize();
313 _initUidRef(cell_metis_uid);
316 cell_metis_uid.setUsed(
false);
319 metis_xadj.
reserve(nb_own_cell+1);
322 metis_adjncy.
reserve(nb_max_face_neighbour_cell);
333 if (!cellUsedWithConstraints(item))
336 metis_xadj.
add(metis_adjncy.
size());
338 getNeighbourCellsUidWithConstraints(item, neighbour_cells, &edge_weights);
340 metis_adjncy.
add(
static_cast<idxtype
>(neighbour_cells[z]));
342 metis_xadj.
add(metis_adjncy.
size());
346 idxtype numflags = 0;
347 idxtype nparts =
static_cast<int>(nb_part);
354 if (!s.
null() && (s.
upper() ==
"PARTICLES"))
360 bool cond = (nb_weight == 1);
376 metis_ubvec.
fill(tolerance_target);
381 cond = converter.isBalancable<realtype>(cells_weights, metis_ubvec, nb_part);
383 info() <<
"We have to increase imbalance :";
384 info() << metis_ubvec;
385 for (
auto& tol: metis_ubvec ) {
395 const bool do_print_weight =
false;
396 if (do_print_weight){
398 Integer iteration =
mesh->subDomain()->commonVariables().globalIteration();
405 for (
int i = 0; i < metis_xadj.
size()-1 ; ++i) {
406 dumpy <<
" Weight uid=" << i;
407 for(
int j=0 ; j < nb_weight ; ++j ){
408 dumpy <<
" w=" << *(metis_vwgt.
begin()+i*nb_weight+j)
409 <<
"(" << cells_weights[i*nb_weight+j] <<
")";
417 converter.computeContrib(edge_weights);
424 if (!initial_partition) {
425 cells_size = cellsSizeWithConstraints();
426 converter.computeContrib(cells_size, (
Real)itr);
431 idxtype metis_options[4];
432 metis_options[0] = 0;
440 metis_options[0] = 1;
441 metis_options[1] = 0;
442 metis_options[2] = m_random_seed;
443 metis_options[3] = 1;
449 idxtype metis_edgecut = 0;
454 std::chrono::high_resolution_clock clock;
455 auto start_time = clock.now();
467 bool redistribute =
true;
468 bool redistribute_by_wrapper =
false;
470 if (call_strategy == MetisCallStrategy::all_processors || call_strategy == MetisCallStrategy::two_scattered_processors) {
471 redistribute =
false;
474 if (call_strategy == MetisCallStrategy::two_processors_two_nodes || call_strategy == MetisCallStrategy::two_scattered_processors) {
475 redistribute_by_wrapper =
true;
492 bool gd_allow_only_one_rank =
false;
493 if (is_shared_memory || (nb_empty_part!=0 && initial_partition))
494 gd_allow_only_one_rank =
true;
499 if ((nb_empty_part+1)==nb_rank){
500 info() <<
"Initialize GraphDistributor with max rank=1";
501 gd.initWithMaxRank(1);
503 else if (call_strategy == MetisCallStrategy::two_gathered_processors && (nb_rank > 2)) {
505 info() <<
"Initialize GraphDistributor with max rank=2";
506 gd.initWithMaxRank(2);
509 info() <<
"Initialize GraphDistributor with one rank per node"
510 <<
" (allow_one?=" << gd_allow_only_one_rank <<
")";
516 info() <<
"Using GraphDistributor?=" << redistribute;
518 metis_xadj = gd.convert<idxtype>(metis_xadj, &metis_part,
true);
519 metis_vwgt = gd.convert<idxtype>(metis_vwgt);
520 metis_adjncy = gd.convert<idxtype>(metis_adjncy);
521 metis_ewgt = gd.convert<idxtype>(metis_ewgt);
522 if (!initial_partition)
523 metis_vsize = gd.convert<idxtype>(metis_vsize);
524 metis_pm = gd.subParallelMng();
526 metis_vtkdist.
resize(gd.size()+1);
527 if (gd.contribute()) {
529 Integer nbVertices = metis_part.size();
531 metis_vtkdist[0] = 0;
532 for (
int i = 1 ; i < metis_vtkdist.
size() ; ++i) {
533 metis_vtkdist[i] = metis_vtkdist[i-1] + buffer[i-1];
536 metis_options[3] = PARMETIS_PSR_UNCOUPLED ;
542 MPI_Comm metis_mpicomm =
static_cast<MPI_Comm
>(metis_pm->
communicator());
544 bool do_call_metis =
true;
546 do_call_metis = gd.contribute();
553 _writeGraph(metis_pm, metis_vtkdist, metis_xadj, metis_adjncy, metis_vwgt, metis_ewgt, name.
toString());
559 float fill_value = (float)(1.0/(
double)nparts);
560 tpwgt.
fill(fill_value);
562 int retval = METIS_OK;
566 force_partition |= initial_partition;
567 force_full_repartition |= force_partition;
568 if (force_full_repartition){
570 if (force_partition) {
571 info() <<
"Metis: use a complete partitioning.";
573 redistribute_by_wrapper,
574 metis_vtkdist.
data(),
579 &wgtflag,&numflags,&nb_weight,
580 &nparts, tpwgt.
data(),
582 metis_options,&metis_edgecut,
586 info() <<
"Metis: use a complete REpartitioning.";
588 redistribute_by_wrapper,
589 metis_vtkdist.
data(),
595 &wgtflag,&numflags,&nb_weight,
596 &nparts,tpwgt.
data(),
598 &itr, metis_options,&metis_edgecut,
606 info() <<
"Metis: use a diffusive REpartitioning";
607 retval = ParMETIS_V3_RefineKway(metis_vtkdist.
data(),
612 &wgtflag,&numflags,&nb_weight,
613 &nparts,tpwgt.
data(),
615 metis_options,&metis_edgecut,
619 if (retval != METIS_OK) {
620 ARCANE_FATAL(
"Error while computing ParMetis partitioning r='{0}'",retval);
624 metis_part = gd.convertBack<idxtype>(metis_part, nb_own_cell);
627 auto end_time = clock.now();
628 std::chrono::microseconds duration = std::chrono::duration_cast<std::chrono::microseconds>(end_time - start_time);
629 info() <<
"Time to partition using parmetis = " << duration.count() <<
" us ";
635 switch(
options()->emptyPartitionStrategy()){
636 case MetisEmptyPartitionStrategy::DoNothing:
638 case MetisEmptyPartitionStrategy::TakeFromBiggestPartitionV1:
641 case MetisEmptyPartitionStrategy::TakeFromBiggestPartitionV2:
654 if (!cellUsedWithConstraints(item))
659 changeCellOwner(item, cells_new_owner, new_owner);
666 cells_new_owner.synchronize();
699 for (
int i =0 ; i < nb_own_cell ; i++) {
702 pm->
computeMinMaxSum(elem_by_part, min_part, max_part, sum_part, min_roc, max_proc);
705 Int32 max_part_id = -1;
706 Int32 max_part_nbr = -1;
708 for(
int i = 0; i < nb_part ; i++) {
709 if (sum_part[i] == 0) {
712 if(max_part_nbr < sum_part[i]) {
713 max_part_nbr = sum_part[i];
717 info() <<
"Parmetis: number empty parts " << nb_hole;
721 if(my_rank == max_proc[max_part_id]) {
723 for(
int i = 0; i < nb_part ; i++) {
725 if (sum_part[i] == 0 && offset < nb_own_cell) {
726 while(offset < nb_own_cell && metis_part[offset] != max_part_id) {
730 if(offset == nb_own_cell)
733 metis_part[offset] = i;
754 const Int32 nb_own_cell = metis_part.
size();
762 for (
int i =0 ; i < nb_own_cell ; i++) {
763 elem_by_part[metis_part[i]]++;
765 pm->
computeMinMaxSum(elem_by_part, min_part, max_part, sum_part, min_rank, max_rank);
770 Int32 max_part_id = -1;
771 Int64 max_part_nbr = -1;
773 Int64 total_nb_cell = 0;
774 for(
int i = 0; i < nb_part ; i++) {
775 Int64 current_sum_part = sum_part[i];
776 total_nb_cell += current_sum_part;
777 if (current_sum_part == 0) {
778 empty_part_ranks.
add(i);
781 if (max_part_nbr < current_sum_part) {
782 max_part_nbr = current_sum_part;
787 ARCANE_FATAL(
"Bad value max_part_id ({0})",max_part_id);
788 info() <<
"Parmetis: check empty parts: (" << algo_iteration <<
") nb_empty_parts=" << nb_hole
789 <<
" nb_max_part=" << max_part_nbr
790 <<
" max_part_rank=" << max_part_id
791 <<
" max_proc_max_part_id=" << max_rank[max_part_id]
792 <<
" empty_part_ranks=" << empty_part_ranks
793 <<
" total_nb_cell=" << total_nb_cell;
800 if (my_rank == max_rank[max_part_id]) {
802 Int32 max_remove_cell = nb_own_cell - 1;
804 for(
int i = 0; i < nb_part ; i++) {
806 if (sum_part[i] == 0 && offset < nb_own_cell) {
807 while (offset < max_remove_cell && metis_part[offset] != max_part_id) {
811 if (offset == max_remove_cell)
814 metis_part[offset] = i;
834 bool do_continue =
true;
835 Int32 last_nb_hole = 0;
838 info() <<
"Parmetis: nb_empty_partition=" << nb_hole <<
" last_nb_partition=" << last_nb_hole;
845 if (last_nb_hole>0 && last_nb_hole<=nb_hole){
846 pwarning() <<
"Can not remove all empty partitions. This is probably because you try"
847 <<
" to cut in too many partitions";
850 last_nb_hole = nb_hole;
873 bool have_ewgt =
false;
877 info() <<
"COMM_SIZE=" << commsize <<
" RANK=" << commrank;
885 info() <<
"MetisVtkDist=" << metis_vtkdist;
886 info() <<
"MetisXAdj =" << metis_xadj;
887 info() <<
"MetisAdjncy =" << metis_adjncy;
888 info() <<
"MetisVWgt =" << metis_vwgt;
889 info() <<
"MetisEWgt =" << metis_ewgt;
891#define METIS_ERROR ARCANE_FATAL("_writeGraph")
895 filename += commrank;
898 if (metis_vtkdist.
size() != commsize + 1)
901 nvtx = metis_vtkdist[commrank+1] - metis_vtkdist[commrank];
902 if (metis_xadj.
size() != nvtx + 1) {
903 std::cerr <<
"_writeGraph : nvtx+1 = " << nvtx <<
" != " << metis_xadj.
size() << std::endl;
908 nwgt = metis_vwgt.
size()/nvtx;
909 if (nwgt != 0 && metis_vwgt.
size() % nwgt != 0)
912 have_ewgt = (metis_ewgt.
size() != 0) ;
913 if (have_ewgt && metis_ewgt.
size() != metis_adjncy.
size())
919 Int64 localEdges = metis_xadj[nvtx];
922 file <<
"2" << std::endl;
923 file << commsize <<
"\t" << commrank << std::endl;
924 file << metis_vtkdist[commsize] <<
"\t" << globalEdges << std::endl;
925 file << nvtx <<
"\t" << localEdges << std::endl;
926 file <<
"0\t" <<
"0" << have_ewgt << nwgt << std::endl;
937 for (idxtype vertnum = 0 ; vertnum < nvtx ; vertnum++) {
939 for (
int dim = 0 ; dim < nwgt ; ++ dim)
940 file << metis_vwgt[vertnum*nwgt+dim] <<
'\t';
941 file << metis_xadj[vertnum + 1] - metis_xadj[vertnum];
942 for (idxtype edgenum = metis_xadj[vertnum] ;
943 edgenum < metis_xadj[vertnum + 1] ; ++edgenum) {
945 file <<
'\t' << metis_ewgt[edgenum];
946 file <<
'\t' << metis_adjncy[edgenum];
967#if ARCANE_DEFAULT_PARTITIONER == METIS_DEFAULT_PARTITIONER
#define ARCANE_THROW(exception_class,...)
Macro pour envoyer une exception avec formattage.
#define ARCANE_FATAL(...)
Macro envoyant une exception FatalErrorException.
#define ARCANE_SERVICE_INTERFACE(ainterface)
Macro pour déclarer une interface lors de l'enregistrement d'un service.
Integer size() const
Nombre d'éléments du vecteur.
ArcaneMetisMeshPartitionerObject(const Arcane::ServiceBuildInfo &sbi)
Constructeur.
CaseOptionsMetisMeshPartitioner * options() const
Options du jeu de données du service.
Exception lorsqu'un argument est invalide.
Conversion d'un tableau d'un type vers un autre type.
Vue modifiable d'un tableau d'un type T.
constexpr Integer size() const noexcept
Retourne la taille du tableau.
void fill(const DataType &data)
Remplissage du tableau.
ConstArrayView< T > constView() const
Vue constante sur ce tableau.
void resize(Int64 s)
Change le nombre d'éléments du tableau à s.
void reserve(Int64 new_capacity)
Réserve le mémoire pour new_capacity éléments.
const T * data() const
Accès à la racine du tableau hors toute protection.
iterator begin()
Itérateur sur le premier élément du tableau.
void add(ConstReferenceType val)
Ajoute l'élément val à la fin du tableau.
Int32 nbEdge() const
Nombre d'arêtes de la maille.
Int32 nbFace() const
Nombre de faces de la maille.
Int32 globalIteration() const
Numéro de l'itération courante.
Vue constante d'un tableau de type T.
Classe permettant d'activer/désactiver temporairement les exceptions flottantes du processeur.
Redistribute graph data to another "communicator".
void initWithOneRankPerNode(bool allow_only_one_rank)
Automatic distribution : one partitioning process per node.
Interface d'un partitionneur de maillage.
Interface d'un partitionneur de maillage.
Interface du gestionnaire de parallélisme pour un sous-domaine.
virtual bool isThreadImplementation() const =0
Indique si l'implémentation utilise les threads.
virtual void computeMinMaxSum(char val, char &min_val, char &max_val, char &sum_val, Int32 &min_rank, Int32 &max_rank)=0
Calcule en une opération la somme, le min, le max d'une valeur.
virtual Int32 commRank() const =0
Rang de cette instance dans le communicateur.
virtual Int32 commSize() const =0
Nombre d'instance dans le communicateur.
virtual void allGather(ConstArrayView< char > send_buf, ArrayView< char > recv_buf)=0
Effectue un regroupement sur tous les processeurs. Il s'agit d'une opération collective....
virtual Parallel::Communicator communicator() const =0
Communicateur MPI associé à ce gestionnaire.
virtual char reduce(eReduceType rt, char v)=0
Effectue la réduction de type rt sur le réel v et retourne la valeur.
virtual void barrier()=0
Effectue une barière.
Interface du gestionnaire d'un sous-domaine.
virtual const CommonVariables & commonVariables() const =0
Informations sur les variables standards.
virtual void flush()=0
Flush tous les flots.
@ PNoDump
Indique que la variable ne doit pas être sauvegardée.
bool hasFlags(Int32 flags) const
Retourne si les flags flags sont positionnées pour l'entité
Real imbalance() const override
Déséquilibre de temps de calcul.
Real maxImbalance() const override
Déséquilibre maximal autorisé
virtual void changeOwnersFromCells()
Positionne les nouveaux propriétaires des noeuds, arêtes et faces à partir des mailles.
Partitioneur de maillage utilisant la bibliothèque PARMetis.
void partitionMesh(bool initial_partition) override
void _removeEmptyPartsV2(Int32 nb_part, ArrayView< idxtype > metis_part)
Comble les partitions vides (version 2).
Int32 _removeEmptyPartsV2Helper(Int32 nb_part, ArrayView< idxtype > metis_part, Int32 algo_iteration)
Applique une itération de l'algorithme de suppression des partitions vides.
void build() override
Construction de niveau build du service.
void _partitionMesh(bool initial_partition, Int32 nb_part)
int _writeGraph(IParallelMng *pm, Span< const idxtype > metis_vtkdist, Span< const idxtype > metis_xadj, Span< const idxtype > metis_adjncy, Span< const idxtype > metis_vwgt, Span< const idxtype > metis_ewgt, const String &Name) const
void _removeEmptyPartsV1(Int32 nb_part, Int32 nb_own_cell, ArrayView< idxtype > metis_part)
Comble les partitions vides (version 1).
Wrapper autour des appels de Parmetis.
int callPartKway(const bool print_digest, const bool gather, idx_t *vtxdist, idx_t *xadj, idx_t *adjncy, idx_t *vwgt, idx_t *adjwgt, idx_t *wgtflag, idx_t *numflag, idx_t *ncon, idx_t *nparts, real_t *tpwgts, real_t *ubvec, idx_t *options, idx_t *edgecut, idx_t *part)
Simple wrapper autour de la routine ParMetis "ParMETIS_V3_PartKway".
int callAdaptiveRepart(const bool print_digest, const bool gather, idx_t *vtxdist, idx_t *xadj, idx_t *adjncy, idx_t *vwgt, idx_t *vsize, idx_t *adjwgt, idx_t *wgtflag, idx_t *numflag, idx_t *ncon, idx_t *nparts, real_t *tpwgts, real_t *ubvec, real_t *ipc2redist, idx_t *options, idx_t *edgecut, idx_t *part)
Simple wrapper autour de la routine ParMetis "ParMETIS_V3_AdaptiveRepart".
Conversion d'un tableau de flottants vers un tableau d'entiers/longs. \abstract Cette classe gere le ...
Structure contenant les informations pour créer un service.
Propriétés de création d'un service.
Vecteur 1D de données avec sémantique par référence.
constexpr __host__ __device__ SizeType size() const noexcept
Retourne la taille du tableau.
Vue d'un tableau d'éléments de type T.
Constructeur de chaîne de caractère unicode.
String toString() const
Retourne la chaîne de caractères construite.
Chaîne de caractères unicode.
bool null() const
Retourne true si la chaîne est nulle.
String upper() const
Transforme tous les caractères de la chaîne en majuscules.
const char * localstr() const
Retourne la conversion de l'instance dans l'encodage UTF-8.
Postionne le nom de l'action en cours d'exécution.
TraceMessage info() const
Flot pour un message d'information.
ITraceMng * traceMng() const
Gestionnaire de trace.
TraceMessage pwarning() const
Vecteur 1D de données avec sémantique par valeur (style STL).
Paramètres nécessaires à la construction d'une variable.
ItemGroupT< Cell > CellGroup
Groupe de mailles.
#define ARCANE_REGISTER_SERVICE(aclass, a_service_property,...)
Macro pour enregistrer un service.
MeshVariableScalarRefT< Cell, Integer > VariableCellInteger
Grandeur au centre des mailles de type entier.
ItemVariableScalarRefT< Int32 > VariableItemInt32
Grandeur de type entier 32 bits.
Integer toInteger(Real r)
Converti un Int64 en un Integer.
Int32 toInt32(Int64 v)
Converti un Int64 en un Int32.
Integer toInteger(Real r)
Converti un Real en Integer.
@ ReduceSum
Somme des valeurs.
@ ReduceMax
Maximum des valeurs.
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
UniqueArray< Int64 > Int64UniqueArray
Tableau dynamique à une dimension d'entiers 64 bits.
std::int64_t Int64
Type entier signé sur 64 bits.
Int32 Integer
Type représentant un entier.
@ ST_SubDomain
Le service s'utilise au niveau du sous-domaine.
@ IK_Cell
Entité de maillage de genre maille.
double Real
Type représentant un réel.
std::int32_t Int32
Type entier signé sur 32 bits.