faceCentredCubic.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  Copyright (C) 2018 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "faceCentredCubic.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 defineTypeNameAndDebug(faceCentredCubic, 0);
37 addToRunTimeSelectionTable(initialPointsMethod, faceCentredCubic, dictionary);
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/4.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(4*nk);
119 
120  for (label k = 0; k < nk; k++)
121  {
122  point p
123  (
124  x0 + i*delta.x(),
125  y0 + j*delta.y(),
126  z0 + k*delta.z()
127  );
128 
129  if (randomiseInitialGrid_)
130  {
131  p.x() += pert*(rndGen().sample01<scalar>() - 0.5);
132  p.y() += pert*(rndGen().sample01<scalar>() - 0.5);
133  p.z() += pert*(rndGen().sample01<scalar>() - 0.5);
134  }
135 
136  if (Pstream::parRun())
137  {
138  if (decomposition().positionOnThisProcessor(p))
139  {
140  // Add this point in parallel only if this position is
141  // on this processor.
142  points[pI++] = p;
143  }
144  }
145  else
146  {
147  points[pI++] = p;
148  }
149 
150  p = point
151  (
152  x0 + i*delta.x(),
153  y0 + (j + 0.5)*delta.y(),
154  z0 + (k + 0.5)*delta.z()
155  );
156 
157  if (randomiseInitialGrid_)
158  {
159  p.x() += pert*(rndGen().sample01<scalar>() - 0.5);
160  p.y() += pert*(rndGen().sample01<scalar>() - 0.5);
161  p.z() += pert*(rndGen().sample01<scalar>() - 0.5);
162  }
163 
164  if (Pstream::parRun())
165  {
166  if (decomposition().positionOnThisProcessor(p))
167  {
168  // Add this point in parallel only if this position is
169  // on this processor.
170  points[pI++] = p;
171  }
172  }
173  else
174  {
175  points[pI++] = p;
176  }
177 
178  p = point
179  (
180  x0 + (i + 0.5)*delta.x(),
181  y0 + j*delta.y(),
182  z0 + (k + 0.5)*delta.z()
183  );
184 
185  if (randomiseInitialGrid_)
186  {
187  p.x() += pert*(rndGen().sample01<scalar>() - 0.5);
188  p.y() += pert*(rndGen().sample01<scalar>() - 0.5);
189  p.z() += pert*(rndGen().sample01<scalar>() - 0.5);
190  }
191 
192  if (Pstream::parRun())
193  {
194  if (decomposition().positionOnThisProcessor(p))
195  {
196  // Add this point in parallel only if this position is
197  // on this processor.
198  points[pI++] = p;
199  }
200  }
201  else
202  {
203  points[pI++] = p;
204  }
205 
206  p = point
207  (
208  x0 + (i + 0.5)*delta.x(),
209  y0 + (j + 0.5)*delta.y(),
210  z0 + k*delta.z()
211  );
212 
213  if (randomiseInitialGrid_)
214  {
215  p.x() += pert*(rndGen().sample01<scalar>() - 0.5);
216  p.y() += pert*(rndGen().sample01<scalar>() - 0.5);
217  p.z() += pert*(rndGen().sample01<scalar>() - 0.5);
218  }
219 
220  if (Pstream::parRun())
221  {
222  if (decomposition().positionOnThisProcessor(p))
223  {
224  // Add this point in parallel only if this position is
225  // on this processor.
226  points[pI++] = p;
227  }
228  }
229  else
230  {
231  points[pI++] = p;
232  }
233  }
234 
235  points.setSize(pI);
236 
237  Field<bool> insidePoints =
238  geometryToConformTo().wellInside
239  (
240  points,
241  minimumSurfaceDistanceCoeffSqr_
242  *sqr(cellShapeControls().cellSize(points))
243  );
244 
245  forAll(insidePoints, i)
246  {
247  if (insidePoints[i])
248  {
249  const point& p(points[i]);
250 
251  initialPoints.append(Vb::Point(p.x(), p.y(), p.z()));
252  }
253  }
254  }
255  }
256 
257  return initialPoints.shrink();
258 }
259 
260 
261 // ************************************************************************* //
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
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::faceCentredCubic::faceCentredCubic
faceCentredCubic(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.
faceCentredCubic.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::faceCentredCubic::initialPoints
virtual List< Vb::Point > initialPoints() const
Return the initial points for the conformalVoronoiMesh.
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)