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-2024 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 
38 namespace Foam
39 {
40 namespace fv
41 {
42  defineTypeNameAndDebug(jouleHeatingSource, 0);
43  addToRunTimeSelectionTable(option, jouleHeatingSource, dictionary);
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
51 Foam::fv::jouleHeatingSource::transformSigma
52 (
53  const volVectorField& sigmaLocal
54 ) const
55 {
56  auto tsigma = volSymmTensorField::New
57  (
58  IOobject::scopedName(typeName, "sigma"),
60  mesh_,
61  dimensionedSymmTensor(sigmaLocal.dimensions(), Zero),
63  );
64  auto& sigma = tsigma.ref();
65 
66  // This check should be unnecessary
67  if (!csysPtr_)
68  {
70  << "Coordinate system undefined"
71  << abort(FatalError);
72  }
73 
74  const auto& csys = *csysPtr_;
75 
76  if (csys.uniform())
77  {
78  sigma.primitiveFieldRef() =
79  csys.transformPrincipal(sigmaLocal);
80  }
81  else
82  {
83  sigma.primitiveFieldRef() =
84  csys.transformPrincipal(mesh_.cellCentres(), sigmaLocal);
85  }
86 
87  sigma.correctBoundaryConditions();
88 
89  return tsigma;
90 }
91 
92 
93 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
94 
96 (
97  const word& sourceName,
98  const word& modelType,
99  const dictionary& dict,
100  const fvMesh& mesh
101 )
102 :
103  fv::option(sourceName, modelType, dict, mesh),
104  TName_("T"),
105  V_
106  (
107  IOobject
108  (
109  IOobject::scopedName(typeName, "V"),
110  mesh.time().timeName(),
111  mesh,
112  IOobject::MUST_READ,
113  IOobject::AUTO_WRITE,
114  IOobject::REGISTER
115  ),
116  mesh
117  ),
118  anisotropicElectricalConductivity_(false),
119  scalarSigmaVsTPtr_(nullptr),
120  vectorSigmaVsTPtr_(nullptr),
121  csysPtr_(nullptr),
122  curTimeIndex_(-1)
123 {
124  // Set the field name to that of the energy
125  // field from which the temperature is obtained
126 
128 
129  fieldNames_.resize(1, thermo.he().name());
130 
132 
134 }
135 
136 
137 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
138 
140 (
141  const volScalarField& rho,
142  fvMatrix<scalar>& eqn,
143  const label fieldi
144 )
145 {
146  DebugInfo<< name() << ": applying source to " << eqn.psi().name() << endl;
147 
148  if (curTimeIndex_ != mesh_.time().timeIndex())
149  {
150  if (anisotropicElectricalConductivity_)
151  {
152  // Update sigma as a function of T if required
153  const volVectorField& sigmaLocal = updateSigma(vectorSigmaVsTPtr_);
154 
155  tmp<volSymmTensorField> sigma = transformSigma(sigmaLocal);
156 
157  // Solve the electrical potential equation
159  VEqn.relax();
160  VEqn.solve();
161  }
162  else
163  {
164  // Update sigma as a function of T if required
165  const volScalarField& sigma = updateSigma(scalarSigmaVsTPtr_);
166 
167  // Solve the electrical potential equation
169  VEqn.relax();
170  VEqn.solve();
171  }
172 
173  curTimeIndex_ = mesh_.time().timeIndex();
174  }
175 
176  // Add the Joule heating contribution
177  const word sigmaName
178  (
179  IOobject::scopedName(typeName, "sigma")
180  );
181 
182  const volVectorField gradV(fvc::grad(V_));
183  if (anisotropicElectricalConductivity_)
184  {
185  const auto& sigmaLocal = mesh_.lookupObject<volVectorField>(sigmaName);
186 
187  tmp<volSymmTensorField> sigma = transformSigma(sigmaLocal);
188 
189  eqn += (sigma & gradV) & gradV;
190  }
191  else
192  {
193  const auto& sigma = mesh_.lookupObject<volScalarField>(sigmaName);
194 
195  eqn += (sigma*gradV) & gradV;
196  }
197 }
198 
199 
201 {
202  if (fv::option::read(dict))
203  {
204  coeffs_.readIfPresent("T", TName_);
205 
206  anisotropicElectricalConductivity_ =
207  coeffs_.get<bool>("anisotropicElectricalConductivity");
208 
209  if (anisotropicElectricalConductivity_)
210  {
211  Info<< " Using vector electrical conductivity" << endl;
212 
213  initialiseSigma(coeffs_, vectorSigmaVsTPtr_);
214 
215  csysPtr_ =
216  coordinateSystem::New
217  (
218  mesh_,
219  coeffs_,
220  coordinateSystem::typeName
221  );
222  }
223  else
224  {
225  Info<< " Using scalar electrical conductivity" << endl;
226 
227  initialiseSigma(coeffs_, scalarSigmaVsTPtr_);
228 
229  csysPtr_.clear(); // Do not need coordinate system
230  }
231 
232  return true;
233  }
234 
235  return false;
236 }
237 
238 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
virtual bool read(const dictionary &dict)
Read source dictionary.
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:59
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: fvOption.H:157
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:153
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
jouleHeatingSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Calculate the matrix for the laplacian of the field.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Macros for easy insertion into run-time selection tables.
word timeName
Definition: getTimeIndex.H:3
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to compressible momentum equation.
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Calculate the gradient of the given field.
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvOptionIO.C:48
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:64
const vectorField & cellCentres() const
errorManip< error > abort(error &err)
Definition: errorManip.H:139
#define DebugInfo
Report an information message using Foam::Info.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
static tmp< GeometricField< symmTensor, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< symmTensor >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
A special matrix type and solver, designed for finite volume solutions of scalar equations.
messageStream Info
Information stream (stdout output on master, null elsewhere)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Do not request registration (bool: false)
Namespace for OpenFOAM.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:123
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127