fieldSmootherTemplates.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 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
29 
30 template <class Type>
32 (
33  const label nIter,
34  const bitSet& isPatchMasterPoint,
35  const bitSet& isPatchMasterEdge,
36  const indirectPrimitivePatch& adaptPatch,
37  const scalarField& fieldMinMag,
39 ) const
40 {
41  const edgeList& edges = adaptPatch.edges();
42  const labelList& meshPoints = adaptPatch.meshPoints();
43 
44  scalarField edgeWeights(edges.size());
45  scalarField invSumWeight(meshPoints.size());
47  (
48  mesh_,
49  isPatchMasterEdge,
50  meshPoints,
51  edges,
52  edgeWeights,
53  invSumWeight
54  );
55 
56  // Get smoothly varying patch field.
57  Info<< typeName << " : Smoothing field ..." << endl;
58 
59  for (label iter = 0; iter < nIter; iter++)
60  {
61  Field<Type> average(adaptPatch.nPoints());
63  (
64  mesh_,
65  isPatchMasterEdge,
66  meshPoints,
67  edges,
68  edgeWeights,
69  field,
70  average
71  );
72  average *= invSumWeight;
73 
74  // Transfer to field
75  forAll(field, pointI)
76  {
77  //full smoothing neighbours + point value
78  average[pointI] = 0.5*(field[pointI]+average[pointI]);
79 
80  // perform monotonic smoothing
81  if
82  (
83  mag(average[pointI]) < mag(field[pointI])
84  && mag(average[pointI]) >= mag(fieldMinMag[pointI])
85  )
86  {
87  field[pointI] = average[pointI];
88  }
89  }
90 
91  // Do residual calculation every so often.
92  if ((iter % 10) == 0)
93  {
94  scalar resid = meshRefinement::gAverage
95  (
96  isPatchMasterPoint,
97  mag(field-average)()
98  );
99  Info<< " Iteration " << iter << " residual " << resid << endl;
100  }
101  }
102 }
103 
104 
105 // ************************************************************************* //
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)
Definition: meshRefinementTemplates.C:264
Foam::meshRefinement::gAverage
static T gAverage(const bitSet &isMasterElem, const UList< T > &values)
Definition: meshRefinementTemplates.C:55
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
Foam::PrimitivePatch::edges
const edgeList & edges() const
Definition: PrimitivePatch.C:176
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:597
Foam::fac::average
tmp< GeometricField< Type, faPatchField, areaMesh > > average(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facAverage.C:40
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::meshRefinement::calculateEdgeWeights
static void calculateEdgeWeights(const polyMesh &mesh, const bitSet &isMasterEdge, const labelList &meshPoints, const edgeList &edges, scalarField &edgeWeights, scalarField &invSumWeight)
Definition: meshRefinement.C:2076
Foam::Field
Generic templated field type.
Definition: Field.H:59
Foam::Info
messageStream Info
field
rDeltaTY field()
Foam::PrimitivePatch::nPoints
label nPoints() const
Definition: PrimitivePatch.H:312
Foam::fieldSmoother::minSmoothField
void minSmoothField(const label nIter, const bitSet &isPatchMasterPoint, const bitSet &isPatchMasterEdge, const indirectPrimitivePatch &adaptPatch, const scalarField &fieldMin, Field< Type > &field) const
Definition: fieldSmootherTemplates.C:25
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Definition: PrimitivePatch.C:323
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:321
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:75