All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
arcEdge.C
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-2019 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 #include "arcEdge.H"
27 #include "unitConversion.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace blockEdges
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::blockEdges::arcEdge::calc(const point& pA)
45 {
46  const tensor dps(pA - p0_, pA - p1_, (pA - p0_)^(pA - p1_));
47  const tensor pMs((p0_ + pA)/2, (p1_ + pA)/2, pA);
48 
49  if (det(dps) < vSmall)
50  {
52  << "Arc point " << pA << " lies on the line of edge "
53  << edge(start_, end_) << abort(FatalError);
54  }
55 
56  const point c =
57  inv(dps)
58  & vector(dps.x() & pMs.x(), dps.y() & pMs.y(), dps.z() & pMs.z());
59 
60  const vector r0 = p0_ - c, r1 = p1_ - c;
61  const scalar cosT = (r0 & r1)/(mag(r0)*mag(r1));
62  scalar t = acos(max(-1, min(cosT, 1)));
63  if (((r0 ^ r1) & dps.z()) > 0)
64  {
66  }
67 
68  centre_ = c;
69  axis_ = - normalised(dps.z());
70  theta_ = t;
71  length_ = 0;
72 }
73 
74 
75 void Foam::blockEdges::arcEdge::calc(const scalar theta, const vector& axis)
76 {
77  if (0 >= theta || theta >= 360)
78  {
80  << "Arc angle for edge " << edge(start_, end_)
81  << " must take a value between 0 and 360 degrees"
82  << abort(FatalError);
83  }
84 
85  const vector dp = p1_ - p0_;
86 
87  const vector pM = (p0_ + p1_)/2;
88  const vector rM = normalised(dp ^ axis);
89 
90  const scalar l = dp & axis;
91 
92  const vector chord = dp - l*axis;
93  const scalar magChord = mag(chord);
94 
95  centre_ = pM - l*axis/2 - rM*magChord/2/tan(degToRad(theta)/2);
96  axis_ = axis;
97  theta_ = degToRad(theta);
98  length_ = l;
99 }
100 
101 
102 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
103 
105 (
106  const dictionary& dict,
107  const label index,
108  const searchableSurfaces& geometry,
109  const pointField& points,
110  Istream& is
111 )
112 :
113  blockEdge(dict, index, points, is),
114  p0_(points_[start_]),
115  p1_(points_[end_]),
116  centre_(),
117  axis_(),
118  theta_(),
119  length_()
120 {
121  token firstToken(is);
122  is.putBack(firstToken);
123 
124  if (firstToken == token::BEGIN_LIST)
125  {
126  const point pA(is);
127  calc(pA);
128  }
129  else
130  {
131  const scalar theta = readScalar(is);
132  const vector axis(is);
133  calc(theta, normalised(axis));
134  }
135 }
136 
137 
138 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
139 
141 {
142  if (lambda < - small || lambda > 1 + small)
143  {
145  << "Parameter out of range, lambda = " << lambda
146  << abort(FatalError);
147  }
148 
149  const scalar t = theta_*lambda;
150  const vector r1 = p0_ - centre_, r2 = axis_ ^ r1;
151 
152  return centre_ + r1*cos(t) + r2*sin(t) + axis_*length_*lambda;
153 }
154 
155 
157 {
158  const vector r1 = p0_ - centre_;
159 
160  // Length of a helical segment
161  return degToRad(theta_*sqrt(magSqr(r1) + sqr(length_)));
162 }
163 
164 
165 // ************************************************************************* //
dimensionedScalar acos(const dimensionedScalar &ds)
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
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
point position(const scalar) const
Return the point position corresponding to the curve parameter.
Definition: arcEdge.C:140
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Vector< Cmpt > x() const
Definition: TensorI.H:289
A token holds items read from Istream.
Definition: token.H:72
void putBack(const token &)
Put back token.
Definition: Istream.C:30
dimensionedScalar sqrt(const dimensionedScalar &ds)
Vector< Cmpt > y() const
Definition: TensorI.H:296
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
dimensionedScalar det(const dimensionedSphericalTensor &dt)
const dimensionedScalar & c
Speed of light in a vacuum.
Macros for easy insertion into run-time selection tables.
An arcEdge between two points on a circle. The arc is defined either by a third point that the arc pa...
Definition: arcEdge.H:57
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
Vector< Cmpt > z() const
Definition: TensorI.H:303
arcEdge(const dictionary &dict, const label index, const searchableSurfaces &geometry, const pointField &points, Istream &)
Construct from Istream setting pointsList.
Definition: arcEdge.C:105
Form normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:413
dimensionedScalar cos(const dimensionedScalar &ds)
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Container for searchableSurfaces.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:75
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensionedScalar sin(const dimensionedScalar &ds)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
defineTypeNameAndDebug(arcEdge, 0)
addToRunTimeSelectionTable(blockEdge, arcEdge, Istream)
Define a curved edge that is parameterized for 0<lambda<1 between the start and end point...
Definition: blockEdge.H:56
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar length() const
Return the length of the curve.
Definition: arcEdge.C:156
dimensionedScalar tan(const dimensionedScalar &ds)
Namespace for OpenFOAM.