compressibleVoFphase.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) 2023-2025 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 "compressibleVoFphase.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
31 (
32  const word& name,
33  const fvMesh& mesh,
34  const volScalarField& T
35 )
36 :
37  VoFphase(name, mesh),
38  thermo_(nullptr),
39  vDot_
40  (
41  IOobject
42  (
43  IOobject::groupName("vDot", name),
44  mesh.time().name(),
45  mesh,
46  IOobject::READ_IF_PRESENT,
47  IOobject::AUTO_WRITE
48  ),
49  mesh,
51  )
52 {
53  {
55  Tp.write();
56  }
57 
58  thermo_ = rhoFluidThermo::New(mesh, name);
59  thermo_->validate(name, "e");
60 }
61 
62 
63 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
64 
66 {
68  return autoPtr<VoFphase>(nullptr);
69 }
70 
71 
73 (
74  const volScalarField& p,
75  const volScalarField& T
76 )
77 {
78  thermo_->he() = thermo_->he(p, T);
79  thermo_->correct();
80 }
81 
82 
83 // ************************************************************************* //
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
tmp< GeometricField< Type, GeoMesh, Field > > T() const
Return transpose (only if it is a tensor field)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word groupName(Name name, const word &group)
Single incompressible VoFphase derived from the phase-fraction. Used as part of the multiPhaseMixture...
Definition: VoFphase.H:53
const word & name() const
Definition: VoFphase.H:98
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
compressibleVoFphase(const word &name, const fvMesh &mesh, const volScalarField &T)
Construct from components.
void correct(const volScalarField &p, const volScalarField &T)
virtual autoPtr< VoFphase > clone() const
Return clone.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
virtual bool write(const bool write=true) const
Write using setting from DB.
static autoPtr< rhoFluidThermo > New(const fvMesh &, const word &phaseName=word::null)
Standard selection based on fvMesh.
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
const dimensionSet dimless
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTime
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
volScalarField & p