momentumTransportModel.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-2023 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 "momentumTransportModel.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "wallFvPatch.H"
30 #include "wallDist.H"
31 #include "nearWallDist.H"
32 #include "fvcFlux.H"
33 #include "fvmDiv.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
40 }
41 
42 
43 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
44 
47 (
48  const objectRegistry& obr,
49  const word& group,
50  bool registerObject
51 )
52 {
53  typeIOobject<IOdictionary> momentumTransport
54  (
55  IOobject::groupName(typeName, group),
56  obr.time().constant(),
57  obr,
61  );
62 
63  if (momentumTransport.headerOk())
64  {
65  return momentumTransport;
66  }
67  else
68  {
69  typeIOobject<IOdictionary> turbulenceProperties
70  (
71  IOobject::groupName("turbulenceProperties", group),
72  obr.time().constant(),
73  obr,
77  );
78 
79  if (turbulenceProperties.headerOk())
80  {
81  return turbulenceProperties;
82  }
83  else
84  {
85  return momentumTransport;
86  }
87  }
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
92 
94 (
95  const volVectorField& U,
96  const surfaceScalarField& alphaRhoPhi,
97  const surfaceScalarField& phi,
98  const viscosity& viscosity
99 )
100 :
101  IOdictionary(readModelDict(U.db(), alphaRhoPhi.group(), true)),
102 
103  runTime_(U.time()),
104  mesh_(U.mesh()),
105 
106  U_(U),
107  alphaRhoPhi_(alphaRhoPhi),
108  phi_(phi),
109  viscosity_(viscosity)
110 {
111  // Ensure name of IOdictionary is typeName
113 }
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 {
120  return phi_;
121 }
122 
123 
125 {
126  return wallDist::New(mesh_).y();
127 }
128 
129 
131 {
132  return nearWallDist::New(mesh_).y();
133 }
134 
135 
137 {
138  return regIOobject::read();
139 }
140 
141 
144 (
145  const tmp<volTensorField>& devTauCorr,
147 ) const
148 {
149  return fvm::divc
150  (
152  (
154  (
155  groupName("devTauCorrFlux"),
156  fvc::flux(devTauCorr)
157  )
158  ),
159  U
160  );
161 }
162 
163 
165 {}
166 
167 
169 {}
170 
171 
173 {}
174 
175 
176 // ************************************************************************* //
static wallDist & New(const word &name, const fvMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Generic GeometricBoundaryField class.
Generic GeometricField class.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:119
bool & registerObject()
Register object created from this IOobject with registry if true.
Definition: IOobject.H:346
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:134
word group() const
Return group (extension part of name)
Definition: IOobject.C:324
static word groupName(Name name, const word &group)
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
Abstract base class for turbulence models (RAS, LES and laminar).
const volScalarField::Boundary & yb() const
Return the near wall distance.
virtual void validate()
Validate the fields after construction.
static typeIOobject< IOdictionary > readModelDict(const objectRegistry &obr, const word &group, bool registerObject=false)
virtual bool read()=0
Read model coefficients if they have changed.
momentumTransportModel(const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
virtual tmp< fvVectorMatrix > divDevTauCorr(const tmp< volTensorField > &devTauCorr, volVectorField &U) const
Return the explicit stress correction matrix.
virtual void correct()=0
Solve the momentum transport model equations.
virtual void predict()=0
Predict the momentum transport coefficients if possible.
const volScalarField & y() const
Return the wall distance.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux field.
const surfaceScalarField & alphaRhoPhi() const
Access function to phase flux field.
const volScalarField::Boundary & y() const
Access to the near-wall distances.
Definition: nearWallDist.H:99
Registry of regIOobjects.
const Time & time() const
Return time.
virtual void rename(const word &newName)
Rename.
Definition: regIOobject.C:417
virtual bool read()
Read object.
A class for managing temporary objects.
Definition: tmp.H:55
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
bool headerOk()
Read header (uses typeGlobalFile to find file) and check.
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
const volScalarField & y() const
Return reference to cached distance-to-wall field.
Definition: wallDist.H:151
A class for handling words, derived from string.
Definition: word.H:62
Calculate the face-flux of the given field.
Calculate the matrix for the divergence of the given field and flux.
U
Definition: pEqn.H:72
const char *const group
Group name for atomic constants.
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > flux(const VolField< Type > &vf)
Return the face-flux field obtained from the given volVectorField.
tmp< fvMatrix< Type > > divc(const tmp< SurfaceField< Type >> &tflux, const VolField< Type > &vf)
Return the explicit div matrix.
Definition: fvmDiv.C:105
Namespace for OpenFOAM.
defineTypeNameAndDebug(combustionModel, 0)
Foam::surfaceFields.