jouleHeatingSource.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2016-2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "jouleHeatingSource.H"
29#include "fvMatrices.H"
30#include "fvmLaplacian.H"
31#include "fvcGrad.H"
33#include "basicThermo.H"
35
36// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
37
38namespace Foam
39{
40namespace fv
41{
44}
45}
46
47const Foam::word Foam::fv::jouleHeatingSource::sigmaName(typeName + ":sigma");
48
49
50// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51
53Foam::fv::jouleHeatingSource::transformSigma
54(
55 const volVectorField& sigmaLocal
56) const
57{
59 (
60 IOobject
61 (
62 sigmaName,
64 mesh_,
67 false
68 ),
69 mesh_,
70 dimensionedSymmTensor(sigmaLocal.dimensions(), Zero),
72 );
73 auto& sigma = tsigma.ref();
74
75 // This check should be unnecessary
76 if (!csysPtr_)
77 {
79 << "Coordinate system undefined"
80 << abort(FatalError);
81 }
82
83 const auto& csys = *csysPtr_;
84
85 if (csys.uniform())
86 {
87 sigma.primitiveFieldRef() =
88 csys.transformPrincipal(sigmaLocal);
89 }
90 else
91 {
92 sigma.primitiveFieldRef() =
93 csys.transformPrincipal(mesh_.cellCentres(), sigmaLocal);
94 }
95
96 sigma.correctBoundaryConditions();
97
98 return tsigma;
99}
100
101
102// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
103
105(
106 const word& sourceName,
107 const word& modelType,
108 const dictionary& dict,
109 const fvMesh& mesh
110)
111:
112 fv::option(sourceName, modelType, dict, mesh),
113 TName_("T"),
114 V_
115 (
117 (
118 typeName + ":V",
119 mesh.time().timeName(),
120 mesh,
121 IOobject::MUST_READ,
122 IOobject::AUTO_WRITE
123 ),
124 mesh
125 ),
126 anisotropicElectricalConductivity_(false),
127 scalarSigmaVsTPtr_(nullptr),
128 vectorSigmaVsTPtr_(nullptr),
129 csysPtr_(nullptr),
130 curTimeIndex_(-1)
131{
132 // Set the field name to that of the energy
133 // field from which the temperature is obtained
134
136
137 fieldNames_.resize(1, thermo.he().name());
138
140
141 read(dict);
142}
143
144
145// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
146
148(
149 const volScalarField& rho,
150 fvMatrix<scalar>& eqn,
151 const label fieldi
152)
153{
154 DebugInfo<< name() << ": applying source to " << eqn.psi().name() << endl;
155
156 if (curTimeIndex_ != mesh_.time().timeIndex())
157 {
158 if (anisotropicElectricalConductivity_)
159 {
160 // Update sigma as a function of T if required
161 const volVectorField& sigmaLocal = updateSigma(vectorSigmaVsTPtr_);
162
163 tmp<volSymmTensorField> sigma = transformSigma(sigmaLocal);
164
165 // Solve the electrical potential equation
166 fvScalarMatrix VEqn(fvm::laplacian(sigma, V_));
167 VEqn.relax();
168 VEqn.solve();
169 }
170 else
171 {
172 // Update sigma as a function of T if required
173 const volScalarField& sigma = updateSigma(scalarSigmaVsTPtr_);
174
175 // Solve the electrical potential equation
176 fvScalarMatrix VEqn(fvm::laplacian(sigma, V_));
177 VEqn.relax();
178 VEqn.solve();
179 }
180
181 curTimeIndex_ = mesh_.time().timeIndex();
182 }
183
184 // Add the Joule heating contribution
185
186 const volVectorField gradV(fvc::grad(V_));
187 if (anisotropicElectricalConductivity_)
188 {
189 const auto& sigmaLocal = mesh_.lookupObject<volVectorField>(sigmaName);
190
191 tmp<volSymmTensorField> sigma = transformSigma(sigmaLocal);
192
193 eqn += (sigma & gradV) & gradV;
194 }
195 else
196 {
197 const auto& sigma = mesh_.lookupObject<volScalarField>(sigmaName);
198
199 eqn += (sigma*gradV) & gradV;
200 }
201}
202
203
205{
207 {
208 coeffs_.readIfPresent("T", TName_);
209
210 anisotropicElectricalConductivity_ =
211 coeffs_.get<bool>("anisotropicElectricalConductivity");
212
213 if (anisotropicElectricalConductivity_)
214 {
215 Info<< " Using vector electrical conductivity" << endl;
216
217 initialiseSigma(coeffs_, vectorSigmaVsTPtr_);
218
219 csysPtr_ =
221 (
222 mesh_,
223 coeffs_,
224 coordinateSystem::typeName_()
225 );
226 }
227 else
228 {
229 Info<< " Using scalar electrical conductivity" << endl;
230
231 initialiseSigma(coeffs_, scalarSigmaVsTPtr_);
232
233 csysPtr_.clear(); // Do not need coordinate system
234 }
235
236 return true;
237 }
238
239 return false;
240}
241
242// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
virtual bool read()
Re-read model coefficients if they have changed.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:66
static const word dictName
Definition: basicThermo.H:256
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Evolves an electrical potential equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1092
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:412
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Evolves an electrical potential equation.
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to compressible momentum equation.
virtual bool read(const dictionary &dict)
Read source dictionary.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:127
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: fvOption.H:148
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:48
const Type & lookupObject(const word &name, const bool recursive=false) const
const vectorField & cellCentres() const
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Calculate the gradient of the given field.
Calculate the matrix for the laplacian of the field.
word timeName
Definition: getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
labelList fv(nPoints)
dictionary dict
static const char *const typeName
The type name used in ensight case files.