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"
37 
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 defineTypeNameAndDebug(hexCellLooper, 0);
45 addToRunTimeSelectionTable(cellLooper, hexCellLooper, word);
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 // Starting from cut edge start walking.
52 bool Foam::hexCellLooper::walkHex
53 (
54  const label celli,
55  const label startFacei,
56  const label startEdgeI,
57 
58  labelList& loop,
59  scalarField& loopWeights
60 ) const
61 {
62  label facei = startFacei;
63 
64  label edgeI = startEdgeI;
65 
66  label cutI = 0;
67 
68  do
69  {
70  if (debug & 2)
71  {
72  Pout<< " walkHex : inserting cut onto edge:" << edgeI
73  << " vertices:" << mesh().edges()[edgeI] << endl;
74  }
75 
76  // Store cut through edge. For now cut edges halfway.
77  loop[cutI] = edgeToEVert(edgeI);
78  loopWeights[cutI] = 0.5;
79  cutI++;
80 
81  facei = meshTools::otherFace(mesh(), celli, facei, edgeI);
82 
83  const edge& e = mesh().edges()[edgeI];
84 
85  // Walk two edges further
86  edgeI = meshTools::walkFace(mesh(), facei, edgeI, e.end(), 2);
87 
88  if (edgeI == startEdgeI)
89  {
90  break;
91  }
92  }
93  while (true);
94 
95  // Checks.
96  if (cutI > 4)
97  {
98  Pout<< "hexCellLooper::walkHex" << "Problem : cell:" << celli
99  << " collected loop:";
100  writeCuts(Pout, loop, loopWeights);
101  Pout<< "loopWeights:" << loopWeights << endl;
102 
103  return false;
104  }
105 
106  return true;
107 }
108 
109 
110 
111 void Foam::hexCellLooper::makeFace
112 (
113  const labelList& loop,
114  const scalarField& loopWeights,
115 
116  labelList& faceVerts,
117  pointField& facePoints
118 ) const
119 {
120  facePoints.setSize(loop.size());
121  faceVerts.setSize(loop.size());
122 
123  forAll(loop, cutI)
124  {
125  label cut = loop[cutI];
126 
127  if (isEdge(cut))
128  {
129  label edgeI = getEdge(cut);
130 
131  const edge& e = mesh().edges()[edgeI];
132 
133  const point& v0 = mesh().points()[e.start()];
134  const point& v1 = mesh().points()[e.end()];
135 
136  facePoints[cutI] =
137  loopWeights[cutI]*v1 + (1.0-loopWeights[cutI])*v0;
138  }
139  else
140  {
141  label vertI = getVertex(cut);
142 
143  facePoints[cutI] = mesh().points()[vertI];
144  }
145 
146  faceVerts[cutI] = cutI;
147  }
148 }
149 
150 
151 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
152 
153 // Construct from components
154 Foam::hexCellLooper::hexCellLooper(const polyMesh& mesh)
155 :
157  hex_(cellModel::ref(cellModel::HEX))
158 {}
159 
160 
161 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
162 
164 {}
165 
166 
167 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
168 
170 (
171  const vector& refDir,
172  const label celli,
173  const boolList& vertIsCut,
174  const boolList& edgeIsCut,
175  const scalarField& edgeWeight,
176 
177  labelList& loop,
178  scalarField& loopWeights
179 ) const
180 {
181  bool success = false;
182 
183  if (mesh().cellShapes()[celli].model() == hex_)
184  {
185  // Get starting edge. Note: should be compatible with way refDir is
186  // determined.
187  label edgeI = meshTools::cutDirToEdge(mesh(), celli, refDir);
188 
189  // Get any face using edge
190  label face0;
191  label face1;
192  meshTools::getEdgeFaces(mesh(), celli, edgeI, face0, face1);
193 
194  // Walk circumference of hex, cutting edges only
195  loop.setSize(4);
196  loopWeights.setSize(4);
197 
198  success = walkHex(celli, face0, edgeI, loop, loopWeights);
199  }
200  else
201  {
203  (
204  refDir,
205  celli,
206  vertIsCut,
207  edgeIsCut,
208  edgeWeight,
209 
210  loop,
211  loopWeights
212  );
213  }
214 
215  if (debug)
216  {
217  if (loop.empty())
218  {
220  << "could not cut cell " << celli << endl;
221 
222  fileName cutsFile("hexCellLooper_" + name(celli) + ".obj");
223 
224  Pout<< "hexCellLooper : writing cell to " << cutsFile << endl;
225 
226  OFstream cutsStream(cutsFile);
227 
229  (
230  cutsStream,
231  mesh().cells(),
232  mesh().faces(),
233  mesh().points(),
234  labelList(1, celli)
235  );
236 
237  return false;
238  }
239 
240  // Check for duplicate cuts.
241  labelHashSet loopSet(loop.size());
242 
243  forAll(loop, elemI)
244  {
245  label elem = loop[elemI];
246 
247  if (loopSet.found(elem))
248  {
250  << abort(FatalError);
251  }
252  loopSet.insert(elem);
253  }
254 
255 
256  face faceVerts(loop.size());
257  pointField facePoints(loop.size());
258 
259  makeFace(loop, loopWeights, faceVerts, facePoints);
260 
261  if ((faceVerts.mag(facePoints) < SMALL) || (loop.size() < 3))
262  {
264  << " on points:" << facePoints << endl
265  << UIndirectList<point>(facePoints, faceVerts)
266  << abort(FatalError);
267  }
268  }
269  return success;
270 }
271 
272 
273 // No shortcuts for cutting with arbitrary plane.
275 (
276  const plane& cutPlane,
277  const label celli,
278  const boolList& vertIsCut,
279  const boolList& edgeIsCut,
280  const scalarField& edgeWeight,
281 
282  labelList& loop,
283  scalarField& loopWeights
284 ) const
285 {
286  return
288  (
289  cutPlane,
290  celli,
291  vertIsCut,
292  edgeIsCut,
293  edgeWeight,
294 
295  loop,
296  loopWeights
297  );
298 }
299 
300 
301 // ************************************************************************* //
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:74
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:1038
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
meshTools.H
ref
rDeltaT ref()
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
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:228
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::primitiveMesh::edges
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Definition: primitiveMeshEdges.C:505
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
polyMesh.H
Foam::HashSet< label, Hash< label > >
Foam::hexCellLooper::~hexCellLooper
virtual ~hexCellLooper()
Destructor.
Definition: hexCellLooper.C:163
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:290
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:176
Foam::label
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:62
Foam::Field< scalar >
plane.H
cellModel.H
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
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:65
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:137
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:99
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
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::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:109
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
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:170
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:241
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
hexCellLooper.H
Foam::edgeVertex::mesh
const polyMesh & mesh() const
Definition: edgeVertex.H:100