StationaryPhaseModel.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) 2015-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::StationaryPhaseModel
26 
27 Description
28  Class which represents a stationary (and therefore probably solid) phase.
29  Generates, but does not store, zero velocity and flux field and turbulent
30  qauantities. Throws an error when non-const access is requested to the
31  motion fields or when the momentum equation is requested. Usage must
32  must protect against such calls.
33 
34 See also
35  MovingPhaseModel
36 
37 SourceFiles
38  StationaryPhaseModel.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef StationaryPhaseModel_H
43 #define StationaryPhaseModel_H
44 
45 #include "phaseModel.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class StationaryPhaseModel Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 template<class BasePhaseModel>
58 :
59  public BasePhaseModel
60 {
61 private:
62 
63  // Private member functions
64 
65  //- Create a zero geometric field
66  template<class Type, template<class> class PatchField, class GeoMesh>
68  (
69  const word& name,
70  const dimensionSet& dims,
71  const bool cache = false
72  ) const;
73 
74  //- Create a zero vol field
75  template<class Type>
77  (
78  const word& name,
79  const dimensionSet& dims,
80  const bool cache = false
81  ) const;
82 
83  //- Create a zero surface field
84  template<class Type>
86  (
87  const word& name,
88  const dimensionSet& dims,
89  const bool cache = false
90  ) const;
91 
92 
93 public:
94 
95  // Constructors
96 
98  (
99  const phaseSystem& fluid,
100  const word& phaseName,
101  const label index
102  );
103 
104 
105  //- Destructor
106  virtual ~StationaryPhaseModel();
107 
108 
109  // Member Functions
110 
111  // Momentum
112 
113  //- Return whether the phase is stationary
114  virtual bool stationary() const;
115 
116  //- Return the momentum equation
117  virtual tmp<fvVectorMatrix> UEqn();
118 
119  //- Return the momentum equation for the face-based algorithm
120  virtual tmp<fvVectorMatrix> UfEqn();
121 
122  //- Return the velocity
123  virtual tmp<volVectorField> U() const;
124 
125  //- Access the velocity
126  virtual volVectorField& URef();
127 
128  //- Return the volumetric flux
129  virtual tmp<surfaceScalarField> phi() const;
130 
131  //- Access the volumetric flux
132  virtual surfaceScalarField& phiRef();
133 
134  //- Return the volumetric flux of the phase
135  virtual tmp<surfaceScalarField> alphaPhi() const;
136 
137  //- Access the volumetric flux of the phase
138  virtual surfaceScalarField& alphaPhiRef();
139 
140  //- Return the mass flux of the phase
141  virtual tmp<surfaceScalarField> alphaRhoPhi() const;
142 
143  //- Access the mass flux of the phase
145 
146  //- Return the substantive acceleration
147  virtual tmp<volVectorField> DUDt() const;
148 
149  //- Return the substantive acceleration on the faces
150  virtual tmp<surfaceScalarField> DUDtf() const;
151 
152  //- Return the continuity error
153  virtual tmp<volScalarField> continuityError() const;
154 
155  //- Return the continuity error due to the flow field
157 
158  //- Return the continuity error due to any sources
160 
161  //- Return the phase kinetic energy
162  virtual tmp<volScalarField> K() const;
163 
164 
165  // Compressibility (variable density)
166 
167  //- Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
168  virtual tmp<volScalarField> divU() const;
169 
170  //- Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
171  virtual void divU(tmp<volScalarField> divU);
172 
173 
174  // Turbulence
175 
176  //- Return the turbulent dynamic viscosity
177  virtual tmp<volScalarField> mut() const;
178 
179  //- Return the effective dynamic viscosity
180  virtual tmp<volScalarField> muEff() const;
181 
182  //- Return the turbulent kinematic viscosity
183  virtual tmp<volScalarField> nut() const;
184 
185  //- Return the effective kinematic viscosity
186  virtual tmp<volScalarField> nuEff() const;
187 
188  //- Effective thermal turbulent diffusivity for temperature
189  // of mixture for patch [J/m/s/K]
190  using BasePhaseModel::kappaEff;
191 
192  //- Return the effective thermal conductivity
193  virtual tmp<volScalarField> kappaEff() const;
194 
195  //- Return the effective thermal conductivity on a patch
196  virtual tmp<scalarField> kappaEff(const label patchi) const;
197 
198  //- Effective thermal turbulent diffusivity of mixture [kg/m/s]
200 
201  //- Return the effective thermal diffusivity
202  virtual tmp<volScalarField> alphaEff() const;
203 
204  //- Return the effective thermal conductivity on a patch
205  virtual tmp<scalarField> alphaEff(const label patchi) const;
206 
207  //- Return the turbulent kinetic energy
208  virtual tmp<volScalarField> k() const;
209 
210  //- Return the phase-pressure'
211  // (derivative of phase-pressure w.r.t. phase-fraction)
212  virtual tmp<volScalarField> pPrime() const;
213 };
214 
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 } // End namespace Foam
219 
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 
222 #ifdef NoRepository
223  #include "StationaryPhaseModel.C"
224 #endif
225 
226 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 
228 #endif
229 
230 // ************************************************************************* //
virtual tmp< volScalarField > mut() const
Return the turbulent dynamic viscosity.
virtual tmp< surfaceScalarField > alphaPhi() const
Return the volumetric flux of the phase.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
multiphaseSystem & fluid
Definition: createFields.H:11
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux.
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
virtual tmp< volScalarField > continuityErrorSources() const
Return the continuity error due to any sources.
virtual tmp< volScalarField > divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual tmp< volScalarField > nuEff() const
Return the effective kinematic viscosity.
virtual tmp< volScalarField > kappaEff() const
Return the effective thermal conductivity.
virtual tmp< volVectorField > U() const
Return the velocity.
Dimension set for the base types.
Definition: dimensionSet.H:120
Class which represents a stationary (and therefore probably solid) phase. Generates, but does not store, zero velocity and flux field and turbulent qauantities. Throws an error when non-const access is requested to the motion fields or when the momentum equation is requested. Usage must must protect against such calls.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:67
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition: zeroField.H:50
A class for handling words, derived from string.
Definition: word.H:59
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
virtual tmp< volScalarField > continuityErrorFlow() const
Return the continuity error due to the flow field.
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
label patchi
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
virtual volVectorField & URef()
Access the velocity.
virtual tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > alphaEff() const
Return the effective thermal diffusivity.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
virtual ~StationaryPhaseModel()
Destructor.
virtual tmp< volScalarField > nut() const
Return the turbulent kinematic viscosity.
Namespace for OpenFOAM.
virtual bool stationary() const
Return whether the phase is stationary.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
StationaryPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)