uniformGrid.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-2015 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 "uniformGrid.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 defineTypeNameAndDebug(uniformGrid, 0);
37 addToRunTimeSelectionTable(initialPointsMethod, uniformGrid, dictionary);
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 (
43  const dictionary& initialPointsDict,
44  const Time& runTime,
45  Random& rndGen,
46  const conformationSurfaces& geometryToConformTo,
47  const cellShapeControl& cellShapeControls,
48  const autoPtr<backgroundMeshDecomposition>& decomposition
49 )
50 :
51  initialPointsMethod
52  (
53  typeName,
54  initialPointsDict,
55  runTime,
56  rndGen,
57  geometryToConformTo,
58  cellShapeControls,
59  decomposition
60  ),
61  initialCellSize_(readScalar(detailsDict().lookup("initialCellSize"))),
62  randomiseInitialGrid_(detailsDict().lookup("randomiseInitialGrid")),
63  randomPerturbationCoeff_
64  (
65  readScalar(detailsDict().lookup("randomPerturbationCoeff"))
66  )
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
72 List<Vb::Point> uniformGrid::initialPoints() const
73 {
74  boundBox bb;
75 
76  // Pick up the bounds of this processor, or the whole geometry, depending
77  // on whether this is a parallel run.
78  if (Pstream::parRun())
79  {
80  bb = decomposition().procBounds();
81  }
82  else
83  {
85  }
86 
87  scalar x0 = bb.min().x();
88  scalar xR = bb.max().x() - x0;
89  label ni = label(xR/initialCellSize_);
90 
91  scalar y0 = bb.min().y();
92  scalar yR = bb.max().y() - y0;
93  label nj = label(yR/initialCellSize_);
94 
95  scalar z0 = bb.min().z();
96  scalar zR = bb.max().z() - z0;
97  label nk = label(zR/initialCellSize_);
98 
99  vector delta(xR/ni, yR/nj, zR/nk);
100 
101  delta *= pow((1.0),-(1.0/3.0));
102 
103  scalar pert = randomPerturbationCoeff_*cmptMin(delta);
104 
105  // Initialise points list
106  DynamicList<Vb::Point> initialPoints(scalar(ni)*nj*nk/10);
107 
108  for (label i = 0; i < ni; i++)
109  {
110  for (label j = 0; j < nj; j++)
111  {
112  // Generating, testing and adding points one line at a time to
113  // reduce the memory requirement for cases with bounding boxes that
114  // are very large in comparison to the volume to be filled
115 
116  label pI = 0;
117 
118  pointField points(nk);
119 
120  for (label k = 0; k < nk; k++)
121  {
122  point p
123  (
124  x0 + (i + 0.5)*delta.x(),
125  y0 + (j + 0.5)*delta.y(),
126  z0 + (k + 0.5)*delta.z()
127  );
128 
129  if (randomiseInitialGrid_)
130  {
131  p.x() += pert*(rndGen().scalar01() - 0.5);
132  p.y() += pert*(rndGen().scalar01() - 0.5);
133  p.z() += pert*(rndGen().scalar01() - 0.5);
134  }
135 
136  if
137  (
139  && !decomposition().positionOnThisProcessor(p)
140  )
141  {
142  // Skip this point if, in parallel, this position is not on
143  // this processor.
144  continue;
145  }
146 
147  points[pI++] = p;
148  }
149 
150  points.setSize(pI);
151 
152  Field<bool> insidePoints =
154  (
155  points,
157  *sqr
158  (
159  cellShapeControls().cellSize(points)
160  )
161  );
162 
163  forAll(insidePoints, i)
164  {
165  if (insidePoints[i])
166  {
167  const point& p(points[i]);
168 
169  initialPoints.append(Vb::Point(p.x(), p.y(), p.z()));
170  }
171  }
172  }
173  }
174 
175  return initialPoints.shrink();
176 }
177 
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 } // End namespace Foam
182 
183 // ************************************************************************* //
scalar delta
#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 treeBoundBox & globalBounds() const
Return the global bounds.
Field< bool > wellInside(const pointField &samplePts, const scalarField &testDistSqr) const
Check if point is inside surfaces to conform to by at least.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar y0(const dimensionedScalar &ds)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
label k
Boltzmann constant.
virtual List< Vb::Point > initialPoints() const
Return the initial points for the conformalVoronoiMesh.
const cellShapeControl & cellShapeControls() const
Macros for easy insertion into run-time selection tables.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
stressControl lookup("compactNormalStress") >> compactNormalStress
const pointField & points
const backgroundMeshDecomposition & decomposition() const
cachedRandom rndGen(label(0), -1)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
const Cmpt & x() const
Definition: VectorI.H:75
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
const conformationSurfaces & geometryToConformTo() const
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:394
vector point
Point is a vector.
Definition: point.H:41
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
volScalarField & p
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
const treeBoundBox & procBounds() const
Return the boundBox of this processor.
scalar scalar01()
Scalar [0..1] (so including 0,1)
Definition: Random.C:67
scalar minimumSurfaceDistanceCoeffSqr_
Only allow the placement of initial points that are within the.
Namespace for OpenFOAM.
virtual List< Vb::Point > initialPoints() const
Return the initial points for the conformalVoronoiMesh.
uniformGrid(const dictionary &initialPointsDict, const Time &runTime, Random &rndGen, const conformationSurfaces &geometryToConformTo, const cellShapeControl &cellShapeControls, const autoPtr< backgroundMeshDecomposition > &decomposition)
Construct from components.