surfaceRedistributePar.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-2017 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 Application
25  surfaceRedistributePar
26 
27 Description
28  (Re)distribution of triSurface. Either takes an undecomposed surface
29  or an already decomposed surface and redistributes it so that each
30  processor has all triangles that overlap its mesh.
31 
32 Note
33  - best decomposition option is hierarchGeomDecomp since
34  guarantees square decompositions.
35  - triangles might be present on multiple processors.
36  - merging uses geometric tolerance so take care with writing precision.
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #include "argList.H"
41 #include "Time.H"
42 #include "polyMesh.H"
44 #include "mapDistribute.H"
45 #include "localIOdictionary.H"
46 
47 using namespace Foam;
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 // Print on master all the per-processor surface stats.
52 void writeProcStats
53 (
54  const triSurface& s,
56 )
57 {
58  // Determine surface bounding boxes, faces, points
60  {
61  surfBb[Pstream::myProcNo()] = treeBoundBox(s.points());
62  Pstream::gatherList(surfBb);
63  Pstream::scatterList(surfBb);
64  }
65 
70 
71  labelList nFaces(Pstream::nProcs());
72  nFaces[Pstream::myProcNo()] = s.size();
73  Pstream::gatherList(nFaces);
74  Pstream::scatterList(nFaces);
75 
76  forAll(surfBb, proci)
77  {
78  const List<treeBoundBox>& bbs = meshBb[proci];
79 
80  Info<< "processor" << proci << nl
81  << "\tMesh bounds : " << bbs[0] << nl;
82  for (label i = 1; i < bbs.size(); i++)
83  {
84  Info<< "\t " << bbs[i]<< nl;
85  }
86  Info<< "\tSurface bounding box : " << surfBb[proci] << nl
87  << "\tTriangles : " << nFaces[proci] << nl
88  << "\tVertices : " << nPoints[proci]
89  << endl;
90  }
91  Info<< endl;
92 }
93 
94 
95 
96 int main(int argc, char *argv[])
97 {
99  (
100  "redistribute a triSurface"
101  );
102 
103  argList::validArgs.append("surface file");
104  argList::validArgs.append("distribution type");
106  (
107  "keepNonMapped",
108  "preserve surface outside of mesh bounds"
109  );
110 
111  #include "setRootCase.H"
112  #include "createTime.H"
113  runTime.functionObjects().off();
114 
115  const fileName surfFileName = args[1];
116  const word distType = args[2];
117 
118  Info<< "Reading surface from " << surfFileName << nl << nl
119  << "Using distribution method "
121  << " " << distType << nl << endl;
122 
123  const bool keepNonMapped = args.options().found("keepNonMapped");
124 
125  if (keepNonMapped)
126  {
127  Info<< "Preserving surface outside of mesh bounds." << nl << endl;
128  }
129  else
130  {
131  Info<< "Removing surface outside of mesh bounds." << nl << endl;
132  }
133 
134 
135  if (!Pstream::parRun())
136  {
138  << "Please run this program on the decomposed case."
139  << " It will read surface " << surfFileName
140  << " and decompose it such that it overlaps the mesh bounding box."
141  << exit(FatalError);
142  }
143 
144 
145  #include "createPolyMesh.H"
146 
147  Random rndGen(653213);
148 
149  // Determine mesh bounding boxes:
151  {
153  (
154  1,
156  (
157  boundBox(mesh.points(), false)
158  ).extend(rndGen, 1e-3)
159  );
162  }
163 
164  IOobject io
165  (
166  surfFileName, // name
167  //runTime.findInstance("triSurface", surfFileName), // instance
168  runTime.constant(), // instance
169  "triSurface", // local
170  runTime, // registry
173  );
174 
175  // Look for file (using searchableSurface rules)
176  const fileName actualPath(typeFilePath<searchableSurface>(io));
177  fileName localPath(actualPath);
178  localPath.replace(runTime.rootPath() + '/', "");
179 
180  if (actualPath == io.objectPath())
181  {
182  Info<< "Loading local (decomposed) surface " << localPath << nl <<endl;
183  }
184  else
185  {
186  Info<< "Loading undecomposed surface " << localPath << nl << endl;
187  }
188 
189 
190  // Create dummy dictionary for bounding boxes if does not exist.
191  if (!isFile(actualPath / "Dict"))
192  {
194  dict.add("bounds", meshBb[Pstream::myProcNo()]);
195  dict.add("distributionType", distType);
196  dict.add("mergeDistance", SMALL);
197 
198  localIOdictionary ioDict
199  (
200  IOobject
201  (
202  io.name() + "Dict",
203  io.instance(),
204  io.local(),
205  io.db(),
208  false
209  ),
210  dict
211  );
212 
213  Info<< "Writing dummy bounds dictionary to " << ioDict.name()
214  << nl << endl;
215 
216  // Force writing in ascii
217  ioDict.regIOobject::writeObject
218  (
221  ioDict.time().writeCompression(),
222  true
223  );
224  }
225 
226 
227  // Load surface
229  Info<< "Loaded surface" << nl << endl;
230 
231 
232  // Generate a test field
233  {
234  const triSurface& s = static_cast<const triSurface&>(surfMesh);
235 
237  (
239  (
240  IOobject
241  (
242  "faceCentres", // name
243  surfMesh.searchableSurface::time().timeName(), // instance
244  surfMesh,
247  ),
248  surfMesh,
249  dimLength,
250  s.faceCentres()
251  )
252  );
253 
254  // Steal pointer and store object on surfMesh
255  fcPtr.ptr()->store();
256  }
257 
258 
259  // Write per-processor stats
260  Info<< "Before redistribution:" << endl;
261  writeProcStats(surfMesh, meshBb);
262 
263 
264  // Do redistribution
265  Info<< "Redistributing surface" << nl << endl;
267  autoPtr<mapDistribute> pointMap;
268  surfMesh.distribute
269  (
271  keepNonMapped,
272  faceMap,
273  pointMap
274  );
275  faceMap.clear();
276  pointMap.clear();
277 
278  Info<< endl;
279 
280 
281  // Write per-processor stats
282  Info<< "After redistribution:" << endl;
283  writeProcStats(surfMesh, meshBb);
284 
285 
286  Info<< "Writing surface." << nl << endl;
287  surfMesh.objectRegistry::write();
288 
289  Info<< "End\n" << endl;
290 
291  return 0;
292 }
293 
294 
295 // ************************************************************************* //
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
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
A class for handling file names.
Definition: fileName.H:69
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const double e
Elementary charge.
Definition: doubleFloat.H:78
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
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:418
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:154
const Field< PointType > & faceCentres() const
Return face centres for patch.
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
static const NamedEnum< distributionType, 3 > distributionTypeNames_
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:814
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1011
localIOdictionary is derived from IOdictionary but excludes parallel master reading.
void clear()
Delete object (if the pointer is valid) and set pointer to.
Definition: autoPtrI.H:126
const Foam::HashTable< string > & options() const
Return options.
Definition: argListI.H:90
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
cachedRandom rndGen(label(0), -1)
const Field< PointType > & points() const
Return reference to global points.
Simple random number generator.
Definition: Random.H:49
A surface mesh consisting of general polygon faces.
Definition: surfMesh.H:55
bool isFile(const fileName &, const bool checkGzip=true, const bool followLink=true)
Does the name exist as a FILE in the file system?
Definition: POSIX.C:551
IOoject and searching on distributed triSurface. All processor hold (possibly overlapping) part of th...
static const char nl
Definition: Ostream.H:262
treeBoundBox extend(Random &, const scalar s) const
Return slightly wider bounding box.
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:394
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:400
static const versionNumber currentVersion
Current version number.
Definition: IOstream.H:206
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
messageStream Info
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:85
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:126
Triangulated surface description with patch information.
Definition: triSurface.H:65
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Namespace for OpenFOAM.