maxDeltaxyz.H
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-2018 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 Class
25  Foam::LESModels::maxDeltaxyz
26 
27 Description
28  Delta calculated by taking the maximum distance between the cell centre
29  and any face centre. For a regular hex cell, the computed delta will
30  equate to half of the cell width; accordingly, the deltaCoeff model
31  coefficient should be set to 2 for this case.
32 
33 SourceFiles
34  maxDeltaxyz.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef maxDeltaxyz_H
39 #define maxDeltaxyz_H
40 
41 #include "LESdelta.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 namespace LESModels
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class maxDeltaxyz Declaration
52 \*---------------------------------------------------------------------------*/
53 
54 class maxDeltaxyz
55 :
56  public LESdelta
57 {
58  // Private data
59 
60  //- Model coefficient
61  scalar deltaCoeff_;
62 
63 
64  // Private Member Functions
65 
66  //- Disallow default bitwise copy construct and assignment
67  maxDeltaxyz(const maxDeltaxyz&);
68  void operator=(const maxDeltaxyz&);
69 
70  // Calculate the delta values
71  void calcDelta();
72 
73 
74 public:
75 
76  //- Runtime type information
77  TypeName("maxDeltaxyz");
78 
79 
80  // Constructors
81 
82  //- Construct from name, turbulenceModel and dictionary
84  (
85  const word& name,
87  const dictionary&
88  );
89 
90 
91  //- Destructor
92  virtual ~maxDeltaxyz()
93  {}
94 
95 
96  // Member Functions
97 
98  //- Read the LESdelta dictionary
99  virtual void read(const dictionary&);
100 
101  // Correct values
102  virtual void correct();
103 };
104 
105 
106 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
107 
108 } // End namespace LESModels
109 } // End namespace Foam
110 
111 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
112 
113 
114 #endif
115 
116 // ************************************************************************* //
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
Definition: LESdelta.H:131
Abstract base class for LES deltas.
Definition: LESdelta.H:50
Abstract base class for turbulence models (RAS, LES and laminar).
A class for handling words, derived from string.
Definition: word.H:59
TypeName("maxDeltaxyz")
Runtime type information.
Delta calculated by taking the maximum distance between the cell centre and any face centre...
Definition: maxDeltaxyz.H:53
virtual ~maxDeltaxyz()
Destructor.
Definition: maxDeltaxyz.H:91
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Definition: maxDeltaxyz.C:118
Namespace for OpenFOAM.