hexCellLooper.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "hexCellLooper.H"
30 #include "cellFeatures.H"
31 #include "polyMesh.H"
32 #include "cellModel.H"
33 #include "plane.H"
34 #include "ListOps.H"
35 #include "meshTools.H"
36 #include "OFstream.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(hexCellLooper, 0);
44  addToRunTimeSelectionTable(cellLooper, hexCellLooper, word);
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 // Starting from cut edge start walking.
51 bool Foam::hexCellLooper::walkHex
52 (
53  const label celli,
54  const label startFacei,
55  const label startEdgeI,
56 
57  labelList& loop,
58  scalarField& loopWeights
59 ) const
60 {
61  label facei = startFacei;
62 
63  label edgeI = startEdgeI;
64 
65  label cutI = 0;
66 
67  do
68  {
69  if (debug & 2)
70  {
71  Pout<< " walkHex : inserting cut onto edge:" << edgeI
72  << " vertices:" << mesh().edges()[edgeI] << endl;
73  }
74 
75  // Store cut through edge. For now cut edges halfway.
76  loop[cutI] = edgeToEVert(edgeI);
77  loopWeights[cutI] = 0.5;
78  cutI++;
79 
80  facei = meshTools::otherFace(mesh(), celli, facei, edgeI);
81 
82  const edge& e = mesh().edges()[edgeI];
83 
84  // Walk two edges further
85  edgeI = meshTools::walkFace(mesh(), facei, edgeI, e.end(), 2);
86 
87  if (edgeI == startEdgeI)
88  {
89  break;
90  }
91  }
92  while (true);
93 
94  // Checks.
95  if (cutI > 4)
96  {
97  Pout<< "hexCellLooper::walkHex" << "Problem : cell:" << celli
98  << " collected loop:";
99  writeCuts(Pout, loop, loopWeights);
100  Pout<< "loopWeights:" << loopWeights << endl;
101 
102  return false;
103  }
104 
105  return true;
106 }
107 
108 
109 
110 void Foam::hexCellLooper::makeFace
111 (
112  const labelList& loop,
113  const scalarField& loopWeights,
114 
115  labelList& faceVerts,
116  pointField& facePoints
117 ) const
118 {
119  facePoints.setSize(loop.size());
120  faceVerts.setSize(loop.size());
121 
122  forAll(loop, cutI)
123  {
124  label cut = loop[cutI];
125 
126  if (isEdge(cut))
127  {
128  label edgeI = getEdge(cut);
129 
130  const edge& e = mesh().edges()[edgeI];
131 
132  const point& v0 = mesh().points()[e.start()];
133  const point& v1 = mesh().points()[e.end()];
134 
135  facePoints[cutI] =
136  loopWeights[cutI]*v1 + (1.0-loopWeights[cutI])*v0;
137  }
138  else
139  {
140  label vertI = getVertex(cut);
141 
142  facePoints[cutI] = mesh().points()[vertI];
143  }
144 
145  faceVerts[cutI] = cutI;
146  }
147 }
148 
149 
150 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
151 
152 Foam::hexCellLooper::hexCellLooper(const polyMesh& mesh)
153 :
155  hex_(cellModel::ref(cellModel::HEX))
156 {}
157 
158 
159 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
160 
162 (
163  const vector& refDir,
164  const label celli,
165  const boolList& vertIsCut,
166  const boolList& edgeIsCut,
167  const scalarField& edgeWeight,
168 
169  labelList& loop,
170  scalarField& loopWeights
171 ) const
172 {
173  bool success = false;
174 
175  if (mesh().cellShapes()[celli].model() == hex_)
176  {
177  // Get starting edge. Note: should be compatible with way refDir is
178  // determined.
179  label edgeI = meshTools::cutDirToEdge(mesh(), celli, refDir);
180 
181  // Get any face using edge
182  label face0;
183  label face1;
184  meshTools::getEdgeFaces(mesh(), celli, edgeI, face0, face1);
185 
186  // Walk circumference of hex, cutting edges only
187  loop.setSize(4);
188  loopWeights.setSize(4);
189 
190  success = walkHex(celli, face0, edgeI, loop, loopWeights);
191  }
192  else
193  {
195  (
196  refDir,
197  celli,
198  vertIsCut,
199  edgeIsCut,
200  edgeWeight,
201 
202  loop,
203  loopWeights
204  );
205  }
206 
207  if (debug)
208  {
209  if (loop.empty())
210  {
212  << "could not cut cell " << celli << endl;
213 
214  fileName cutsFile("hexCellLooper_" + name(celli) + ".obj");
215 
216  Pout<< "hexCellLooper : writing cell to " << cutsFile << endl;
217 
218  OFstream cutsStream(cutsFile);
219 
221  (
222  cutsStream,
223  mesh().cells(),
224  mesh().faces(),
225  mesh().points(),
226  labelList(1, celli)
227  );
228 
229  return false;
230  }
231 
232  // Check for duplicate cuts.
233  labelHashSet loopSet(loop.size());
234 
235  forAll(loop, elemI)
236  {
237  label elem = loop[elemI];
238 
239  if (loopSet.found(elem))
240  {
242  << abort(FatalError);
243  }
244  loopSet.insert(elem);
245  }
246 
247 
248  face faceVerts(loop.size());
249  pointField facePoints(loop.size());
250 
251  makeFace(loop, loopWeights, faceVerts, facePoints);
252 
253  if ((faceVerts.mag(facePoints) < SMALL) || (loop.size() < 3))
254  {
256  << " on points:" << facePoints << endl
257  << UIndirectList<point>(facePoints, faceVerts)
258  << abort(FatalError);
259  }
260  }
261  return success;
262 }
263 
264 
265 // No shortcuts for cutting with arbitrary plane.
267 (
268  const plane& cutPlane,
269  const label celli,
270  const boolList& vertIsCut,
271  const boolList& edgeIsCut,
272  const scalarField& edgeWeight,
273 
274  labelList& loop,
275  scalarField& loopWeights
276 ) const
277 {
278  return
280  (
281  cutPlane,
282  celli,
283  vertIsCut,
284  edgeIsCut,
285  edgeWeight,
286 
287  loop,
288  loopWeights
289  );
290 }
291 
292 
293 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
meshTools.H
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::meshTools::otherFace
label otherFace(const primitiveMesh &mesh, const label celli, const label facei, const label edgeI)
Return face on cell using edgeI but not facei. Throws error.
Definition: meshTools.C:555
success
bool success
Definition: LISASMDCalcMethod1.H:16
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Foam::geomCellLooper::cut
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.
Definition: geomCellLooper.C:217
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::primitiveMesh::edges
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Definition: primitiveMeshEdges.C:505
ref
rDeltaT ref()
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
polyMesh.H
Foam::HashSet< label, Hash< label > >
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
OFstream.H
Foam::meshTools::walkFace
label walkFace(const primitiveMesh &mesh, 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:603
Foam::plane
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:89
Foam::edgeVertex::edgeToEVert
static label edgeToEVert(const primitiveMesh &mesh, const label edgeI)
Convert edgeI to eVert.
Definition: edgeVertex.H:177
Foam::Field< scalar >
plane.H
cellModel.H
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::meshTools::cutDirToEdge
label cutDirToEdge(const primitiveMesh &mesh, const label celli, const vector &cutDir)
Reverse of edgeToCutDir: given direction find edge bundle and.
Definition: meshTools.C:810
Foam::meshTools::getEdgeFaces
void getEdgeFaces(const primitiveMesh &mesh, 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:479
Foam::geomCellLooper
Implementation of cellLooper. Does pure geometric cut through cell.
Definition: geomCellLooper.H:66
cut
Patchify triangles based on orientation w.r.t other (triangulated or triangulatable) surfaces.
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
cellFeatures.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::Vector< scalar >
Foam::List< bool >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
cellShapes
cellShapeList cellShapes
Definition: createBlockMesh.H:3
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
ListOps.H
Various functions to operate on Lists.
Foam::cellModel
Maps a geometry to a set of cell primitives.
Definition: cellModel.H:72
Foam::hexCellLooper::cut
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.
Definition: hexCellLooper.C:162
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::edgeVertex::writeCuts
Ostream & writeCuts(Ostream &os, const labelList &, const scalarField &) const
Write cut descriptions to Ostream.
Definition: edgeVertex.C:240
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
hexCellLooper.H
Foam::edgeVertex::mesh
const polyMesh & mesh() const
Definition: edgeVertex.H:101