phaseChange.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) 2021-2024 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 "phaseChange.H"
27 #include "fluidThermo.H"
28 #include "multicomponentThermo.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
43 
45 {
46  for (label i = 0; i < 2; ++ i)
47  {
48  if (isA<fluidThermo>(thermos()[i]))
49  {
50  return refCast<const fluidThermo>(thermos()[i]).p();
51  }
52  }
53 
54  return mesh().lookupObject<volScalarField>("p");
55 }
56 
57 
59 (
61 )
62 {
65  (
66  tvif().name(),
67  tvif().mesh(),
68  tvif().dimensions(),
70  );
71 
72  tvf->internalFieldRef() = tvif();
73  tvf->correctBoundaryConditions();
74 
75  tvif.clear();
76 
77  return tvf;
78 }
79 
80 
82 (
83  const tmp<volScalarField>& tvf
84 )
85 {
87 
88  tvf.clear();
89 
90  return tvif;
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
97 (
98  const word& name,
99  const word& modelType,
100  const fvMesh& mesh,
101  const dictionary& dict,
102  const Pair<bool>& fluidThermosRequired,
103  const Pair<bool>& specieThermosRequired
104 )
105 :
106  massTransfer(name, modelType, mesh, dict),
107  thermos_(mesh, phaseNames()),
108  fluidThermos_(thermos_),
109  specieThermos_(thermos_),
110  heNames_(thermos_.first().he().name(), thermos_.second().he().name())
111 {
112  forAll(fluidThermos_.valid(), i)
113  {
114  if (!fluidThermos_.valid()[i] && fluidThermosRequired[i])
115  {
117  << "Model " << name << " of type " << modelType
118  << " requires a fluid thermo for phase "
119  << phaseNames()[i] << exit(FatalError);
120  }
121  }
122 
123  forAll(specieThermos_.valid(), i)
124  {
125  if (!specieThermos_.valid()[i] && specieThermosRequired[i])
126  {
128  << "Model " << name << " of type " << modelType
129  << " requires a multicomponent thermo for phase "
130  << phaseNames()[i] << exit(FatalError);
131  }
132  }
133 }
134 
135 
136 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
137 
140 {
141  const volScalarField::Internal mDot(this->mDot());
142 
143  return pos0(mDot)*thermos().first().T() + neg(mDot)*thermos().second().T();
144 }
145 
146 
149 {
150  const volScalarField& kappa1 = thermos().first().kappa();
151  const volScalarField& kappa2 = thermos().second().kappa();
152 
153  return vfToVif(kappa2/(kappa1 + kappa2));
154 }
155 
156 
157 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
T & first()
Return the first element of the list.
Definition: UListI.H:114
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
dimensioned< Type > T() const
Return transpose.
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:53
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:47
Base class for mass transfers between phases.
Definition: massTransfer.H:55
const Pair< word > & phaseNames() const
Return the names of the phases.
Definition: massTransferI.H:52
Base class for phase change models.
Definition: phaseChange.H:59
phaseChange(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict, const Pair< bool > &fluidThermosRequired, const Pair< bool > &specieThermosRequired)
Construct from explicit source name and mesh.
Definition: phaseChange.C:97
virtual tmp< DimensionedField< scalar, volMesh > > Lfraction() const
Return the fraction of the latent heat that is transferred into.
Definition: phaseChange.C:148
static tmp< DimensionedField< scalar, volMesh > > vfToVif(const tmp< volScalarField > &tvf)
Remove the boundary field from the given geometric field.
Definition: phaseChange.C:82
const volScalarField & p() const
Access the pressure field.
Definition: phaseChange.C:44
virtual tmp< DimensionedField< scalar, volMesh > > Tchange() const
Return the temperature at which the phases are considered to be.
Definition: phaseChange.C:139
const ThermoRefPair< basicThermo > & thermos() const
Return the thermo references.
Definition: phaseChangeI.H:31
static tmp< volScalarField > vifToVf(const tmp< DimensionedField< scalar, volMesh >> &tvif)
Add a boundary field to the given internal field.
Definition: phaseChange.C:59
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
A class for managing temporary objects.
Definition: tmp.H:55
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:205
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedScalar pos0(const dimensionedScalar &ds)
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
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
dimensionedScalar neg(const dimensionedScalar &ds)
error FatalError
thermo he()
labelList fv(nPoints)
dictionary dict