extendedFeatureEdgeMeshI.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-2016 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 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
27 
29 {
30  return convexStart_;
31 }
32 
33 
35 {
36  return concaveStart_;
37 }
38 
39 
41 {
42  return mixedStart_;
43 }
44 
45 
47 {
48  return nonFeatureStart_;
49 }
50 
51 
53 {
54  return externalStart_;
55 }
56 
57 
59 {
60  return internalStart_;
61 }
62 
63 
65 {
66  return flatStart_;
67 }
68 
69 
71 {
72  return openStart_;
73 }
74 
75 
77 {
78  return multipleStart_;
79 }
80 
81 
83 {
84  return ptI < nonFeatureStart_;
85 }
86 
87 
89 {
90  return normals_;
91 }
92 
95 {
96  return normalVolumeTypes_;
97 }
98 
100 const
101 {
102  return edgeDirections_;
103 }
104 
105 inline const Foam::labelListList&
107 {
108  return normalDirections_;
109 }
110 
111 
113 (
114  label edgeI,
115  label ptI
116 ) const
117 {
118  const edge& e = edges()[edgeI];
119 
120  if (ptI == e.start())
121  {
122  return edgeDirections()[edgeI];
123  }
124  else if (ptI == e.end())
125  {
126  return -edgeDirections()[edgeI];
127  }
128  else
129  {
131  << "Requested ptI " << ptI << " is not a point on the requested "
132  << "edgeI " << edgeI << ". edgeI start and end: "
133  << e.start() << " " << e.end()
134  << exit(FatalError);
135 
136  return Zero;
137  }
138 }
139 
140 
142 const
143 {
144  return edgeNormals_;
145 }
146 
147 
149 (
150  const labelList& edgeNormIs
151 ) const
152 {
153  vectorField norms(edgeNormIs.size());
154 
155  forAll(edgeNormIs, i)
156  {
157  norms[i] = normals_[edgeNormIs[i]];
158  }
159 
160  return norms;
161 }
162 
163 
165 const
166 {
167  return edgeNormals(edgeNormals_[edgeI]);
168 }
169 
170 
171 inline const Foam::labelListList&
173 {
174  return featurePointNormals_;
175 }
176 
177 
179 (
180  label ptI
181 ) const
182 {
183  if (!featurePoint(ptI))
184  {
186  << "Requesting the normals of a non-feature point. "
187  << "Returned zero length vectorField."
188  << endl;
189 
190  return vectorField(0);
191  }
192 
193  labelList featPtNormIs(featurePointNormals_[ptI]);
194 
195  vectorField norms(featPtNormIs.size());
196 
197  forAll(featPtNormIs, i)
198  {
199  norms[i] = normals_[featPtNormIs[i]];
200  }
201 
202  return norms;
203 }
204 
205 
206 inline const Foam::labelListList&
208 {
209  return featurePointEdges_;
210 }
211 
212 
214 {
215  return regionEdges_;
216 }
217 
218 
221 {
222  if (ptI < concaveStart_)
223  {
224  return CONVEX;
225  }
226  else if (ptI < mixedStart_)
227  {
228  return CONCAVE;
229  }
230  else if (ptI < nonFeatureStart_)
231  {
232  return MIXED;
233  }
234  else
235  {
236  return NONFEATURE;
237  }
238 }
239 
240 
243 {
244  if (edgeI < internalStart_)
245  {
246  return EXTERNAL;
247  }
248  else if (edgeI < flatStart_)
249  {
250  return INTERNAL;
251  }
252  else if (edgeI < openStart_)
253  {
254  return FLAT;
255  }
256  else if (edgeI < multipleStart_)
257  {
258  return OPEN;
259  }
260  else
261  {
262  return MULTIPLE;
263  }
264 }
265 
266 
268 (
269  label edgeI
270 ) const
271 {
272  const labelList& eNormals = edgeNormals_[edgeI];
273 
274  DynamicList<label> edgeBaffles(eNormals.size());
275 
276  forAll(eNormals, enI)
277  {
278  const label normI = eNormals[enI];
279 
280  if (normalVolumeTypes_[normI])
281  {
282  edgeBaffles.append(normI);
283  }
284  }
285 
286  return PackedList<2>(edgeBaffles);
287 }
288 
289 
290 // ************************************************************************* //
label convexStart() const
Return the index of the start of the convex feature points.
label concaveStart() const
Return the index of the start of the concave feature points.
label nonFeatureStart() const
Return the index of the start of the non-feature points.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
vectorField normals_
Normals of the features, to be referred to by index by both feature.
const labelListList & featurePointEdges() const
Return the edge labels for a given feature point. Edges are.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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
label mixedStart_
Index of the start of the mixed type feature points.
label openStart_
Index of the start of the open feature edges.
label externalStart() const
Return the index of the start of the external feature edges.
bool featurePoint(label ptI) const
Return whether or not the point index is a feature point.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
label concaveStart_
Index of the start of the concave feature points.
List< sideVolumeType > normalVolumeTypes_
Type per normal: which side of normal to mesh.
PackedList< 2 > edgeBaffles(label edgeI) const
Return the baffle faces of a specified edge.
label internalStart_
Index of the start of the internal feature edges.
labelListList normalDirections_
Starting directions for the edges.
const labelListList & featurePointNormals() const
Return the indices of the normals that are adjacent to the.
labelList regionEdges_
Feature edges which are on the boundary between regions.
static label convexStart_
Index of the start of the convex feature points - static as 0.
label multipleStart() const
Return the index of the start of the multiply-connected feature.
label flatStart_
Index of the start of the flat feature edges.
pointStatus getPointStatus(label ptI) const
Return the pointStatus of a specified point.
static label externalStart_
Index of the start of the external feature edges - static as 0.
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:91
vector edgeDirection(label edgeI, label ptI) const
Return the direction of edgeI, pointing away from ptI.
const labelList & regionEdges() const
Return the feature edges which are on the boundary between.
PackedList< nBits > & append(const unsigned int val)
Append a value at the end of the list.
Definition: PackedListI.H:1031
label openStart() const
Return the index of the start of the open feature edges.
const List< sideVolumeType > & normalVolumeTypes() const
Return.
label internalStart() const
Return the index of the start of the internal feature edges.
edgeStatus getEdgeStatus(label edgeI) const
Return the edgeStatus of a specified edge.
label multipleStart_
Index of the start of the multiply-connected feature edges.
#define WarningInFunction
Report a warning using Foam::Warning.
label flatStart() const
Return the index of the start of the flat feature edges.
const vectorField & normals() const
Return the normals of the surfaces adjacent to the feature edges.
const labelListList & normalDirections() const
const edgeList & edges() const
Return edges.
Definition: edgeMeshI.H:45
const labelListList & edgeNormals() const
Return the indices of the normals that are adjacent to the.
labelListList featurePointEdges_
Indices of feature edges attached to feature points. The edges are.
Field< vector > vectorField
Specialisation of Field<T> for vector.
label mixedStart() const
Return the index of the start of the mixed type feature points.
labelListList featurePointNormals_
Indices of the normals that are adjacent to the feature points.
const vectorField & edgeDirections() const
Return the edgeDirection vectors.
labelListList edgeNormals_
Indices of the normals that are adjacent to the feature edges.
label nonFeatureStart_
Index of the start of the non-feature points.
vectorField edgeDirections_
Flat and open edges require the direction of the edge.