bodyCentredCubic.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2015 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "bodyCentredCubic.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 defineTypeNameAndDebug(bodyCentredCubic, 0);
36 addToRunTimeSelectionTable(initialPointsMethod, bodyCentredCubic, dictionary);
37 }
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43 (
44  const dictionary& initialPointsDict,
45  const Time& runTime,
46  Random& rndGen,
47  const conformationSurfaces& geometryToConformTo,
48  const cellShapeControl& cellShapeControls,
49  const autoPtr<backgroundMeshDecomposition>& decomposition
50 )
51 :
52  initialPointsMethod
53  (
54  typeName,
55  initialPointsDict,
56  runTime,
57  rndGen,
58  geometryToConformTo,
59  cellShapeControls,
60  decomposition
61  ),
62  initialCellSize_(detailsDict().get<scalar>("initialCellSize")),
63  randomiseInitialGrid_(detailsDict().get<Switch>("randomiseInitialGrid")),
64  randomPerturbationCoeff_
65  (
66  detailsDict().get<scalar>("randomPerturbationCoeff")
67  )
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
74 {
75  boundBox bb;
76 
77  // Pick up the bounds of this processor, or the whole geometry, depending
78  // on whether this is a parallel run.
79  if (Pstream::parRun())
80  {
81  bb = decomposition().procBounds();
82  }
83  else
84  {
85  bb = geometryToConformTo().globalBounds();
86  }
87 
88  scalar x0 = bb.min().x();
89  scalar xR = bb.max().x() - x0;
90  label ni = label(xR/initialCellSize_);
91 
92  scalar y0 = bb.min().y();
93  scalar yR = bb.max().y() - y0;
94  label nj = label(yR/initialCellSize_);
95 
96  scalar z0 = bb.min().z();
97  scalar zR = bb.max().z() - z0;
98  label nk = label(zR/initialCellSize_);
99 
100  vector delta(xR/ni, yR/nj, zR/nk);
101 
102  delta *= pow((1.0/2.0),-(1.0/3.0));
103 
104  scalar pert = randomPerturbationCoeff_*cmptMin(delta);
105 
106  DynamicList<Vb::Point> initialPoints(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(2*nk);
119 
120  for (label k = 0; k < nk; k++)
121  {
122  point pA
123  (
124  x0 + i*delta.x(),
125  y0 + j*delta.y(),
126  z0 + k*delta.z()
127  );
128 
129  point pB = pA + 0.5*delta;
130 
131  if (randomiseInitialGrid_)
132  {
133  pA.x() += pert*(rndGen().sample01<scalar>() - 0.5);
134  pA.y() += pert*(rndGen().sample01<scalar>() - 0.5);
135  pA.z() += pert*(rndGen().sample01<scalar>() - 0.5);
136  }
137 
138  if (Pstream::parRun())
139  {
140  if (decomposition().positionOnThisProcessor(pA))
141  {
142  // Add this point in parallel only if this position is
143  // on this processor.
144  points[pI++] = pA;
145  }
146  }
147  else
148  {
149  points[pI++] = pA;
150  }
151 
152  if (randomiseInitialGrid_)
153  {
154  pB.x() += pert*(rndGen().sample01<scalar>() - 0.5);
155  pB.y() += pert*(rndGen().sample01<scalar>() - 0.5);
156  pB.z() += pert*(rndGen().sample01<scalar>() - 0.5);
157  }
158 
159  if (Pstream::parRun())
160  {
161  if (decomposition().positionOnThisProcessor(pB))
162  {
163  // Add this point in parallel only if this position is
164  // on this processor.
165  points[pI++] = pB;
166  }
167  }
168  else
169  {
170  points[pI++] = pB;
171  }
172  }
173 
174  points.setSize(pI);
175 
176  Field<bool> insidePoints =
177  geometryToConformTo().wellInside
178  (
179  points,
180  minimumSurfaceDistanceCoeffSqr_
181  *sqr(cellShapeControls().cellSize(points))
182  );
183 
184  forAll(insidePoints, i)
185  {
186  if (insidePoints[i])
187  {
188  const point& p(points[i]);
189 
190  initialPoints.append(Vb::Point(p.x(), p.y(), p.z()));
191  }
192  }
193  }
194  }
195 
196  return initialPoints.shrink();
197 }
198 
199 
200 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::bodyCentredCubic::initialPoints
virtual List< Vb::Point > initialPoints() const
Return the initial points for the conformalVoronoiMesh.
p
volScalarField & p
Definition: createFieldRefs.H:8
CGAL::indexedVertex::Point
Vb::Point Point
Definition: indexedVertex.H:135
Foam::Random::sample01
Type sample01()
Return a sample whose components lie in the range [0,1].
Definition: RandomTemplates.C:36
Foam::bodyCentredCubic::bodyCentredCubic
bodyCentredCubic(const dictionary &initialPointsDict, const Time &runTime, Random &rndGen, const conformationSurfaces &geometryToConformTo, const cellShapeControl &cellShapeControls, const autoPtr< backgroundMeshDecomposition > &decomposition)
Construct from components.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::cmptMin
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:302
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::y0
dimensionedScalar y0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:281
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
bodyCentredCubic.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
points
const pointField & points
Definition: gmvOutputHeader.H:1
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
insidePoints
insidePoints((1 2 3))
rndGen
Random rndGen
Definition: createFields.H:23
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)