maxDeltaxyz.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 "maxDeltaxyz.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace LESModels
34 {
35  defineTypeNameAndDebug(maxDeltaxyz, 0);
36  addToRunTimeSelectionTable(LESdelta, maxDeltaxyz, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 void Foam::LESModels::maxDeltaxyz::calcDelta()
44 {
45  const fvMesh& mesh = momentumTransportModel_.mesh();
46 
47  label nD = mesh.nGeometricD();
48 
49  const cellList& cells = mesh.cells();
50  scalarField hmax(cells.size());
51 
52  forAll(cells,celli)
53  {
54  scalar deltaMaxTmp = 0.0;
55  const labelList& cFaces = mesh.cells()[celli];
56  const point& centrevector = mesh.cellCentres()[celli];
57 
58  forAll(cFaces, cFacei)
59  {
60  label facei = cFaces[cFacei];
61  const point& facevector = mesh.faceCentres()[facei];
62  scalar tmp = mag(facevector - centrevector);
63  if (tmp > deltaMaxTmp)
64  {
65  deltaMaxTmp = tmp;
66  }
67  }
68 
69  hmax[celli] = deltaCoeff_*deltaMaxTmp;
70  }
71 
72  if (nD == 3)
73  {
74  delta_.primitiveFieldRef() = hmax;
75  }
76  else if (nD == 2)
77  {
79  << "Case is 2D, LES is not strictly applicable\n"
80  << endl;
81 
82  delta_.primitiveFieldRef() = hmax;
83  }
84  else
85  {
87  << "Case is not 3D or 2D, LES is not applicable"
88  << exit(FatalError);
89  }
90 }
91 
92 
93 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
94 
96 (
97  const word& name,
98  const momentumTransportModel& turbulence,
99  const dictionary& dict
100 )
101 :
102  LESdelta(name, turbulence),
103  deltaCoeff_
104  (
105  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<scalar>
106  (
107  "deltaCoeff",
108  1
109  )
110  )
111 {
112  calcDelta();
113 }
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 {
120  const dictionary& coeffsDict(dict.optionalSubDict(type() + "Coeffs"));
121 
122  coeffsDict.readIfPresent<scalar>("deltaCoeff", deltaCoeff_);
123 
124  calcDelta();
125 }
126 
127 
129 {
130  if (momentumTransportModel_.mesh().changing())
131  {
132  calcDelta();
133  }
134 }
135 
136 
137 // ************************************************************************* //
#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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
Abstract base class for LES deltas.
Definition: LESdelta.H:50
addToRunTimeSelectionTable(LESdelta, cubeRootVolDelta, dictionary)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Macros for easy insertion into run-time selection tables.
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
List< label > labelList
A List of labels.
Definition: labelList.H:56
defineTypeNameAndDebug(cubeRootVolDelta, 0)
Abstract base class for turbulence models (RAS, LES and laminar).
vector point
Point is a vector.
Definition: point.H:41
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
maxDeltaxyz(const word &name, const momentumTransportModel &turbulence, const dictionary &)
Construct from name, momentumTransportModel and dictionary.
Definition: maxDeltaxyz.C:96
dimensioned< scalar > mag(const dimensioned< Type > &)
List< cell > cellList
list of cells
Definition: cellList.H:42
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Definition: maxDeltaxyz.C:118
Namespace for OpenFOAM.