circleSet.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 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "circleSet.H"
30#include "sampledSet.H"
31#include "meshSearch.H"
32#include "DynamicList.H"
33#include "polyMesh.H"
35#include "word.H"
36#include "unitConversion.H"
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
44}
45
46
47// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48
49void Foam::circleSet::calcSamples
50(
51 DynamicList<point>& samplingPts,
52 DynamicList<label>& samplingCells,
53 DynamicList<label>& samplingFaces,
54 DynamicList<label>& samplingSegments,
55 DynamicList<scalar>& samplingCurveDist
56) const
57{
58 // Set start point
59 label celli = searchEngine().findCell(startPoint_);
60 if (celli != -1)
61 {
62 samplingPts.append(startPoint_);
63 samplingCells.append(celli);
64 samplingFaces.append(-1);
65 samplingSegments.append(0);
66 samplingCurveDist.append(0.0);
67 }
68 else
69 {
71 << "Unable to find cell at point id " << 0
72 << " at location " << startPoint_ << endl;
73 }
74
75 // Add remaining points
76 const scalar alpha = degToRad(dTheta_);
77 const scalar sinAlpha = sin(alpha);
78 const scalar cosAlpha = cos(alpha);
79
80 // First axis
81 vector axis1 = startPoint_ - origin_;
82 const scalar radius = mag(axis1);
83
84 if (mag(axis1 & circleAxis_) > SMALL)
85 {
87 << "Vector defined by (startPoint - origin) not orthogonal to "
88 << "circleAxis:" << nl
89 << " startPoint - origin = " << axis1 << nl
90 << " circleAxis = " << circleAxis_ << nl
91 << endl;
92 }
93
94 axis1 /= mag(axis1);
95
96 scalar theta = dTheta_;
97 label nPoint = 1;
98 while (theta < 360)
99 {
100 axis1 = normalised(axis1*cosAlpha + (axis1^circleAxis_)*sinAlpha);
101 point pt = origin_ + radius*axis1;
102
103 label celli = searchEngine().findCell(pt);
104
105 if (celli != -1)
106 {
107 samplingPts.append(pt);
108 samplingCells.append(celli);
109 samplingFaces.append(-1);
110 samplingSegments.append(nPoint);
111 samplingCurveDist.append(radius*degToRad(theta));
112
113 ++nPoint;
114 }
115 else
116 {
118 << "Unable to find cell at point id " << nPoint
119 << " at location " << pt << endl;
120 }
121 theta += dTheta_;
122 }
123}
124
125
126void Foam::circleSet::genSamples()
127{
128 // Storage for sample points
129 DynamicList<point> samplingPts;
130 DynamicList<label> samplingCells;
131 DynamicList<label> samplingFaces;
132 DynamicList<label> samplingSegments;
133 DynamicList<scalar> samplingCurveDist;
134
135 calcSamples
136 (
137 samplingPts,
138 samplingCells,
139 samplingFaces,
140 samplingSegments,
141 samplingCurveDist
142 );
143
144 samplingPts.shrink();
145 samplingCells.shrink();
146 samplingFaces.shrink();
147 samplingSegments.shrink();
148 samplingCurveDist.shrink();
149
150 // Move into *this
151 setSamples
152 (
153 std::move(samplingPts),
154 std::move(samplingCells),
155 std::move(samplingFaces),
156 std::move(samplingSegments),
157 std::move(samplingCurveDist)
158 );
159
160 if (debug)
161 {
162 write(Info);
163 }
164}
165
166
167// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
168
170(
171 const word& name,
172 const polyMesh& mesh,
173 const meshSearch& searchEngine,
174 const word& axis,
175 const point& origin,
176 const vector& circleAxis,
177 const point& startPoint,
178 const scalar dTheta
179)
180:
181 sampledSet(name, mesh, searchEngine, axis),
182 origin_(origin),
183 circleAxis_(circleAxis),
184 startPoint_(startPoint),
185 dTheta_(dTheta)
186{
187 genSamples();
188}
189
190
192(
193 const word& name,
194 const polyMesh& mesh,
195 const meshSearch& searchEngine,
196 const dictionary& dict
197)
198:
199 sampledSet(name, mesh, searchEngine, dict),
200 origin_(dict.get<point>("origin")),
201 circleAxis_(normalised(dict.get<vector>("circleAxis"))),
202 startPoint_(dict.get<point>("startPoint")),
203 dTheta_(dict.get<scalar>("dTheta"))
204{
205 genSamples();
206}
207
208
209// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
210
212{
213 return startPoint_;
214}
215
216
217// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
Samples along a circular path.
Definition: circleSet.H:103
virtual point getRefPoint(const List< point > &pts) const
Get reference point.
Definition: circleSet.C:211
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:61
label findCell(const point &location, const label seedCelli=-1, const bool useTreeSearch=true) const
Find cell containing location.
Definition: meshSearch.C:780
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:86
const meshSearch & searchEngine() const noexcept
Definition: sampledSet.H:324
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
vector point
Point is a vector.
Definition: point.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
runTime write()
volScalarField & alpha
dictionary dict
Unit conversion functions.