extendedEdgeMeshI.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 
30 inline Foam::label Foam::extendedEdgeMesh::convexStart() const
31 {
32  return convexStart_;
33 }
34 
35 
36 inline Foam::label Foam::extendedEdgeMesh::concaveStart() const
37 {
38  return concaveStart_;
39 }
40 
41 
42 inline Foam::label Foam::extendedEdgeMesh::mixedStart() const
43 {
44  return mixedStart_;
45 }
46 
47 
48 inline Foam::label Foam::extendedEdgeMesh::nonFeatureStart() const
49 {
50  return nonFeatureStart_;
51 }
52 
53 
54 inline Foam::label Foam::extendedEdgeMesh::externalStart() const
55 {
56  return externalStart_;
57 }
58 
59 
60 inline Foam::label Foam::extendedEdgeMesh::internalStart() const
61 {
62  return internalStart_;
63 }
64 
65 
66 inline Foam::label Foam::extendedEdgeMesh::flatStart() const
67 {
68  return flatStart_;
69 }
70 
71 
72 inline Foam::label Foam::extendedEdgeMesh::openStart() const
73 {
74  return openStart_;
75 }
76 
77 
78 inline Foam::label Foam::extendedEdgeMesh::multipleStart() const
79 {
80  return multipleStart_;
81 }
82 
83 
84 inline bool Foam::extendedEdgeMesh::featurePoint(label ptI) const
85 {
86  return ptI < nonFeatureStart_;
87 }
88 
89 
91 {
92  return normals_;
93 }
94 
95 
98 {
99  return normalVolumeTypes_;
100 }
101 
102 
104 const
105 {
106  return edgeDirections_;
107 }
108 
109 
110 inline const Foam::labelListList&
112 {
113  return normalDirections_;
114 }
115 
116 
118 (
119  label edgeI,
120  label ptI
121 ) const
122 {
123  const edge& e = edges()[edgeI];
124 
125  if (ptI == e.start())
126  {
127  return edgeDirections()[edgeI];
128  }
129  else if (ptI == e.end())
130  {
131  return -edgeDirections()[edgeI];
132  }
133  else
134  {
136  << "Requested ptI " << ptI << " is not a point on the requested "
137  << "edgeI " << edgeI << ". edgeI start and end: "
138  << e.start() << " " << e.end()
140 
141  return Zero;
142  }
143 }
144 
145 
147 const
148 {
149  return edgeNormals_;
150 }
151 
152 
154 (
155  const labelList& edgeNormIs
156 ) const
157 {
158  vectorField norms(edgeNormIs.size());
159 
160  forAll(edgeNormIs, i)
161  {
162  norms[i] = normals_[edgeNormIs[i]];
163  }
164 
165  return norms;
166 }
167 
168 
170 const
171 {
172  return edgeNormals(edgeNormals_[edgeI]);
173 }
174 
175 
176 inline const Foam::labelListList&
178 {
179  return featurePointNormals_;
180 }
181 
182 
184 (
185  label ptI
186 ) const
187 {
188  if (!featurePoint(ptI))
189  {
191  << "Requesting the normals of a non-feature point. "
192  << "Returned zero length vectorField."
193  << endl;
194 
195  return vectorField(0);
196  }
197 
198  labelList featPtNormIs(featurePointNormals_[ptI]);
199 
200  vectorField norms(featPtNormIs.size());
201 
202  forAll(featPtNormIs, i)
203  {
204  norms[i] = normals_[featPtNormIs[i]];
205  }
206 
207  return norms;
208 }
209 
210 
211 inline const Foam::labelListList&
213 {
214  return featurePointEdges_;
215 }
216 
217 
219 {
220  return regionEdges_;
221 }
222 
223 
226 {
227  if (ptI < concaveStart_)
228  {
229  return CONVEX;
230  }
231  else if (ptI < mixedStart_)
232  {
233  return CONCAVE;
234  }
235  else if (ptI < nonFeatureStart_)
236  {
237  return MIXED;
238  }
239 
240  return NONFEATURE;
241 }
242 
243 
245 Foam::extendedEdgeMesh::getEdgeStatus(label edgeI) const
246 {
247  if (edgeI < internalStart_)
248  {
249  return EXTERNAL;
250  }
251  else if (edgeI < flatStart_)
252  {
253  return INTERNAL;
254  }
255  else if (edgeI < openStart_)
256  {
257  return FLAT;
258  }
259  else if (edgeI < multipleStart_)
260  {
261  return OPEN;
262  }
263 
264  return MULTIPLE;
265 }
266 
267 
269 (
270  label edgeI
271 ) const
272 {
273  const labelList& eNormals = edgeNormals_[edgeI];
274 
275  DynamicList<label> edgeBaffles(eNormals.size());
276 
277  for (const label normi : eNormals)
278  {
279  if (normalVolumeTypes_[normi])
280  {
281  edgeBaffles.append(normi);
282  }
283  }
284 
285  return PackedList<2>(edgeBaffles);
286 }
287 
288 
289 // ************************************************************************* //
label nonFeatureStart() const
Return the index of the start of the non-feature points.
pointStatus getPointStatus(label ptI) const
Return the pointStatus of a specified point.
label convexStart() const
Return the index of the start of the convex feature points.
const labelListList & featurePointNormals() const
Return the indices of the normals that are adjacent to the.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const labelListList & edgeNormals() const
Return the indices of the normals that are adjacent to the.
const labelListList & normalDirections() const
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Mixed uniform/non-uniform (eg, after reduction)
Definition: ListPolicy.H:108
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:517
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
const List< sideVolumeType > & normalVolumeTypes() const
Return.
label multipleStart() const
Return the index of the start of the multiply-connected feature.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
bool featurePoint(label ptI) const
Return whether or not the point index is a feature point.
const vectorField & normals() const
Return the normals of the surfaces adjacent to the feature edges.
label openStart() const
Return the index of the start of the open feature edges.
const labelList & regionEdges() const
Return the feature edges which are on the boundary between.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
edgeStatus getEdgeStatus(label edgeI) const
Return the edgeStatus of a specified edge.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
label flatStart() const
Return the index of the start of the flat feature edges.
static label convexStart_
Index of the start of the convex feature points - static as 0.
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
label internalStart() const
Return the index of the start of the internal feature edges.
label mixedStart() const
Return the index of the start of the mixed type feature points.
label externalStart() const
Return the index of the start of the external feature edges.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
PackedList< 2 > edgeBaffles(label edgeI) const
Return the baffle faces of a specified edge.
const labelListList & featurePointEdges() const
Return the edge labels for a given feature point. Edges are.
const vectorField & edgeDirections() const
Return the edgeDirection vectors.
#define WarningInFunction
Report a warning using Foam::Warning.
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
vector edgeDirection(label edgeI, label ptI) const
Return the direction of edgeI, pointing away from ptI.
label concaveStart() const
Return the index of the start of the concave feature points.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127