fieldSmoother.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) 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 "fieldSmoother.H"
27 #include "pointFields.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  defineTypeNameAndDebug(fieldSmoother, 0);
34 }
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 :
40  mesh_(mesh)
41 {}
42 
43 
44 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
45 
47 {}
48 
49 
50 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
51 
53 (
54  const label nIter,
55  const PackedBoolList& isMeshMasterPoint,
56  const PackedBoolList& isMeshMasterEdge,
57  const labelList& fixedPoints,
58  pointVectorField& normals
59 ) const
60 {
61  // Get smoothly varying internal normals field.
62  Info<< typeName
63  << " : Smoothing normals in interior ..." << endl;
64 
65  const edgeList& edges = mesh_.edges();
66 
67  // Points that do not change.
68  PackedBoolList isFixedPoint(mesh_.nPoints());
69 
70  // Internal points that are fixed
71  forAll(fixedPoints, i)
72  {
73  label meshPointI = fixedPoints[i];
74  isFixedPoint.set(meshPointI, 1);
75  }
76 
77  // Make sure that points that are coupled to meshPoints but not on a patch
78  // are fixed as well
79  syncTools::syncPointList(mesh_, isFixedPoint, maxEqOp<unsigned int>(), 0);
80 
81 
82  // Correspondence between local edges/points and mesh edges/points
83  const labelList meshPoints(identity(mesh_.nPoints()));
84 
85  // Calculate inverse sum of weights
86 
87  scalarField edgeWeights(mesh_.nEdges());
88  scalarField invSumWeight(meshPoints.size());
90  (
91  mesh_,
92  isMeshMasterEdge,
93  meshPoints,
94  edges,
95  edgeWeights,
96  invSumWeight
97  );
98 
100  for (label iter = 0; iter < nIter; iter++)
101  {
103  (
104  mesh_,
105  isMeshMasterEdge,
106  meshPoints,
107  edges,
108  edgeWeights,
109  normals,
110  average
111  );
112  average *= invSumWeight;
113 
114  // Do residual calculation every so often.
115  if ((iter % 10) == 0)
116  {
117  scalar resid = meshRefinement::gAverage
118  (
119  isMeshMasterPoint,
120  mag(normals-average)()
121  );
122  Info<< " Iteration " << iter << " residual " << resid << endl;
123  }
124 
125  // Transfer to normals vector field
126  forAll(average, pointI)
127  {
128  if (isFixedPoint.get(pointI) == 0)
129  {
130  //full smoothing neighbours + point value
131  average[pointI] = 0.5*(normals[pointI]+average[pointI]);
132  normals[pointI] = average[pointI];
133  normals[pointI] /= mag(normals[pointI]) + VSMALL;
134  }
135  }
136  }
137 }
138 
139 
141 (
142  const label nIter,
143  const PackedBoolList& isPatchMasterPoint,
144  const PackedBoolList& isPatchMasterEdge,
145  const indirectPrimitivePatch& adaptPatch,
146  pointField& normals
147 ) const
148 {
149  const edgeList& edges = adaptPatch.edges();
150  const labelList& meshPoints = adaptPatch.meshPoints();
151 
152  // Get smoothly varying internal normals field.
153  Info<< typeName << " : Smoothing normals ..." << endl;
154 
155  scalarField edgeWeights(edges.size());
156  scalarField invSumWeight(meshPoints.size());
158  (
159  mesh_,
160  isPatchMasterEdge,
161  meshPoints,
162  edges,
163  edgeWeights,
164  invSumWeight
165  );
166 
167 
169  for (label iter = 0; iter < nIter; iter++)
170  {
172  (
173  mesh_,
174  isPatchMasterEdge,
175  meshPoints,
176  edges,
177  edgeWeights,
178  normals,
179  average
180  );
181  average *= invSumWeight;
182 
183  // Do residual calculation every so often.
184  if ((iter % 10) == 0)
185  {
186  scalar resid = meshRefinement::gAverage
187  (
188  isPatchMasterPoint,
189  mag(normals-average)()
190  );
191  Info<< " Iteration " << iter << " residual " << resid << endl;
192  }
193 
194  // Transfer to normals vector field
195  forAll(average, pointI)
196  {
197  // full smoothing neighbours + point value
198  average[pointI] = 0.5*(normals[pointI]+average[pointI]);
199  normals[pointI] = average[pointI];
200  normals[pointI] /= mag(normals[pointI]) + VSMALL;
201  }
202  }
203 }
204 
205 
207 (
208  const label nIter,
209  const PackedBoolList& isMeshMasterPoint,
210  const PackedBoolList& isMeshMasterEdge,
211  const PackedBoolList& isToBeSmoothed,
212  vectorField& displacement
213 ) const
214 {
215  const edgeList& edges = mesh_.edges();
216 
217  // Correspondence between local edges/points and mesh edges/points
218  const labelList meshPoints(identity(mesh_.nPoints()));
219 
220  // Calculate inverse sum of weights
221  scalarField edgeWeights(mesh_.nEdges());
222  scalarField invSumWeight(meshPoints.size());
224  (
225  mesh_,
226  isMeshMasterEdge,
227  meshPoints,
228  edges,
229  edgeWeights,
230  invSumWeight
231  );
232 
233  // Get smoothly varying patch field.
234  Info<< typeName << " : Smoothing displacement ..." << endl;
235 
236  const scalar lambda = 0.33;
237  const scalar mu = -0.34;
238 
240 
241  for (label iter = 0; iter < nIter; iter++)
242  {
244  (
245  mesh_,
246  isMeshMasterEdge,
247  meshPoints,
248  edges,
249  edgeWeights,
250  displacement,
251  average
252  );
253  average *= invSumWeight;
254 
255  forAll(displacement, i)
256  {
257  if (isToBeSmoothed[i])
258  {
259  displacement[i] = (1-lambda)*displacement[i]+lambda*average[i];
260  }
261  }
262 
264  (
265  mesh_,
266  isMeshMasterEdge,
267  meshPoints,
268  edges,
269  edgeWeights,
270  displacement,
271  average
272  );
273  average *= invSumWeight;
274 
275 
276  forAll(displacement, i)
277  {
278  if (isToBeSmoothed[i])
279  {
280  displacement[i] = (1-mu)*displacement[i]+mu*average[i];
281  }
282  }
283 
284 
285  // Do residual calculation every so often.
286  if ((iter % 10) == 0)
287  {
288  scalar resid = meshRefinement::gAverage
289  (
290  isMeshMasterPoint,
291  mag(displacement-average)()
292  );
293  Info<< " Iteration " << iter << " residual " << resid << endl;
294  }
295  }
296 }
297 
298 
299 // ************************************************************************* //
Foam::PackedBoolList
A bit-packed bool list.
Definition: PackedBoolList.H:63
Foam::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatchTemplate.C:212
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFields.H:13
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::PackedBoolList::set
void set(const PackedList< 1 > &)
Set specified bits.
Definition: PackedBoolList.C:169
Foam::maxEqOp
Definition: ops.H:77
Foam::fieldSmoother::~fieldSmoother
virtual ~fieldSmoother()
Destructor.
Definition: fieldSmoother.C:46
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::meshRefinement::calculateEdgeWeights
static void calculateEdgeWeights(const polyMesh &mesh, const PackedBoolList &isMasterEdge, const labelList &meshPoints, const edgeList &edges, scalarField &edgeWeights, scalarField &invSumWeight)
Helper: calculate edge weights (1/length)
Definition: meshRefinement.C:1835
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::fieldSmoother::smoothNormals
void smoothNormals(const label nIter, const PackedBoolList &isMeshMasterPoint, const PackedBoolList &isMeshMasterEdge, const labelList &fixedPoints, pointVectorField &normals) const
Smooth interior normals.
Definition: fieldSmoother.C:53
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::fieldSmoother::smoothLambdaMuDisplacement
void smoothLambdaMuDisplacement(const label nIter, const PackedBoolList &isMeshMasterPoint, const PackedBoolList &isMeshMasterEdge, const PackedBoolList &isSmoothable, vectorField &displacement) const
Smooth and then un-smooth a displacement.
Definition: fieldSmoother.C:207
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::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Info
messageStream Info
Foam::identity
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
fieldSmoother.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::meshRefinement::weightedSum
static void weightedSum(const polyMesh &mesh, const PackedBoolList &isMasterEdge, const labelList &meshPoints, const edgeList &edges, const scalarField &edgeWeights, const Field< Type > &data, Field< Type > &sum)
Helper: weighted sum (over all subset of mesh points) by.
Definition: meshRefinementTemplates.C:280
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::fieldSmoother::smoothPatchNormals
void smoothPatchNormals(const label nIter, const PackedBoolList &isPatchMasterPoint, const PackedBoolList &isPatchMasterEdge, const indirectPrimitivePatch &adaptPatch, pointField &normals) const
Smooth patch normals.
Definition: fieldSmoother.C:141
Foam::PackedList::get
unsigned int get(const label) const
Get value at index I.
Definition: PackedListI.H:964
Foam::meshRefinement::gAverage
static T gAverage(const PackedBoolList &isMasterElem, const UList< T > &values)
Helper: calculate average.
Definition: meshRefinementTemplates.C:59
Foam::fieldSmoother::fieldSmoother
fieldSmoother(const fieldSmoother &)
Disallow default bitwise copy construct.
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatchTemplate.C:392
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
pointFields.H
lambda
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
Foam::syncTools::syncPointList
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
Definition: syncToolsTemplates.C:984
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:335
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88