surfaceRedistributePar.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 "distributionMap.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  // Determine mesh bounding boxes:
149  {
151  (
152  1,
153  treeBoundBox(boundBox(mesh.points(), false)).extend(1e-3)
154  );
157  }
158 
160  (
161  surfFileName,
162  runTime.constant(),
164  runTime,
167  );
168 
169  // Look for file (using searchableSurface rules)
170  const fileName actualPath(io.filePath());
171  fileName localPath(actualPath);
172  localPath.replace(runTime.rootPath() + '/', "");
173 
174  if (actualPath == io.objectPath(false))
175  {
176  Info<< "Loading local (decomposed) surface " << localPath << nl <<endl;
177  }
178  else
179  {
180  Info<< "Loading undecomposed surface " << localPath << nl << endl;
181  }
182 
183 
184  // Create dummy dictionary for bounding boxes if does not exist.
185  if (!isFile(actualPath / "Dict"))
186  {
188  dict.add("bounds", meshBb[Pstream::myProcNo()]);
189  dict.add("distributionType", distType);
190  dict.add("mergeDistance", small);
191 
192  localIOdictionary ioDict
193  (
194  IOobject
195  (
196  io.name() + "Dict",
197  io.instance(),
198  io.local(),
199  io.db(),
202  false
203  ),
204  dict
205  );
206 
207  Info<< "Writing dummy bounds dictionary to " << ioDict.name()
208  << nl << endl;
209 
210  // Force writing in ascii
211  ioDict.regIOobject::writeObject
212  (
215  ioDict.time().writeCompression(),
216  true
217  );
218  }
219 
220 
221  // Load surface
223  Info<< "Loaded surface" << nl << endl;
224 
225 
226  // Generate a test field
227  {
228  const triSurface& s = static_cast<const triSurface&>(surfMesh);
229 
231  (
233  (
234  IOobject
235  (
236  "faceCentres", // name
237  surfMesh.searchableSurface::time().timeName(), // instance
238  surfMesh,
241  ),
242  surfMesh,
243  dimLength,
244  s.faceCentres()
245  )
246  );
247 
248  // Steal pointer and store object on surfMesh
249  fcPtr.ptr()->store();
250  }
251 
252 
253  // Write per-processor stats
254  Info<< "Before redistribution:" << endl;
255  writeProcStats(surfMesh, meshBb);
256 
257 
258  // Do redistribution
259  Info<< "Redistributing surface" << nl << endl;
261  autoPtr<distributionMap> pointMap;
262  surfMesh.distribute
263  (
265  keepNonMapped,
266  faceMap,
267  pointMap
268  );
269  faceMap.clear();
270  pointMap.clear();
271 
272  Info<< endl;
273 
274 
275  // Write per-processor stats
276  Info<< "After redistribution:" << endl;
277  writeProcStats(surfMesh, meshBb);
278 
279 
280  Info<< "Writing surface." << nl << endl;
281  surfMesh.objectRegistry::write();
282 
283  Info<< "End\n" << endl;
284 
285  return 0;
286 }
287 
288 
289 // ************************************************************************* //
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:434
A class for handling file names.
Definition: fileName.H:79
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:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const Field< PointType > & faceCentres() const
Return face centres for patch.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
bool isFile(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist as a file in the file system?
Definition: POSIX.C:555
Templated form of IOobject providing type information for file reading and header type checking...
Definition: IOobject.H:537
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
fvMesh & mesh
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:1153
const dimensionSet dimLength
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
localIOdictionary derived from IOdictionary with global set false to disable 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:96
static const word & geometryDir()
Return the geometry directory name.
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
const Field< PointType > & points() const
Return reference to global points.
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(1e-3))
A surface mesh consisting of general polygon faces.
Definition: surfMesh.H:55
IOoject and searching on distributed triSurface. All processor hold (possibly overlapping) part of th...
static const char nl
Definition: Ostream.H:260
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static const versionNumber currentVersion
Current version number.
Definition: IOstream.H:203
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
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
Triangulated surface description with patch information.
Definition: triSurface.H:66
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
treeBoundBox extend(const scalar s) const
Return asymetrically extended bounding box, with guaranteed.
Namespace for OpenFOAM.