All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
faceAreaIntersectI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
27 
28 inline void Foam::faceAreaIntersect::setTriPoints
29 (
30  const point& a,
31  const point& b,
32  const point& c,
33  label& count,
34  FixedList<triPoints, 10>& tris
35 ) const
36 {
37  triPoints& tp = tris[count++];
38  tp[0] = a;
39  tp[1] = b;
40  tp[2] = c;
41 }
42 
43 
44 inline Foam::faceAreaIntersect::triPoints Foam::faceAreaIntersect::getTriPoints
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 
70 inline Foam::point Foam::faceAreaIntersect::planeIntersection
71 (
72  const FixedList<scalar, 3>& d,
73  const triPoints& t,
74  const label negI,
75  const label posI
76 ) const
77 {
78  return (d[posI]*t[negI] - d[negI]*t[posI])/(-d[negI] + d[posI]);
79 }
80 
81 
82 inline Foam::scalar Foam::faceAreaIntersect::triArea(const triPoints& t) const
83 {
84  return mag(0.5*((t[1] - t[0])^(t[2] - t[0])));
85 }
86 
87 
88 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
89 
91 {
92  return tol;
93 }
94 
95 
96 // ************************************************************************* //
static scalar & tolerance()
Fraction of local length scale to use as intersection tolerance.
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
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
const dimensionedScalar & c
Speed of light in a vacuum.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
const dimensionedScalar & b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
vector point
Point is a vector.
Definition: point.H:41
FixedList< point, 3 > triPoints
dimensioned< scalar > mag(const dimensioned< Type > &)