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) 2019-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 "fam.H"
30#include "faScalarMatrix.H"
32
33// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace fa
38{
41}
42}
43
44
45// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46
48(
49 const word& sourceName,
50 const word& modelType,
51 const dictionary& dict,
52 const fvPatch& patch
53)
54:
55 fa::faceSetOption(sourceName, modelType, dict, patch),
56 TName_(dict.getOrDefault<word>("T", "T")),
57 V_
58 (
60 (
61 typeName + ":V_" + regionName_,
62 mesh().time().timeName(),
63 mesh(),
64 IOobject::MUST_READ,
65 IOobject::AUTO_WRITE
66 ),
67 regionMesh()
68 ),
69 scalarSigmaVsTPtr_(nullptr),
70 tensorSigmaVsTPtr_(nullptr),
71 curTimeIndex_(-1),
72 nIter_(1),
73 anisotropicElectricalConductivity_(false)
74{
75 fieldNames_.resize(1, TName_);
76
78
79 if (anisotropicElectricalConductivity_)
80 {
81 Info<< " Using tensor electrical conductivity" << endl;
82
83 initialiseSigma(coeffs_, tensorSigmaVsTPtr_);
84 }
85 else
86 {
87 Info<< " Using scalar electrical conductivity" << endl;
88
89 initialiseSigma(coeffs_, scalarSigmaVsTPtr_);
90 }
91
92 read(dict);
93}
94
95
96// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97
99(
100 const areaScalarField& h,
101 const areaScalarField& rho,
102 faMatrix<scalar>& eqn,
103 const label fieldi
104)
105{
106 if (isActive())
107 {
108 DebugInfo<< name() << ": applying source to " << eqn.psi().name()
109 << endl;
110
111 if (curTimeIndex_ != mesh().time().timeIndex())
112 {
113 for (label i = 0; i < nIter_; ++i)
114 {
115 if (anisotropicElectricalConductivity_)
116 {
117 // Update sigma as a function of T if required
118 const areaTensorField& sigma =
119 updateSigma(tensorSigmaVsTPtr_);
120
121 // Solve the electrical potential equation
122 faScalarMatrix VEqn(fam::laplacian(h*sigma, V_));
123 VEqn.relax();
124 VEqn.solve();
125 }
126 else
127 {
128 // Update sigma as a function of T if required
129 const areaScalarField& sigma =
130 updateSigma(scalarSigmaVsTPtr_);
131
132 // Solve the electrical potential equation
133 faScalarMatrix VEqn(fam::laplacian(h*sigma, V_));
134 VEqn.relax();
135 VEqn.solve();
136 }
137 }
138
139 curTimeIndex_ = mesh().time().timeIndex();
140 }
141
142 // Add the Joule heating contribution
143 areaVectorField gradV("gradV", fac::grad(V_));
144
145 if (anisotropicElectricalConductivity_)
146 {
147 const auto& sigma =
148 mesh_.lookupObject<areaTensorField>
149 (
150 typeName + ":sigma_" + regionName_
151 );
152
153 eqn += (h*sigma & gradV) & gradV;
154 }
155 else
156 {
157 const auto& sigma =
158 mesh_.lookupObject<areaScalarField>
159 (
160 typeName + ":sigma_" + regionName_
161 );
162
163 eqn += (h*sigma*gradV) & gradV;
164
165 if (mesh().time().outputTime() && debug)
166 {
167 areaScalarField qgradV("gradVSource", (gradV & gradV));
168 qgradV.write();
169 }
170 }
171 }
172}
173
174
176{
178 {
179 dict.readIfPresent("T", TName_);
180
181 dict.readIfPresent("nIter", nIter_);
182
183 anisotropicElectricalConductivity_ =
184 dict.get<bool>("anisotropicElectricalConductivity");
185
186 return true;
187 }
188
189 return false;
190}
191
192
193// ************************************************************************* //
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.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A special matrix type and solver, designed for finite area solutions of scalar equations....
Definition: faMatrix.H:76
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: faMatrix.C:491
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: faMatrixSolve.C:55
const GeometricField< Type, faPatchField, areaMesh > & psi() const
Definition: faMatrix.H:270
Intermediate abstract class for handling face-set options for the derived faOptions.
Evolves an electrical potential equation.
virtual void addSup(const areaScalarField &h, const areaScalarField &rho, faMatrix< 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 area options (i.e. faOption).
Definition: faOption.H:134
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: faOption.H:164
dictionary coeffs_
Dictionary containing source coefficients.
Definition: faOption.H:161
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: faOption.C:45
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
virtual bool write(const bool valid=true) const
Write using setting from DB.
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
Namespace of functions to calculate implicit derivatives returning a matrix. Time derivatives are cal...
word timeName
Definition: getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh > > grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:56
tmp< faMatrix< Type > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famLaplacian.C:49
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
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
label timeIndex
Definition: getTimeIndex.H:30
dictionary dict
volScalarField & h