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 -------------------------------------------------------------------------------
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 const Foam::word Foam::fv::jouleHeatingSource::sigmaName(typeName + ":sigma");
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
53 Foam::fv::jouleHeatingSource::transformSigma
54 (
55  const volVectorField& sigmaLocal
56 ) const
57 {
58  auto tsigma = tmp<volSymmTensorField>::New
59  (
60  IOobject
61  (
62  sigmaName,
63  mesh_.time().timeName(),
64  mesh_,
67  false
68  ),
69  mesh_,
70  dimensionedSymmTensor(sigmaLocal.dimensions(), Zero),
71  zeroGradientFvPatchField<symmTensor>::typeName
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  (
116  IOobject
117  (
118  typeName + ":V",
119  mesh.time().timeName(),
120  mesh,
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 
135  const auto& thermo = mesh_.lookupObject<basicThermo>(basicThermo::dictName);
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
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
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 {
206  if (fv::option::read(dict))
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 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
basicThermo.H
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::dimensionedSymmTensor
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
Definition: dimensionedSymmTensor.H:50
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::fv::option::resetApplied
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:48
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::basicThermo
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:63
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
zeroGradientFvPatchField.H
rho
rho
Definition: readInitialConditions.H:88
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::fv::option
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:126
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Foam::fv::jouleHeatingSource::addSup
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to compressible momentum equation.
Definition: jouleHeatingSource.C:148
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1183
fvmLaplacian.H
Calculate the matrix for the laplacian of the field.
Foam::coordinateSystem::New
static autoPtr< coordinateSystem > New(word modelType, const objectRegistry &obr, const dictionary &dict)
Definition: coordinateSystemNew.C:84
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:399
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
jouleHeatingSource.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::fv::option::read
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvOptionIO.C:55
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
fv
labelList fv(nPoints)
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::fv::jouleHeatingSource::read
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: jouleHeatingSource.C:204
fvcGrad.H
Calculate the gradient of the given field.
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:319
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::fv::jouleHeatingSource::jouleHeatingSource
jouleHeatingSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
Definition: jouleHeatingSource.C:105
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:60