Loading...
Searching...
No Matches
maxDeltaxyz.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-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "maxDeltaxyz.H"
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace LESModels
37{
40}
41}
42
43
44// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45
47{
49
50 label nD = mesh.nGeometricD();
51
52 const cellList& cells = mesh.cells();
53 const vectorField& cellC = mesh.cellCentres();
54 const vectorField& faceC = mesh.faceCentres();
55 const vectorField faceN(mesh.faceAreas()/mag(mesh.faceAreas()));
56 scalarField hmax(cells.size());
57
58 forAll(cells, celli)
59 {
60 scalar deltaMaxTmp = 0.0;
61 const labelList& cFaces = cells[celli];
62 const point& cc = cellC[celli];
63
64 forAll(cFaces, cFacei)
65 {
66 label facei = cFaces[cFacei];
67 const point& fc = faceC[facei];
68 const vector& n = faceN[facei];
69
70 scalar tmp = mag(n & (fc - cc));
71 if (tmp > deltaMaxTmp)
72 {
73 deltaMaxTmp = tmp;
74 }
75 }
76
77 hmax[celli] = deltaCoeff_*deltaMaxTmp;
78 }
79
80 if (nD == 3)
81 {
82 delta_.primitiveFieldRef() = hmax;
83 }
84 else if (nD == 2)
85 {
87 << "Case is 2D, LES is not strictly applicable" << nl
88 << endl;
89
90 delta_.primitiveFieldRef() = hmax;
91 }
92 else
93 {
95 << "Case is not 3D or 2D, LES is not applicable"
96 << exit(FatalError);
97 }
98
99 // Handle coupled boundaries
100 delta_.correctBoundaryConditions();
101}
102
103
104// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
105
106Foam::LESModels::maxDeltaxyz::maxDeltaxyz
107(
108 const word& name,
110 const dictionary& dict
111)
112:
114 deltaCoeff_
115 (
116 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
117 (
118 "deltaCoeff",
119 2
120 )
121 )
123 calcDelta();
124}
125
126
127// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128
130{
131 const dictionary& coeffsDict(dict.optionalSubDict(type() + "Coeffs"));
133 coeffsDict.readIfPresent<scalar>("deltaCoeff", deltaCoeff_);
134
135 calcDelta();
136}
137
138
140{
141 if (turbulenceModel_.mesh().changing())
142 {
143 calcDelta();
144 }
145}
146
147
148// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void calcDelta()
Calculate the delta values.
Definition maxDeltaxyz.C:39
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Abstract base class for LES deltas.
Definition LESdelta.H:50
LESdelta(const LESdelta &)=delete
No copy construct.
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
Definition LESdelta.H:146
volScalarField delta_
Definition LESdelta.H:57
const turbulenceModel & turbulenceModel_
Definition LESdelta.H:55
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A class for managing temporary objects.
Definition tmp.H:75
Abstract base class for turbulence models (RAS, LES and laminar).
const fvMesh & mesh() const
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const cellShapeList & cells
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for LES SGS models.
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299