hexCellLooper.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-2024 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 \*---------------------------------------------------------------------------*/
25 
26 #include "hexCellLooper.H"
27 #include "cellFeatures.H"
28 #include "polyMesh.H"
29 #include "cellModeller.H"
30 #include "plane.H"
31 #include "ListOps.H"
32 #include "meshTools.H"
33 #include "OFstream.H"
34 
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 // Starting from cut edge start walking.
48 bool Foam::hexCellLooper::walkHex
49 (
50  const label celli,
51  const label startFacei,
52  const label startEdgeI,
53 
54  labelList& loop,
55  scalarField& loopWeights
56 ) const
57 {
58  label facei = startFacei;
59 
60  label edgeI = startEdgeI;
61 
62  label cutI = 0;
63 
64  do
65  {
66  if (debug & 2)
67  {
68  Pout<< " walkHex : inserting cut onto edge:" << edgeI
69  << " vertices:" << mesh().edges()[edgeI] << endl;
70  }
71 
72  // Store cut through edge. For now cut edges halfway.
73  loop[cutI] = edgeToEVert(edgeI);
74  loopWeights[cutI] = 0.5;
75  cutI++;
76 
77  facei = meshTools::otherFace(mesh(), celli, facei, edgeI);
78 
79  const edge& e = mesh().edges()[edgeI];
80 
81  // Walk two edges further
82  edgeI = meshTools::walkFace(mesh(), facei, edgeI, e.end(), 2);
83 
84  if (edgeI == startEdgeI)
85  {
86  break;
87  }
88  }
89  while (true);
90 
91  // Checks.
92  if (cutI > 4)
93  {
94  Pout<< "hexCellLooper::walkHex" << "Problem : cell:" << celli
95  << " collected loop:";
96  writeCuts(Pout, loop, loopWeights);
97  Pout<< "loopWeights:" << loopWeights << endl;
98 
99  return false;
100  }
101  else
102  {
103  return true;
104  }
105 }
106 
107 
108 
109 void Foam::hexCellLooper::makeFace
110 (
111  const labelList& loop,
112  const scalarField& loopWeights,
113 
114  labelList& faceVerts,
115  pointField& facePoints
116 ) const
117 {
118  facePoints.setSize(loop.size());
119  faceVerts.setSize(loop.size());
120 
121  forAll(loop, cutI)
122  {
123  label cut = loop[cutI];
124 
125  if (isEdge(cut))
126  {
127  label edgeI = getEdge(cut);
128 
129  const edge& e = mesh().edges()[edgeI];
130 
131  const point& v0 = mesh().points()[e.start()];
132  const point& v1 = mesh().points()[e.end()];
133 
134  facePoints[cutI] =
135  loopWeights[cutI]*v1 + (1.0-loopWeights[cutI])*v0;
136  }
137  else
138  {
139  label vertI = getVertex(cut);
140 
141  facePoints[cutI] = mesh().points()[vertI];
142  }
143 
144  faceVerts[cutI] = cutI;
145  }
146 }
147 
148 
149 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
150 
151 // Construct from components
153 :
154  geomCellLooper(mesh),
155  hex_(*(cellModeller::lookup("hex")))
156 {}
157 
158 
159 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
160 
162 {}
163 
164 
165 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
166 
168 (
169  const vector& refDir,
170  const label celli,
171  const boolList& vertIsCut,
172  const boolList& edgeIsCut,
173  const scalarField& edgeWeight,
174 
175  labelList& loop,
176  scalarField& loopWeights
177 ) const
178 {
179  bool success = false;
180 
181  if (mesh().cellShapes()[celli].model() == hex_)
182  {
183  // Get starting edge. Note: should be compatible with way refDir is
184  // determined.
185  label edgeI = meshTools::cutDirToEdge(mesh(), celli, refDir);
186 
187  // Get any face using edge
188  label face0;
189  label face1;
190  meshTools::getEdgeFaces(mesh(), celli, edgeI, face0, face1);
191 
192  // Walk circumference of hex, cutting edges only
193  loop.setSize(4);
194  loopWeights.setSize(4);
195 
196  success = walkHex(celli, face0, edgeI, loop, loopWeights);
197  }
198  else
199  {
201  (
202  refDir,
203  celli,
204  vertIsCut,
205  edgeIsCut,
206  edgeWeight,
207 
208  loop,
209  loopWeights
210  );
211  }
212 
213  if (debug)
214  {
215  if (loop.empty())
216  {
218  << "could not cut cell " << celli << endl;
219 
220  fileName cutsFile("hexCellLooper_" + name(celli) + ".obj");
221 
222  Pout<< "hexCellLooper : writing cell to " << cutsFile << endl;
223 
224  OFstream cutsStream(cutsFile);
225 
227  (
228  cutsStream,
229  mesh().cells(),
230  mesh().faces(),
231  mesh().points(),
232  labelList(1, celli)
233  );
234 
235  return false;
236  }
237 
238  // Check for duplicate cuts.
239  labelHashSet loopSet(loop.size());
240 
241  forAll(loop, elemI)
242  {
243  label elem = loop[elemI];
244 
245  if (loopSet.found(elem))
246  {
248  << abort(FatalError);
249  }
250  loopSet.insert(elem);
251  }
252 
253 
254  face faceVerts(loop.size());
255  pointField facePoints(loop.size());
256 
257  makeFace(loop, loopWeights, faceVerts, facePoints);
258 
259  if ((faceVerts.mag(facePoints) < small) || (loop.size() < 3))
260  {
262  << " on points:" << facePoints << endl
263  << UIndirectList<point>(facePoints, faceVerts)()
264  << abort(FatalError);
265  }
266  }
267  return success;
268 }
269 
270 
271 // No shortcuts for cutting with arbitrary plane.
273 (
274  const plane& cutPlane,
275  const label celli,
276  const boolList& vertIsCut,
277  const boolList& edgeIsCut,
278  const scalarField& edgeWeight,
279 
280  labelList& loop,
281  scalarField& loopWeights
282 ) const
283 {
284  return
286  (
287  cutPlane,
288  celli,
289  vertIsCut,
290  edgeIsCut,
291  edgeWeight,
292 
293  loop,
294  loopWeights
295  );
296 }
297 
298 
299 // ************************************************************************* //
bool success
Various functions to operate on Lists.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
Output to file stream.
Definition: OFstream.H:86
A List with indirect addressing.
Definition: UIndirectList.H:60
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
A static collection of cell models, and a means of looking them up.
Definition: cellModeller.H:53
Ostream & writeCuts(Ostream &os, const labelList &, const scalarField &) const
Write cut descriptions to Ostream.
Definition: edgeVertex.C:239
static label edgeToEVert(const primitiveMesh &mesh, const label edgeI)
Convert edgeI to eVert.
Definition: edgeVertex.H:169
const polyMesh & mesh() const
Definition: edgeVertex.H:93
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
scalar mag(const pointField &) const
Return scalar magnitude.
Definition: faceI.H:81
A class for handling file names.
Definition: fileName.H:82
Implementation of cellLooper. Does pure geometric cut through cell.
virtual bool cut(const vector &refDir, const label celli, const boolList &vertIsCut, const boolList &edgeIsCut, const scalarField &edgeWeight, labelList &loop, scalarField &loopWeights) const
Create cut along circumference of celli. Gets current mesh cuts.
Implementation of cellLooper.
Definition: hexCellLooper.H:63
hexCellLooper(const polyMesh &mesh)
Construct from components.
virtual ~hexCellLooper()
Destructor.
virtual bool cut(const vector &refDir, const label celli, const boolList &vertIsCut, const boolList &edgeIsCut, const scalarField &edgeWeight, labelList &loop, scalarField &loopWeights) const
Create cut along circumference of celli. Gets current mesh cuts.
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:61
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const cellShapeList & cellShapes
volScalarField scalarField(fieldObject, mesh)
const pointField & points
const cellShapeList & cells
#define WarningInFunction
Report a warning using Foam::Warning.
void getEdgeFaces(const primitiveMesh &, const label celli, const label edgeI, label &face0, label &face1)
Get faces on cell using edgeI. Throws error if no two found.
Definition: meshTools.C:458
label walkFace(const primitiveMesh &, const label facei, const label startEdgeI, const label startVertI, const label nEdges)
Returns label of edge nEdges away from startEdge (in the direction.
Definition: meshTools.C:587
label cutDirToEdge(const primitiveMesh &, const label celli, const vector &cutDir)
Reverse of edgeToCutDir: given direction find edge bundle and.
Definition: meshTools.C:795
label otherFace(const primitiveMesh &, const label celli, const label facei, const label edgeI)
Return face on cell using edgeI but not facei. Throws error.
Definition: meshTools.C:536
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
List< label > labelList
A List of labels.
Definition: labelList.H:56
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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
errorManip< error > abort(error &err)
Definition: errorManip.H:131
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
vector point
Point is a vector.
Definition: point.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError