solidEquilibriumEnergySource.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) 2019-2022 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 
27 #include "fvmDdt.H"
28 #include "fvmLaplacian.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
37  defineTypeNameAndDebug(solidEquilibriumEnergySource, 0);
39  (
40  fvModel,
41  solidEquilibriumEnergySource,
42  dictionary
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void Foam::fv::solidEquilibriumEnergySource::readCoeffs()
51 {
52  phaseName_ = coeffs().lookupOrDefault<word>("phase", word::null);
53 
54  solidPhaseName_ = coeffs().lookup<word>("solidPhase");
55 }
56 
57 
59 Foam::fv::solidEquilibriumEnergySource::solidAlpha() const
60 {
61  const word alphaName = IOobject::groupName("alpha", solidPhaseName_);
62 
63  if (!mesh().foundObject<volScalarField>(alphaName))
64  {
65  volScalarField* alphaPtr =
66  new volScalarField
67  (
68  IOobject
69  (
70  alphaName,
71  mesh().time().constant(),
72  mesh(),
75  ),
76  mesh()
77  );
78 
79  alphaPtr->store();
80  }
81 
82  return mesh().lookupObject<volScalarField>(alphaName);
83 }
84 
85 
86 const Foam::solidThermo&
87 Foam::fv::solidEquilibriumEnergySource::solidThermo() const
88 {
89  const word thermoName =
90  IOobject::groupName(physicalProperties::typeName, solidPhaseName_);
91 
92  if (!mesh().foundObject<Foam::solidThermo>(thermoName))
93  {
94  Foam::solidThermo* thermoPtr =
95  solidThermo::New(mesh(), solidPhaseName_).ptr();
96 
97  thermoPtr->properties().store();
98  }
99 
100  return mesh().lookupObject<Foam::solidThermo>(thermoName);
101 }
102 
103 
104 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
105 
107 (
108  const word& name,
109  const word& modelType,
110  const dictionary& dict,
111  const fvMesh& mesh
112 )
113 :
114  fvModel(name, modelType, dict, mesh),
115  phaseName_(word::null),
116  solidPhaseName_(word::null)
117 {
118  read(dict);
119  solidAlpha();
120  solidThermo();
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
125 
127 {}
128 
129 
130 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
131 
133 {
134  const basicThermo& thermo =
135  mesh().lookupObject<basicThermo>
136  (
137  IOobject::groupName(physicalProperties::typeName, phaseName_)
138  );
139 
140  return wordList(1, thermo.he().name());
141 }
142 
143 
145 (
146  const volScalarField& rho,
147  fvMatrix<scalar>& eqn,
148  const word& fieldName
149 ) const
150 {
151  const volScalarField alphahe(solidThermo().alphahe());
152 
153  const volScalarField& A = solidAlpha();
154  const volScalarField B(1 - A);
155 
156  eqn -=
157  A/B*fvm::ddt(solidThermo().rho(), eqn.psi());
158  - 1/B*fvm::laplacian
159  (
160  A*alphahe,
161  eqn.psi(),
162  "laplacian(" + alphahe.name() + "," + eqn.psi().name() + ")"
163  );
164 }
165 
166 
168 (
169  const volScalarField& alpha,
170  const volScalarField& rho,
171  fvMatrix<scalar>& eqn,
172  const word& fieldName
173 ) const
174 {
175  const volScalarField alphahe(alpha*solidThermo().alphahe());
176 
177  const volScalarField& A = solidAlpha();
178  const volScalarField B(1 - A);
179 
180  eqn -=
181  A/B*fvm::ddt(alpha, solidThermo().rho(), eqn.psi());
182  - 1/B*fvm::laplacian
183  (
184  A*alphahe,
185  eqn.psi(),
186  "laplacian(" + alphahe.name() + "," + eqn.psi().name() + ")"
187  );
188 }
189 
190 
192 {
193  return true;
194 }
195 
196 
198 (
199  const polyTopoChangeMap&
200 )
201 {}
202 
203 
205 {}
206 
207 
209 (
210  const polyDistributionMap&
211 )
212 {}
213 
214 
216 {
217  if (fvModel::read(dict))
218  {
219  readCoeffs();
220  return true;
221  }
222  else
223  {
224  return false;
225  }
226 }
227 
228 
229 // ************************************************************************* //
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
const word & name() const
Return name.
Definition: IOobject.H:315
fluidReactionThermo & thermo
Definition: createFields.H:28
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:77
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual bool read(const dictionary &dict)
Read dictionary.
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:290
static autoPtr< solidThermo > New(const fvMesh &, const word &phaseName=word::null)
Standard selection based on fvMesh.
Definition: solidThermo.C:112
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:164
Calculate the matrix for the laplacian of the field.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
Finite volume model abstract base class.
Definition: fvModel.H:57
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
virtual const IOdictionary & properties() const =0
Properties dictionary.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
virtual bool movePoints()
Update for mesh motion.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
static word groupName(Name name, const word &group)
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calculate the matrix for the first temporal derivative.
static const word null
An empty word.
Definition: word.H:77
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:40
Base-class for solid thermodynamic properties.
Definition: solidThermo.H:53
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
List< word > wordList
A List of words.
Definition: fileName.H:54
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
virtual void addSup(const volScalarField &, fvMatrix< scalar > &, const word &fieldName) const
Explicit and implicit sources for compressible equations.
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
solidEquilibriumEnergySource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Namespace for OpenFOAM.