stitchMesh.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) 2023-2025 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  Utility to stitch or conform pairs of patches,
29  converting the patch faces either into internal faces
30  or conformal faces or another patch.
31 
32 Usage
33  \b stitchMesh (<list of patch pairs>)
34 
35  E.g. to stitch patches \c top1 to \c top2 and \c bottom1 to \c bottom2
36  stitchMesh "((top1 top2) (bottom1 bottom2))"
37 
38  Options:
39  - \par -noOverwrite \n
40  Do not replace the old mesh with the new one, writing the new one
41  into a separate time directory
42 
43  - \par -region <name>
44  Specify an alternative mesh region.
45 
46  - \par -fields
47  Update vol and point fields
48 
49  - \par -tol
50  Merge tolerance relative to local edge length (default 1e-4)
51 
52 See also
53  Foam::mergePatchPairs
54 
55 \*---------------------------------------------------------------------------*/
56 
57 #include "argList.H"
58 #include "fvMesh.H"
59 #include "timeSelector.H"
60 #include "ReadFields.H"
61 #include "volFields.H"
62 #include "pointFields.H"
63 #include "mergePatchPairs.H"
64 
65 using namespace Foam;
66 
67 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69 int main(int argc, char *argv[])
70 {
71  timeSelector::addOptions(true, false);
72 
74  (
75  "Stitch the faces on the specified patch pairs\n"
76  "converting the overlapping patch faces into internal faces.\n"
77  );
78 
80  #include "addNoOverwriteOption.H"
81  #include "addMeshOption.H"
82  #include "addRegionOption.H"
84  (
85  "fields",
86  "update fields"
87  );
89  (
90  "tol",
91  "scalar",
92  "merge tolerance relative to local edge length (default 1e-4)"
93  );
94 
95  argList::validArgs.append("patchPairs");
96 
97  #include "setRootCase.H"
99 
100  // Select time if specified
102 
104 
105  const scalar snapTol = args.optionLookupOrDefault("tol", 1e-4);
106  #include "setNoOverwrite.H"
107  const bool fields = args.optionFound("fields");
108 
109  const List<Pair<word>> patchPairNames((IStringStream(args[1])()));
110 
111  const word oldInstance = mesh.pointsInstance();
112 
113  // Search for list of objects for this time
114  IOobjectList objects(mesh, runTime.name());
115 
116  if (fields) Info<< "Reading geometric fields" << nl << endl;
117 
118  #include "readVolFields.H"
119  #include "readPointFields.H"
120 
121  Info<< endl;
122 
123  if (!overwrite)
124  {
125  runTime++;
126  }
127 
128  // Stitch all the patch-pairs
129  mergePatchPairs(mesh, patchPairNames, snapTol);
130 
131  if (overwrite)
132  {
133  mesh.setInstance(oldInstance);
134  }
135 
136  Info<< "Writing mesh to " << mesh.facesInstance() << endl;
137  mesh.write();
138 
139  Info<< "End\n" << endl;
140 
141  return 0;
142 }
143 
144 
145 // ************************************************************************* //
Field reading functions for post-processing utilities.
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
Input from memory buffer stream.
Definition: IStringStream.H:52
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:128
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
static void noParallel()
Remove the parallel options.
Definition: argList.C:175
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
T optionLookupOrDefault(const word &opt, const T &deflt) const
Read a value from the named option if present.
Definition: argListI.H:294
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1785
const word & name() const
Return reference to name.
Definition: fvMesh.H:434
Class to stitch mesh by merging patch-pairs.
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:994
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:988
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:91
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
static instantList selectIfPresent(Time &runTime, const argList &args)
If any time option provided return the set of times (as select0)
Definition: timeSelector.C:283
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
int main(int argc, char *argv[])
Definition: financialFoam.C:44
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:234
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
static const char nl
Definition: Ostream.H:267
objects
const bool overwrite
Definition: setNoOverwrite.H:1
Foam::argList args(argc, argv)