patchSeedSet.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) 2012-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 "patchSeedSet.H"
27 #include "polyMesh.H"
29 #include "treeBoundBox.H"
30 #include "treeDataFace.H"
31 #include "Time.H"
32 #include "meshTools.H"
33 //#include "Random.H"
34 // For 'facePoint' helper function only
35 #include "mappedPatchBase.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(patchSeedSet, 0);
42  addToRunTimeSelectionTable(sampledSet, patchSeedSet, word);
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::patchSeedSet::calcSamples
49 (
50  DynamicList<point>& samplingPts,
51  DynamicList<label>& samplingCells,
52  DynamicList<label>& samplingFaces,
53  DynamicList<label>& samplingSegments,
54  DynamicList<scalar>& samplingCurveDist
55 )
56 {
57  if (debug)
58  {
59  Info<< "patchSeedSet : sampling on patches :" << endl;
60  }
61 
62  // Construct search tree for all patch faces.
63  label sz = 0;
64  forAllConstIter(labelHashSet, patchSet_, iter)
65  {
66  const polyPatch& pp = mesh().boundaryMesh()[iter.key()];
67 
68  sz += pp.size();
69 
70  if (debug)
71  {
72  Info<< " " << pp.name() << " size " << pp.size() << endl;
73  }
74  }
75 
76  labelList patchFaces(sz);
77  sz = 0;
78  forAllConstIter(labelHashSet, patchSet_, iter)
79  {
80  const polyPatch& pp = mesh().boundaryMesh()[iter.key()];
81  forAll(pp, i)
82  {
83  patchFaces[sz++] = pp.start()+i;
84  }
85  }
86 
87 
88  label totalSize = returnReduce(sz, sumOp<label>());
89 
90 
91  // Shuffle and truncate if in random mode
92  if (maxPoints_ < totalSize)
93  {
94  // Check what fraction of maxPoints_ I need to generate locally.
95  label myMaxPoints = label(scalar(sz)/totalSize*maxPoints_);
96 
97  rndGenPtr_.reset(new Random(123456));
98  Random& rndGen = rndGenPtr_();
99 
100  labelList subset = identity(sz);
101  for (label iter = 0; iter < 4; iter++)
102  {
103  forAll(subset, i)
104  {
105  label j = rndGen.integer(0, subset.size()-1);
106  Swap(subset[i], subset[j]);
107  }
108  }
109  // Truncate
110  subset.setSize(myMaxPoints);
111 
112  // Subset patchFaces
113  patchFaces = UIndirectList<label>(patchFaces, subset)();
114 
115  if (debug)
116  {
117  Pout<< "In random mode : selected " << patchFaces.size()
118  << " faces out of " << sz << endl;
119  }
120  }
121 
122 
123  // Get points on patchFaces.
124  globalIndex globalSampleNumbers(patchFaces.size());
125 
126  samplingPts.setCapacity(patchFaces.size());
127  samplingCells.setCapacity(patchFaces.size());
128  samplingFaces.setCapacity(patchFaces.size());
129  samplingSegments.setCapacity(patchFaces.size());
130  samplingCurveDist.setCapacity(patchFaces.size());
131 
132  // For calculation of min-decomp tet base points
133  (void)mesh().tetBasePtIs();
134 
135  forAll(patchFaces, i)
136  {
137  label facei = patchFaces[i];
139  (
140  mesh(),
141  facei,
143  );
144  label celli = mesh().faceOwner()[facei];
145 
146  if (info.hit())
147  {
148  // Move the point into the cell
149  const point& cc = mesh().cellCentres()[celli];
150  samplingPts.append
151  (
152  info.hitPoint() + 1e-1*(cc-info.hitPoint())
153  );
154  }
155  else
156  {
157  samplingPts.append(info.rawPoint());
158  }
159  samplingCells.append(celli);
160  samplingFaces.append(facei);
161  samplingSegments.append(0);
162  samplingCurveDist.append(globalSampleNumbers.toGlobal(i));
163  }
164 }
165 
166 
167 void Foam::patchSeedSet::genSamples()
168 {
169  // Storage for sample points
170  DynamicList<point> samplingPts;
171  DynamicList<label> samplingCells;
172  DynamicList<label> samplingFaces;
173  DynamicList<label> samplingSegments;
174  DynamicList<scalar> samplingCurveDist;
175 
176  calcSamples
177  (
178  samplingPts,
179  samplingCells,
180  samplingFaces,
181  samplingSegments,
182  samplingCurveDist
183  );
184 
185  samplingPts.shrink();
186  samplingCells.shrink();
187  samplingFaces.shrink();
188  samplingSegments.shrink();
189  samplingCurveDist.shrink();
190 
191  setSamples
192  (
193  samplingPts,
194  samplingCells,
195  samplingFaces,
196  samplingSegments,
197  samplingCurveDist
198  );
199 }
200 
201 
202 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
203 
205 (
206  const word& name,
207  const polyMesh& mesh,
208  const meshSearch& searchEngine,
209  const dictionary& dict
210 )
211 :
212  sampledSet(name, mesh, searchEngine, dict),
213  patchSet_
214  (
215  mesh.boundaryMesh().patchSet
216  (
217  wordReList(dict.lookup("patches"))
218  )
219  ),
220  //searchDist_(readScalar(dict.lookup("maxDistance"))),
221  //offsetDist_(readScalar(dict.lookup("offsetDist"))),
222  maxPoints_(readLabel(dict.lookup("maxPoints")))
223 {
224  genSamples();
225 
226  if (debug)
227  {
228  write(Info);
229  }
230 }
231 
232 
233 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
234 
236 {}
237 
238 
239 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
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:428
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
const double e
Elementary charge.
Definition: doubleFloat.H:78
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:841
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
Macros for easy insertion into run-time selection tables.
patchSeedSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
Definition: patchSeedSet.C:205
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
dynamicFvMesh & mesh
void Swap(T &a, T &b)
Definition: Swap.H:43
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
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:184
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1049
List< label > labelList
A List of labels.
Definition: labelList.H:56
const vectorField & cellCentres() const
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
label readLabel(Istream &is)
Definition: label.H:64
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
virtual ~patchSeedSet()
Destructor.
Definition: patchSeedSet.C:235
defineTypeNameAndDebug(combustionModel, 0)
static pointIndexHit facePoint(const polyMesh &, const label facei, const polyMesh::cellDecomposition)
Get a point on the face given a face decomposition method:
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
vector point
Point is a vector.
Definition: point.H:41
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
ListType subset(const UList< T > &select, const T &value, const ListType &)
Extract elements of List when select is a certain value.
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576