Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / src / turbulenceModels / LES / LESdeltas / smoothDelta / smoothDelta.C
blob28fef61d433739bd07a6d93e6c33ebb9785fb3fe
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation, either version 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "smoothDelta.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "FaceCellWave.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 defineTypeNameAndDebug(smoothDelta, 0);
38 addToRunTimeSelectionTable(LESdelta, smoothDelta, dictionary);
40 scalar smoothDelta::deltaData::maxDeltaRatio = 1.2;
43 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
45 // Fill changedFaces (with face labels) and changedFacesInfo (with delta)
46 // This is the initial set of faces from which to start the waves.
47 // Since there might be lots of places with delta jumps we can follow various
48 // strategies for this initial 'seed'.
49 // - start from single cell/face and let FaceCellWave pick up all others
50 //   from there. might be quite a few waves before everything settles.
51 // - start from all faces. Lots of initial transfers.
52 // We do something inbetween:
53 // - start from all faces where there is a jump. Since we cannot easily
54 //   determine this across coupled patches (cyclic, processor) introduce
55 //   all faces of these and let FaceCellWave sort it out.
56 void smoothDelta::setChangedFaces
58     const polyMesh& mesh,
59     const volScalarField& delta,
60     dynamicLabelList& changedFaces,
61     DynamicList<deltaData>& changedFacesInfo
64     for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
65     {
66         scalar ownDelta = delta[mesh.faceOwner()[faceI]];
68         scalar neiDelta = delta[mesh.faceNeighbour()[faceI]];
70         // Check if owner delta much larger than neighbour delta or vice versa
72         if (ownDelta > deltaData::maxDeltaRatio * neiDelta)
73         {
74             changedFaces.append(faceI);
75             changedFacesInfo.append(deltaData(ownDelta));
76         }
77         else if (neiDelta > deltaData::maxDeltaRatio * ownDelta)
78         {
79             changedFaces.append(faceI);
80             changedFacesInfo.append(deltaData(neiDelta));
81         }
82     }
84     // Insert all faces of coupled patches no matter what. Let FaceCellWave
85     // sort it out.
86     forAll(mesh.boundaryMesh(), patchI)
87     {
88         const polyPatch& patch = mesh.boundaryMesh()[patchI];
90         if (patch.coupled())
91         {
92             forAll(patch, patchFaceI)
93             {
94                 label meshFaceI = patch.start() + patchFaceI;
96                 scalar ownDelta = delta[mesh.faceOwner()[meshFaceI]];
98                 changedFaces.append(meshFaceI);
99                 changedFacesInfo.append(deltaData(ownDelta));
100             }
101         }
102     }
104     changedFaces.shrink();
105     changedFacesInfo.shrink();
109 void smoothDelta::calcDelta()
111     deltaData::maxDeltaRatio = maxDeltaRatio_;
112     const volScalarField& geometricDelta = geometricDelta_();
114     // Fill changed faces with info
115     dynamicLabelList changedFaces(mesh_.nFaces()/100 + 100);
116     DynamicList<deltaData> changedFacesInfo(changedFaces.size());
118     setChangedFaces(mesh_, geometricDelta, changedFaces, changedFacesInfo);
120     // Set initial field on cells.
121     List<deltaData> cellDeltaData(mesh_.nCells());
123     forAll(geometricDelta, cellI)
124     {
125         cellDeltaData[cellI] = geometricDelta[cellI];
126     }
128     // Set initial field on faces.
129     List<deltaData> faceDeltaData(mesh_.nFaces());
132     // Propagate information over whole domain.
133     FaceCellWave<deltaData> deltaCalc
134     (
135         mesh_,
136         changedFaces,
137         changedFacesInfo,
138         faceDeltaData,
139         cellDeltaData,
140         mesh_.nCells()  // max iterations
141     );
143     forAll(delta_, cellI)
144     {
145         delta_[cellI] = cellDeltaData[cellI].delta();
146     }
150 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
152 smoothDelta::smoothDelta
154     const word& name,
155     const fvMesh& mesh,
156     const dictionary& dd
159     LESdelta(name, mesh),
160     geometricDelta_
161     (
162         LESdelta::New("geometricDelta", mesh, dd.subDict(type() + "Coeffs"))
163     ),
164     maxDeltaRatio_
165     (
166         readScalar(dd.subDict(type() + "Coeffs").lookup("maxDeltaRatio"))
167     )
169     calcDelta();
173 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
175 void smoothDelta::read(const dictionary& d)
177     const dictionary& dd(d.subDict(type() + "Coeffs"));
179     geometricDelta_().read(dd);
180     dd.lookup("maxDeltaRatio") >> maxDeltaRatio_;
181     calcDelta();
185 void smoothDelta::correct()
187     geometricDelta_().correct();
189     if (mesh_.changing())
190     {
191         calcDelta();
192     }
196 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198 } // End namespace Foam
200 // ************************************************************************* //