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 | Copyright (C) 2011-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 "conformalVoronoiMesh.H"
27 
28 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 
30 template<class Triangulation>
32 {
33  if (!Pstream::parRun())
34  {
35  return false;
36  }
37 
38  Info<< nl << "Redistributing points" << endl;
39 
40  timeCheck("Before distribute");
41 
42  label iteration = 0;
43 
44  scalar previousLoadUnbalance = 0;
45 
46  while (true)
47  {
48  scalar maxLoadUnbalance = mesh.calculateLoadUnbalance();
49 
50  if
51  (
52  maxLoadUnbalance <= foamyHexMeshControls().maxLoadUnbalance()
53  || maxLoadUnbalance <= previousLoadUnbalance
54  )
55  {
56  // If this is the first iteration, return false, if it was a
57  // subsequent one, return true;
58  return iteration != 0;
59  }
60 
61  previousLoadUnbalance = maxLoadUnbalance;
62 
63  Info<< " Total number of vertices before redistribution "
64  << returnReduce(label(mesh.number_of_vertices()), sumOp<label>())
65  << endl;
66 
67  const fvMesh& bMesh = decomposition_().mesh();
68 
69  volScalarField cellWeights
70  (
71  IOobject
72  (
73  "cellWeights",
74  bMesh.time().timeName(),
75  bMesh,
78  ),
79  bMesh,
80  dimensionedScalar("weight", dimless, 1e-2),
81  zeroGradientFvPatchScalarField::typeName
82  );
83 
84  meshSearch cellSearch(bMesh, polyMesh::FACE_PLANES);
85 
86  labelList cellVertices(bMesh.nCells(), label(0));
87 
88  for
89  (
90  typename Triangulation::Finite_vertices_iterator vit
91  = mesh.finite_vertices_begin();
92  vit != mesh.finite_vertices_end();
93  ++vit
94  )
95  {
96  // Only store real vertices that are not feature vertices
97  if (vit->real() && !vit->featurePoint())
98  {
99  pointFromPoint v = topoint(vit->point());
100 
101  label cellI = cellSearch.findCell(v);
102 
103  if (cellI == -1)
104  {
105 // Pout<< "findCell conformalVoronoiMesh::distribute "
106 // << "findCell "
107 // << vit->type() << " "
108 // << vit->index() << " "
109 // << v << " "
110 // << cellI
111 // << " find nearest cellI ";
112 
113  cellI = cellSearch.findNearestCell(v);
114  }
115 
116  cellVertices[cellI]++;
117  }
118  }
119 
120  forAll(cellVertices, cI)
121  {
122  // Give a small but finite weight for empty cells. Some
123  // decomposition methods have difficulty with integer overflows in
124  // the sum of the normalised weight field.
125  cellWeights.internalField()[cI] = max
126  (
127  cellVertices[cI],
128  1e-2
129  );
130  }
131 
132  autoPtr<mapDistributePolyMesh> mapDist = decomposition_().distribute
133  (
134  cellWeights
135  );
136 
138 
139  distribute();
140 
141  timeCheck("After distribute");
142 
143  iteration++;
144  }
145 
146  return true;
147 }
148 
149 
150 // ************************************************************************* //
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::conformalVoronoiMesh::timeCheck
static void timeCheck(const Time &runTime, const string &description=string::null, const bool check=true)
Write the elapsedCpuTime and memory usage, with an optional.
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::conformalVoronoiMesh::decomposition_
autoPtr< backgroundMeshDecomposition > decomposition_
Background mesh decomposition, only available in parallel.
Definition: conformalVoronoiMesh.H:157
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::topoint
pointFromPoint topoint(const Point &P)
Definition: pointConversion.H:70
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::label
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
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::bMesh
PrimitivePatch< face, List, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points)
Definition: bMesh.H:45
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::conformalVoronoiMesh::cellShapeControl_
cellShapeControl cellShapeControl_
The cell shape control object.
Definition: conformalVoronoiMesh.H:160
Foam::cellShapeControlMesh::distribute
void distribute(const backgroundMeshDecomposition &decomposition)
Foam::cellShapeControl::shapeControlMesh
cellShapeControlMesh & shapeControlMesh()
Definition: cellShapeControlI.H:29
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::conformalVoronoiMesh::foamyHexMeshControls
const cvControls & foamyHexMeshControls() const
Return the foamyHexMeshControls object.
Definition: conformalVoronoiMeshI.H:578
Foam::Vector< scalar >
Foam::polyMesh::FACE_PLANES
@ FACE_PLANES
Definition: polyMesh.H:100
Foam::conformalVoronoiMesh::distributeBackground
bool distributeBackground(const Triangulation &mesh)
In parallel redistribute the backgroundMeshDecomposition and.
Foam::conformalVoronoiMesh::distribute
void distribute()
conformalVoronoiMesh.H