stitchMesh.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 Application
25  stitchMesh
26 
27 Description
28  'Stitches' a mesh.
29 
30  Takes a mesh and two patches and merges the faces on the two patches
31  (if geometrically possible) so the faces become internal.
32 
33  Can do
34  - 'perfect' match: faces and points on patches align exactly. Order might
35  be different though.
36  - 'integral' match: where the surfaces on both patches exactly
37  match but the individual faces not
38  - 'partial' match: where the non-overlapping part of the surface remains
39  in the respective patch.
40 
41  Note : Is just a front-end to perfectInterface/slidingInterface.
42 
43  Comparable to running a meshModifier of the form
44  (if masterPatch is called "M" and slavePatch "S"):
45 
46  \verbatim
47  couple
48  {
49  type slidingInterface;
50  masterFaceZoneName MSMasterZone
51  slaveFaceZoneName MSSlaveZone
52  cutPointZoneName MSCutPointZone
53  cutFaceZoneName MSCutFaceZone
54  masterPatchName M;
55  slavePatchName S;
56  typeOfMatch partial or integral
57  }
58  \endverbatim
59 
60 
61 \*---------------------------------------------------------------------------*/
62 
63 #include "fvCFD.H"
64 #include "polyTopoChanger.H"
65 #include "mapPolyMesh.H"
66 #include "ListOps.H"
67 #include "slidingInterface.H"
68 #include "perfectInterface.H"
69 #include "IOobjectList.H"
70 #include "ReadFields.H"
71 
72 
73 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74 
75 label addPointZone(const polyMesh& mesh, const word& name)
76 {
77  label zoneID = mesh.pointZones().findZoneID(name);
78 
79  if (zoneID != -1)
80  {
81  Info<< "Reusing existing pointZone "
82  << mesh.pointZones()[zoneID].name()
83  << " at index " << zoneID << endl;
84  }
85  else
86  {
87  pointZoneMesh& pointZones = const_cast<polyMesh&>(mesh).pointZones();
88  zoneID = pointZones.size();
89  Info<< "Adding pointZone " << name << " at index " << zoneID << endl;
90 
91  pointZones.setSize(zoneID+1);
92  pointZones.set
93  (
94  zoneID,
95  new pointZone
96  (
97  name,
98  labelList(0),
99  zoneID,
100  pointZones
101  )
102  );
103  }
104  return zoneID;
105 }
106 
107 
108 label addFaceZone(const polyMesh& mesh, const word& name)
109 {
110  label zoneID = mesh.faceZones().findZoneID(name);
111 
112  if (zoneID != -1)
113  {
114  Info<< "Reusing existing faceZone " << mesh.faceZones()[zoneID].name()
115  << " at index " << zoneID << endl;
116  }
117  else
118  {
119  faceZoneMesh& faceZones = const_cast<polyMesh&>(mesh).faceZones();
120  zoneID = faceZones.size();
121  Info<< "Adding faceZone " << name << " at index " << zoneID << endl;
122 
123  faceZones.setSize(zoneID+1);
124  faceZones.set
125  (
126  zoneID,
127  new faceZone
128  (
129  name,
130  labelList(0),
131  boolList(),
132  zoneID,
133  faceZones
134  )
135  );
136  }
137  return zoneID;
138 }
139 
140 
141 label addCellZone(const polyMesh& mesh, const word& name)
142 {
143  label zoneID = mesh.cellZones().findZoneID(name);
144 
145  if (zoneID != -1)
146  {
147  Info<< "Reusing existing cellZone " << mesh.cellZones()[zoneID].name()
148  << " at index " << zoneID << endl;
149  }
150  else
151  {
152  cellZoneMesh& cellZones = const_cast<polyMesh&>(mesh).cellZones();
153  zoneID = cellZones.size();
154  Info<< "Adding cellZone " << name << " at index " << zoneID << endl;
155 
156  cellZones.setSize(zoneID+1);
157  cellZones.set
158  (
159  zoneID,
160  new cellZone
161  (
162  name,
163  labelList(0),
164  zoneID,
165  cellZones
166  )
167  );
168  }
169  return zoneID;
170 }
171 
172 
173 // Checks whether patch present
174 void checkPatch(const polyBoundaryMesh& bMesh, const word& name)
175 {
176  const label patchi = bMesh.findPatchID(name);
177 
178  if (patchi == -1)
179  {
181  << "Cannot find patch " << name << endl
182  << "It should be present and of non-zero size" << endl
183  << "Valid patches are " << bMesh.names()
184  << exit(FatalError);
185  }
186 
187  if (bMesh[patchi].empty())
188  {
190  << "Patch " << name << " is present but zero size"
191  << exit(FatalError);
192  }
193 }
194 
195 
196 
197 int main(int argc, char *argv[])
198 {
199  argList::addNote
200  (
201  "Merge the faces on the specified patches (if geometrically possible)\n"
202  "so the faces become internal.\n"
203  "Integral matching is used when the options -partial and -perfect are "
204  "omitted.\n"
205  );
206 
207  argList::noParallel();
208  #include "addOverwriteOption.H"
209  #include "addRegionOption.H"
210 
211  argList::validArgs.append("masterPatch");
212  argList::validArgs.append("slavePatch");
213 
214  argList::addBoolOption
215  (
216  "partial",
217  "couple partially overlapping patches (optional)"
218  );
219  argList::addBoolOption
220  (
221  "perfect",
222  "couple perfectly aligned patches (optional)"
223  );
224  argList::addOption
225  (
226  "toleranceDict",
227  "file",
228  "dictionary file with tolerances"
229  );
230 
231  #include "setRootCase.H"
232  #include "createTime.H"
233  runTime.functionObjects().off();
234  #include "createNamedMesh.H"
235 
236  const word oldInstance = mesh.pointsInstance();
237 
238  const word masterPatchName = args[1];
239  const word slavePatchName = args[2];
240 
241  const bool partialCover = args.optionFound("partial");
242  const bool perfectCover = args.optionFound("perfect");
243  const bool overwrite = args.optionFound("overwrite");
244 
245  if (partialCover && perfectCover)
246  {
248  << "Cannot supply both partial and perfect." << endl
249  << "Use perfect match option if the patches perfectly align"
250  << " (both vertex positions and face centres)" << endl
251  << exit(FatalError);
252  }
253 
254 
255  const word mergePatchName(masterPatchName + slavePatchName);
256  const word cutZoneName(mergePatchName + "CutFaceZone");
257 
258  slidingInterface::typeOfMatch tom = slidingInterface::INTEGRAL;
259 
260  if (partialCover)
261  {
262  Info<< "Coupling partially overlapping patches "
263  << masterPatchName << " and " << slavePatchName << nl
264  << "Resulting internal faces will be in faceZone " << cutZoneName
265  << nl
266  << "Any uncovered faces will remain in their patch"
267  << endl;
268 
269  tom = slidingInterface::PARTIAL;
270  }
271  else if (perfectCover)
272  {
273  Info<< "Coupling perfectly aligned patches "
274  << masterPatchName << " and " << slavePatchName << nl
275  << "Resulting (internal) faces will be in faceZone " << cutZoneName
276  << nl << nl
277  << "Note: both patches need to align perfectly." << nl
278  << "Both the vertex"
279  << " positions and the face centres need to align to within" << nl
280  << "a tolerance given by the minimum edge length on the patch"
281  << endl;
282  }
283  else
284  {
285  Info<< "Coupling patches " << masterPatchName << " and "
286  << slavePatchName << nl
287  << "Resulting (internal) faces will be in faceZone " << cutZoneName
288  << nl << nl
289  << "Note: the overall area covered by both patches should be"
290  << " identical (\"integral\" interface)." << endl
291  << "If this is not the case use the -partial option" << nl << endl;
292  }
293 
294  // set up the tolerances for the sliding mesh
295  dictionary slidingTolerances;
296  if (args.options().found("toleranceDict"))
297  {
298  IOdictionary toleranceFile
299  (
300  IOobject
301  (
302  args.options()["toleranceDict"],
303  runTime.constant(),
304  mesh,
305  IOobject::MUST_READ_IF_MODIFIED,
306  IOobject::NO_WRITE
307  )
308  );
309  slidingTolerances += toleranceFile;
310  }
311 
312  // Check for non-empty master and slave patches
313  checkPatch(mesh.boundaryMesh(), masterPatchName);
314  checkPatch(mesh.boundaryMesh(), slavePatchName);
315 
316  // Create and add face zones and mesh modifiers
317 
318  // Master patch
319  const polyPatch& masterPatch = mesh.boundaryMesh()[masterPatchName];
320 
321  // Make list of masterPatch faces
322  labelList isf(masterPatch.size());
323 
324  forAll(isf, i)
325  {
326  isf[i] = masterPatch.start() + i;
327  }
328 
329  polyTopoChanger stitcher(mesh);
330  stitcher.setSize(1);
331 
332  mesh.pointZones().clearAddressing();
333  mesh.faceZones().clearAddressing();
334  mesh.cellZones().clearAddressing();
335 
336  if (perfectCover)
337  {
338  // Add empty zone for resulting internal faces
339  label cutZoneID = addFaceZone(mesh, cutZoneName);
340 
341  mesh.faceZones()[cutZoneID].resetAddressing
342  (
343  isf,
344  boolList(masterPatch.size(), false)
345  );
346 
347  // Add the perfect interface mesh modifier
348  stitcher.set
349  (
350  0,
351  new perfectInterface
352  (
353  "couple",
354  0,
355  stitcher,
356  cutZoneName,
357  masterPatchName,
358  slavePatchName
359  )
360  );
361  }
362  else
363  {
364  label pointZoneID = addPointZone(mesh, mergePatchName + "CutPointZone");
365  mesh.pointZones()[pointZoneID] = labelList(0);
366 
367  label masterZoneID = addFaceZone(mesh, mergePatchName + "MasterZone");
368 
369  mesh.faceZones()[masterZoneID].resetAddressing
370  (
371  isf,
372  boolList(masterPatch.size(), false)
373  );
374 
375  // Slave patch
376  const polyPatch& slavePatch = mesh.boundaryMesh()[slavePatchName];
377 
378  labelList osf(slavePatch.size());
379 
380  forAll(osf, i)
381  {
382  osf[i] = slavePatch.start() + i;
383  }
384 
385  label slaveZoneID = addFaceZone(mesh, mergePatchName + "SlaveZone");
386  mesh.faceZones()[slaveZoneID].resetAddressing
387  (
388  osf,
389  boolList(slavePatch.size(), false)
390  );
391 
392  // Add empty zone for cut faces
393  label cutZoneID = addFaceZone(mesh, cutZoneName);
394  mesh.faceZones()[cutZoneID].resetAddressing
395  (
396  labelList(0),
397  boolList(0, false)
398  );
399 
400 
401  // Add the sliding interface mesh modifier
402  stitcher.set
403  (
404  0,
405  new slidingInterface
406  (
407  "couple",
408  0,
409  stitcher,
410  mergePatchName + "MasterZone",
411  mergePatchName + "SlaveZone",
412  mergePatchName + "CutPointZone",
413  cutZoneName,
414  masterPatchName,
415  slavePatchName,
416  tom, // integral or partial
417  true // couple/decouple mode
418  )
419  );
420  static_cast<slidingInterface&>(stitcher[0]).setTolerances
421  (
422  slidingTolerances,
423  true
424  );
425  }
426 
427 
428  // Search for list of objects for this time
429  IOobjectList objects(mesh, runTime.timeName());
430 
431  // Read all current fvFields so they will get mapped
432  Info<< "Reading all current volfields" << endl;
433  PtrList<volScalarField> volScalarFields;
434  ReadFields(mesh, objects, volScalarFields);
435 
436  PtrList<volVectorField> volVectorFields;
437  ReadFields(mesh, objects, volVectorFields);
438 
439  PtrList<volSphericalTensorField> volSphericalTensorFields;
440  ReadFields(mesh, objects, volSphericalTensorFields);
441 
442  PtrList<volSymmTensorField> volSymmTensorFields;
443  ReadFields(mesh, objects, volSymmTensorFields);
444 
445  PtrList<volTensorField> volTensorFields;
446  ReadFields(mesh, objects, volTensorFields);
447 
448  //- Uncomment if you want to interpolate surface fields (usually bad idea)
449  //Info<< "Reading all current surfaceFields" << endl;
450  //PtrList<surfaceScalarField> surfaceScalarFields;
451  //ReadFields(mesh, objects, surfaceScalarFields);
452  //
453  //PtrList<surfaceVectorField> surfaceVectorFields;
454  //ReadFields(mesh, objects, surfaceVectorFields);
455  //
456  //PtrList<surfaceTensorField> surfaceTensorFields;
457  //ReadFields(mesh, objects, surfaceTensorFields);
458 
459  if (!overwrite)
460  {
461  runTime++;
462  }
463 
464  // Execute all polyMeshModifiers
465  autoPtr<mapPolyMesh> morphMap = stitcher.changeMesh(true);
466 
467  mesh.movePoints(morphMap->preMotionPoints());
468 
469  // Write mesh
470  if (overwrite)
471  {
472  mesh.setInstance(oldInstance);
473  stitcher.instance() = oldInstance;
474  }
475  Info<< nl << "Writing polyMesh to time " << runTime.timeName() << endl;
476 
477  IOstream::defaultPrecision(max(10u, IOstream::defaultPrecision()));
478 
479  // Bypass runTime write (since only writes at writeTime)
480  if
481  (
482  !runTime.objectRegistry::writeObject
483  (
484  runTime.writeFormat(),
485  IOstream::currentVersion,
486  runTime.writeCompression()
487  )
488  )
489  {
491  << "Failed writing polyMesh."
492  << exit(FatalError);
493  }
494 
495  mesh.faceZones().write();
496  mesh.pointZones().write();
497  mesh.cellZones().write();
498 
499  // Write fields
500  runTime.write();
501 
502  Info<< "End\n" << endl;
503 
504  return 0;
505 }
506 
507 
508 // ************************************************************************* //
DynamicID< pointZoneMesh > pointZoneID
Foam::pointZoneID.
Definition: ZoneIDs.H:47
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
wordList ReadFields(const Mesh &mesh, const IOobjectList &objects, PtrList< GeoField > &fields, const bool syncPar=true)
Read all fields of the specified type.
#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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Field reading functions for post-processing utilities.
const Foam::HashTable< string > & options() const
Return options.
Definition: argListI.H:90
Various functions to operate on Lists.
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
dynamicFvMesh & mesh
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
List< label > labelList
A List of labels.
Definition: labelList.H:56
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
static const char nl
Definition: Ostream.H:262
messageStream Info
Foam::argList args(argc, argv)
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with the type cellZone.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29