Loading...
Searching...
No Matches
electricPotential.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) 2021-2024 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 "electricPotential.H"
29#include "fvc.H"
30#include "fvm.H"
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace functionObjects
38{
41}
42}
43
44
45// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46
47Foam::volScalarField& Foam::functionObjects::electricPotential::getOrReadField
48(
49 const word& fieldName
50) const
51{
52 auto* ptr = mesh_.getObjectPtr<volScalarField>(fieldName);
53
54 if (!ptr)
55 {
56 ptr = new volScalarField
57 (
59 (
60 fieldName,
62 mesh_.thisDb(),
66 ),
67 mesh_
68 );
70 }
71
72 return *ptr;
73}
74
75
76Foam::tmp<Foam::volScalarField>
77Foam::functionObjects::electricPotential::sigma() const
78{
79 const IOobject sigmaIO
80 (
81 mesh_.thisDb().newIOobject(IOobject::scopedName(typeName, "sigma"))
82 );
83
84 if (phases_.size())
85 {
86 tmp<volScalarField> tsigma = phases_[0]*sigmas_[0];
87
88 for (label i = 1; i < phases_.size(); ++i)
89 {
90 tsigma.ref() += phases_[i]*sigmas_[i];
91 }
92
94 (
95 sigmaIO,
96 tsigma,
98 );
99 }
100
102 (
103 sigmaIO,
104 mesh_,
105 sigma_,
107 );
108}
109
110
111Foam::tmp<Foam::volScalarField>
112Foam::functionObjects::electricPotential::epsilonm() const
113{
114 // Vacuum permittivity (aka the electric constant) [A^2 s^4/(kg m^3)]
116 (
118 8.8541878128e-12 // CODATA value
119 );
120
121 const IOobject epsilonrIO
122 (
123 mesh_.thisDb().newIOobject(IOobject::scopedName(typeName, "epsilonr"))
124 );
125
126 if (phases_.size())
127 {
128 tmp<volScalarField> tepsilonr = phases_[0]*epsilonrs_[0];
129
130 for (label i = 1; i < phases_.size(); ++i)
131 {
132 tepsilonr.ref() += phases_[i]*epsilonrs_[i];
133 }
134
136 (
137 epsilonrIO,
138 epsilon0*tepsilonr,
140 );
141 }
142
144 (
145 epsilonrIO,
146 mesh_,
147 epsilon0*epsilonr_,
149 );
150}
151
152
153// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
154
155Foam::functionObjects::electricPotential::electricPotential
156(
157 const word& name,
158 const Time& runTime,
159 const dictionary& dict
160)
161:
163 phasesDict_(dict.subOrEmptyDict("phases")),
164 phaseNames_(),
165 phases_(),
166 sigmas_(),
167 sigma_
168 (
170 (
172 dict.getCheckOrDefault<scalar>
173 (
174 "sigma",
175 scalar(1),
176 scalarMinMax::ge(SMALL)
177 )
178 )
179 ),
180 epsilonrs_(),
181 epsilonr_
182 (
184 (
185 dimless,
186 dict.getCheckOrDefault<scalar>
187 (
188 "epsilonr",
189 scalar(1),
190 scalarMinMax::ge(SMALL)
191 )
192 )
193 ),
194 Vname_
195 (
196 dict.getOrDefault<word>
197 (
198 "V",
199 IOobject::scopedName(typeName, "V")
200 )
201 ),
202 Ename_
203 (
204 dict.getOrDefault<word>
205 (
206 "E",
207 IOobject::scopedName(typeName, "E")
208 )
209 ),
210 fvOptions_(mesh_),
211 tol_(1),
212 nCorr_(1),
213 writeDerivedFields_(false),
214 electricField_(false)
215{
216 read(dict);
217
218 // Force creation of transported field so any BCs using it can
219 // look it up
220 volScalarField& eV = getOrReadField(Vname_);
222
223 if (electricField_)
224 {
225 auto* ptr = getObjectPtr<volVectorField>(Ename_);
226
227 if (!ptr)
228 {
229 ptr = new volVectorField
230 (
232 (
233 Ename_,
234 mesh_.time().timeName(),
235 mesh_.thisDb(),
239 ),
240 -fvc::grad(eV)
241 );
244 }
245}
246
247
248// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
249
251{
253 {
254 return false;
255 }
256
257 Log << type() << " read: " << name() << endl;
258
259 dict.readIfPresent("sigma", sigma_);
260 dict.readIfPresent("epsilonr", epsilonr_);
261 dict.readIfPresent("nCorr", nCorr_);
262 dict.readIfPresent("tolerance", tol_);
263 dict.readIfPresent("writeDerivedFields", writeDerivedFields_);
264 dict.readIfPresent("electricField", electricField_);
265
266 // If flow is multiphase
267 if (!phasesDict_.empty())
268 {
269 phaseNames_.setSize(phasesDict_.size());
270 phases_.setSize(phasesDict_.size());
271 sigmas_.setSize(phasesDict_.size());
272
273 if (writeDerivedFields_)
274 {
275 epsilonrs_.setSize(phasesDict_.size());
276 }
277
278 label phasei = 0;
279 for (const entry& dEntry : phasesDict_)
280 {
281 const word& key = dEntry.keyword();
282
283 if (!dEntry.isDict())
284 {
285 FatalIOErrorInFunction(phasesDict_)
286 << "Entry " << key << " is not a dictionary" << nl
287 << exit(FatalIOError);
288 }
289
290 const dictionary& subDict = dEntry.dict();
291
292 phaseNames_[phasei] = key;
293
294 sigmas_.set
295 (
296 phasei,
298 (
300 subDict.getCheck<scalar>
301 (
302 "sigma",
303 scalarMinMax::ge(SMALL)
304 )
305 )
306 );
307
308 if (writeDerivedFields_)
309 {
310 epsilonrs_.set
311 (
312 phasei,
314 (
315 dimless,
316 subDict.getCheck<scalar>
317 (
318 "epsilonr",
319 scalarMinMax::ge(SMALL)
320 )
321 )
322 );
323 }
324
325 ++phasei;
326 }
327
328 forAll(phaseNames_, i)
329 {
330 phases_.set
331 (
332 i,
333 mesh_.getObjectPtr<volScalarField>(phaseNames_[i])
334 );
335 }
336 }
337
338 if (const dictionary* dictptr = dict.findDict("fvOptions"))
339 {
340 fvOptions_.reset(*dictptr);
341 }
342
343 return true;
344}
345
346
348{
349 Log << type() << " execute: " << name() << endl;
350
351 // Convergence monitor parameters
352 bool converged = false;
353 label iter = 0;
354
355 tmp<volScalarField> tsigma = this->sigma();
356 const auto& sigma = tsigma();
357
358 volScalarField& eV = getOrReadField(Vname_);
359
360 for (int i = 1; i <= nCorr_; ++i)
361 {
362 fvScalarMatrix eVEqn
363 (
364 - fvm::laplacian(sigma, eV)
365 );
366
367 eVEqn.relax();
368
369 fvOptions_.constrain(eVEqn);
370
371 ++iter;
372 converged = (eVEqn.solve().initialResidual() < tol_);
373 if (converged) break;
374 }
375
376 if (electricField_)
377 {
378 auto& E = lookupObjectRef<volVectorField>(Ename_);
379 E == -fvc::grad(eV);
380 }
381
382 if (converged)
383 {
384 Log << type() << ": " << name() << ": "
385 << eV.name() << " is converged." << nl
386 << tab << "initial-residual tolerance: " << tol_ << nl
387 << tab << "outer iteration: " << iter << nl;
388 }
390 Log << endl;
391
392 return true;
393}
394
395
397{
398 Log << type() << " write: " << name() << nl
399 << tab << Vname_
400 << endl;
401
402 volScalarField& eV = getOrReadField(Vname_);
403
404 if (electricField_)
405 {
406 const auto& E = lookupObject<volVectorField>(Ename_);
407
408 Log << tab << E.name() << endl;
409
410 E.write();
411 }
412
413 if (writeDerivedFields_)
414 {
415 // Write the current density field
416 tmp<volScalarField> tsigma = this->sigma();
417
418 auto eJ = volVectorField::New
419 (
422 -tsigma*fvc::grad(eV),
424 );
425
426 Log << tab << eJ().name() << endl;
427
428 eJ->write();
429
430
431 // Write the free-charge density field
432 tmp<volScalarField> tepsilonm = this->epsilonm();
433
434 auto erho = volScalarField::New
435 (
438 fvc::div(tepsilonm*(-fvc::grad(eV))),
440 );
441
442 Log << tab << erho().name() << endl;
443
444 erho->write();
445 }
446
447 eV.write();
448
449 return true;
450}
451
452
453// ************************************************************************* //
#define Log
Definition PDRblock.C:28
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< vector >::calculatedType())
void correctBoundaryConditions()
Correct boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
static MinMax< scalar > ge(const scalar &minVal)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T getCheck(const word &keyword, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T with additional checking FatalIOError if not found, or if the number of tokens is...
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
Abstract base-class for Time/database function objects.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
word scopedName(const word &name) const
Return a scoped (prefixed) name.
Computes the steady-state equation of charge conservation to obtain the electric potential by strictl...
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function-object results.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
ObjectType * getObjectPtr(const word &fieldName) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
ObjectType & lookupObjectRef(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition fvMatrix.C:1094
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition fvMesh.H:376
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
bool store()
Register object with its registry and transfer ownership to the registry.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
engineTime & runTime
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
auto & name
label phasei
Definition pEqn.H:27
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
constexpr auto key(const Type &t) noexcept
Helper function to return the enum value.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
fvMatrix< scalar > fvScalarMatrix
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimCurrent(0, 0, 0, 0, 0, 1, 0)
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar pow4(const dimensionedScalar &ds)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
dictionary dict
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299