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 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  option,
41  solidEquilibriumEnergySource,
42  dictionary
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
51 Foam::fv::solidEquilibriumEnergySource::alpha() const
52 {
53  const word alphaName = IOobject::groupName("alpha", phaseName_);
54 
55  if (!mesh_.foundObject<volScalarField>(alphaName))
56  {
57  volScalarField* alphaPtr =
58  new volScalarField
59  (
60  IOobject
61  (
62  alphaName,
63  mesh_.time().constant(),
64  mesh_,
67  ),
68  mesh_
69  );
70 
71  alphaPtr->store();
72  }
73 
74  return mesh_.lookupObject<volScalarField>(alphaName);
75 }
76 
77 
78 const Foam::solidThermo& Foam::fv::solidEquilibriumEnergySource::thermo() const
79 {
80  const word thermoName =
82 
83  if (!mesh_.foundObject<solidThermo>(thermoName))
84  {
85  solidThermo* thermoPtr = solidThermo::New(mesh_, phaseName_).ptr();
86 
87  thermoPtr->store();
88  }
89 
90  return mesh_.lookupObject<solidThermo>(thermoName);
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
97 (
98  const word& name,
99  const word& modelType,
100  const dictionary& dict,
101  const fvMesh& mesh
102 )
103 :
104  option(name, modelType, dict, mesh),
105  phaseName_(dict.lookup<word>("phase"))
106 {
107  read(dict);
108  alpha();
109  thermo();
110 }
111 
112 
113 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
114 
116 {}
117 
118 
119 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
120 
122 (
123  const volScalarField& rho,
124  fvMatrix<scalar>& eqn,
125  const label fieldi
126 )
127 {
128  const volScalarField alphahe(thermo().alphahe());
129 
130  const volScalarField& A = this->alpha();
131  const volScalarField B(1 - A);
132 
133  eqn -=
134  A/B*fvm::ddt(thermo().rho(), eqn.psi());
135  - 1/B*fvm::laplacian
136  (
137  A*alphahe,
138  eqn.psi(),
139  "laplacian(" + alphahe.name() + "," + eqn.psi().name() + ")"
140  );
141 }
142 
143 
145 (
146  const volScalarField& alpha,
147  const volScalarField& rho,
148  fvMatrix<scalar>& eqn,
149  const label fieldi
150 )
151 {
152  const volScalarField alphahe(alpha*thermo().alphahe());
153 
154  const volScalarField& A = this->alpha();
155  const volScalarField B(1 - A);
156 
157  eqn -=
158  A/B*fvm::ddt(alpha, thermo().rho(), eqn.psi());
159  - 1/B*fvm::laplacian
160  (
161  A*alphahe,
162  eqn.psi(),
163  "laplacian(" + alphahe.name() + "," + eqn.psi().name() + ")"
164  );
165 }
166 
167 
169 {
170  if (option::read(dict))
171  {
172  fieldNames_ = wordList(1, coeffs_.lookup<word>("field"));
173 
174  applied_.setSize(fieldNames_.size(), false);
175 
176  return true;
177  }
178  else
179  {
180  return false;
181  }
182 }
183 
184 
185 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
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 word & name() const
Return name.
Definition: IOobject.H:303
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
virtual bool read(const dictionary &dict)
Read dictionary.
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:282
static autoPtr< solidThermo > New(const fvMesh &, const word &phaseName=word::null)
Return a pointer to a new solidThermo created from.
Definition: solidThermo.C:92
virtual void addSup(const volScalarField &, fvMatrix< scalar > &, const label)
Explicit and implicit sources for compressible equations.
Calculate the matrix for the laplacian of the field.
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
rhoReactionThermo & thermo
Definition: createFields.H:28
Macros for easy insertion into run-time selection tables.
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
bool read(const char *, int32_t &)
Definition: int32IO.C:85
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
static word groupName(Name name, const word &group)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calculate the matrix for the first temporal derivative.
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvOptionIO.C:53
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
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:48
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:78
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
solidEquilibriumEnergySource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812
Finite volume options abstract base class. Provides a base set of controls, e.g.: ...
Definition: fvOption.H:66