splitMesh.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-2023 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  splitMesh
26 
27 Description
28  Splits mesh by making internal faces external. Uses attachDetach.
29 
30  Generates a meshModifier of the form:
31 
32  Splitter
33  {
34  type attachDetach;
35  faceZoneName membraneFaces;
36  masterPatchName masterPatch;
37  slavePatchName slavePatch;
38  triggerTimes runTime.value();
39  }
40 
41  so will detach at the current time and split all faces in membraneFaces
42  into masterPatch and slavePatch (which have to be present but of 0 size)
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #include "argList.H"
47 #include "polyMesh.H"
48 #include "Time.H"
49 #include "polyTopoChange.H"
50 #include "polyTopoChangeMap.H"
51 #include "faceSet.H"
52 #include "attachDetach.H"
53 #include "attachPolyTopoChanger.H"
54 #include "regionSide.H"
55 #include "primitiveFacePatch.H"
56 
57 using namespace Foam;
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61 // Find edge between points v0 and v1.
62 label findEdge(const primitiveMesh& mesh, const label v0, const label v1)
63 {
64  const labelList& pEdges = mesh.pointEdges()[v0];
65 
66  forAll(pEdges, pEdgeI)
67  {
68  label edgeI = pEdges[pEdgeI];
69 
70  const edge& e = mesh.edges()[edgeI];
71 
72  if (e.otherVertex(v0) == v1)
73  {
74  return edgeI;
75  }
76  }
77 
79  << "Cannot find edge between mesh points " << v0 << " and " << v1
80  << abort(FatalError);
81 
82  return -1;
83 }
84 
85 
86 // Checks whether patch present
87 void checkPatch(const polyBoundaryMesh& bMesh, const word& name)
88 {
90 
91  if (patchi == -1)
92  {
94  << "Cannot find patch " << name << nl
95  << "It should be present but of zero size" << endl
96  << "Valid patches are " << bMesh.names()
97  << exit(FatalError);
98  }
99 
100  if (bMesh[patchi].size())
101  {
103  << "Patch " << name << " is present but non-zero size"
104  << exit(FatalError);
105  }
106 }
107 
108 
109 
110 int main(int argc, char *argv[])
111 {
113  #include "addOverwriteOption.H"
114 
115  argList::validArgs.append("faceSet");
116  argList::validArgs.append("masterPatch");
117  argList::validArgs.append("slavePatch");
118 
119  #include "setRootCase.H"
121  #include "createPolyMesh.H"
122  const word oldInstance = mesh.pointsInstance();
123 
124  const word setName = args[1];
125  const word masterPatch = args[2];
126  const word slavePatch = args[3];
127  const bool overwrite = args.optionFound("overwrite");
128 
129  // List of faces to split
130  faceSet facesSet(mesh, setName);
131 
132  Info<< "Read " << facesSet.size() << " faces to split" << endl << endl;
133 
134 
135  // Convert into labelList and check
136 
137  labelList faces(facesSet.toc());
138 
139  forAll(faces, i)
140  {
141  if (!mesh.isInternalFace(faces[i]))
142  {
144  << "Face " << faces[i] << " in faceSet " << setName
145  << " is not an internal face."
146  << exit(FatalError);
147  }
148  }
149 
150 
151  // Check for empty master and slave patches
152  checkPatch(mesh.boundaryMesh(), masterPatch);
153  checkPatch(mesh.boundaryMesh(), slavePatch);
154 
155 
156  //
157  // Find 'side' of all faces on splitregion. Uses regionSide which needs
158  // set of edges on side of this region. Use PrimitivePatch to find these.
159  //
160 
161  // Addressing on faces only in mesh vertices.
162  primitiveFacePatch fPatch
163  (
164  faceList
165  (
167  (
168  mesh.faces(),
169  faces
170  )
171  ),
172  mesh.points()
173  );
174 
175  const labelList& meshPoints = fPatch.meshPoints();
176 
177  // Mark all fence edges : edges on boundary of fPatch but not on boundary
178  // of polyMesh
179  labelHashSet fenceEdges(fPatch.size());
180 
181  const labelListList& allEdgeFaces = fPatch.edgeFaces();
182 
183  forAll(allEdgeFaces, patchEdgeI)
184  {
185  if (allEdgeFaces[patchEdgeI].size() == 1)
186  {
187  const edge& e = fPatch.edges()[patchEdgeI];
188 
189  label edgeI =
190  findEdge
191  (
192  mesh,
193  meshPoints[e.start()],
194  meshPoints[e.end()]
195  );
196 
197  fenceEdges.insert(edgeI);
198  }
199  }
200 
201  // Find sides reachable from 0th face of faceSet
202  label startFacei = faces[0];
203 
204  regionSide regionInfo
205  (
206  mesh,
207  facesSet,
208  fenceEdges,
209  mesh.faceOwner()[startFacei],
210  startFacei
211  );
212 
213  // Determine flip state for all faces in faceSet
214  boolList zoneFlip(faces.size());
215 
216  forAll(faces, i)
217  {
218  zoneFlip[i] = !regionInfo.sideOwner().found(faces[i]);
219  }
220 
221 
222  // Create and add face zones and mesh modifiers
223  List<pointZone*> pz(0);
224  List<faceZone*> fz(1);
225  List<cellZone*> cz(0);
226 
227  fz[0] =
228  new faceZone
229  (
230  "membraneFaces",
231  faces,
232  zoneFlip,
233  0,
234  mesh.faceZones()
235  );
236 
237  Info<< "Adding point and face zones" << endl;
238  mesh.addZones(pz, fz, cz);
239 
240  attachPolyTopoChanger splitter(mesh);
241  splitter.setSize(1);
242 
243  // Add the sliding interface mesh modifier to start working at current
244  // time
245  splitter.set
246  (
247  0,
248  new attachDetach
249  (
250  "Splitter",
251  0,
252  splitter,
253  "membraneFaces",
254  masterPatch,
255  slavePatch,
256  scalarField(1, runTime.value())
257  )
258  );
259 
260  Info<< nl << "Constructed topologyModifier:" << endl;
261  splitter[0].writeDict(Info);
262 
263  if (!overwrite)
264  {
265  runTime++;
266  }
267 
268  splitter.attach();
269 
270  if (overwrite)
271  {
272  mesh.setInstance(oldInstance);
273  }
274  else
275  {
276  mesh.setInstance(runTime.name());
277  }
278 
279  Info<< "Writing mesh to " << runTime.name() << endl;
280  if (!mesh.write())
281  {
283  << "Failed writing polyMesh."
284  << exit(FatalError);
285  }
286 
287  Info<< "End\n" << endl;
288  return 0;
289 }
290 
291 
292 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A list of faces which address into the list of points.
A List with indirect addressing.
Definition: UIndirectList.H:60
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
Attach/detach boundary mesh modifier. This modifier takes a set of internal faces and converts them i...
Definition: attachDetach.H:62
This class is derived from polyMesh and serves as a tool for statically connecting pieces of a mesh b...
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A list of face labels.
Definition: faceSet.H:51
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:68
Foam::polyBoundaryMesh.
label findPatchID(const word &patchName) const
Find patch index given a name.
wordList names() const
Return a list of patch names.
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
const labelListList & pointEdges() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
virtual const faceList & faces() const =0
Return faces.
virtual const labelList & faceOwner() const =0
Face face-owner addressing.
virtual const pointField & points() const =0
Return mesh points.
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
Determines the 'side' for every face and connected to a singly-connected (through edges) region of fa...
Definition: regionSide.H:62
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
const polyBoundaryMesh & bMesh
label patchi
label findEdge(const edgeList &edges, const labelList &candidates, const label v0, const label v1)
Return edge among candidates that uses the two vertices.
Definition: meshTools.C:337
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const doubleScalar e
Definition: doubleScalar.H:105
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
error FatalError
static const char nl
Definition: Ostream.H:260
Foam::argList args(argc, argv)