Newer
Older
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Description
Splits mesh by making internal faces external. Uses attachDetach.
Generates a meshModifier of the form:
Splitter
{
type attachDetach;
faceZoneName membraneFaces;
masterPatchName masterPatch;
slavePatchName slavePatch;
triggerTimes runTime.value();
}
so will detach at the current time and split all faces in membraneFaces
into masterPatch and slavePatch (which have to be present but of 0 size)
\*---------------------------------------------------------------------------*/
#include "argList.H"
#include "polyMesh.H"
#include "Time.H"
#include "polyTopoChange.H"
#include "mapPolyMesh.H"
#include "faceSet.H"
#include "attachDetach.H"
#include "attachPolyTopoChanger.H"
#include "regionSide.H"
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
using namespace Foam;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// Find edge between points v0 and v1.
label findEdge(const primitiveMesh& mesh, const label v0, const label v1)
{
const labelList& pEdges = mesh.pointEdges()[v0];
forAll(pEdges, pEdgeI)
{
label edgeI = pEdges[pEdgeI];
const edge& e = mesh.edges()[edgeI];
if (e.otherVertex(v0) == v1)
{
return edgeI;
}
}
FatalErrorIn
(
"findEdge(const primitiveMesh&, const label, const label)"
) << "Cannot find edge between mesh points " << v0 << " and " << v1
<< abort(FatalError);
return -1;
}
// Checks whether patch present
void checkPatch(const polyBoundaryMesh& bMesh, const word& name)
{
label patchI = bMesh.findPatchID(name);
if (patchI == -1)
{
FatalErrorIn("checkPatch(const polyBoundaryMesh&, const word&)")
<< "Cannot find patch " << name << endl
<< "It should be present but of zero size" << endl
<< "Valid patches are " << bMesh.names()
<< exit(FatalError);
}
if (bMesh[patchI].size())
{
FatalErrorIn("checkPatch(const polyBoundaryMesh&, const word&)")
<< "Patch " << name << " is present but non-zero size"
<< exit(FatalError);
}
}
// Main program:
int main(int argc, char *argv[])
{
Foam::argList::noParallel();
Foam::argList::validArgs.append("faceSet");
Foam::argList::validArgs.append("masterPatch");
Foam::argList::validArgs.append("slavePatch");
Foam::argList::validOptions.insert("overwrite", "");
# include "setRootCase.H"
# include "createTime.H"
# include "createPolyMesh.H"
word setName(args.additionalArgs()[0]);
word masterPatch(args.additionalArgs()[1]);
word slavePatch(args.additionalArgs()[2]);
bool overwrite = args.optionFound("overwrite");
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
// List of faces to split
faceSet facesSet(mesh, setName);
Info<< "Read " << facesSet.size() << " faces to split" << endl << endl;
// Convert into labelList and check
labelList faces(facesSet.toc());
forAll(faces, i)
{
if (!mesh.isInternalFace(faces[i]))
{
FatalErrorIn(args.executable())
<< "Face " << faces[i] << " in faceSet " << setName
<< " is not an internal face."
<< exit(FatalError);
}
}
// Check for empty master and slave patches
checkPatch(mesh.boundaryMesh(), masterPatch);
checkPatch(mesh.boundaryMesh(), slavePatch);
//
// Find 'side' of all faces on splitregion. Uses regionSide which needs
// set of edges on side of this region. Use PrimitivePatch to find these.
//
// Addressing on faces only in mesh vertices.
UIndirectList<face>
(
mesh.faces(),
faces
)
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
const labelList& meshPoints = fPatch.meshPoints();
// Mark all fence edges : edges on boundary of fPatch but not on boundary
// of polyMesh
labelHashSet fenceEdges(fPatch.size());
const labelListList& allEdgeFaces = fPatch.edgeFaces();
forAll(allEdgeFaces, patchEdgeI)
{
if (allEdgeFaces[patchEdgeI].size() == 1)
{
const edge& e = fPatch.edges()[patchEdgeI];
label edgeI =
findEdge
(
mesh,
meshPoints[e.start()],
meshPoints[e.end()]
);
fenceEdges.insert(edgeI);
}
}
// Find sides reachable from 0th face of faceSet
label startFaceI = faces[0];
regionSide regionInfo
(
mesh,
facesSet,
fenceEdges,
mesh.faceOwner()[startFaceI],
startFaceI
);
// Determine flip state for all faces in faceSet
boolList zoneFlip(faces.size());
forAll(faces, i)
{
zoneFlip[i] = !regionInfo.sideOwner().found(faces[i]);
}
// Create and add face zones and mesh modifiers
List<pointZone*> pz(0);
List<faceZone*> fz(1);
List<cellZone*> cz(0);
fz[0] =
new faceZone
(
"membraneFaces",
faces,
zoneFlip,
0,
mesh.faceZones()
);
Info << "Adding point and face zones" << endl;
mesh.addZones(pz, fz, cz);
attachPolyTopoChanger splitter(mesh);
splitter.setSize(1);
// Add the sliding interface mesh modifier to start working at current
// time
splitter.set
(
0,
new attachDetach
(
"Splitter",
0,
splitter,
"membraneFaces",
masterPatch,
slavePatch,
scalarField(1, runTime.value())
)
);
Info<< nl << "Constructed topologyModifier:" << endl;
splitter[0].writeDict(Info);
splitter.attach();
if (overwrite)
{
mesh.setInstance(oldInstance);
}
else
{
mesh.setInstance(runTime.timeName());
}
Info<< "Writing mesh to " << runTime.timeName() << endl;
if (!mesh.write())
{
FatalErrorIn(args.executable())
<< "Failed writing polyMesh."
<< exit(FatalError);
}
Info<< nl << "end" << endl;
return 0;
}
// ************************************************************************* //