boxUniform.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 "boxUniform.H"
27 #include "sampledSet.H"
28 #include "meshSearch.H"
29 #include "DynamicList.H"
30 #include "polyMesh.H"
32 #include "word.H"
33 #include "transform.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace sampledSets
40 {
41  defineTypeNameAndDebug(boxUniform, 0);
42  addToRunTimeSelectionTable(sampledSet, boxUniform, word);
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::sampledSets::boxUniform::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  for (label k = 0; k < nPoints_.z(); ++ k)
59  {
60  for (label j = 0; j < nPoints_.y(); ++ j)
61  {
62  for (label i = 0; i < nPoints_.x(); ++ i)
63  {
64  const vector t =
65  cmptDivide(vector(i, j, k), vector(nPoints_) - vector::one);
66 
67  const point pt =
68  cmptMultiply(vector::one - t, box_.min())
69  + cmptMultiply(t, box_.max());
70 
71  const label celli = searchEngine().findCell(pt);
72 
73  if (celli != -1)
74  {
75  samplingPts.append(pt);
76  samplingCells.append(celli);
77  samplingFaces.append(-1);
78  samplingSegments.append(0);
79  samplingCurveDist.append(scalar(i));
80  }
81  }
82  }
83  }
84 }
85 
86 
87 void Foam::sampledSets::boxUniform::genSamples()
88 {
89  // Storage for sample points
90  DynamicList<point> samplingPts;
91  DynamicList<label> samplingCells;
92  DynamicList<label> samplingFaces;
93  DynamicList<label> samplingSegments;
94  DynamicList<scalar> samplingCurveDist;
95 
96  calcSamples
97  (
98  samplingPts,
99  samplingCells,
100  samplingFaces,
101  samplingSegments,
102  samplingCurveDist
103  );
104 
105  samplingPts.shrink();
106  samplingCells.shrink();
107  samplingFaces.shrink();
108  samplingSegments.shrink();
109  samplingCurveDist.shrink();
110 
111  setSamples
112  (
113  samplingPts,
114  samplingCells,
115  samplingFaces,
116  samplingSegments,
117  samplingCurveDist
118  );
119 }
120 
121 
122 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
123 
125 (
126  const word& name,
127  const polyMesh& mesh,
128  const meshSearch& searchEngine,
129  const dictionary& dict
130 )
131 :
132  sampledSet(name, mesh, searchEngine, dict),
133  box_(dict.lookup("box")),
134  nPoints_(dict.lookup("nPoints"))
135 {
136  genSamples();
137 
138  if (debug)
139  {
140  write(Info);
141  }
142 }
143 
144 
145 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
146 
148 {}
149 
150 
151 // ************************************************************************* //
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:156
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
label k
Boltzmann constant.
Macros for easy insertion into run-time selection tables.
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
static const Form min
Definition: VectorSpace.H:116
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
3D tensor transformation operations.
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
addToRunTimeSelectionTable(sampledSet, arcUniform, word)
boxUniform(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
Definition: boxUniform.C:125
vector point
Point is a vector.
Definition: point.H:41
virtual ~boxUniform()
Destructor.
Definition: boxUniform.C:147
virtual bool write()
Sample and write.
Definition: sampledSets.C:228
messageStream Info
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:844