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-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 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  //- Access the velocity
98  virtual const volVectorField& URef() const;
99 
100  //- Return the volumetric flux
101  virtual tmp<surfaceScalarField> phi() const;
102 
103  //- Access the volumetric flux
104  virtual surfaceScalarField& phiRef();
105 
106  //- Access the volumetric flux
107  virtual const surfaceScalarField& phiRef() const;
108 
109  //- Return the face velocity
110  // Required for moving mesh cases
111  virtual const autoPtr<surfaceVectorField>& Uf() const;
112 
113  //- Access the face velocity
114  // Required for moving mesh cases
115  virtual surfaceVectorField& UfRef();
116 
117  //- Access the face velocity
118  // Required for moving mesh cases
119  virtual const surfaceVectorField& UfRef() const;
120 
121  //- Return the volumetric flux of the phase
122  virtual tmp<surfaceScalarField> alphaPhi() const;
123 
124  //- Access the volumetric flux of the phase
125  virtual surfaceScalarField& alphaPhiRef();
126 
127  //- Access the volumetric flux of the phase
128  virtual const surfaceScalarField& alphaPhiRef() const;
129 
130  //- Return the mass flux of the phase
131  virtual tmp<surfaceScalarField> alphaRhoPhi() const;
132 
133  //- Access the mass flux of the phase
135 
136  //- Access the mass flux of the phase
137  virtual const surfaceScalarField& alphaRhoPhiRef() const;
138 
139  //- Return the substantive acceleration
140  virtual tmp<volVectorField> DUDt() const;
141 
142  //- Return the substantive acceleration on the faces
143  virtual tmp<surfaceScalarField> DUDtf() const;
144 
145  //- Return the continuity error
146  virtual tmp<volScalarField> continuityError() const;
147 
148  //- Return the phase kinetic energy
149  virtual tmp<volScalarField> K() const;
150 
151 
152  // Compressibility (variable density)
153 
154  //- Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
155  virtual const autoPtr<volScalarField>& divU() const;
156 
157  //- Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
158  virtual void divU(tmp<volScalarField> divU);
159 
160 
161  // Momentum transport
162 
163  //- Return the effective thermal conductivity on a patch
164  virtual tmp<scalarField> kappaEff(const label patchi) const;
165 
166  //- Return the turbulent kinetic energy
167  virtual tmp<volScalarField> k() const;
168 
169  //- Return the phase-pressure'
170  // (derivative of phase-pressure w.r.t. phase-fraction)
171  virtual tmp<volScalarField> pPrime() const;
172 
173 
174  // Thermophysical transport
175 
176  //- Return the source term for the energy equation
177  virtual tmp<fvScalarMatrix> divq(volScalarField& he) const;
178 };
179 
180 
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 
183 } // End namespace Foam
184 
185 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186 
187 #ifdef NoRepository
188  #include "StationaryPhaseModel.C"
189 #endif
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 #endif
194 
195 // ************************************************************************* //
Generic GeometricField class.
Class which represents a stationary (and therefore probably solid) phase. Generates,...
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
virtual tmp< volVectorField > U() const
Return the velocity.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
virtual const autoPtr< volScalarField > & divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
virtual tmp< surfaceScalarField > alphaPhi() const
Return the volumetric flux of the phase.
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
virtual const autoPtr< surfaceVectorField > & Uf() const
Return the face velocity.
virtual surfaceVectorField & UfRef()
Access the face velocity.
virtual ~StationaryPhaseModel()
Destructor.
StationaryPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux.
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
virtual volVectorField & URef()
Access the velocity.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
virtual bool stationary() const
Return whether the phase is stationary.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
virtual tmp< scalarField > kappaEff(const label patchi) const
Return the effective thermal conductivity on a patch.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
label patchi
Namespace for OpenFOAM.
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
thermo he()