conformalVoronoiMeshTemplates.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) 2011-2016 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 "conformalVoronoiMesh.H"
29 
30 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
31 
32 template<class Triangulation>
33 bool Foam::conformalVoronoiMesh::distributeBackground(const Triangulation& mesh)
34 {
35  if (!Pstream::parRun())
36  {
37  return false;
38  }
39 
40  Info<< nl << "Redistributing points" << endl;
41 
42  timeCheck("Before distribute");
43 
44  label iteration = 0;
45 
46  scalar previousLoadUnbalance = 0;
47 
48  while (true)
49  {
50  scalar maxLoadUnbalance = mesh.calculateLoadUnbalance();
51 
52  if
53  (
54  maxLoadUnbalance <= foamyHexMeshControls().maxLoadUnbalance()
55  || maxLoadUnbalance <= previousLoadUnbalance
56  )
57  {
58  // If this is the first iteration, return false, if it was a
59  // subsequent one, return true;
60  return iteration != 0;
61  }
62 
63  previousLoadUnbalance = maxLoadUnbalance;
64 
65  Info<< " Total number of vertices before redistribution "
66  << returnReduce(label(mesh.number_of_vertices()), sumOp<label>())
67  << endl;
68 
69  const fvMesh& bMesh = decomposition_().mesh();
70 
71  volScalarField cellWeights
72  (
73  IOobject
74  (
75  "cellWeights",
76  bMesh.time().timeName(),
77  bMesh,
80  ),
81  bMesh,
84  );
85 
86  meshSearch cellSearch(bMesh, polyMesh::FACE_PLANES);
87 
88  labelList cellVertices(bMesh.nCells(), Zero);
89 
90  for
91  (
92  typename Triangulation::Finite_vertices_iterator vit
93  = mesh.finite_vertices_begin();
94  vit != mesh.finite_vertices_end();
95  ++vit
96  )
97  {
98  // Only store real vertices that are not feature vertices
99  if (vit->real() && !vit->featurePoint())
100  {
101  pointFromPoint v = topoint(vit->point());
102 
103  label celli = cellSearch.findCell(v);
104 
105  if (celli == -1)
106  {
107 // Pout<< "findCell conformalVoronoiMesh::distribute "
108 // << "findCell "
109 // << vit->type() << " "
110 // << vit->index() << " "
111 // << v << " "
112 // << celli
113 // << " find nearest celli ";
114 
115  celli = cellSearch.findNearestCell(v);
116  }
117 
118  cellVertices[celli]++;
119  }
120  }
121 
122  scalarField& cwi = cellWeights.primitiveFieldRef();
123 
124  forAll(cellVertices, cI)
125  {
126  // Give a small but finite weight for empty cells. Some
127  // decomposition methods have difficulty with integer overflows in
128  // the sum of the normalised weight field.
129  cwi[cI] = max(cellVertices[cI], 1e-2);
130  }
131 
132  autoPtr<mapDistributePolyMesh> mapDist = decomposition_().distribute
133  (
134  cellWeights
135  );
136 
137  cellShapeControl_.shapeControlMesh().distribute(decomposition_());
138 
139  distribute();
140 
141  timeCheck("After distribute");
142 
143  iteration++;
144  }
145 
146  return true;
147 }
148 
149 
150 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:218
pointFromPoint topoint(const Point &P)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1004
static void timeCheck(const Time &runTime, const string &description=string::null, const bool check=true)
Write the elapsedCpuTime and memory usage, with an optional.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static const word null
An empty word.
Definition: word.H:84
void distribute(const backgroundMeshDecomposition &decomposition)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Nothing to be read.
messageStream Info
Information stream (stdout output on master, null elsewhere)
cellShapeControlMesh & shapeControlMesh()
const cvControls & foamyHexMeshControls() const
Return the foamyHexMeshControls object.
List< label > labelList
A List of labels.
Definition: List.H:62
PrimitivePatch< List< face >, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points) ...
Definition: bMesh.H:39
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133