zoltanRenumber.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::zoltanRenumber
26 
27 Description
28  Renumber using Zoltan
29 
30  Zoltan install:
31  - in your ~/.bashrc:
32  export ZOLTAN_ARCH_DIR=\
33  $WM_THIRD_PARTY_DIR/platforms/linux64Gcc/Zoltan_XXX
34  - unpack into $WM_THIRD_PARTY_DIR
35  - cd Zoltan_XXX
36  - mkdir build
37  - cd build
38  - export CCFLAGS="-fPIC"
39  - export CXXFLAGS="-fPIC"
40  - export CFLAGS="-fPIC"
41  - export LDFLAGS="-shared"
42  - ../configure \
43  --prefix=$ZOLTAN_ARCH_DIR \
44  --with-ccflags=-fPIC --with-cxxflags=-fPIC --with-ldflags=-shared
45 
46 SourceFiles
47  zoltanRenumber.C
48 
49 \*---------------------------------------------------------------------------*/
50 
51 #include "zoltanRenumber.H"
53 #include "IFstream.H"
54 #include "labelIOList.H"
55 #include "polyMesh.H"
56 #include "globalMeshData.H"
57 #include "globalIndex.H"
58 #include "uint.H"
59 
60 #include "zoltan.h"
61 #include <mpi.h>
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 namespace Foam
66 {
67  defineTypeNameAndDebug(zoltanRenumber, 0);
68 
70  (
71  renumberMethod,
72  zoltanRenumber,
73  dictionary
74  );
75 }
76 
77 
78 static int get_number_of_vertices(void *data, int *ierr)
79 {
80  const Foam::polyMesh& mesh = *static_cast<const Foam::polyMesh*>(data);
81  *ierr = ZOLTAN_OK;
82  return mesh.nCells();
83 }
84 
85 
86 static void get_vertex_list(void *data, int sizeGID, int sizeLID,
87  ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID,
88  int wgt_dim, float *obj_wgts, int *ierr)
89 {
90  const Foam::polyMesh& mesh = *static_cast<const Foam::polyMesh*>(data);
91 
92  *ierr = ZOLTAN_OK;
93 
94  /* In this example, return the IDs of our vertices, but no weights.
95  * Zoltan will assume equally weighted vertices.
96  */
97 
98  wgt_dim = 0;
99  obj_wgts = NULL;
100 
101  for (Foam::label i=0; i<mesh.nCells(); i++)
102  {
103  globalID[i] = i; // should be global
104  localID[i] = i;
105  }
106 }
107 
109 static void get_num_edges_list(void *data, int sizeGID, int sizeLID,
110  int num_obj,
111  ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID,
112  int *numEdges, int *ierr)
113 {
114  const Foam::polyMesh& mesh = *static_cast<const Foam::polyMesh*>(data);
115 
116  if ((sizeGID != 1) || (sizeLID != 1) || (num_obj != mesh.nCells()))
117  {
118  *ierr = ZOLTAN_FATAL;
119  return;
120  }
121 
122  for (Foam::label i=0; i < num_obj ;i++)
123  {
124  Foam::label celli = localID[i];
125  const Foam::cell& cFaces = mesh.cells()[celli];
126  forAll(cFaces, cFacei)
127  {
128  Foam::label n = 0;
129  if (mesh.isInternalFace(cFaces[cFacei]))
130  {
131  n++;
132  }
133  numEdges[i] = n;
134  }
135  }
136 
137  *ierr = ZOLTAN_OK;
138 }
139 
141 static void get_edge_list(void *data, int sizeGID, int sizeLID,
142  int num_obj, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID,
143  int *num_edges,
144  ZOLTAN_ID_PTR nborGID, int *nborProc,
145  int wgt_dim, float *ewgts, int *ierr)
146 {
147  const Foam::polyMesh& mesh = *static_cast<const Foam::polyMesh*>(data);
148 
149  if
150  (
151  (sizeGID != 1)
152  || (sizeLID != 1)
153  || (num_obj != mesh.nCells())
154  || (wgt_dim != 1)
155  )
156  {
157  *ierr = ZOLTAN_FATAL;
158  return;
159  }
160 
161  ZOLTAN_ID_TYPE* nextNbor = nborGID;
162  int* nextProc = nborProc;
163  float* nextWgt = ewgts;
164 
165  for (Foam::label i=0; i < num_obj; i++)
166  {
167  Foam::label celli = localID[i];
168 
169  const Foam::cell& cFaces = mesh.cells()[celli];
170  forAll(cFaces, cFacei)
171  {
172  Foam::label n = 0;
173 
174  Foam::label facei = cFaces[cFacei];
175  if (mesh.isInternalFace(facei))
176  {
177  Foam::label nbr = mesh.faceOwner()[facei];
178  if (nbr == celli)
179  {
180  nbr = mesh.faceNeighbour()[facei];
181  }
182 
183  // Note: global index
184  *nextNbor++ = nbr;
185  *nextProc++ = 0;
186  *nextWgt++ = 1.0;
187 
188  n++;
189  }
190  if (n != num_edges[i])
191  {
192  *ierr = ZOLTAN_FATAL;
193  return;
194  }
195  }
196  }
197  *ierr = ZOLTAN_OK;
198 }
199 
201 static int get_mesh_dim(void *data, int *ierr)
202 {
203  const Foam::polyMesh& mesh = *static_cast<const Foam::polyMesh*>(data);
204 
205  return mesh.nSolutionD();
206 }
207 
208 
209 static void get_geom_list
210 (
211  void *data,
212  int num_gid_entries,
213  int num_lid_entries,
214  int num_obj,
215  ZOLTAN_ID_PTR global_ids,
216  ZOLTAN_ID_PTR local_ids,
217  int num_dim,
218  double *geom_vec,
219  int *ierr
220 )
221 {
222  const Foam::polyMesh& mesh = *static_cast<const Foam::polyMesh*>(data);
223 
224  if
225  (
226  (num_gid_entries != 1)
227  || (num_lid_entries != 1)
228  || (num_obj != mesh.nCells())
229  || (num_dim != mesh.nSolutionD())
230  )
231  {
232  *ierr = ZOLTAN_FATAL;
233  return;
234  }
235 
236  double* p = geom_vec;
237 
238 
239  const Foam::Vector<Foam::label>& sol = mesh.solutionD();
240 
241  const Foam::pointField& cc = mesh.cellCentres();
242 
243  for (Foam::label celli = 0; celli < num_obj; celli++)
244  {
245  const Foam::point& pt = cc[celli];
246 
247  for (Foam::direction cmpt = 0; cmpt < Foam::vector::nComponents; cmpt++)
248  {
249  if (sol[cmpt] == 1)
250  {
251  *p++ = pt[cmpt];
252  }
253  }
254  }
255  *ierr = ZOLTAN_OK;
256 }
257 
258 
259 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
261 Foam::zoltanRenumber::zoltanRenumber(const dictionary& renumberDict)
262 :
263  renumberMethod(renumberDict),
264  coeffsDict_(renumberDict.subDict(typeName+"Coeffs"))
265 {}
266 
267 
268 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
269 
271 (
272  const polyMesh& pMesh,
273  const pointField& points
274 ) const
275 {
276  stringList args(1);
277  args[0] = "zoltanRenumber";
278 
279  int argc = args.size();
280  char* argv[argc];
281  for (label i = 0; i < argc; i++)
282  {
283  argv[i] = strdup(args[i].c_str());
284  }
285 
286  float ver;
287  int rc = Zoltan_Initialize(argc, argv, &ver);
288 
289  Foam::Pout<< "Initialised to " << ver << Foam::endl;
290 
291  if (rc != ZOLTAN_OK)
292  {
294  << "Failed initialising Zoltan" << exit(FatalError);
295  }
296 
297  struct Zoltan_Struct *zz = Zoltan_Create(MPI_COMM_WORLD);
298 
299  polyMesh& mesh = const_cast<polyMesh&>(pMesh);
300 
301 
302  forAllConstIter(IDLList<entry>, coeffsDict_, iter)
303  {
304  if (!iter().isDict())
305  {
306  const word& key = iter().keyword();
307  const word value(iter().stream());
308 
309  Info<< typeName << " : setting parameter " << key
310  << " to " << value << endl;
311 
312  Zoltan_Set_Param(zz, key.c_str(), value.c_str());
313  }
314  }
315 
316 
317  Zoltan_Set_Num_Obj_Fn(zz, get_number_of_vertices, &mesh);
318  Zoltan_Set_Obj_List_Fn(zz, get_vertex_list, &mesh);
319 
320  // Callbacks for geometry
321  Zoltan_Set_Num_Geom_Fn(zz, get_mesh_dim, &mesh);
322  Zoltan_Set_Geom_Multi_Fn(zz, get_geom_list, &mesh);
323 
324  // Callbacks for connectivity
325  Zoltan_Set_Num_Edges_Multi_Fn(zz, get_num_edges_list, &mesh);
326  Zoltan_Set_Edge_List_Multi_Fn(zz, get_edge_list, &mesh);
327 
328 
329 
330  //Note: !global indices
331  List<ZOLTAN_ID_TYPE> wantedCells(mesh.nCells());
332 
333  globalIndex globalCells(mesh.nCells());
334  forAll(wantedCells, i)
335  {
336  //wantedCells[i] = i;
337  wantedCells[i] = globalCells.toGlobal(i);
338  }
339 
340  List<ZOLTAN_ID_TYPE> oldToNew(mesh.nCells());
341 
342  int err = Zoltan_Order
343  (
344  zz,
345  1, //int num_gid_entries,
346  mesh.globalData().nTotalCells(), //int num_obj,
347  wantedCells.begin(),
348  oldToNew.begin()
349  );
350 
351  if (err != ZOLTAN_OK)
352  {
354  << "Failed Zoltan_Order" << exit(FatalError);
355  }
356 
357 
358  for (label i = 0; i < argc; i++)
359  {
360  free(argv[i]);
361  }
362 
363 
364  labelList order(oldToNew.size());
365  forAll(order, i)
366  {
367  order[i] = oldToNew[i];
368  }
369  return order;
370 }
371 
372 
373 // ************************************************************************* //
static int get_mesh_dim(void *data, int *ierr)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
uint8_t direction
Definition: direction.H:46
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:801
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static void get_geom_list(void *data, int num_gid_entries, int num_lid_entries, int num_obj, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int num_dim, double *geom_vec, int *ierr)
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:812
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
const cellList & cells() const
Macros for easy insertion into run-time selection tables.
label nCells() const
Abstract base class for renumbering.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
static void get_edge_list(void *data, int sizeGID, int sizeLID, int num_obj, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID, int *num_edges, ZOLTAN_ID_PTR nborGID, int *nborProc, int wgt_dim, float *ewgts, int *ierr)
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:96
dynamicFvMesh & mesh
System uinteger.
A class for handling words, derived from string.
Definition: word.H:59
const vectorField & cellCentres() const
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1140
Templated 3D Vector derived from VectorSpace adding construction from 3 components, element access using x(), y() and z() member functions and the inner-product (dot-product) and cross product operators.
Definition: Vector.H:57
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
static int get_number_of_vertices(void *data, int *ierr)
static void get_vertex_list(void *data, int sizeGID, int sizeLID, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID, int wgt_dim, float *obj_wgts, int *ierr)
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
virtual labelList renumber(const pointField &) const
Return the order in which cells need to be visited, i.e.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
messageStream Info
label n
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
Foam::argList args(argc, argv)
label nTotalCells() const
Return total number of cells in decomposed mesh.
Namespace for OpenFOAM.
static void get_num_edges_list(void *data, int sizeGID, int sizeLID, int num_obj, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID, int *numEdges, int *ierr)