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-2020 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 
45  (
46  option,
47  jouleHeatingSource,
48  dictionary
49  );
50 }
51 }
52 
53 const Foam::word Foam::fv::jouleHeatingSource::sigmaName(typeName + ":sigma");
54 
55 
56 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
57 
59 Foam::fv::jouleHeatingSource::transformSigma
60 (
61  const volVectorField& sigmaLocal
62 ) const
63 {
64  auto tsigma = tmp<volSymmTensorField>::New
65  (
66  IOobject
67  (
68  sigmaName,
69  mesh_.time().timeName(),
70  mesh_,
73  false
74  ),
75  mesh_,
76  dimensionedSymmTensor(sigmaLocal.dimensions(), Zero),
77  zeroGradientFvPatchField<symmTensor>::typeName
78  );
79  auto& sigma = tsigma.ref();
80 
81  // This check should be unnecessary
82  if (!csysPtr_)
83  {
85  << "Coordinate system undefined"
86  << abort(FatalError);
87  }
88 
89  const auto& csys = *csysPtr_;
90 
91  if (csys.uniform())
92  {
93  sigma.primitiveFieldRef() =
94  csys.transformPrincipal(sigmaLocal);
95  }
96  else
97  {
98  sigma.primitiveFieldRef() =
99  csys.transformPrincipal(mesh_.cellCentres(), sigmaLocal);
100  }
101 
102  sigma.correctBoundaryConditions();
103 
104  return tsigma;
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
109 
110 Foam::fv::jouleHeatingSource::jouleHeatingSource
111 (
112  const word& sourceName,
113  const word& modelType,
114  const dictionary& dict,
115  const fvMesh& mesh
116 )
117 :
118  option(sourceName, modelType, dict, mesh),
119  TName_("T"),
120  V_
121  (
122  IOobject
123  (
124  typeName + ":V",
125  mesh.time().timeName(),
126  mesh,
129  ),
130  mesh
131  ),
132  anisotropicElectricalConductivity_(false),
133  scalarSigmaVsTPtr_(nullptr),
134  vectorSigmaVsTPtr_(nullptr),
135  csysPtr_(nullptr),
136  curTimeIndex_(-1)
137 {
138  // Set the field name to that of the energy field from which the temperature
139  // is obtained
140 
141  const basicThermo& thermo =
142  mesh_.lookupObject<basicThermo>(basicThermo::dictName);
143 
144  fieldNames_.setSize(1, thermo.he().name());
145 
146  applied_.setSize(fieldNames_.size(), false);
147 
148  read(dict);
149 }
150 
151 
152 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
153 
155 (
156  const volScalarField& rho,
157  fvMatrix<scalar>& eqn,
158  const label fieldi
159 )
160 {
161  DebugInfo<< name() << ": applying source to " << eqn.psi().name() << endl;
162 
163  if (curTimeIndex_ != mesh_.time().timeIndex())
164  {
165  if (anisotropicElectricalConductivity_)
166  {
167  // Update sigma as a function of T if required
168  const volVectorField& sigmaLocal = updateSigma(vectorSigmaVsTPtr_);
169 
170  tmp<volSymmTensorField> sigma = transformSigma(sigmaLocal);
171 
172  // Solve the electrical potential equation
174  VEqn.relax();
175  VEqn.solve();
176  }
177  else
178  {
179  // Update sigma as a function of T if required
180  const volScalarField& sigma = updateSigma(scalarSigmaVsTPtr_);
181 
182  // Solve the electrical potential equation
184  VEqn.relax();
185  VEqn.solve();
186  }
187 
188  curTimeIndex_ = mesh_.time().timeIndex();
189  }
190 
191  // Add the Joule heating contribution
192 
193  const volVectorField gradV(fvc::grad(V_));
194  if (anisotropicElectricalConductivity_)
195  {
196  const volVectorField& sigmaLocal =
197  mesh_.lookupObject<volVectorField>(sigmaName);
198 
199  tmp<volSymmTensorField> sigma = transformSigma(sigmaLocal);
200 
201  eqn += (sigma & gradV) & gradV;
202  }
203  else
204  {
205  const volScalarField& sigma =
206  mesh_.lookupObject<volScalarField>(sigmaName);
207 
208  eqn += (sigma*gradV) & gradV;
209  }
210 }
211 
212 
214 {
215  if (option::read(dict))
216  {
217  coeffs_.readIfPresent("T", TName_);
218 
219  anisotropicElectricalConductivity_ =
220  coeffs_.get<bool>("anisotropicElectricalConductivity");
221 
222  if (anisotropicElectricalConductivity_)
223  {
224  Info<< " Using vector electrical conductivity" << endl;
225 
226  initialiseSigma(coeffs_, vectorSigmaVsTPtr_);
227 
228  csysPtr_ =
230  (
231  mesh_,
232  coeffs_,
233  coordinateSystem::typeName_()
234  );
235  }
236  else
237  {
238  Info<< " Using scalar electrical conductivity" << endl;
239 
240  initialiseSigma(coeffs_, scalarSigmaVsTPtr_);
241 
242  csysPtr_.clear(); // Do not need coordinate system
243  }
244 
245  return true;
246  }
247 
248  return false;
249 }
250 
251 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
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:62
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:59
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:107
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:785
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:350
zeroGradientFvPatchField.H
rho
rho
Definition: readInitialConditions.H:88
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:285
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:82
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::fv::option
Finite volume options abstract base class. Provides a base set of controls, e.g.:
Definition: fvOption.H:69
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
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:155
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:560
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:82
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionary.H:458
dict
dictionary dict
Definition: searchingEngine.H:14
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:121
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:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
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:60
fv
labelList fv(nPoints)
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:175
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:354
Foam::fv::jouleHeatingSource::read
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: jouleHeatingSource.C:213
jouleHeatingSource.H
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:76
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:301
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:248
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:123
Foam::IOobject::MUST_READ
Definition: IOobject.H:120