circleRandom.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-2021 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 "circleRandom.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 #include "Random.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace sampledSets
41 {
42  defineTypeNameAndDebug(circleRandom, 0);
43  addToRunTimeSelectionTable(sampledSet, circleRandom, word);
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void Foam::sampledSets::circleRandom::calcSamples
51 (
52  DynamicList<point>& samplingPositions,
53  DynamicList<label>& samplingSegments,
54  DynamicList<label>& samplingCells,
55  DynamicList<label>& samplingFaces
56 ) const
57 {
58  Random rndGen(261782);
59 
60  const vector radial1 = normalised(perpendicular(normal_));
61  const vector radial2 = normalised(normal_ ^ radial1);
62 
63  for (label i = 0; i < nPoints_; ++ i)
64  {
65  // Request all random numbers simultaneously on all processors so that
66  // the generator state stays consistent
67 
68  const scalar r = radius_*rndGen.scalar01();
69  const scalar theta = 2*constant::mathematical::pi*rndGen.scalar01();
70  const scalar c = cos(theta), s = sin(theta);
71 
72  const point pt = centre_ + r*(c*radial1 + s*radial2);
73  const label celli = searchEngine().findCell(pt);
74 
75  if (celli != -1)
76  {
77  samplingPositions.append(pt);
78  samplingSegments.append(i);
79  samplingCells.append(celli);
80  samplingFaces.append(-1);
81  }
82  }
83 }
84 
85 
86 void Foam::sampledSets::circleRandom::genSamples()
87 {
88  DynamicList<point> samplingPositions;
89  DynamicList<label> samplingSegments;
90  DynamicList<label> samplingCells;
91  DynamicList<label> samplingFaces;
92 
93  calcSamples
94  (
95  samplingPositions,
96  samplingSegments,
97  samplingCells,
98  samplingFaces
99  );
100 
101  samplingPositions.shrink();
102  samplingSegments.shrink();
103  samplingCells.shrink();
104  samplingFaces.shrink();
105 
106  setSamples
107  (
108  samplingPositions,
109  samplingSegments,
110  samplingCells,
111  samplingFaces
112  );
113 }
114 
115 
116 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
117 
119 (
120  const word& name,
121  const polyMesh& mesh,
122  const meshSearch& searchEngine,
123  const dictionary& dict
124 )
125 :
126  sampledSet(name, mesh, searchEngine, dict),
127  centre_(dict.lookup("centre")),
128  normal_(normalised(dict.lookup<vector>("normal"))),
129  radius_(dict.lookup<scalar>("radius")),
130  nPoints_(dict.lookup<label>("nPoints"))
131 {
132  genSamples();
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
137 
139 {}
140 
141 
142 // ************************************************************************* //
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
Definition: meshSearch.H:57
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Macros for easy insertion into run-time selection tables.
circleRandom(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
Definition: circleRandom.C:119
Random rndGen(label(0))
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))
Form normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:413
dimensionedScalar cos(const dimensionedScalar &ds)
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)
dimensionedScalar sin(const dimensionedScalar &ds)
Vector< Cmpt > perpendicular(const Vector< Cmpt > &v)
Definition: VectorI.H:166
vector point
Point is a vector.
Definition: point.H:41
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
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:864