smoothDelta.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 "smoothDelta.H"
28 #include "FvFaceCellWave.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace LESModels
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::LESModels::smoothDelta::setChangedFaces
45 (
46  const fvMesh& mesh,
47  const volScalarField& delta,
48  DynamicList<labelPair>& changedFaces,
49  DynamicList<deltaData>& changedFacesInfo
50 )
51 {
52  forAll(mesh.owner(), facei)
53  {
54  const scalar ownDelta = delta[mesh.owner()[facei]];
55  const scalar neiDelta = delta[mesh.neighbour()[facei]];
56 
57  // Check if owner delta much larger than neighbour delta or vice versa
58  if (ownDelta > maxDeltaRatio_*neiDelta)
59  {
60  changedFaces.append(labelPair(-1, facei));
61  changedFacesInfo.append(deltaData(ownDelta));
62  }
63  else if (neiDelta > maxDeltaRatio_*ownDelta)
64  {
65  changedFaces.append(labelPair(-1, facei));
66  changedFacesInfo.append(deltaData(neiDelta));
67  }
68  }
69 
70  // Insert all faces of coupled patches no matter what. Let FvFaceCellWave
71  // sort it out.
72  forAll(mesh.boundary(), patchi)
73  {
74  const fvPatch& patch = mesh.boundary()[patchi];
75 
76  if (patch.coupled())
77  {
78  forAll(patch, patchFacei)
79  {
80  const scalar ownDelta = delta[patch.faceCells()[patchFacei]];
81 
82  changedFaces.append(labelPair(patchi, patchFacei));
83  changedFacesInfo.append(deltaData(ownDelta));
84  }
85  }
86  }
87 
88  changedFaces.shrink();
89  changedFacesInfo.shrink();
90 }
91 
92 
93 void Foam::LESModels::smoothDelta::calcDelta()
94 {
95  const fvMesh& mesh = momentumTransportModel_.mesh();
96 
97  const volScalarField& geometricDelta = geometricDelta_();
98 
99  // Fill changed faces with info
100  DynamicList<labelPair> changedFaces(mesh.nFaces()/100 + 100);
101  DynamicList<deltaData> changedFacesInfo(changedFaces.size());
102  setChangedFaces(mesh, geometricDelta, changedFaces, changedFacesInfo);
103 
104  // Set initial field on cells.
105  List<deltaData> cellDeltaData(mesh.nCells());
106  forAll(geometricDelta, celli)
107  {
108  cellDeltaData[celli] = geometricDelta[celli];
109  }
110 
111  // Set initial field on faces.
112  List<deltaData> internalFaceDeltaData(mesh.nInternalFaces());
113  List<List<deltaData>> patchFaceDeltaData
114  (
115  FvFaceCellWave<deltaData>::template
116  sizesListList<List<List<deltaData>>>
117  (
118  FvFaceCellWave<deltaData>::template
119  listListSizes(mesh.boundary()),
120  deltaData()
121  )
122  );
123 
124  // Propagate information over whole domain.
125  FvFaceCellWave<deltaData, scalar> deltaCalc
126  (
127  mesh,
128  changedFaces,
129  changedFacesInfo,
130  internalFaceDeltaData,
131  patchFaceDeltaData,
132  cellDeltaData,
133  mesh.globalData().nTotalCells() + 1, // max iterations
134  maxDeltaRatio_
135  );
136 
137  forAll(delta_, celli)
138  {
139  delta_[celli] = cellDeltaData[celli].delta();
140  }
141 }
142 
143 
144 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
145 
147 (
148  const word& name,
150  const dictionary& dict
151 )
152 :
154  geometricDelta_
155  (
156  LESdelta::New
157  (
158  "geometricDelta",
159  turbulence,
160  dict.optionalSubDict(type() + "Coeffs")
161  )
162  ),
163  maxDeltaRatio_
164  (
165  dict.optionalSubDict(type() + "Coeffs").lookup<scalar>("maxDeltaRatio")
166  )
167 {
168  calcDelta();
169 }
170 
171 
172 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
173 
175 {
176  const dictionary& coeffsDict(dict.optionalSubDict(type() + "Coeffs"));
177 
178  geometricDelta_().read(coeffsDict);
179  coeffsDict.lookup("maxDeltaRatio") >> maxDeltaRatio_;
180  calcDelta();
181 }
182 
183 
185 {
186  geometricDelta_().correct();
187 
188  if (momentumTransportModel_.mesh().changing())
189  {
190  calcDelta();
191  }
192 }
193 
194 
195 // ************************************************************************* //
scalar delta
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Smoothed delta which takes a given simple geometric delta and applies smoothing to it such that the r...
Definition: smoothDelta.H:55
smoothDelta(const word &name, const momentumTransportModel &turbulence, const dictionary &)
Construct from name, momentumTransportModel and dictionary.
Definition: smoothDelta.C:147
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Definition: smoothDelta.C:174
Abstract base class for LES deltas.
Definition: LESdelta.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
Abstract base class for turbulence models (RAS, LES and laminar).
A class for handling words, derived from string.
Definition: word.H:62
label patchi
defineTypeNameAndDebug(cubeRootVolDelta, 0)
addToRunTimeSelectionTable(LESdelta, cubeRootVolDelta, dictionary)
Namespace for OpenFOAM.
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))