All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
points.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-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 #include "points.H"
27 #include "meshSearch.H"
28 #include "DynamicList.H"
29 #include "polyMesh.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace sampledSets
37 {
39  addToRunTimeSelectionTable(sampledSet, points, word);
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 void Foam::sampledSets::points::calcSamplesUnordered
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  forAll(points_, i)
56  {
57  const point& pt = points_[i];
58  const label celli = searchEngine().findCell(pt);
59 
60  if (celli != -1)
61  {
62  samplingPts.append(pt);
63  samplingCells.append(celli);
64  samplingFaces.append(-1);
65  samplingSegments.append(samplingSegments.size());
66  samplingCurveDist.append(scalar(i));
67  }
68  }
69 }
70 
71 void Foam::sampledSets::points::calcSamplesOrdered
72 (
73  DynamicList<point>& samplingPts,
74  DynamicList<label>& samplingCells,
75  DynamicList<label>& samplingFaces,
76  DynamicList<label>& samplingSegments,
77  DynamicList<scalar>& samplingCurveDist
78 ) const
79 {
80  const label n = points_.size();
81 
82  label sampleSegmentI = 0;
83  label sampleI = 0;
84  scalar sampleDist = 0;
85 
86  while (sampleI < n)
87  {
88  const point pt = points_[sampleI];
89 
90  const label sampleCellI = searchEngine().findCell(pt);
91 
92  if (sampleCellI == -1)
93  {
94  if (++ sampleI < n)
95  {
96  sampleDist += mag(points_[sampleI] - points_[sampleI - 1]);
97  }
98  }
99  else
100  {
101  passiveParticle sampleParticle(mesh(), pt, sampleCellI);
102 
103  do
104  {
105  samplingPts.append(sampleParticle.position());
106  samplingCells.append(sampleParticle.cell());
107  samplingFaces.append(-1);
108  samplingSegments.append(sampleSegmentI);
109  samplingCurveDist.append(sampleDist);
110 
111  if (++ sampleI < n)
112  {
113  const vector s = points_[sampleI] - points_[sampleI - 1];
114  sampleDist += mag(s);
115  sampleParticle.track(s, 0);
116  }
117  }
118  while (sampleI < n && !sampleParticle.onBoundaryFace());
119 
120  ++ sampleSegmentI;
121  }
122  }
123 }
124 
125 
126 void Foam::sampledSets::points::genSamples()
127 {
128  DynamicList<point> samplingPts;
129  DynamicList<label> samplingCells;
130  DynamicList<label> samplingFaces;
131  DynamicList<label> samplingSegments;
132  DynamicList<scalar> samplingCurveDist;
133 
134  if (!ordered_)
135  {
136  calcSamplesUnordered
137  (
138  samplingPts,
139  samplingCells,
140  samplingFaces,
141  samplingSegments,
142  samplingCurveDist
143  );
144  }
145  else
146  {
147  calcSamplesOrdered
148  (
149  samplingPts,
150  samplingCells,
151  samplingFaces,
152  samplingSegments,
153  samplingCurveDist
154  );
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 dictionary& dict
182 )
183 :
184  sampledSet(name, mesh, searchEngine, dict),
185  points_(dict.lookup("points")),
186  ordered_(dict.lookup("ordered"))
187 {
188  genSamples();
189 
190  if (debug)
191  {
192  write(Info);
193  }
194 }
195 
196 
197 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
198 
200 {}
201 
202 
203 // ************************************************************************* //
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
Definition: meshSearch.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:158
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
virtual ~points()
Destructor.
Definition: points.C:199
Macros for easy insertion into run-time selection tables.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dynamicFvMesh & mesh
const pointField & points
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:61
label findCell(const point &location, const label seedCelli=-1, const bool useTreeSearch=true) const
Find cell containing location.
Definition: meshSearch.C:750
A class for handling words, derived from string.
Definition: word.H:59
addToRunTimeSelectionTable(sampledSet, arcUniform, word)
vector point
Point is a vector.
Definition: point.H:41
virtual bool write()
Sample and write.
Definition: sampledSets.C:241
points(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
Definition: points.C:177
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
defineTypeNameAndDebug(arcUniform, 0)
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812