TrioCFD 1.9.8
TrioCFD documentation
Loading...
Searching...
No Matches
Partitionneur_Parmetis.cpp
1/****************************************************************************
2* Copyright (c) 2026, CEA
3* All rights reserved.
4*
5* Redistribution and use in source and binary forms, with or without modification, are permitted provided that the following conditions are met:
6* 1. Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer.
7* 2. Redistributions in binary form must reproduce the above copyright notice, this list of conditions and the following disclaimer in the documentation and/or other materials provided with the distribution.
8* 3. Neither the name of the copyright holder nor the names of its contributors may be used to endorse or promote products derived from this software without specific prior written permission.
9*
10* THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED.
11* IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS;
12* OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
13*
14*****************************************************************************/
15#include <Partitionneur_Parmetis.h>
16#include <Domaine.h>
17#include <Static_Int_Lists.h>
18#include <parmetis++.h>
19#include <Param.h>
20
21#include <communications.h>
22#include <Poly_geom_base.h>
23#include <Matrix_tools.h>
24#include <Matrice_Morse.h>
25#include <Array_tools.h>
26#include <Comm_Group_MPI.h>
27
28#include <Domain_Graph.h>
29#include <MD_Vector_tools.h>
30
31inline void not_implemented(const Nom& chaine)
32{
33 Cerr << chaine << " is not implemented yet to the PARMETIS API." << finl;
35}
36
37// XD partitionneur_parmetis partitionneur_deriv parmetis INHERITS_BRACE ParMETIS is the MPI-parallel counterpart of
38// XD_CONT METIS — same general-purpose mesh partitioning, run in parallel.
39Implemente_instanciable_sans_constructeur(Partitionneur_Parmetis,"Partitionneur_Parmetis",Partitionneur_base);
40
41Partitionneur_Parmetis::Partitionneur_Parmetis()
42{
43}
44
46{
47 Cerr << "Partitionneur_Parmetis::printOn invalid\n" << finl;
48 exit();
49 return os;
50}
51
53{
54 param.ajouter("nb_parts",&nb_parties_,Param::REQUIRED);
55 param.ajouter_flag("use_weights",&use_weights_); // XD attr use_weights rien use_weights OPT Weight element-element
56 // XD_CONT links in the graph to keep periodic opposite elements on the same processor. Same meaning as in metis;
57 // XD_CONT costs extra memory and time but is otherwise non-mandatory (a correction pass runs unconditionally).
58}
59
61{
62 if (nb_parties_ < 1 || nb_parties_ > 100000)
63 {
64 Process::exit("Partitionneur_Parmetis::validate_params: The following condition must be satisfied : 1 <= nb_parts <= 100000");
65 }
66}
67
72
73void Partitionneur_Parmetis::associer_domaine(const Domaine& domaine)
74{
75 ref_domaine_ = domaine;
76}
77
78/*! @brief Calcule le graphe de connectivite pour parmetis, appelle le partitionneur et remplit elem_part (pour chaque element, numero de la partie qui lui
79 *
80 * est attribuee).
81 *
82 */
83void Partitionneur_Parmetis::construire_partition(IntVect& elem_part, int& nb_parts_tot) const
84{
85#ifndef PARMETIS_
86 Cerr << "PARMETIS is not compiled with this version. Use another partition tool like Tranche." << finl;
88#else
89
90 if (!ref_domaine_)
91 {
92 Cerr << "Error in Partitionneur_Parmetis::construire_partition\n";
93 Cerr << " The domain has not been associated" << finl;
94 exit();
95 }
96 if (nb_parties_ <= 0)
97 {
98 Cerr << "Error in Partitionneur_Parmetis::construire_partition\n";
99 Cerr << " The parts number has not been initialized" << finl;
100 exit();
101 }
102
103 // Cas particulier: si nb_parts == 1, METIS ne veut rien faire...
104 //ToDo: I don't know is that's the case with ParMetis
105 if (nb_parties_ == 1)
106 {
107
108 int nb_elem = ref_domaine_->nb_elem();
109 elem_part.resize_array(nb_elem);
110 elem_part = 0;
111 return;
112 }
113
114 if (ref_domaine_->nb_elem() == 0)
115 return;
116
117 Cerr << "Partitionneur_Parmetis::construire_partition" << finl;
118 Cerr << " Construction of graph connectivity..." << finl;
119 Static_Int_Lists graph_elements_perio;
120 //const Domaine& dom = ref_domaine_.valeur();
121 Domain_Graph graph;
122 graph.construire_graph_elem_elem(ref_domaine_.valeur(),
123 use_weights_,
124 graph_elements_perio);
125
126
127 std::vector<idx_t> partition(graph.nvtxs);
128
129 idx_t int_parts = nb_parties_;
130 idx_t edgecut = 0; // valeur renvoyee par parmetis (nombre total de faces de joint)
131 Cerr << " Call for PARMETIS" << finl;
132 idx_t options[3];
133 options[0] = 1; //personnalized options
134 options[1] = 111111111; // Mode verbose maximal;
135 options[2] = 0; //random seed
136
137 idx_t ncon=1;
138 real_t ubvec = 1.05f; //recommanded value
139 idx_t numflag = 0; //numerotation C
140 std::vector<real_t> tpwgts(ncon*int_parts, (real_t)(1.0/nb_parties_)); //we want the weight to be equally distributed on each sub_somain
141 MPI_Comm comm = Comm_Group_MPI::get_trio_u_world();
142 int status = ParMETIS_V3_PartKway(graph.vtxdist.addr(), graph.xadj.addr(), graph.adjncy.addr(),
143 graph.vwgts.addr(), graph.ewgts.addr(), &graph.weightflag,
144 &numflag, &ncon, &int_parts, tpwgts.data(), &ubvec, options ,
145 &edgecut, partition.data(), &comm);
146 if (status != METIS_OK)
147 {
148 Cerr << "Call to PARMETIS failed." << finl;
149 if (status == METIS_ERROR) Cerr << "It seems there is a PARMETIS internal error." << finl;
150 Cerr << "Contact TRUST support." << finl;
151 exit();
152 }
153
154 Cerr << "Partitioning quality : edgecut = " << edgecut << finl;
155 Cerr << "-> It is roughly the total number of edges (faces) which will be shared by the processors." << finl;
156 Cerr << "-> The lesser this number is, the lesser the total volume of communication between processors." << finl;
157 Cerr << "-> You can increase nb_essais option (default 1) to try to reduce (but at a higher CPU cost) this number." << finl;
158 Cerr << "===============" << finl;
159
160 MD_Vector_tools::creer_tableau_distribue(ref_domaine_->md_vector_elements(), elem_part);
161 const int n = ref_domaine_->nb_elem();
162 for (int i = 0; i < n; i++)
163 elem_part[i] = static_cast<int>(partition[i]); // partition[i] is a a proc number...
164
165 // Correction de la partition pour la periodicite. (***)
166 if (graph_elements_perio.get_nb_lists() > 0)
167 {
168 Cerr << "Correction of the partition for the periodicity" << finl;
169 corriger_bords_avec_liste(ref_domaine_.valeur(),
170 0,
171 elem_part);
172 Cerr << " If this number is high, we can improve the splitting with the option use_weights\n"
173 << " but it takes more memory)" << finl;
174 }
175
176 Cerr << "Correction elem0 on processor 0" << finl;
177 corriger_elem0_sur_proc0(elem_part);
178 elem_part.echange_espace_virtuel();
179#endif
180}
181
Build the graph of the domain that the METIS/PARMETIS/PTSCOTCH libraries need.
ArrOfIdx_ vwgts
void construire_graph_elem_elem(const Domaine_32_64< _SIZE_ > &dom, bool use_weights, Static_Int_Lists_32_64< _SIZE_ > &graph_elements_perio)
ArrOfIdx_ ewgts
ArrOfIdx_ adjncy
ArrOfIdx_ vtxdist
ArrOfIdx_ xadj
static void creer_tableau_distribue(const MD_Vector &, Array_base &, RESIZE_OPTIONS opt=RESIZE_OPTIONS::COPY_INIT)
transforme v en un tableau parallele ayant la structure md.
Une chaine de caractere (Nom) en majuscules.
Definition Motcle.h:26
class Nom Une chaine de caractere pour nommer les objets de TRUST
Definition Nom.h:31
virtual void set_param(Param &) const
Definition Objet_U.h:135
friend class Entree
Definition Objet_U.h:76
virtual Sortie & printOn(Sortie &) const
Ecriture de l'objet sur un flot de sortie Methode a surcharger.
Definition Objet_U.cpp:282
Helper class to factorize the readOn method of Objet_U classes.
Definition Param.h:112
void ajouter_flag(const char *keyword, const bool *value)
Register a boolean flag whose mere presence switches it to true.
Definition Param.cpp:474
void ajouter(const char *keyword, const int *value, Param::Nature nat=Param::OPTIONAL)
Register an integer parameter.
Definition Param.cpp:364
@ REQUIRED
Definition Param.h:115
Partition d'un domaine en nb_parties parties equilibrees en utilisant la librairie PARMETIS.
void associer_domaine(const Domaine &domaine) override
void construire_partition(IntVect &elem_part, int &nb_parts_tot) const override
Calcule le graphe de connectivite pour parmetis, appelle le partitionneur et remplit elem_part (pour ...
void validate_params() const override
Called in the readOn of Objet_U_With_Params, after reading the params.
int lire_motcle_non_standard(const Motcle &, Entree &) override
Lecture des parametres de type non simple d'un objet_U a partir d'un flot d'entree.
int lire_motcle_non_standard(const Motcle &, Entree &) override
static void corriger_bords_avec_liste(const Domaine_t &dom, const int_t my_offset, BigIntVect_ &elem_part)
static void corriger_elem0_sur_proc0(BigIntVect_ &elem_part)
static void exit(int exit_code=-1)
Routine de sortie de TRUST dans une region Kokkos.
Definition Process.cpp:455
Classe de base des flux de sortie.
Definition Sortie.h:52
int_t get_nb_lists() const
renvoie le nombre de listes stockees
_TYPE_ * addr()
void resize_array(_SIZE_ new_size, RESIZE_OPTIONS opt=RESIZE_OPTIONS::COPY_INIT)
virtual void echange_espace_virtuel(IsExchangeBlocking exchange_type=IsExchangeBlocking::DefaultBlocking, const std::string kernel_name="noname")