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-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 \*---------------------------------------------------------------------------*/
25 
27 #include "fvmDdt.H"
28 #include "fvmLaplacian.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
39  (
40  fvModel,
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 fvMesh& mesh,
111  const dictionary& dict
112 )
113 :
114  fvModel(name, modelType, mesh, dict),
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
152  (
153  "alphahe",
154  solidThermo().kappa()/solidThermo().Cpv()
155  );
156 
157  const volScalarField& A = solidAlpha();
158  const volScalarField B(1 - A);
159 
160  eqn -=
161  A/B*fvm::ddt(solidThermo().rho(), eqn.psi());
162  - 1/B*fvm::laplacian
163  (
164  A*alphahe,
165  eqn.psi(),
166  "laplacian(" + alphahe.name() + "," + eqn.psi().name() + ")"
167  );
168 }
169 
170 
172 (
173  const volScalarField& alpha,
174  const volScalarField& rho,
175  fvMatrix<scalar>& eqn,
176  const word& fieldName
177 ) const
178 {
179  const volScalarField alphahe
180  (
181  "alphahe",
182  solidThermo().kappa()/solidThermo().Cpv()
183  );
184 
185  const volScalarField& A = solidAlpha();
186  const volScalarField B(1 - A);
187 
188  eqn -=
189  A/B*fvm::ddt(alpha, solidThermo().rho(), eqn.psi());
190  - 1/B*fvm::laplacian
191  (
192  A*alphahe,
193  eqn.psi(),
194  "laplacian(" + alphahe.name() + "," + eqn.psi().name() + ")"
195  );
196 }
197 
198 
200 {
201  return true;
202 }
203 
204 
206 (
207  const polyTopoChangeMap&
208 )
209 {}
210 
211 
213 {}
214 
215 
217 (
218  const polyDistributionMap&
219 )
220 {}
221 
222 
224 {
225  if (fvModel::read(dict))
226  {
227  readCoeffs();
228  return true;
229  }
230  else
231  {
232  return false;
233  }
234 }
235 
236 
237 // ************************************************************************* //
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("B", Foam::dimless, 18.678)
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
const word & name() const
Return name.
Definition: IOobject.H:310
static word groupName(Name name, const word &group)
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:78
virtual const IOdictionary & properties() const =0
Properties dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:40
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:187
This fvModel adds the thermal inertia of a solid phase into the energy equation. It assumes that the ...
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual void addSup(const volScalarField &, fvMatrix< scalar > &, const word &fieldName) const
Explicit and implicit sources for compressible equations.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read dictionary.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
solidEquilibriumEnergySource(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:40
Base-class for solid thermodynamic properties.
Definition: solidThermo.H:56
static autoPtr< solidThermo > New(const fvMesh &, const word &phaseName=word::null)
Standard selection based on fvMesh.
Definition: solidThermo.C:88
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the laplacian of the field.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31