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-2020 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 "FaceCellWave.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace LESModels
35 {
36  defineTypeNameAndDebug(smoothDelta, 0);
37  addToRunTimeSelectionTable(LESdelta, smoothDelta, dictionary);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::LESModels::smoothDelta::setChangedFaces
45 (
46  const polyMesh& mesh,
47  const volScalarField& delta,
48  DynamicList<label>& changedFaces,
49  DynamicList<deltaData>& changedFacesInfo
50 )
51 {
52  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
53  {
54  scalar ownDelta = delta[mesh.faceOwner()[facei]];
55 
56  scalar neiDelta = delta[mesh.faceNeighbour()[facei]];
57 
58  // Check if owner delta much larger than neighbour delta or vice versa
59 
60  if (ownDelta > maxDeltaRatio_ * neiDelta)
61  {
62  changedFaces.append(facei);
63  changedFacesInfo.append(deltaData(ownDelta));
64  }
65  else if (neiDelta > maxDeltaRatio_ * ownDelta)
66  {
67  changedFaces.append(facei);
68  changedFacesInfo.append(deltaData(neiDelta));
69  }
70  }
71 
72  // Insert all faces of coupled patches no matter what. Let FaceCellWave
73  // sort it out.
74  forAll(mesh.boundaryMesh(), patchi)
75  {
76  const polyPatch& patch = mesh.boundaryMesh()[patchi];
77 
78  if (patch.coupled())
79  {
80  forAll(patch, patchFacei)
81  {
82  label meshFacei = patch.start() + patchFacei;
83 
84  scalar ownDelta = delta[mesh.faceOwner()[meshFacei]];
85 
86  changedFaces.append(meshFacei);
87  changedFacesInfo.append(deltaData(ownDelta));
88  }
89  }
90  }
91 
92  changedFaces.shrink();
93  changedFacesInfo.shrink();
94 }
95 
96 
97 void Foam::LESModels::smoothDelta::calcDelta()
98 {
99  const fvMesh& mesh = momentumTransportModel_.mesh();
100 
101  const volScalarField& geometricDelta = geometricDelta_();
102 
103  // Fill changed faces with info
104  DynamicList<label> changedFaces(mesh.nFaces()/100 + 100);
105  DynamicList<deltaData> changedFacesInfo(changedFaces.size());
106 
107  setChangedFaces(mesh, geometricDelta, changedFaces, changedFacesInfo);
108 
109  // Set initial field on cells.
110  List<deltaData> cellDeltaData(mesh.nCells());
111 
112  forAll(geometricDelta, celli)
113  {
114  cellDeltaData[celli] = geometricDelta[celli];
115  }
116 
117  // Set initial field on faces.
118  List<deltaData> faceDeltaData(mesh.nFaces());
119 
120 
121  // Propagate information over whole domain.
122  FaceCellWave<deltaData, scalar> deltaCalc
123  (
124  mesh,
125  changedFaces,
126  changedFacesInfo,
127  faceDeltaData,
128  cellDeltaData,
129  mesh.globalData().nTotalCells()+1, // max iterations
130  maxDeltaRatio_
131  );
132 
133  forAll(delta_, celli)
134  {
135  delta_[celli] = cellDeltaData[celli].delta();
136  }
137 }
138 
139 
140 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
141 
143 (
144  const word& name,
145  const momentumTransportModel& turbulence,
146  const dictionary& dict
147 )
148 :
149  LESdelta(name, turbulence),
150  geometricDelta_
151  (
153  (
154  "geometricDelta",
155  turbulence,
156  dict.optionalSubDict(type() + "Coeffs")
157  )
158  ),
159  maxDeltaRatio_
160  (
161  dict.optionalSubDict(type() + "Coeffs").lookup<scalar>("maxDeltaRatio")
162  )
163 {
164  calcDelta();
165 }
166 
167 
168 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
169 
171 {
172  const dictionary& coeffsDict(dict.optionalSubDict(type() + "Coeffs"));
173 
174  geometricDelta_().read(coeffsDict);
175  coeffsDict.lookup("maxDeltaRatio") >> maxDeltaRatio_;
176  calcDelta();
177 }
178 
179 
181 {
182  geometricDelta_().correct();
183 
184  if (momentumTransportModel_.mesh().changing())
185  {
186  calcDelta();
187  }
188 }
189 
190 
191 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
static autoPtr< LESdelta > New(const word &name, const momentumTransportModel &turbulence, const dictionary &dict)
Return a reference to the selected LES delta.
Definition: LESdelta.C:66
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Definition: smoothDelta.C:170
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Abstract base class for LES deltas.
Definition: LESdelta.H:50
addToRunTimeSelectionTable(LESdelta, cubeRootVolDelta, dictionary)
bool read(Istream &, const bool keepHeader=false)
Read dictionary from Istream, optionally keeping the header.
Definition: dictionaryIO.C:104
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:1040
A class for handling words, derived from string.
Definition: word.H:59
defineTypeNameAndDebug(cubeRootVolDelta, 0)
Abstract base class for turbulence models (RAS, LES and laminar).
label patchi
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
smoothDelta(const word &name, const momentumTransportModel &turbulence, const dictionary &)
Construct from name, momentumTransportModel and dictionary.
Definition: smoothDelta.C:143
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844