solidEqulibriumEnergySource.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(solidEqulibriumEnergySource, 0);
39  (
40  option,
41  solidEqulibriumEnergySource,
42  dictionary
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 const Foam::volScalarField& Foam::fv::solidEqulibriumEnergySource::alpha() const
51 {
52  const word alphaName = IOobject::groupName("alpha", phaseName_);
53 
54  if (!mesh_.foundObject<volScalarField>(alphaName))
55  {
56  volScalarField* alphaPtr =
57  new volScalarField
58  (
59  IOobject
60  (
61  alphaName,
62  mesh_.time().constant(),
63  mesh_,
66  ),
67  mesh_
68  );
69 
70  alphaPtr->store();
71  }
72 
73  return mesh_.lookupObject<volScalarField>(alphaName);
74 }
75 
76 
77 const Foam::solidThermo& Foam::fv::solidEqulibriumEnergySource::thermo() const
78 {
79  const word thermoName =
81 
82  if (!mesh_.foundObject<solidThermo>(thermoName))
83  {
84  solidThermo* thermoPtr = solidThermo::New(mesh_, phaseName_).ptr();
85 
86  thermoPtr->store();
87  }
88 
89  return mesh_.lookupObject<solidThermo>(thermoName);
90 }
91 
92 
93 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
94 
96 (
97  const word& name,
98  const word& modelType,
99  const dictionary& dict,
100  const fvMesh& mesh
101 )
102 :
103  option(name, modelType, dict, mesh),
104  phaseName_(dict.lookupType<word>("phase"))
105 {
106  read(dict);
107  alpha();
108  thermo();
109 }
110 
111 
112 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
113 
115 {}
116 
117 
118 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
119 
121 (
122  const volScalarField& rho,
123  fvMatrix<scalar>& eqn,
124  const label fieldi
125 )
126 {
127  const volScalarField alphahe(thermo().alphahe());
128 
129  const volScalarField& A = this->alpha();
130  const volScalarField B(1 - A);
131 
132  eqn -=
133  A/B*fvm::ddt(thermo().rho(), eqn.psi());
134  - 1/B*fvm::laplacian
135  (
136  A*alphahe,
137  eqn.psi(),
138  "laplacian(" + alphahe.name() + "," + eqn.psi().name() + ")"
139  );
140 }
141 
142 
144 (
145  const volScalarField& alpha,
146  const volScalarField& rho,
147  fvMatrix<scalar>& eqn,
148  const label fieldi
149 )
150 {
151  const volScalarField alphahe(alpha*thermo().alphahe());
152 
153  const volScalarField& A = this->alpha();
154  const volScalarField B(1 - A);
155 
156  eqn -=
157  A/B*fvm::ddt(alpha, thermo().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  if (option::read(dict))
170  {
171  fieldNames_ = wordList(1, coeffs_.lookupType<word>("field"));
172 
173  applied_.setSize(fieldNames_.size(), false);
174 
175  return true;
176  }
177  else
178  {
179  return false;
180  }
181 }
182 
183 
184 // ************************************************************************* //
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:295
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
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
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:52
T lookupType(const word &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
virtual bool read(const dictionary &dict)
Read dictionary.
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
solidEqulibriumEnergySource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:48
virtual void addSup(const volScalarField &, fvMatrix< scalar > &, const label)
Explicit and implicit sources for compressible equations.
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))
Namespace for OpenFOAM.
Finite volume options abstract base class. Provides a base set of controls, e.g.: ...
Definition: fvOption.H:66