faceAreaIntersectI.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-2012 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
27 
29 (
30  const point& a,
31  const point& b,
32  const point& c,
33  label& count,
35 ) const
36 {
37  triPoints& tp = tris[count++];
38  tp[0] = a;
39  tp[1] = b;
40  tp[2] = c;
41 }
42 
43 
45 (
46  const pointField& points,
47  const face& f,
48  const bool reverse
49 ) const
50 {
51  triPoints result;
52 
53  if (reverse)
54  {
55  result[2] = points[f[0]];
56  result[1] = points[f[1]];
57  result[0] = points[f[2]];
58  }
59  else
60  {
61  result[0] = points[f[0]];
62  result[1] = points[f[1]];
63  result[2] = points[f[2]];
64  }
65 
66  return result;
67 }
68 
69 
71 (
72  const face& f,
73  DynamicList<face>& faces
74 ) const
75 {
76  if (f.size() > 2)
77  {
78  const label v0 = 0;
79 
80  labelList indices(3);
81 
82  for (label i = 1; i < f.size() - 1; i++)
83  {
84  indices[0] = f[v0];
85  indices[1] = f[i];
86  indices[2] = f[i + 1];
87  faces.append(face(indices));
88  }
89  }
90 }
91 
92 
94 (
95  const FixedList<scalar, 3>& d,
96  const triPoints& t,
97  const label negI,
98  const label posI
99 ) const
100 {
101  return (d[posI]*t[negI] - d[negI]*t[posI])/(-d[negI] + d[posI]);
102 }
103 
104 
105 inline Foam::scalar Foam::faceAreaIntersect::triArea(const triPoints& t) const
106 {
107  return mag(0.5*((t[1] - t[0])^(t[2] - t[0])));
108 }
109 
110 
111 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
112 
114 {
115  return tol;
116 }
117 
118 
119 // ************************************************************************* //
Foam::faceAreaIntersect::triangleFan
void triangleFan(const face &f, DynamicList< face > &faces) const
Decompose face into triangle fan.
Definition: faceAreaIntersectI.H:71
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::faceAreaIntersect::setTriPoints
void setTriPoints(const point &a, const point &b, const point &c, label &count, FixedList< triPoints, 10 > &tris) const
Set triPoints into tri list.
Definition: faceAreaIntersectI.H:29
tp
volScalarField tp(IOobject("tp", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedScalar("tp", dimensionSet(0, 0, 1, 0, 0, 0, 0), 0))
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::faceAreaIntersect::planeIntersection
point planeIntersection(const FixedList< scalar, 3 > &d, const triPoints &t, const label negI, const label posI) const
Return point of intersection between plane and triangle edge.
Definition: faceAreaIntersectI.H:94
Foam::triPoints
Triangle storage. Null constructable (unfortunately triangle<point, point> is not)
Definition: triPoints.H:49
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::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::faceAreaIntersect::getTriPoints
triPoints getTriPoints(const pointField &points, const face &f, const bool reverse) const
Get triPoints from face.
Definition: faceAreaIntersectI.H:45
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::faceAreaIntersect::triArea
scalar triArea(const triPoints &t) const
Return triangle area.
Definition: faceAreaIntersectI.H:105
f
labelList f(nPoints)
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::faceAreaIntersect::tolerance
static scalar & tolerance()
Fraction of local length scale to use as intersection tolerance.
Definition: faceAreaIntersectI.H:113
Foam::FixedList
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:53
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::reverse
void reverse(UList< T > &, const label n)
Definition: UListI.H:322