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-------------------------------------------------------------------------------
10License
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
33namespace Foam
34{
36}
37
38// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39
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// ************************************************************************* //
A list of faces which address into the list of points.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:590
bool test(const label pos) const
Test value at specified position, never auto-vivify entries.
Definition: bitSetI.H:521
Utility functions for mesh motion solvers.
Definition: fieldSmoother.H:51
void smoothLambdaMuDisplacement(const label nIter, const bitSet &isMeshMasterPoint, const bitSet &isMeshMasterEdge, const bitSet &isSmoothable, vectorField &displacement) const
Smooth and then un-smooth a displacement.
virtual ~fieldSmoother()
Destructor.
Definition: fieldSmoother.C:48
void smoothNormals(const label nIter, const bitSet &isMeshMasterPoint, const bitSet &isMeshMasterEdge, const labelList &fixedPoints, pointVectorField &normals) const
Smooth interior normals.
Definition: fieldSmoother.C:55
void smoothPatchNormals(const label nIter, const bitSet &isPatchMasterPoint, const bitSet &isPatchMasterEdge, const indirectPrimitivePatch &adaptPatch, pointField &normals) const
Smooth patch normals.
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)
static T gAverage(const bitSet &isMasterElem, const UList< T > &values)
Helper: calculate average.
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.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & mu
dynamicFvMesh & mesh
Namespace for OpenFOAM.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333