geomCellLooper.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2014 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 Class
25  Foam::geomCellLooper
26 
27 Description
28  Implementation of cellLooper. Does pure geometric cut through cell.
29 
30  Handles all cell shapes in the same way: cut edges with plane through
31  cell centre and normal in direction of provided direction. Snaps cuts
32  close to edge endpoints (close = snapTol * minEdgeLen) to vertices.
33 
34  Currently determines cuts through edges (and edge endpoints close to plane)
35  in random order and then sorts them acc. to angle. Could be converted to
36  use walk but problem is that face can be cut multiple times (since does
37  not need to be convex). Another problem is that edges parallel to plane
38  might not be cut. So these are handled by looking at the distance from
39  edge endpoints to the plane.
40 
41 SourceFiles
42  geomCellLooper.C
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #ifndef geomCellLooper_H
47 #define geomCellLooper_H
48 
49 #include "cellLooper.H"
50 #include "typeInfo.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 
57 // Forward declaration of classes
58 class plane;
59 
60 /*---------------------------------------------------------------------------*\
61  Class geomCellLooper Declaration
62 \*---------------------------------------------------------------------------*/
63 
64 class geomCellLooper
65 :
66  public cellLooper
67 {
68 
69  // Static
70 
71  //- Tolerance for point equal test. Fraction of edge length.
72  static const scalar pointEqualTol_;
73 
74  //- Tolerance for cut through edges to get snapped to edge end point.
75  // Fraction of length of minimum connected edge length.
76  static scalar snapTol_;
77 
78 
79  // Private Member Functions
80 
81  //- Min length of attached edges
82  scalar minEdgeLen(const label vertI) const;
83 
84  //- Return true and set weight if edge is cut
85  bool cutEdge
86  (
87  const plane& cutPlane,
88  const label edgeI,
89  scalar& weight
90  ) const;
91 
92  //- Snaps cut through edge by cut through vertex (if weight closer than
93  // tol to 0 or 1). Returns vertex label snapped to or -1.
95  (
96  const scalar tol,
97  const label edgeI,
98  const scalar weight
99  ) const;
100 
101  //- Gets two (random) vectors perpendicular to n and each other to be
102  // used as base.
103  void getBase
104  (
105  const vector& n,
106  vector& e0,
107  vector& e1
108  ) const;
109 
110  //- Return true if the cut edge at loop[index] is inbetween the cuts
111  // through the edge end points.
112  bool edgeEndsCut(const labelList&, const label index) const;
113 
114 
115  //- Disallow default bitwise copy construct
117 
118  //- Disallow default bitwise assignment
119  void operator=(const geomCellLooper&);
120 
121 
122 public:
123 
124  //- Runtime type information
125  TypeName("geomCellLooper");
126 
127 
128  // Static Functions
129 
130  static scalar snapTol()
131  {
132  return snapTol_;
133  }
134 
135  static void setSnapTol(const scalar tol)
136  {
137  snapTol_ = tol;
138  }
139 
140 
141 
142 
143  // Constructors
144 
145  //- Construct from components
146  geomCellLooper(const polyMesh& mesh);
147 
148 
149  //- Destructor
150  virtual ~geomCellLooper();
151 
152 
153  // Member Functions
154 
155 
156 
157  //- Create cut along circumference of cellI. Gets current mesh cuts.
158  // Cut along circumference is expressed as loop of cuts plus weights
159  // for cuts along edges (only valid for edge cuts).
160  // Return true if successful cut.
161  virtual bool cut
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  //- Same but now also base point of cut provided (instead of always
174  // cell centre)
175  virtual bool cut
176  (
177  const plane& cutPlane,
178  const label cellI,
179  const boolList& vertIsCut,
180  const boolList& edgeIsCut,
181  const scalarField& edgeWeight,
182 
183  labelList& loop,
184  scalarField& loopWeights
185  ) const;
186 };
187 
188 
189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190 
191 } // End namespace Foam
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
195 #endif
196 
197 // ************************************************************************* //
typeInfo.H
Foam::geomCellLooper::snapTol_
static scalar snapTol_
Tolerance for cut through edges to get snapped to edge end point.
Definition: geomCellLooper.H:75
cellLooper.H
Foam::geomCellLooper::snapToVert
label snapToVert(const scalar tol, const label edgeI, const scalar weight) const
Snaps cut through edge by cut through vertex (if weight closer than.
Definition: geomCellLooper.C:111
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::plane
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
Foam::geomCellLooper::snapTol
static scalar snapTol()
Definition: geomCellLooper.H:129
n
label n
Definition: TABSMDCalcMethod2.H:31
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:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::geomCellLooper::getBase
void getBase(const vector &n, vector &e0, vector &e1) const
Gets two (random) vectors perpendicular to n and each other to be.
Definition: geomCellLooper.C:134
Foam::cellLooper
Abstract base class. Concrete implementations know how to cut a cell (i.e. determine a loop around th...
Definition: cellLooper.H:69
Foam::geomCellLooper
Implementation of cellLooper. Does pure geometric cut through cell.
Definition: geomCellLooper.H:63
Foam::geomCellLooper::TypeName
TypeName("geomCellLooper")
Runtime type information.
Foam::geomCellLooper::cutEdge
bool cutEdge(const plane &cutPlane, const label edgeI, scalar &weight) const
Return true and set weight if edge is cut.
Definition: geomCellLooper.C:81
Foam::geomCellLooper::geomCellLooper
geomCellLooper(const geomCellLooper &)
Disallow default bitwise copy construct.
Foam::geomCellLooper::pointEqualTol_
static const scalar pointEqualTol_
Tolerance for point equal test. Fraction of edge length.
Definition: geomCellLooper.H:71
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::geomCellLooper::setSnapTol
static void setSnapTol(const scalar tol)
Definition: geomCellLooper.H:134
Foam::geomCellLooper::~geomCellLooper
virtual ~geomCellLooper()
Destructor.
Definition: geomCellLooper.C:222
Foam::geomCellLooper::minEdgeLen
scalar minEdgeLen(const label vertI) const
Min length of attached edges.
Definition: geomCellLooper.C:62
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::geomCellLooper::operator=
void operator=(const geomCellLooper &)
Disallow default bitwise assignment.
Foam::geomCellLooper::edgeEndsCut
bool edgeEndsCut(const labelList &, const label index) const
Return true if the cut edge at loop[index] is inbetween the cuts.
Definition: geomCellLooper.C:179
Foam::edgeVertex::mesh
const polyMesh & mesh() const
Definition: edgeVertex.H:98
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:229