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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 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 "fieldSmoother.H"
29 #include "pointFields.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(fieldSmoother, 0);
36 }
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
40 Foam::fieldSmoother::fieldSmoother(const polyMesh& mesh)
41 :
42  mesh_(mesh)
43 {}
44 
45 
46 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
47 
49 {}
50 
51 
52 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
53 
55 (
56  const label nIter,
57  const bitSet& isMeshMasterPoint,
58  const bitSet& isMeshMasterEdge,
59  const labelList& fixedPoints,
60  pointVectorField& normals
61 ) const
62 {
63  // Get smoothly varying internal normals field.
64  Info<< typeName
65  << " : Smoothing normals in interior ..." << endl;
66 
67  const edgeList& edges = mesh_.edges();
68 
69  // Points that do not change.
70  bitSet isFixedPoint(mesh_.nPoints());
71 
72  // Internal points that are fixed
73  forAll(fixedPoints, i)
74  {
75  label meshPointI = fixedPoints[i];
76  isFixedPoint.set(meshPointI);
77  }
78 
79  // Make sure that points that are coupled to meshPoints but not on a patch
80  // are fixed as well
81  syncTools::syncPointList(mesh_, isFixedPoint, maxEqOp<unsigned int>(), 0);
82 
83 
84  // Correspondence between local edges/points and mesh edges/points
85  const labelList meshPoints(identity(mesh_.nPoints()));
86 
87  // Calculate inverse sum of weights
88 
89  scalarField edgeWeights(mesh_.nEdges());
90  scalarField invSumWeight(meshPoints.size());
92  (
93  mesh_,
94  isMeshMasterEdge,
95  meshPoints,
96  edges,
97  edgeWeights,
98  invSumWeight
99  );
100 
102  for (label iter = 0; iter < nIter; iter++)
103  {
105  (
106  mesh_,
107  isMeshMasterEdge,
108  meshPoints,
109  edges,
110  edgeWeights,
111  normals,
112  average
113  );
114  average *= invSumWeight;
115 
116  // Do residual calculation every so often.
117  if ((iter % 10) == 0)
118  {
119  scalar resid = meshRefinement::gAverage
120  (
121  isMeshMasterPoint,
122  mag(normals-average)()
123  );
124  Info<< " Iteration " << iter << " residual " << resid << endl;
125  }
126 
127  // Transfer to normals vector field
128  forAll(average, pointI)
129  {
130  if (!isFixedPoint.test(pointI))
131  {
132  //full smoothing neighbours + point value
133  average[pointI] = 0.5*(normals[pointI]+average[pointI]);
134  normals[pointI] = normalised(average[pointI]);
135  }
136  }
137  }
138 }
139 
140 
142 (
143  const label nIter,
144  const bitSet& isPatchMasterPoint,
145  const bitSet& isPatchMasterEdge,
146  const indirectPrimitivePatch& adaptPatch,
147  pointField& normals
148 ) const
149 {
150  const edgeList& edges = adaptPatch.edges();
151  const labelList& meshPoints = adaptPatch.meshPoints();
152 
153  // Get smoothly varying internal normals field.
154  Info<< typeName << " : Smoothing normals ..." << endl;
155 
156  scalarField edgeWeights(edges.size());
157  scalarField invSumWeight(meshPoints.size());
159  (
160  mesh_,
161  isPatchMasterEdge,
162  meshPoints,
163  edges,
164  edgeWeights,
165  invSumWeight
166  );
167 
168 
170  for (label iter = 0; iter < nIter; iter++)
171  {
173  (
174  mesh_,
175  isPatchMasterEdge,
176  meshPoints,
177  edges,
178  edgeWeights,
179  normals,
180  average
181  );
182  average *= invSumWeight;
183 
184  // Do residual calculation every so often.
185  if ((iter % 10) == 0)
186  {
187  scalar resid = meshRefinement::gAverage
188  (
189  isPatchMasterPoint,
190  mag(normals-average)()
191  );
192  Info<< " Iteration " << iter << " residual " << resid << endl;
193  }
194 
195  // Transfer to normals vector field
196  forAll(average, pointI)
197  {
198  // full smoothing neighbours + point value
199  average[pointI] = 0.5*(normals[pointI]+average[pointI]);
200  normals[pointI] = normalised(average[pointI]);
201  }
202  }
203 }
204 
205 
207 (
208  const label nIter,
209  const bitSet& isMeshMasterPoint,
210  const bitSet& isMeshMasterEdge,
211  const bitSet& 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.test(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.test(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::meshRefinement::weightedSum
static void weightedSum(const polyMesh &mesh, const bitSet &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:271
Foam::meshRefinement::gAverage
static T gAverage(const bitSet &isMasterElem, const UList< T > &values)
Helper: calculate average.
Definition: meshRefinementTemplates.C:62
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatch.C:183
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
Foam::fieldSmoother::~fieldSmoother
virtual ~fieldSmoother()
Destructor.
Definition: fieldSmoother.C:48
Foam::fieldSmoother::smoothPatchNormals
void smoothPatchNormals(const label nIter, const bitSet &isPatchMasterPoint, const bitSet &isPatchMasterEdge, const indirectPrimitivePatch &adaptPatch, pointField &normals) const
Smooth patch normals.
Definition: fieldSmoother.C:142
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::fieldSmoother::smoothLambdaMuDisplacement
void smoothLambdaMuDisplacement(const label nIter, const bitSet &isMeshMasterPoint, const bitSet &isMeshMasterEdge, const bitSet &isSmoothable, vectorField &displacement) const
Smooth and then un-smooth a displacement.
Definition: fieldSmoother.C:207
Foam::bitSet::test
bool test(const label pos) const
Test value at specified position, never auto-vivify entries.
Definition: bitSetI.H:520
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::syncTools::syncPointList
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
Definition: syncToolsTemplates.C:721
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::meshRefinement::calculateEdgeWeights
static void calculateEdgeWeights(const polyMesh &mesh, const bitSet &isMasterEdge, const labelList &meshPoints, const edgeList &edges, scalarField &edgeWeights, scalarField &invSumWeight)
Helper: calculate edge weights (1/length)
Definition: meshRefinement.C:2081
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
fieldSmoother.H
lambda
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::maxEqOp
Definition: ops.H:80
Foam::fieldSmoother::smoothNormals
void smoothNormals(const label nIter, const bitSet &isMeshMasterPoint, const bitSet &isMeshMasterEdge, const labelList &fixedPoints, pointVectorField &normals) const
Smooth interior normals.
Definition: fieldSmoother.C:55
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:487
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::List::set
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:341
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatch.C:330
Foam::GeometricField< vector, pointPatchField, pointMesh >
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
pointFields.H
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:328
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79