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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2016 OpenCFD Ltd.
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 {
37  addToRunTimeSelectionTable(LESdelta, smoothDelta, dictionary);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
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 
98 {
99  const fvMesh& mesh = turbulenceModel_.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.
123  (
124  mesh,
125  changedFaces,
126  changedFacesInfo,
127  faceDeltaData,
128  cellDeltaData,
129  mesh.globalData().nTotalCells()+1, // max iterations
131  );
132 
133  forAll(delta_, cellI)
134  {
135  delta_[cellI] = cellDeltaData[cellI].delta();
136  }
137 
138  // Handle coupled boundaries
140 }
141 
142 
143 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
144 
146 (
147  const word& name,
149  const dictionary& dict
150 )
151 :
153  geometricDelta_
154  (
156  (
157  "geometricDelta",
158  turbulence,
159  dict.subDict(type() + "Coeffs")
160  )
161  ),
162  maxDeltaRatio_
163  (
164  readScalar(dict.subDict(type() + "Coeffs").lookup("maxDeltaRatio"))
165  )
166 {
167  calcDelta();
168 }
169 
170 
171 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
172 
174 {
175  const dictionary& coeffsDict(dict.subDict(type() + "Coeffs"));
176 
177  geometricDelta_().read(coeffsDict);
178  coeffsDict.lookup("maxDeltaRatio") >> maxDeltaRatio_;
179  calcDelta();
180 }
181 
182 
184 {
185  geometricDelta_().correct();
186 
187  if (turbulenceModel_.mesh().changing())
188  {
189  calcDelta();
190  }
191 }
192 
193 
194 // ************************************************************************* //
Foam::LESModels::smoothDelta::read
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Definition: smoothDelta.C:173
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::LESModels::smoothDelta::deltaData
Public member class used by mesh-wave to propagate the delta-ratio.
Definition: smoothDelta.H:59
Foam::DynamicList< label >
Foam::polyPatch::coupled
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:322
Foam::LESdelta::delta_
volScalarField delta_
Definition: LESdelta.H:59
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::LESModels::smoothDelta::smoothDelta
smoothDelta(const smoothDelta &)
Disallow default bitwise copy construct and assignment.
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
smoothDelta.H
smoothDelta
Smoothed delta which takes a given simple geometric delta and applies smoothing to it such that the r...
Foam::LESModels::smoothDelta::geometricDelta_
autoPtr< LESdelta > geometricDelta_
Definition: smoothDelta.H:216
Foam::label
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
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::LESModels::smoothDelta::maxDeltaRatio_
scalar maxDeltaRatio_
Definition: smoothDelta.H:217
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::LESdelta::turbulenceModel_
const turbulenceModel & turbulenceModel_
Definition: LESdelta.H:57
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::DynamicList::shrink
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:258
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:60
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::LESModels::smoothDelta::calcDelta
void calcDelta()
Definition: smoothDelta.C:97
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
Foam::FaceCellWave
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:75
Foam::LESdelta
Abstract base class for LES deltas.
Definition: LESdelta.H:50
Foam::LESModels::defineTypeNameAndDebug
defineTypeNameAndDebug(cubeRootVolDelta, 0)
Foam::LESdelta::New
static autoPtr< LESdelta > New(const word &name, const turbulenceModel &turbulence, const dictionary &dict)
Return a reference to the selected LES delta.
Definition: LESdelta.C:66
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::turbulenceModel::mesh
const fvMesh & mesh() const
Definition: turbulenceModel.H:127
FaceCellWave.H
Foam::LESModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
Foam::LESModels::smoothDelta::correct
virtual void correct()
Definition: smoothDelta.C:183
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress
Foam::LESModels::smoothDelta::setChangedFaces
void setChangedFaces(const polyMesh &mesh, const volScalarField &delta, DynamicList< label > &changedFaces, DynamicList< deltaData > &changedFacesInfo)
Fill changedFaces (with face labels) and changedFacesInfo.
Definition: smoothDelta.C:45