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 -------------------------------------------------------------------------------
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 "electricPotential.H"
29 #include "fvc.H"
30 #include "fvm.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace functionObjects
38 {
39  defineTypeNameAndDebug(electricPotential, 0);
40  addToRunTimeSelectionTable(functionObject, electricPotential, dictionary);
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 Foam::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  (
58  IOobject
59  (
60  fieldName,
61  mesh_.time().timeName(),
62  mesh_.thisDb(),
66  ),
67  mesh_
68  );
69  regIOobject::store(ptr);
70  }
71 
72  return *ptr;
73 }
74 
75 
77 Foam::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 
112 Foam::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 
155 Foam::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  nCorr_(1),
212  writeDerivedFields_(false),
213  electricField_(false)
214 {
215  read(dict);
216 
217  // Force creation of transported field so any BCs using it can
218  // look it up
219  volScalarField& eV = getOrReadField(Vname_);
221 
222  if (electricField_)
223  {
224  auto* ptr = getObjectPtr<volVectorField>(Ename_);
225 
226  if (!ptr)
227  {
228  ptr = new volVectorField
229  (
230  IOobject
231  (
232  Ename_,
233  mesh_.time().timeName(),
234  mesh_.thisDb(),
238  ),
239  -fvc::grad(eV)
240  );
241  regIOobject::store(ptr);
242  }
243  }
244 }
245 
246 
247 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
248 
250 {
252  {
253  return false;
254  }
255 
256  Log << type() << " read: " << name() << endl;
257 
258  dict.readIfPresent("sigma", sigma_);
259  dict.readIfPresent("epsilonr", epsilonr_);
260  dict.readIfPresent("nCorr", nCorr_);
261  dict.readIfPresent("writeDerivedFields", writeDerivedFields_);
262  dict.readIfPresent("electricField", electricField_);
263 
264  // If flow is multiphase
265  if (!phasesDict_.empty())
266  {
267  phaseNames_.setSize(phasesDict_.size());
268  phases_.setSize(phasesDict_.size());
269  sigmas_.setSize(phasesDict_.size());
270 
271  if (writeDerivedFields_)
272  {
273  epsilonrs_.setSize(phasesDict_.size());
274  }
275 
276  label phasei = 0;
277  for (const entry& dEntry : phasesDict_)
278  {
279  const word& key = dEntry.keyword();
280 
281  if (!dEntry.isDict())
282  {
283  FatalIOErrorInFunction(phasesDict_)
284  << "Entry " << key << " is not a dictionary" << nl
285  << exit(FatalIOError);
286  }
287 
288  const dictionary& subDict = dEntry.dict();
289 
290  phaseNames_[phasei] = key;
291 
292  sigmas_.set
293  (
294  phasei,
296  (
298  subDict.getCheck<scalar>
299  (
300  "sigma",
301  scalarMinMax::ge(SMALL)
302  )
303  )
304  );
305 
306  if (writeDerivedFields_)
307  {
308  epsilonrs_.set
309  (
310  phasei,
312  (
313  dimless,
314  subDict.getCheck<scalar>
315  (
316  "epsilonr",
317  scalarMinMax::ge(SMALL)
318  )
319  )
320  );
321  }
322 
323  ++phasei;
324  }
325 
326  forAll(phaseNames_, i)
327  {
328  phases_.set
329  (
330  i,
331  mesh_.getObjectPtr<volScalarField>(phaseNames_[i])
332  );
333  }
334  }
335 
336  if (const dictionary* dictptr = dict.findDict("fvOptions"))
337  {
338  fvOptions_.reset(*dictptr);
339  }
340 
341  return true;
342 }
343 
344 
346 {
347  Log << type() << " execute: " << name() << endl;
348 
349  tmp<volScalarField> tsigma = this->sigma();
350  const auto& sigma = tsigma();
351 
352  volScalarField& eV = getOrReadField(Vname_);
353 
354  for (int i = 1; i <= nCorr_; ++i)
355  {
356  fvScalarMatrix eVEqn
357  (
358  - fvm::laplacian(sigma, eV)
359  );
360 
361  eVEqn.relax();
362 
363  fvOptions_.constrain(eVEqn);
364 
365  eVEqn.solve();
366  }
367 
368  if (electricField_)
369  {
370  auto& E = lookupObjectRef<volVectorField>(Ename_);
371  E == -fvc::grad(eV);
372  }
374  Log << endl;
375 
376  return true;
377 }
378 
379 
381 {
382  Log << type() << " write: " << name() << nl
383  << tab << Vname_
384  << endl;
385 
386  volScalarField& eV = getOrReadField(Vname_);
387 
388  if (electricField_)
389  {
390  const auto& E = lookupObject<volVectorField>(Ename_);
391 
392  Log << tab << E.name() << endl;
393 
394  E.write();
395  }
396 
397  if (writeDerivedFields_)
398  {
399  // Write the current density field
400  tmp<volScalarField> tsigma = this->sigma();
401 
402  auto eJ = volVectorField::New
403  (
404  IOobject::scopedName(typeName, "J"),
406  -tsigma*fvc::grad(eV),
408  );
409 
410  Log << tab << eJ().name() << endl;
411 
412  eJ->write();
413 
414 
415  // Write the free-charge density field
416  tmp<volScalarField> tepsilonm = this->epsilonm();
417 
418  auto erho = volScalarField::New
419  (
420  IOobject::scopedName(typeName, "rho"),
422  fvc::div(tepsilonm*(-fvc::grad(eV))),
424  );
425 
426  Log << tab << erho().name() << endl;
427 
428  erho->write();
429  }
430 
431  eV.write();
432 
433  return true;
434 }
435 
436 
437 // ************************************************************************* //
virtual bool write()
Write the function object output.
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
defineTypeNameAndDebug(ObukhovLength, 0)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
label phasei
Definition: pEqn.H:27
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
bool store()
Register object with its registry and transfer ownership to the registry.
Definition: regIOobjectI.H:36
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
virtual bool read(const dictionary &dict)
Read the function object data.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
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...
static MinMax< scalar > ge(const scalar &minVal)
A semi-infinite range from minVal to the type max.
Definition: MinMaxI.H:24
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
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
virtual bool execute()
Calculate the function object.
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())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
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
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1095
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:637
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:103
const dimensionSet dimCurrent(0, 0, 0, 0, 0, 1, 0)
Definition: dimensionSets.H:54
dimensionedScalar pow4(const dimensionedScalar &ds)
Nothing to be read.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:204
#define Log
Definition: PDRblock.C:28
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
void correctBoundaryConditions()
Correct boundary field.
virtual bool read(const dictionary &dict)
Read optional controls.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
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
Request registration (bool: true)
const fvMesh & mesh_
Reference to the fvMesh.
Do not request registration (bool: false)
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...