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 | Copyright (C) 2011-2016 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 "circleSet.H"
27 #include "sampledSet.H"
28 #include "meshSearch.H"
29 #include "DynamicList.H"
30 #include "polyMesh.H"
32 #include "word.H"
33 #include "mathematicalConstants.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(circleSet, 0);
40  addToRunTimeSelectionTable(sampledSet, circleSet, word);
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 void Foam::circleSet::calcSamples
47 (
48  DynamicList<point>& samplingPts,
49  DynamicList<label>& samplingCells,
50  DynamicList<label>& samplingFaces,
51  DynamicList<label>& samplingSegments,
52  DynamicList<scalar>& samplingCurveDist
53 ) const
54 {
55  static const string funcName =
56  (
57  "void circleSet::calcSamples"
58  "("
59  "DynamicList<point>&, "
60  "DynamicList<label>&, "
61  "DynamicList<label>&, "
62  "DynamicList<label>&, "
63  "DynamicList<scalar>&"
64  ") const"
65  );
66 
67  // set start point
68  label celli = searchEngine().findCell(startPoint_);
69  if (celli != -1)
70  {
71  samplingPts.append(startPoint_);
72  samplingCells.append(celli);
73  samplingFaces.append(-1);
74  samplingSegments.append(0);
75  samplingCurveDist.append(0.0);
76  }
77  else
78  {
80  << "Unable to find cell at point id " << 0
81  << " at location " << startPoint_ << endl;
82  }
83 
84  // add remaining points
85  const scalar alpha = constant::mathematical::pi/180.0*dTheta_;
86  const scalar sinAlpha = sin(alpha);
87  const scalar cosAlpha = cos(alpha);
88 
89  // first axis
90  vector axis1 = startPoint_ - origin_;
91  const scalar radius = mag(axis1);
92 
93  if (mag(axis1 & circleAxis_) > SMALL)
94  {
96  << "Vector defined by (startPoint - origin) not orthogonal to "
97  << "circleAxis:" << nl
98  << " startPoint - origin = " << axis1 << nl
99  << " circleAxis = " << circleAxis_ << nl
100  << endl;
101  }
102 
103  axis1 /= mag(axis1);
104 
105  scalar theta = dTheta_;
106  label nPoint = 1;
107  while (theta < 360)
108  {
109  axis1 = axis1*cosAlpha + (axis1^circleAxis_)*sinAlpha;
110  axis1 /= mag(axis1);
111  point pt = origin_ + radius*axis1;
112 
113  label celli = searchEngine().findCell(pt);
114 
115  if (celli != -1)
116  {
117  samplingPts.append(pt);
118  samplingCells.append(celli);
119  samplingFaces.append(-1);
120  samplingSegments.append(nPoint);
121  samplingCurveDist.append
122  (
123  radius*constant::mathematical::pi/180.0*theta
124  );
125 
126  nPoint++;
127  }
128  else
129  {
131  << "Unable to find cell at point id " << nPoint
132  << " at location " << pt << endl;
133  }
134  theta += dTheta_;
135  }
136 }
137 
138 
139 void Foam::circleSet::genSamples()
140 {
141  // Storage for sample points
142  DynamicList<point> samplingPts;
143  DynamicList<label> samplingCells;
144  DynamicList<label> samplingFaces;
145  DynamicList<label> samplingSegments;
146  DynamicList<scalar> samplingCurveDist;
147 
148  calcSamples
149  (
150  samplingPts,
151  samplingCells,
152  samplingFaces,
153  samplingSegments,
154  samplingCurveDist
155  );
156 
157  samplingPts.shrink();
158  samplingCells.shrink();
159  samplingFaces.shrink();
160  samplingSegments.shrink();
161  samplingCurveDist.shrink();
162 
163  setSamples
164  (
165  samplingPts,
166  samplingCells,
167  samplingFaces,
168  samplingSegments,
169  samplingCurveDist
170  );
171 }
172 
173 
174 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
175 
177 (
178  const word& name,
179  const polyMesh& mesh,
180  const meshSearch& searchEngine,
181  const word& axis,
182  const point& origin,
183  const vector& circleAxis,
184  const point& startPoint,
185  const scalar dTheta
186 )
187 :
188  sampledSet(name, mesh, searchEngine, axis),
189  origin_(origin),
190  circleAxis_(circleAxis),
191  startPoint_(startPoint),
192  dTheta_(dTheta)
193 {
194  genSamples();
195 
196  if (debug)
197  {
198  write(Info);
199  }
200 }
201 
202 
204 (
205  const word& name,
206  const polyMesh& mesh,
207  const meshSearch& searchEngine,
208  const dictionary& dict
209 )
210 :
211  sampledSet(name, mesh, searchEngine, dict),
212  origin_(dict.lookup("origin")),
213  circleAxis_(dict.lookup("circleAxis")),
214  startPoint_(dict.lookup("startPoint")),
215  dTheta_(readScalar(dict.lookup("dTheta")))
216 {
217  // normalise circleAxis
218  circleAxis_ /= mag(circleAxis_);
219 
220  genSamples();
221 
222  if (debug)
223  {
224  write(Info);
225  }
226 }
227 
228 
229 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
230 
232 {}
233 
234 
235 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
236 
238 {
239  return startPoint_;
240 }
241 
242 
243 // ************************************************************************* //
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
Definition: meshSearch.H:57
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 list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
label findCell(const point &location, const label seedCelli=-1, const bool useTreeSearch=true) const
Find cell containing location.
Definition: meshSearch.C:780
virtual point getRefPoint(const List< point > &) const
Get reference point.
Definition: circleSet.C:237
Macros for easy insertion into run-time selection tables.
dimensionedScalar cos(const dimensionedScalar &ds)
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:64
A class for handling words, derived from string.
Definition: word.H:59
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
dimensionedScalar sin(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:262
defineTypeNameAndDebug(combustionModel, 0)
vector point
Point is a vector.
Definition: point.H:41
circleSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const word &axis, const point &origin, const vector &circleAxis, const point &startPoint, const scalar dTheta)
Construct from components.
Definition: circleSet.C:177
#define WarningInFunction
Report a warning using Foam::Warning.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
virtual ~circleSet()
Definition: circleSet.C:231
runTime write()
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451