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-2021 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  quantities. 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 public:
62 
63  // Constructors
64 
66  (
67  const phaseSystem& fluid,
68  const word& phaseName,
69  const bool referencePhase,
70  const label index
71  );
72 
73 
74  //- Destructor
75  virtual ~StationaryPhaseModel();
76 
77 
78  // Member Functions
79 
80  // Momentum
81 
82  //- Return whether the phase is stationary
83  virtual bool stationary() const;
84 
85  //- Return the momentum equation
86  virtual tmp<fvVectorMatrix> UEqn();
87 
88  //- Return the momentum equation for the face-based algorithm
89  virtual tmp<fvVectorMatrix> UfEqn();
90 
91  //- Return the velocity
92  virtual tmp<volVectorField> U() const;
93 
94  //- Access the velocity
95  virtual volVectorField& URef();
96 
97  //- Return the volumetric flux
98  virtual tmp<surfaceScalarField> phi() const;
99 
100  //- Access the volumetric flux
101  virtual surfaceScalarField& phiRef();
102 
103  //- Return the face velocity
104  // Required for moving mesh cases
105  virtual tmp<surfaceVectorField> Uf() const;
106 
107  //- Access the face velocity
108  // Required for moving mesh cases
109  virtual surfaceVectorField& UfRef();
110 
111  //- Return the volumetric flux of the phase
112  virtual tmp<surfaceScalarField> alphaPhi() const;
113 
114  //- Access the volumetric flux of the phase
115  virtual surfaceScalarField& alphaPhiRef();
116 
117  //- Return the mass flux of the phase
118  virtual tmp<surfaceScalarField> alphaRhoPhi() const;
119 
120  //- Access the mass flux of the phase
122 
123  //- Return the substantive acceleration
124  virtual tmp<volVectorField> DUDt() const;
125 
126  //- Return the substantive acceleration on the faces
127  virtual tmp<surfaceScalarField> DUDtf() const;
128 
129  //- Return the continuity error
130  virtual tmp<volScalarField> continuityError() const;
131 
132  //- Return the phase kinetic energy
133  virtual tmp<volScalarField> K() const;
134 
135 
136  // Compressibility (variable density)
137 
138  //- Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
139  virtual tmp<volScalarField> divU() const;
140 
141  //- Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
142  virtual void divU(tmp<volScalarField> divU);
143 
144 
145  // Momentum transport
146 
147  //- Return the effective thermal conductivity on a patch
148  virtual tmp<scalarField> kappaEff(const label patchi) const;
149 
150  //- Return the turbulent kinetic energy
151  virtual tmp<volScalarField> k() const;
152 
153  //- Return the phase-pressure'
154  // (derivative of phase-pressure w.r.t. phase-fraction)
155  virtual tmp<volScalarField> pPrime() const;
156 
157 
158  // Thermophysical transport
159 
160  //- Return the source term for the energy equation
161  virtual tmp<fvScalarMatrix> divq(volScalarField& he) const;
162 };
163 
164 
165 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 
167 } // End namespace Foam
168 
169 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
170 
171 #ifdef NoRepository
172  #include "StationaryPhaseModel.C"
173 #endif
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 #endif
178 
179 // ************************************************************************* //
virtual tmp< surfaceScalarField > alphaPhi() const
Return the volumetric flux of the phase.
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
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 > divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual tmp< volVectorField > U() const
Return the velocity.
Class which represents a stationary (and therefore probably solid) phase. Generates, but does not store, zero velocity and flux field and turbulent quantities. 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.
virtual tmp< scalarField > kappaEff(const label patchi) const
Return the effective thermal conductivity on a patch.
virtual tmp< surfaceVectorField > Uf() const
Return the face velocity.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:68
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.
phaseSystem & fluid
Definition: createFields.H:11
thermo he()
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
label patchi
StationaryPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
virtual volVectorField & URef()
Access the velocity.
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
virtual surfaceVectorField & UfRef()
Access the face velocity.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
virtual ~StationaryPhaseModel()
Destructor.
Namespace for OpenFOAM.
virtual bool stationary() const
Return whether the phase is stationary.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.