Go to the documentation of this file.
37 Foam::electrostaticDepositionFvPatchScalarField::eVPatch
48 refCast<const electrostaticDepositionFvPatchScalarField>(bf[patchi]);
54 void Foam::electrostaticDepositionFvPatchScalarField::setMaster()
const
69 if (isA<electrostaticDepositionFvPatchScalarField>(bf[patchi]))
71 electrostaticDepositionFvPatchScalarField& eVpf = eVPatch(patchi);
78 eVpf.master() = master;
84 void Foam::electrostaticDepositionFvPatchScalarField::round
92 f = std::round(
f*dcml)/dcml;
97 void Foam::electrostaticDepositionFvPatchScalarField::writeFilmFields()
const
104 const fvMesh&
mesh = eV.mesh();
123 if (isA<electrostaticDepositionFvPatchScalarField>(bf[patchi]))
125 electrostaticDepositionFvPatchScalarField& eVpf = eVPatch(patchi);
127 auto& hp =
h.boundaryFieldRef()[patchi];
146 fixedValueFvPatchScalarField(
p, iF),
175 fixedValueFvPatchScalarField(
p, iF,
dict, false),
176 h_(
"h",
dict,
p.size()),
197 Vi_(
dict.getOrDefault<scalar>(
"Vi", 0)),
198 Vanode_(
dict.getOrDefault<scalar>(
"Vanode", GREAT)),
199 phasesDict_(
dict.subOrEmptyDict(
"phases")),
208 dict.getCheckOrDefault<scalar>
219 if (
dict.found(
"value"))
232 if (!phasesDict_.empty())
234 phaseNames_.
setSize(phasesDict_.size());
235 phases_.setSize(phasesDict_.size());
236 sigmas_.setSize(phasesDict_.size());
239 for (
const entry& dEntry : phasesDict_)
241 const word&
key = dEntry.keyword();
243 if (!dEntry.isDict())
246 <<
"Entry " <<
key <<
" is not a dictionary" <<
nl
250 const dictionary& subDict = dEntry.dict();
260 subDict.getCheck<scalar>
276 db().getObjectPtr<volScalarField>(phaseNames_[i])
286 const electrostaticDepositionFvPatchScalarField& ptf,
288 const DimensionedField<scalar, volMesh>& iF,
289 const fvPatchFieldMapper& mapper
292 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
294 qcum_(ptf.qcum_, mapper),
295 Vfilm_(ptf.Vfilm_, mapper),
296 Ceffptr_(ptf.Ceffptr_.clone(
p.
patch())),
297 rptr_(ptf.rptr_.clone(
p.
patch())),
302 Vanode_(ptf.Vanode_),
303 phasesDict_(ptf.phasesDict_),
304 phaseNames_(ptf.phaseNames_),
305 phases_(ptf.phases_),
319 fixedValueFvPatchScalarField(ptf),
323 Ceffptr_(ptf.Ceffptr_.clone(
patch().
patch())),
329 Vanode_(ptf.Vanode_),
330 phasesDict_(ptf.phasesDict_),
331 phaseNames_(ptf.phaseNames_),
332 phases_(ptf.phases_),
347 fixedValueFvPatchScalarField(ptf, iF),
351 Ceffptr_(ptf.Ceffptr_.clone(
patch().
patch())),
357 Vanode_(ptf.Vanode_),
358 phasesDict_(ptf.phasesDict_),
359 phaseNames_(ptf.phaseNames_),
360 phases_(ptf.phases_),
375 fixedValueFvPatchScalarField::autoMap(m);
383 Ceffptr_->autoMap(m);
399 fixedValueFvPatchScalarField::rmap(ptf, addr);
402 refCast<const electrostaticDepositionFvPatchScalarField>(ptf);
404 h_.rmap(tiptf.h_, addr);
405 qcum_.rmap(tiptf.qcum_, addr);
406 Vfilm_.rmap(tiptf.Vfilm_, addr);
410 Ceffptr_->rmap(tiptf.Ceffptr_(), addr);
415 rptr_->rmap(tiptf.rptr_(), addr);
423 const label patchi =
patch().index();
427 tmp<scalarField> tsigma =
428 phases_[0].boundaryField()[patchi]*sigmas_[0].value();
430 for (label i = 1; i < phases_.size(); ++i)
433 phases_[i].boundaryField()[patchi]*sigmas_[i].value();
455 const scalar t = db().time().timeOutputValue();
456 const scalar dt = db().time().deltaTValue();
457 const label patchi =
patch().index();
463 tmp<scalarField> tjnp = -this->
sigma()*eV.boundaryField()[patchi].snGrad();
465 jnp =
max(jnp, scalar(0));
472 tmp<scalarField> tCoulombicEfficiency = Ceffptr_->value(t);
473 tmp<scalarField> tdh = tCoulombicEfficiency*(jnp - jMin_)*dt;
477 dh =
max(dh, scalar(0));
485 if (qcum_[i] < qMin_)
496 tmp<scalarField> tresistivity = rptr_->value(t);
497 tmp<scalarField> tRfilm = tresistivity*tdh;
498 tmp<scalarField> tdV = jnp*tRfilm;
500 Vfilm_ =
min(Vfilm_, Vanode_);
504 tmp<scalarField> tVbody = tjnp*Rbody_;
511 fixedValueFvPatchScalarField::updateCoeffs();
513 timei_ = db().time().timeIndex();
516 const scalar hMin =
gMin(h_);
517 const scalar hMax =
gMax(h_);
523 <<
", h: min = " << hMin
524 <<
", max = " << hMax
525 <<
", average = " << hAvg <<
nl
531 if (db().time().writeTime())
536 if (
patch().index() == master_)
548 h_.writeEntry(
"h",
os);
552 Ceffptr_->writeData(
os);
557 rptr_->writeData(
os);
560 if (!phasesDict_.empty())
562 phasesDict_.writeEntry(phasesDict_.dictName(),
os);
566 sigma_.writeEntry(
"sigma",
os);
577 writeEntry(
"value",
os);
588 electrostaticDepositionFvPatchScalarField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
A keyword and a list of tokens is an 'entry'.
fvPatchField< scalar > fvPatchScalarField
List< label > labelList
A List of labels.
virtual void write(Ostream &) const
virtual void operator=(const UList< Type > &)
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A class for handling words, derived from Foam::string.
virtual void autoMap(const fvPatchFieldMapper &)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Type gAverage(const FieldField< Field, Type > &f)
static word timeName(const scalar t, const int precision=precision_)
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
static bool master(const label communicator=worldComm)
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
Ostream & endl(Ostream &os)
label min(const labelHashSet &set, label minValue=labelMax)
const dimensionSet dimCurrent(0, 0, 0, 0, 0, 1, 0)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Generic templated field type.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionedScalar h
void setSize(const label n)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
T getCheckOrDefault(const word &keyword, const T &deflt, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
label max(const labelHashSet &set, label maxValue=labelMin)
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))
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Generic dimensioned Type class.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
virtual void rmap(const fvPatchScalarField &, const labelList &)
errorManipArg< error, int > exit(error &err, const int errNo=1)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
tmp< scalarField > sigma() const
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static MinMax< T > ge(const T &minVal)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
static tmp< T > New(Args &&... args)
T getCheck(const word &keyword, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
Ostream & writeEntry(const keyType &key, const T &value)
The electrostaticDeposition is a boundary condition to calculate electric potential (V) on a given bo...
const Time & time() const
virtual void updateCoeffs()
#define FatalIOErrorInFunction(ios)
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Type gMin(const FieldField< Field, Type > &f)
static word scopedName(const std::string &scope, const word &name)
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Generic GeometricField class.
A min/max value pair with additional methods. In addition to conveniently storing values,...
virtual void write(Ostream &) const
Type gMax(const FieldField< Field, Type > &f)
const Boundary & boundaryField() const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
electrostaticDepositionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)