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 "nearWallDist.H"
31 #include "fvcFlux.H"
32 #include "fvmDiv.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
39 }
40 
41 
42 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
43 
46 (
47  const objectRegistry& obr,
48  const word& group,
49  bool registerObject
50 )
51 {
52  typeIOobject<IOdictionary> momentumTransport
53  (
54  IOobject::groupName(typeName, group),
55  obr.time().constant(),
56  obr,
60  );
61 
62  if (momentumTransport.headerOk())
63  {
64  return momentumTransport;
65  }
66  else
67  {
68  typeIOobject<IOdictionary> turbulenceProperties
69  (
70  IOobject::groupName("turbulenceProperties", group),
71  obr.time().constant(),
72  obr,
76  );
77 
78  if (turbulenceProperties.headerOk())
79  {
80  return turbulenceProperties;
81  }
82  else
83  {
84  return momentumTransport;
85  }
86  }
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
91 
93 (
94  const volVectorField& U,
95  const surfaceScalarField& alphaRhoPhi,
96  const surfaceScalarField& phi,
97  const viscosity& viscosity
98 )
99 :
100  IOdictionary(readModelDict(U.db(), alphaRhoPhi.group(), true)),
101 
102  runTime_(U.time()),
103  mesh_(U.mesh()),
104 
105  U_(U),
106  alphaRhoPhi_(alphaRhoPhi),
107  phi_(phi),
108  viscosity_(viscosity)
109 {
110  // Ensure name of IOdictionary is typeName
112 }
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
118 {
119  return phi_;
120 }
121 
122 
124 {
125  return nearWallDist::New(mesh_).y();
126 }
127 
128 
130 {
131  return regIOobject::read();
132 }
133 
134 
137 (
138  const tmp<volTensorField>& devTauCorr,
140 ) const
141 {
142  return fvm::divc
143  (
145  (
147  (
148  groupName("devTauCorrFlux"),
149  fvc::flux(devTauCorr)
150  )
151  ),
152  U
153  );
154 }
155 
156 
158 {}
159 
160 
162 {}
163 
164 
166 {}
167 
168 
169 // ************************************************************************* //
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).
virtual void validate()
Validate the fields after construction.
const volScalarField::Boundary & y() const
Return the near wall distance.
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.
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:105
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
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.