smoothDelta.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) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "smoothDelta.H"
31 #include "FaceCellWave.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace LESModels
38 {
39  defineTypeNameAndDebug(smoothDelta, 0);
40  addToRunTimeSelectionTable(LESdelta, smoothDelta, dictionary);
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 void Foam::LESModels::smoothDelta::setChangedFaces
48 (
49  const polyMesh& mesh,
50  const volScalarField& delta,
51  DynamicList<label>& changedFaces,
52  DynamicList<deltaData>& changedFacesInfo
53 )
54 {
55  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
56  {
57  scalar ownDelta = delta[mesh.faceOwner()[facei]];
58 
59  scalar neiDelta = delta[mesh.faceNeighbour()[facei]];
60 
61  // Check if owner delta much larger than neighbour delta or vice versa
62 
63  if (ownDelta > maxDeltaRatio_ * neiDelta)
64  {
65  changedFaces.append(facei);
66  changedFacesInfo.append(deltaData(ownDelta));
67  }
68  else if (neiDelta > maxDeltaRatio_ * ownDelta)
69  {
70  changedFaces.append(facei);
71  changedFacesInfo.append(deltaData(neiDelta));
72  }
73  }
74 
75  // Insert all faces of coupled patches no matter what. Let FaceCellWave
76  // sort it out.
77  forAll(mesh.boundaryMesh(), patchi)
78  {
79  const polyPatch& patch = mesh.boundaryMesh()[patchi];
80 
81  if (patch.coupled())
82  {
83  forAll(patch, patchFacei)
84  {
85  label meshFacei = patch.start() + patchFacei;
86 
87  scalar ownDelta = delta[mesh.faceOwner()[meshFacei]];
88 
89  changedFaces.append(meshFacei);
90  changedFacesInfo.append(deltaData(ownDelta));
91  }
92  }
93  }
94 
95  changedFaces.shrink();
96  changedFacesInfo.shrink();
97 }
98 
99 
100 void Foam::LESModels::smoothDelta::calcDelta()
101 {
102  const fvMesh& mesh = turbulenceModel_.mesh();
103 
104  const volScalarField& geometricDelta = geometricDelta_();
105 
106  // Fill changed faces with info
107  DynamicList<label> changedFaces(mesh.nFaces()/100 + 100);
108  DynamicList<deltaData> changedFacesInfo(changedFaces.size());
109 
110  setChangedFaces(mesh, geometricDelta, changedFaces, changedFacesInfo);
111 
112  // Set initial field on cells.
113  List<deltaData> cellDeltaData(mesh.nCells());
114 
115  forAll(geometricDelta, celli)
116  {
117  cellDeltaData[celli] = geometricDelta[celli];
118  }
119 
120  // Set initial field on faces.
121  List<deltaData> faceDeltaData(mesh.nFaces());
122 
123 
124  // Propagate information over whole domain.
125  FaceCellWave<deltaData, scalar> deltaCalc
126  (
127  mesh,
128  changedFaces,
129  changedFacesInfo,
130  faceDeltaData,
131  cellDeltaData,
132  mesh.globalData().nTotalCells()+1, // max iterations
133  maxDeltaRatio_
134  );
135 
136  forAll(delta_, celli)
137  {
138  delta_[celli] = cellDeltaData[celli].delta();
139  }
140 
141  // Handle coupled boundaries
142  delta_.correctBoundaryConditions();
143 }
144 
145 
146 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
147 
148 Foam::LESModels::smoothDelta::smoothDelta
149 (
150  const word& name,
152  const dictionary& dict
153 )
154 :
156  geometricDelta_
157  (
158  LESdelta::New
159  (
160  IOobject::groupName("geometricDelta", turbulence.U().group()),
161  turbulence,
162  dict.optionalSubDict(type() + "Coeffs")
163  )
164  ),
165  maxDeltaRatio_
166  (
167  dict.optionalSubDict(type() + "Coeffs").get<scalar>("maxDeltaRatio")
168  )
169 {
170  calcDelta();
171 }
172 
173 
174 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
175 
177 {
178  const dictionary& coeffsDict(dict.optionalSubDict(type() + "Coeffs"));
180  geometricDelta_().read(coeffsDict);
181  coeffsDict.readEntry("maxDeltaRatio", maxDeltaRatio_);
182  calcDelta();
183 }
184 
185 
187 {
188  geometricDelta_().correct();
189 
190  if (turbulenceModel_.mesh().changing())
191  {
192  calcDelta();
193  }
194 }
195 
196 
197 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
scalar delta
dictionary dict
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
bool read(Istream &is)
Read dictionary from Istream (discards the header). Reads entries until EOF or when the first token i...
Definition: dictionaryIO.C:134
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Abstract base class for LES deltas.
Definition: LESdelta.H:49
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
constexpr const char *const group
Group name for atomic constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
virtual void read(const dictionary &dict)
Read the LESdelta dictionary.
Definition: smoothDelta.C:169
defineTypeNameAndDebug(cubeRootVolDelta, 0)
U
Definition: pEqn.H:72
const std::string patch
OpenFOAM patch number as a std::string.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Namespace for OpenFOAM.