electrostaticDepositionFvPatchScalarField.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 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 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 
37 Foam::electrostaticDepositionFvPatchScalarField::eVPatch
38 (
39  const label patchi
40 ) const
41 {
42  const auto& eV =
43  db().lookupObject<volScalarField>(this->internalField().name());
44 
45  const volScalarField::Boundary& bf = eV.boundaryField();
46 
47  const auto& eVpf =
48  refCast<const electrostaticDepositionFvPatchScalarField>(bf[patchi]);
49 
50  return const_cast<electrostaticDepositionFvPatchScalarField&>(eVpf);
51 }
52 
53 
54 void Foam::electrostaticDepositionFvPatchScalarField::setMaster() const
55 {
56  if (master_ != -1)
57  {
58  return;
59  }
60 
61  const auto& eV =
62  db().lookupObject<volScalarField>(this->internalField().name());
63 
64  const volScalarField::Boundary& bf = eV.boundaryField();
65 
66  label master = -1;
67  forAll(bf, patchi)
68  {
69  if (isA<electrostaticDepositionFvPatchScalarField>(bf[patchi]))
70  {
71  electrostaticDepositionFvPatchScalarField& eVpf = eVPatch(patchi);
72 
73  if (master == -1)
74  {
75  master = patchi;
76  }
77 
78  eVpf.master() = master;
79  }
80  }
81 }
82 
83 
84 void Foam::electrostaticDepositionFvPatchScalarField::round
85 (
87  const scalar dcml
88 ) const
89 {
90  for (auto& f : fld)
91  {
92  f = std::round(f*dcml)/dcml;
93  }
94 }
95 
96 
97 void Foam::electrostaticDepositionFvPatchScalarField::writeFilmFields() const
98 {
99  const auto& eV =
100  db().lookupObject<volScalarField>(this->internalField().name());
101 
102  const volScalarField::Boundary& bf = eV.boundaryField();
103 
104  const fvMesh& mesh = eV.mesh();
105 
107  (
108  IOobject
109  (
110  IOobject::scopedName("electrostaticDeposition", "h"),
111  mesh.time().timeName(),
112  mesh,
116  ),
117  mesh,
119  );
120 
121  forAll(bf, patchi)
122  {
123  if (isA<electrostaticDepositionFvPatchScalarField>(bf[patchi]))
124  {
125  electrostaticDepositionFvPatchScalarField& eVpf = eVPatch(patchi);
126 
127  auto& hp = h.boundaryFieldRef()[patchi];
128 
129  hp = eVpf.h();
130  }
131  }
132 
133  h.write();
134 }
135 
136 
137 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
138 
141 (
142  const fvPatch& p,
144 )
145 :
146  fixedValueFvPatchScalarField(p, iF),
147  h_(p.size(), 0),
148  qcum_(p.size(), 0),
149  Vfilm_(p.size(), 0),
150  Ceffptr_(nullptr),
151  rptr_(nullptr),
152  jMin_(0),
153  qMin_(0),
154  Rbody_(0),
155  Vi_(0),
156  Vanode_(GREAT),
157  phasesDict_(),
158  phaseNames_(),
159  phases_(),
160  sigmas_(),
161  sigma_(sqr(dimCurrent)*pow3(dimTime)/(dimMass*pow3(dimLength)), scalar(1)),
162  timei_(-1),
163  master_(-1)
164 {}
165 
166 
169 (
170  const fvPatch& p,
172  const dictionary& dict
173 )
174 :
175  fixedValueFvPatchScalarField(p, iF, dict, IOobjectOption::NO_READ),
176  h_("h", dict, p.size()),
177  qcum_("qCumulative", dict, p.size(), IOobjectOption::LAZY_READ),
178  Vfilm_("Vfilm", dict, p.size(), IOobjectOption::LAZY_READ),
179  Ceffptr_
180  (
181  PatchFunction1<scalar>::New(p.patch(), "CoulombicEfficiency", dict)
182  ),
183  rptr_(PatchFunction1<scalar>::New(p.patch(), "resistivity", dict)),
184  jMin_(dict.getCheckOrDefault<scalar>("jMin", 0, scalarMinMax::ge(0))),
185  qMin_(dict.getCheckOrDefault<scalar>("qMin", 0, scalarMinMax::ge(0))),
186  Rbody_(dict.getCheckOrDefault<scalar>("Rbody", 0, scalarMinMax::ge(0))),
187  Vi_(dict.getOrDefault<scalar>("Vi", 0)),
188  Vanode_(dict.getOrDefault<scalar>("Vanode", GREAT)),
189  phasesDict_(dict.subOrEmptyDict("phases")),
190  phaseNames_(),
191  phases_(),
192  sigmas_(),
193  sigma_
194  (
196  (
198  dict.getCheckOrDefault<scalar>
199  (
200  "sigma",
201  scalar(1),
202  scalarMinMax::ge(SMALL)
203  )
204  )
205  ),
206  timei_(-1),
207  master_(-1)
208 {
209  if (!this->readValueEntry(dict))
210  {
211  this->extrapolateInternal(); // Zero-gradient patch values
212  }
213 
214  // If flow is multiphase
215  if (!phasesDict_.empty())
216  {
217  phaseNames_.setSize(phasesDict_.size());
218  phases_.setSize(phasesDict_.size());
219  sigmas_.setSize(phasesDict_.size());
220 
221  label phasei = 0;
222  for (const entry& dEntry : phasesDict_)
223  {
224  const word& key = dEntry.keyword();
225 
226  if (!dEntry.isDict())
227  {
228  FatalIOErrorInFunction(phasesDict_)
229  << "Entry " << key << " is not a dictionary" << nl
230  << exit(FatalIOError);
231  }
232 
233  const dictionary& subDict = dEntry.dict();
234 
235  phaseNames_[phasei] = key;
236 
237  sigmas_.set
238  (
239  phasei,
241  (
243  subDict.getCheck<scalar>
244  (
245  "sigma",
246  scalarMinMax::ge(SMALL)
247  )
248  )
249  );
250 
251  ++phasei;
252  }
253 
254  forAll(phaseNames_, i)
255  {
256  phases_.set
257  (
258  i,
259  db().getObjectPtr<volScalarField>(phaseNames_[i])
260  );
261  }
262  }
263 }
264 
265 
268 (
269  const electrostaticDepositionFvPatchScalarField& ptf,
270  const fvPatch& p,
271  const DimensionedField<scalar, volMesh>& iF,
272  const fvPatchFieldMapper& mapper
273 )
274 :
275  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
276  h_(ptf.h_, mapper),
277  qcum_(ptf.qcum_, mapper),
278  Vfilm_(ptf.Vfilm_, mapper),
279  Ceffptr_(ptf.Ceffptr_.clone(p.patch())),
280  rptr_(ptf.rptr_.clone(p.patch())),
281  jMin_(ptf.jMin_),
282  qMin_(ptf.qMin_),
283  Rbody_(ptf.Rbody_),
284  Vi_(ptf.Vi_),
285  Vanode_(ptf.Vanode_),
286  phasesDict_(ptf.phasesDict_),
287  phaseNames_(ptf.phaseNames_),
288  phases_(ptf.phases_),
289  sigmas_(),
290  sigma_(ptf.sigma_),
291  timei_(ptf.timei_),
292  master_(-1)
293 {}
294 
295 
298 (
300 )
301 :
302  fixedValueFvPatchScalarField(ptf),
303  h_(ptf.h_),
304  qcum_(ptf.qcum_),
305  Vfilm_(ptf.Vfilm_),
306  Ceffptr_(ptf.Ceffptr_.clone(patch().patch())),
307  rptr_(ptf.rptr_.clone(patch().patch())),
308  jMin_(ptf.jMin_),
309  qMin_(ptf.qMin_),
310  Rbody_(ptf.Rbody_),
311  Vi_(ptf.Vi_),
312  Vanode_(ptf.Vanode_),
313  phasesDict_(ptf.phasesDict_),
314  phaseNames_(ptf.phaseNames_),
315  phases_(ptf.phases_),
316  sigmas_(),
317  sigma_(ptf.sigma_),
318  timei_(ptf.timei_),
319  master_(-1)
320 {}
321 
322 
325 (
328 )
329 :
330  fixedValueFvPatchScalarField(ptf, iF),
331  h_(ptf.h_),
332  qcum_(ptf.qcum_),
333  Vfilm_(ptf.Vfilm_),
334  Ceffptr_(ptf.Ceffptr_.clone(patch().patch())),
335  rptr_(ptf.rptr_.clone(patch().patch())),
336  jMin_(ptf.jMin_),
337  qMin_(ptf.qMin_),
338  Rbody_(ptf.Rbody_),
339  Vi_(ptf.Vi_),
340  Vanode_(ptf.Vanode_),
341  phasesDict_(ptf.phasesDict_),
342  phaseNames_(ptf.phaseNames_),
343  phases_(ptf.phases_),
344  sigmas_(),
345  sigma_(ptf.sigma_),
346  timei_(ptf.timei_),
347  master_(-1)
348 {}
349 
350 
351 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
352 
354 (
355  const fvPatchFieldMapper& m
356 )
357 {
358  fixedValueFvPatchScalarField::autoMap(m);
359 
360  h_.autoMap(m);
361  qcum_.autoMap(m);
362  Vfilm_.autoMap(m);
363 
364  if (Ceffptr_)
365  {
366  Ceffptr_->autoMap(m);
367  }
368 
369  if (rptr_)
370  {
371  rptr_->autoMap(m);
372  }
373 }
374 
375 
377 (
378  const fvPatchScalarField& ptf,
379  const labelList& addr
380 )
381 {
382  fixedValueFvPatchScalarField::rmap(ptf, addr);
383 
384  const auto& tiptf =
385  refCast<const electrostaticDepositionFvPatchScalarField>(ptf);
386 
387  h_.rmap(tiptf.h_, addr);
388  qcum_.rmap(tiptf.qcum_, addr);
389  Vfilm_.rmap(tiptf.Vfilm_, addr);
390 
391  if (Ceffptr_)
392  {
393  Ceffptr_->rmap(tiptf.Ceffptr_(), addr);
394  }
395 
396  if (rptr_)
397  {
398  rptr_->rmap(tiptf.rptr_(), addr);
399  }
400 }
401 
402 
405 {
406  const label patchi = patch().index();
407 
408  if (phases_.size())
409  {
410  tmp<scalarField> tsigma =
411  phases_[0].boundaryField()[patchi]*sigmas_[0].value();
412 
413  for (label i = 1; i < phases_.size(); ++i)
414  {
415  tsigma.ref() +=
416  phases_[i].boundaryField()[patchi]*sigmas_[i].value();
417  }
418 
419  return tsigma;
420  }
421 
422  return tmp<scalarField>::New(patch().size(), sigma_.value());
423 }
424 
425 
427 {
428  if (updated())
429  {
430  return;
431  }
432 
433  if (timei_ == db().time().timeIndex())
434  {
435  return;
436  }
437 
438  const scalar t = db().time().timeOutputValue();
439  const scalar dt = db().time().deltaTValue();
440  const label patchi = patch().index();
441 
442  const auto& eV =
443  db().lookupObject<volScalarField>(this->internalField().name());
444 
445  // Current density on film interface
446  tmp<scalarField> tjnp = -this->sigma()*eV.boundaryField()[patchi].snGrad();
447  scalarField& jnp = tjnp.ref();
448  jnp = max(jnp, scalar(0)); // experimental - do not allow any negative jnp
449  // experimental - avoid micro/nano currents/volts
450  // to reduce snowballing effects of lateral gradients on the patch
451  round(jnp);
452 
453 
454  // Calculate film-thickness finite increments
455  tmp<scalarField> tCoulombicEfficiency = Ceffptr_->value(t);
456  tmp<scalarField> tdh = tCoulombicEfficiency*(jnp - jMin_)*dt;
457  scalarField& dh = tdh.ref();
458 
459  // Do not allow any depletion or abrasion of deposition
460  dh = max(dh, scalar(0));
461 
462  // Do not allow any deposition when accumulative specific
463  // charge is less than minimum accumulative specific charge
464  qcum_ += jnp*dt;
465 
466  forAll(dh, i)
467  {
468  if (qcum_[i] < qMin_)
469  {
470  dh[i] = 0;
471  }
472  }
473 
474  // Add finite increments of film thickness to total film thickness
475  h_ += dh;
476 
477 
478  // Calculate incremental electric potential due to film resistance
479  tmp<scalarField> tresistivity = rptr_->value(t);
480  tmp<scalarField> tRfilm = tresistivity*tdh;
481  tmp<scalarField> tdV = jnp*tRfilm;
482  Vfilm_ += tdV;
483  Vfilm_ = min(Vfilm_, Vanode_);
484 
485 
486  // Calculate electric potential due to body resistance
487  tmp<scalarField> tVbody = tjnp*Rbody_;
488 
489 
490  // Add all electric potential contributions
491  operator==(min(Vi_ + Vfilm_ + tVbody, Vanode_));
492 
493 
494  fixedValueFvPatchScalarField::updateCoeffs();
495 
496  timei_ = db().time().timeIndex();
497 
498  {
499  const scalar hMin = gMin(h_);
500  const scalar hMax = gMax(h_);
501  const scalar hAvg = gAverage(h_);
502 
503  if (Pstream::master())
504  {
505  Info<< " patch: " << patch().name()
506  << ", h: min = " << hMin
507  << ", max = " << hMax
508  << ", average = " << hAvg << nl
509  << endl;
510  }
511  }
512 
513  // Write here to avoid any upset to redistributePar-decompose
514  if (db().time().writeTime())
515  {
516  // Write film thickness fields as patch fields of a volScalarField
517  setMaster();
518 
519  if (patch().index() == master_)
520  {
521  writeFilmFields();
522  }
523  }
524 }
525 
526 
528 {
530 
531  h_.writeEntry("h", os);
532 
533  if (Ceffptr_)
534  {
535  Ceffptr_->writeData(os);
536  }
537 
538  if (rptr_)
539  {
540  rptr_->writeData(os);
541  }
542 
543  if (!phasesDict_.empty())
544  {
545  phasesDict_.writeEntry(phasesDict_.dictName(), os);
546  }
547  else
548  {
549  sigma_.writeEntry("sigma", os);
550  }
551 
552  os.writeEntryIfDifferent<scalar>("jMin", 0, jMin_);
553  os.writeEntryIfDifferent<scalar>("qMin", 0, qMin_);
554  os.writeEntryIfDifferent<scalar>("Rbody", 0, Rbody_);
555  os.writeEntryIfDifferent<scalar>("Vi", 0, Vi_);
556  os.writeEntryIfDifferent<scalar>("Vanode", GREAT, Vanode_);
557  qcum_.writeEntry("qCumulative", os);
558  Vfilm_.writeEntry("Vfilm", os);
559 
561 }
562 
563 
564 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
565 
566 namespace Foam
567 {
569  (
571  electrostaticDepositionFvPatchScalarField
572  );
573 }
574 
575 // ************************************************************************* //
Foam::surfaceFields.
tmp< scalarField > sigma() const
Return the isotropic electrical conductivity field of mixture.
dictionary dict
electrostaticDepositionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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
Type gMin(const FieldField< Field, Type > &f)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
static MinMax< scalar > ge(const scalar &minVal)
A semi-infinite range from minVal to the type max.
Definition: MinMaxI.H:24
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
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
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
labelList f(nPoints)
const dimensionedScalar h
Planck constant.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
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:627
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Type gAverage(const FieldField< Field, Type > &f)
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
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
const dimensionSet dimCurrent(0, 0, 0, 0, 0, 1, 0)
Definition: dimensionSets.H:54
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
Nothing to be read.
const std::string patch
OpenFOAM patch number as a std::string.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
The electrostaticDeposition is a boundary condition to calculate electric potential (V) on a given bo...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Do not request registration (bool: false)
Namespace for OpenFOAM.
label timeIndex
Definition: getTimeIndex.H:24
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...