RASModelVariablesI.H
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) 2007-2023 PCOpt/NTUA
9  Copyright (C) 2013-2023 FOSS GP
10  Copyright (C) 2019-2020 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 \*---------------------------------------------------------------------------*/
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace incompressible
34 {
35 
36 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
37 
38 inline const word& RASModelVariables::TMVar1BaseName() const
39 {
40  return TMVar1BaseName_;
41 }
42 
43 
44 inline const word& RASModelVariables::TMVar2BaseName() const
45 {
46  return TMVar2BaseName_;
47 }
48 
49 
50 inline const word& RASModelVariables::nutBaseName() const
51 {
52  return nutBaseName_;
53 }
54 
55 
56 inline bool RASModelVariables::hasTMVar1() const
57 {
58  return bool(TMVar1Ptr_);
59 }
60 
61 
62 inline bool RASModelVariables::hasTMVar2() const
63 {
64  return bool(TMVar2Ptr_);
65 }
66 
67 
68 inline bool RASModelVariables::hasNut() const
69 {
70  return bool(nutPtr_);
71 }
72 
73 
74 inline bool RASModelVariables::hasDist() const
75 {
76  return bool(distPtr_);
77 }
78 
79 
80 inline const volScalarField& RASModelVariables::TMVar1() const
81 {
83  {
84  return TMVar1MeanPtr_.cref();
85  }
86 
87  return TMVar1Ptr_.cref();
88 }
89 
90 
92 {
94  {
95  return TMVar1MeanPtr_.ref();
96  }
97 
98  return TMVar1Ptr_.ref();
99 }
100 
101 
102 inline const volScalarField& RASModelVariables::TMVar2() const
103 {
105  {
106  return TMVar2MeanPtr_.cref();
107  }
108 
109  return TMVar2Ptr_.cref();
110 }
111 
113 {
115  {
116  return TMVar2MeanPtr_.ref();
117  }
118 
119  return TMVar2Ptr_.ref();
120 }
121 
122 inline const volScalarField& RASModelVariables::nutRef() const
123 {
125  {
126  return nutMeanPtr_.cref();
127  }
128 
129  return nutPtr_.cref();
130 }
131 
132 
134 {
136  {
137  return nutMeanPtr_.ref();
138  }
139 
140  return nutPtr_.ref();
141 }
142 
143 
145 {
146  if (hasNut())
147  {
148  return tmp<volScalarField>(nutRef());
149  }
150 
151  return volScalarField::New
152  (
153  "dummylaminarNut",
155  mesh_,
157  );
158 }
159 
161 inline const volScalarField& RASModelVariables::d() const
162 {
163  return distPtr_.cref();
164 }
165 
168 {
169  return distPtr_.ref();
170 }
171 
172 
173 inline tmp<scalarField> RASModelVariables::TMVar1(const label patchi) const
174 {
175  if (hasTMVar1())
176  {
177  return TMVar1().boundaryField()[patchi];
178  }
179 
180  return tmp<scalarField>::New(mesh_.boundary()[patchi].size(), Zero);
181 }
182 
183 
184 inline tmp<scalarField> RASModelVariables::TMVar2(const label patchi) const
185 {
186  if (hasTMVar2())
187  {
188  return TMVar2().boundaryField()[patchi];
189  }
190 
191  return tmp<scalarField>::New(mesh_.boundary()[patchi].size(), Zero);
192 }
193 
194 
195 inline tmp<scalarField> RASModelVariables::nut(const label patchi) const
196 {
197  if (hasNut())
198  {
199  return nutRef().boundaryField()[patchi];
200  }
201 
202  return tmp<scalarField>::New(mesh_.boundary()[patchi].size(), Zero);
203 }
204 
205 
207 RASModelVariables::nutPatchField(const label patchi) const
208 {
209  if (hasNut())
210  {
211  return nutRef().boundaryField()[patchi];
212  }
214  // Using dummy internalField
215  return
217 }
218 
220 inline const volScalarField& RASModelVariables::TMVar1Inst() const
221 {
222  return TMVar1Ptr_.cref();
223 }
224 
227 {
228  return TMVar1Ptr_.ref();
229 }
230 
232 inline const volScalarField& RASModelVariables::TMVar2Inst() const
233 {
234  return TMVar2Ptr_.cref();
235 }
236 
239 {
240  return TMVar2Ptr_.ref();
241 }
242 
244 inline const volScalarField& RASModelVariables::nutRefInst() const
245 {
246  return nutPtr_.cref();
247 }
248 
249 
251 {
252  return nutPtr_.ref();
253 }
254 
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 } // End namespace incompressible
259 } // End namespace Foam
260 
261 // ************************************************************************* //
const volScalarField & TMVar1Inst() const
Return references to instantaneous turbulence fields.
const dimensionSet dimViscosity
const volScalarField & TMVar2() const
const volScalarField & d() const
const volScalarField & nutRef() const
virtual bool hasTMVar1() const
Bools to identify which turbulent fields are present.
tmp< fvPatchScalarField > nutPatchField(const label patchi) const
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
A class for handling words, derived from Foam::string.
Definition: word.H:63
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const volScalarField & TMVar1() const
Return references to turbulence fields.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
const word & TMVar1BaseName() const
Turbulence field names.
const volScalarField & TMVar2Inst() const
const volScalarField & nutRefInst() const
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
A class for managing temporary objects.
Definition: HashPtrTable.H:50
bool useAveragedFields() const
Use averaged fields? For solving the adjoint equations or computing sensitivities based on averaged f...
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Do not request registration (bool: false)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127