arcEdge.C
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-2016 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 #include "arcEdge.H"
29 #include "unitConversion.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace blockEdges
37 {
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 Foam::coordSystem::cylindrical Foam::blockEdges::arcEdge::calcAngle()
47 {
48  const vector a = p2_ - p1_;
49  const vector b = p3_ - p1_;
50 
51  // Find centre of arcEdge
52  const scalar asqr = a & a;
53  const scalar bsqr = b & b;
54  const scalar adotb = a & b;
55 
56  const scalar denom = asqr*bsqr - adotb*adotb;
57 
58  if (mag(denom) < VSMALL)
59  {
61  << denom
62  << abort(FatalError);
63  }
64 
65  const scalar fact = 0.5*(bsqr - adotb)/denom;
66 
67  point centre = 0.5*a + fact*((a ^ b) ^ a);
68 
69  centre += p1_;
70 
71  // Find position vectors w.r.t. the arcEdge centre
72  const vector r1(p1_ - centre);
73  const vector r2(p2_ - centre);
74  const vector r3(p3_ - centre);
75 
76  // Find angle (in degrees)
77  angle_ = radToDeg(acos((r3 & r1)/(mag(r3) * mag(r1))));
78 
79  // Check if the vectors define an exterior or an interior arcEdge
80  if (((r1 ^ r2) & (r1 ^ r3)) < 0.0)
81  {
82  angle_ = 360.0 - angle_;
83  }
84 
85  vector arcAxis;
86 
87  if (angle_ <= 180.0)
88  {
89  arcAxis = r1 ^ r3;
90 
91  if (mag(arcAxis)/(mag(r1)*mag(r3)) < 0.001)
92  {
93  arcAxis = r1 ^ r2;
94  }
95  }
96  else
97  {
98  arcAxis = r3 ^ r1;
99  }
100 
101  radius_ = mag(r3);
102 
103  // The corresponding local cylindrical coordinate system (radians)
104  return coordSystem::cylindrical("arc", centre, arcAxis, r1);
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
109 
110 Foam::blockEdges::arcEdge::arcEdge
111 (
112  const pointField& points,
113  const label start,
114  const label end,
115  const point& pMid
116 )
117 :
119  p1_(points_[start_]),
120  p2_(pMid),
121  p3_(points_[end_]),
122  cs_(calcAngle())
123 {}
124 
125 
126 Foam::blockEdges::arcEdge::arcEdge
127 (
128  const dictionary& dict,
129  const label index,
130  const searchableSurfaces& geometry,
131  const pointField& points,
132  Istream& is
133 )
134 :
135  blockEdge(dict, index, points, is),
136  p1_(points_[start_]),
137  p2_(is),
138  p3_(points_[end_]),
139  cs_(calcAngle())
140 {}
141 
142 
143 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
144 
146 {
147  if (lambda < -SMALL || lambda > 1 + SMALL)
148  {
150  << "Parameter out of range, lambda = " << lambda
151  << abort(FatalError);
152  }
153 
154  if (lambda < SMALL)
155  {
156  return p1_;
157  }
158  else if (lambda > 1 - SMALL)
159  {
160  return p3_;
161  }
162 
163  // The angle is degrees, the coordinate system in radians
164  return cs_.globalPosition(vector(radius_, degToRad(lambda*angle_), 0));
165 }
166 
167 
169 {
170  return degToRad(radius_*angle_);
171 }
172 
173 
174 // ************************************************************************* //
Foam::blockEdges::arcEdge::length
scalar length() const
The length of the curve.
Definition: arcEdge.C:168
unitConversion.H
Unit conversion functions.
Foam::blockEdge
Define a curved edge that is parameterized for 0<lambda<1 between the start and end point.
Definition: blockEdge.H:59
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:62
Foam::blockEdges::defineTypeNameAndDebug
defineTypeNameAndDebug(arcEdge, 0)
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< vector >
Foam::blockEdges::addToRunTimeSelectionTable
addToRunTimeSelectionTable(blockEdge, arcEdge, Istream)
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::radToDeg
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
Definition: unitConversion.H:54
Foam::coordSystem::cylindrical
A cylindrical coordinate system (r-theta-z). The coordinate system angle theta is always in radians.
Definition: cylindricalCS.H:71
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::blockEdges::arcEdge::position
point position(const scalar lambda) const
The point corresponding to the curve parameter [0-1].
Definition: arcEdge.C:145
lambda
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:115
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::Vector< scalar >
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:268
Foam::searchableSurfaces
Container for searchableSurfaces. The collection is specified as a dictionary. For example,...
Definition: searchableSurfaces.H:92
arcEdge.H
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::start
label ListType::const_reference const label start
Definition: ListOps.H:408
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::blockEdges::arcEdge
Defines the arcEdge of a circle in terms of 3 points on its circumference.
Definition: arcEdge.H:54
Foam::point
vector point
Point is a vector.
Definition: point.H:43