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-2023 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_,
66  ),
67  mesh_
68  );
69  mesh_.objectRegistry::store(ptr);
70  }
71 
72  return *ptr;
73 }
74 
75 
77 Foam::functionObjects::electricPotential::sigma() const
78 {
79  const IOobject sigmaIO
80  (
81  IOobject::scopedName(typeName, "sigma"),
82  mesh_.time().timeName(),
83  mesh_.time(),
87  );
88 
89  if (phases_.size())
90  {
91  tmp<volScalarField> tsigma = phases_[0]*sigmas_[0];
92 
93  for (label i = 1; i < phases_.size(); ++i)
94  {
95  tsigma.ref() += phases_[i]*sigmas_[i];
96  }
97 
99  (
100  sigmaIO,
101  tsigma,
103  );
104  }
105 
107  (
108  sigmaIO,
109  mesh_,
110  sigma_,
112  );
113 }
114 
115 
117 Foam::functionObjects::electricPotential::epsilonm() const
118 {
119  // Vacuum permittivity (aka the electric constant) [A^2 s^4/(kg m^3)]
121  (
123  8.8541878128e-12 // CODATA value
124  );
125 
126  const IOobject epsilonrIO
127  (
128  IOobject::scopedName(typeName, "epsilonr"),
129  mesh_.time().timeName(),
130  mesh_.time(),
134  );
135 
136  if (phases_.size())
137  {
138  tmp<volScalarField> tepsilonr = phases_[0]*epsilonrs_[0];
139 
140  for (label i = 1; i < phases_.size(); ++i)
141  {
142  tepsilonr.ref() += phases_[i]*epsilonrs_[i];
143  }
144 
146  (
147  epsilonrIO,
148  epsilon0*tepsilonr,
150  );
151  }
152 
154  (
155  epsilonrIO,
156  mesh_,
157  epsilon0*epsilonr_,
159  );
160 }
161 
162 
163 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
164 
165 Foam::functionObjects::electricPotential::electricPotential
166 (
167  const word& name,
168  const Time& runTime,
169  const dictionary& dict
170 )
171 :
173  phasesDict_(dict.subOrEmptyDict("phases")),
174  phaseNames_(),
175  phases_(),
176  sigmas_(),
177  sigma_
178  (
180  (
182  dict.getCheckOrDefault<scalar>
183  (
184  "sigma",
185  scalar(1),
186  scalarMinMax::ge(SMALL)
187  )
188  )
189  ),
190  epsilonrs_(),
191  epsilonr_
192  (
194  (
195  dimless,
196  dict.getCheckOrDefault<scalar>
197  (
198  "epsilonr",
199  scalar(1),
200  scalarMinMax::ge(SMALL)
201  )
202  )
203  ),
204  Vname_
205  (
206  dict.getOrDefault<word>
207  (
208  "V",
209  IOobject::scopedName(typeName, "V")
210  )
211  ),
212  Ename_
213  (
214  dict.getOrDefault<word>
215  (
216  "E",
217  IOobject::scopedName(typeName, "E")
218  )
219  ),
220  nCorr_(1),
221  writeDerivedFields_(false),
222  electricField_(false)
223 {
224  read(dict);
225 
226  // Force creation of transported field so any BCs using it can
227  // look it up
228  volScalarField& eV = getOrReadField(Vname_);
230 
231  if (electricField_)
232  {
233  auto* ptr = getObjectPtr<volVectorField>(Ename_);
234 
235  if (!ptr)
236  {
237  ptr = new volVectorField
238  (
239  IOobject
240  (
241  Ename_,
242  mesh_.time().timeName(),
243  mesh_,
246  ),
247  -fvc::grad(eV)
248  );
249  mesh_.objectRegistry::store(ptr);
250  }
251  }
252 }
253 
254 
255 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
256 
258 {
260  {
261  return false;
262  }
263 
264  Log << type() << " read: " << name() << endl;
265 
266  dict.readIfPresent("sigma", sigma_);
267  dict.readIfPresent("epsilonr", epsilonr_);
268  dict.readIfPresent("nCorr", nCorr_);
269  dict.readIfPresent("writeDerivedFields", writeDerivedFields_);
270  dict.readIfPresent("electricField", electricField_);
271 
272  // If flow is multiphase
273  if (!phasesDict_.empty())
274  {
275  phaseNames_.setSize(phasesDict_.size());
276  phases_.setSize(phasesDict_.size());
277  sigmas_.setSize(phasesDict_.size());
278 
279  if (writeDerivedFields_)
280  {
281  epsilonrs_.setSize(phasesDict_.size());
282  }
283 
284  label phasei = 0;
285  for (const entry& dEntry : phasesDict_)
286  {
287  const word& key = dEntry.keyword();
288 
289  if (!dEntry.isDict())
290  {
291  FatalIOErrorInFunction(phasesDict_)
292  << "Entry " << key << " is not a dictionary" << nl
293  << exit(FatalIOError);
294  }
295 
296  const dictionary& subDict = dEntry.dict();
297 
298  phaseNames_[phasei] = key;
299 
300  sigmas_.set
301  (
302  phasei,
304  (
306  subDict.getCheck<scalar>
307  (
308  "sigma",
309  scalarMinMax::ge(SMALL)
310  )
311  )
312  );
313 
314  if (writeDerivedFields_)
315  {
316  epsilonrs_.set
317  (
318  phasei,
320  (
321  dimless,
322  subDict.getCheck<scalar>
323  (
324  "epsilonr",
325  scalarMinMax::ge(SMALL)
326  )
327  )
328  );
329  }
330 
331  ++phasei;
332  }
333 
334  forAll(phaseNames_, i)
335  {
336  phases_.set
337  (
338  i,
339  mesh_.getObjectPtr<volScalarField>(phaseNames_[i])
340  );
341  }
342  }
343 
344  return false;
345 }
346 
347 
349 {
350  Log << type() << " execute: " << name() << endl;
351 
352  tmp<volScalarField> tsigma = this->sigma();
353  const auto& sigma = tsigma();
354 
355  volScalarField& eV = getOrReadField(Vname_);
356 
357  for (int i = 1; i <= nCorr_; ++i)
358  {
359  fvScalarMatrix eVEqn
360  (
361  - fvm::laplacian(sigma, eV)
362  );
363 
364  eVEqn.relax();
365 
366  eVEqn.solve();
367  }
368 
369  if (electricField_)
370  {
371  auto& E = lookupObjectRef<volVectorField>(Ename_);
372  E == -fvc::grad(eV);
373  }
375  Log << endl;
376 
377  return true;
378 }
379 
380 
382 {
383  Log << type() << " write: " << name() << nl
384  << tab << Vname_
385  << endl;
386 
387  volScalarField& eV = getOrReadField(Vname_);
388 
389  if (electricField_)
390  {
391  const auto& E = lookupObject<volVectorField>(Ename_);
392 
393  Log << tab << E.name() << endl;
394 
395  E.write();
396  }
397 
398  if (writeDerivedFields_)
399  {
400  // Write the current density field
401  tmp<volScalarField> tsigma = this->sigma();
402 
403  auto eJ = tmp<volVectorField>::New
404  (
405  IOobject
406  (
407  IOobject::scopedName(typeName, "J"),
408  mesh_.time().timeName(),
409  mesh_.time(),
413  ),
414  -tsigma*fvc::grad(eV),
416  );
417 
418  Log << tab << eJ().name() << endl;
419 
420  eJ->write();
421 
422 
423  // Write the free-charge density field
424  tmp<volScalarField> tepsilonm = this->epsilonm();
425 
426  auto erho = tmp<volScalarField>::New
427  (
428  IOobject
429  (
430  IOobject::scopedName(typeName, "rho"),
431  mesh_.time().timeName(),
432  mesh_.time(),
436  ),
437  fvc::div(tepsilonm*(-fvc::grad(eV))),
439  );
440 
441  Log << tab << erho().name() << endl;
442 
443  erho->write();
444  }
445 
446  eV.write();
447 
448  return true;
449 }
450 
451 
452 // ************************************************************************* //
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)
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
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:120
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:49
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:48
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:362
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
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.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
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 INVALID.
Definition: exprTraits.C:52
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:212
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
virtual bool execute()
Calculate the function object.
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:770
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
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:607
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:201
#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:171
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 ...