mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
STYLE: electricPotential: simplify read function
- change 'fieldName' to 'Vname' for better clarity
This commit is contained in:
committed by
Kutalmış Berçin
parent
7269cc1d3b
commit
efe8220a26
@ -201,11 +201,11 @@ Foam::functionObjects::electricPotential::electricPotential
|
|||||||
)
|
)
|
||||||
)
|
)
|
||||||
),
|
),
|
||||||
fieldName_
|
Vname_
|
||||||
(
|
(
|
||||||
dict.getOrDefault<word>
|
dict.getOrDefault<word>
|
||||||
(
|
(
|
||||||
"field",
|
"V",
|
||||||
IOobject::scopedName(typeName, "V")
|
IOobject::scopedName(typeName, "V")
|
||||||
)
|
)
|
||||||
),
|
),
|
||||||
@ -216,7 +216,7 @@ Foam::functionObjects::electricPotential::electricPotential
|
|||||||
|
|
||||||
// Force creation of transported field so any BCs using it can
|
// Force creation of transported field so any BCs using it can
|
||||||
// look it up
|
// look it up
|
||||||
volScalarField& eV = getOrReadField(fieldName_);
|
volScalarField& eV = getOrReadField(Vname_);
|
||||||
eV.correctBoundaryConditions();
|
eV.correctBoundaryConditions();
|
||||||
}
|
}
|
||||||
|
|
||||||
@ -225,88 +225,88 @@ Foam::functionObjects::electricPotential::electricPotential
|
|||||||
|
|
||||||
bool Foam::functionObjects::electricPotential::read(const dictionary& dict)
|
bool Foam::functionObjects::electricPotential::read(const dictionary& dict)
|
||||||
{
|
{
|
||||||
if (fvMeshFunctionObject::read(dict))
|
if (!fvMeshFunctionObject::read(dict))
|
||||||
{
|
{
|
||||||
Log << type() << " read: " << name() << endl;
|
return false;
|
||||||
|
}
|
||||||
|
|
||||||
dict.readIfPresent("sigma", sigma_);
|
Log << type() << " read: " << name() << endl;
|
||||||
dict.readIfPresent("epsilonr", epsilonr_);
|
|
||||||
dict.readIfPresent("nCorr", nCorr_);
|
|
||||||
dict.readIfPresent("writeDerivedFields", writeDerivedFields_);
|
|
||||||
|
|
||||||
// If flow is multiphase
|
dict.readIfPresent("sigma", sigma_);
|
||||||
if (!phasesDict_.empty())
|
dict.readIfPresent("epsilonr", epsilonr_);
|
||||||
|
dict.readIfPresent("nCorr", nCorr_);
|
||||||
|
dict.readIfPresent("writeDerivedFields", writeDerivedFields_);
|
||||||
|
|
||||||
|
// If flow is multiphase
|
||||||
|
if (!phasesDict_.empty())
|
||||||
|
{
|
||||||
|
phaseNames_.setSize(phasesDict_.size());
|
||||||
|
phases_.setSize(phasesDict_.size());
|
||||||
|
sigmas_.setSize(phasesDict_.size());
|
||||||
|
|
||||||
|
if (writeDerivedFields_)
|
||||||
{
|
{
|
||||||
phaseNames_.setSize(phasesDict_.size());
|
epsilonrs_.setSize(phasesDict_.size());
|
||||||
phases_.setSize(phasesDict_.size());
|
}
|
||||||
sigmas_.setSize(phasesDict_.size());
|
|
||||||
|
label phasei = 0;
|
||||||
|
for (const entry& dEntry : phasesDict_)
|
||||||
|
{
|
||||||
|
const word& key = dEntry.keyword();
|
||||||
|
|
||||||
|
if (!dEntry.isDict())
|
||||||
|
{
|
||||||
|
FatalIOErrorInFunction(phasesDict_)
|
||||||
|
<< "Entry " << key << " is not a dictionary" << nl
|
||||||
|
<< exit(FatalIOError);
|
||||||
|
}
|
||||||
|
|
||||||
|
const dictionary& subDict = dEntry.dict();
|
||||||
|
|
||||||
|
phaseNames_[phasei] = key;
|
||||||
|
|
||||||
|
sigmas_.set
|
||||||
|
(
|
||||||
|
phasei,
|
||||||
|
new dimensionedScalar
|
||||||
|
(
|
||||||
|
sqr(dimCurrent)*pow3(dimTime)/(dimMass*pow3(dimLength)),
|
||||||
|
subDict.getCheck<scalar>
|
||||||
|
(
|
||||||
|
"sigma",
|
||||||
|
scalarMinMax::ge(SMALL)
|
||||||
|
)
|
||||||
|
)
|
||||||
|
);
|
||||||
|
|
||||||
if (writeDerivedFields_)
|
if (writeDerivedFields_)
|
||||||
{
|
{
|
||||||
epsilonrs_.setSize(phasesDict_.size());
|
epsilonrs_.set
|
||||||
}
|
|
||||||
|
|
||||||
label phasei = 0;
|
|
||||||
for (const entry& dEntry : phasesDict_)
|
|
||||||
{
|
|
||||||
const word& key = dEntry.keyword();
|
|
||||||
|
|
||||||
if (!dEntry.isDict())
|
|
||||||
{
|
|
||||||
FatalIOErrorInFunction(phasesDict_)
|
|
||||||
<< "Entry " << key << " is not a dictionary" << nl
|
|
||||||
<< exit(FatalIOError);
|
|
||||||
}
|
|
||||||
|
|
||||||
const dictionary& subDict = dEntry.dict();
|
|
||||||
|
|
||||||
phaseNames_[phasei] = key;
|
|
||||||
|
|
||||||
sigmas_.set
|
|
||||||
(
|
(
|
||||||
phasei,
|
phasei,
|
||||||
new dimensionedScalar
|
new dimensionedScalar
|
||||||
(
|
(
|
||||||
sqr(dimCurrent)*pow3(dimTime)/(dimMass*pow3(dimLength)),
|
dimless,
|
||||||
subDict.getCheck<scalar>
|
subDict.getCheck<scalar>
|
||||||
(
|
(
|
||||||
"sigma",
|
"epsilonr",
|
||||||
scalarMinMax::ge(SMALL)
|
scalarMinMax::ge(SMALL)
|
||||||
)
|
)
|
||||||
)
|
)
|
||||||
);
|
);
|
||||||
|
|
||||||
if (writeDerivedFields_)
|
|
||||||
{
|
|
||||||
epsilonrs_.set
|
|
||||||
(
|
|
||||||
phasei,
|
|
||||||
new dimensionedScalar
|
|
||||||
(
|
|
||||||
dimless,
|
|
||||||
subDict.getCheck<scalar>
|
|
||||||
(
|
|
||||||
"epsilonr",
|
|
||||||
scalarMinMax::ge(SMALL)
|
|
||||||
)
|
|
||||||
)
|
|
||||||
);
|
|
||||||
}
|
|
||||||
|
|
||||||
++phasei;
|
|
||||||
}
|
}
|
||||||
|
|
||||||
forAll(phaseNames_, i)
|
++phasei;
|
||||||
{
|
|
||||||
phases_.set
|
|
||||||
(
|
|
||||||
i,
|
|
||||||
mesh_.getObjectPtr<volScalarField>(phaseNames_[i])
|
|
||||||
);
|
|
||||||
}
|
|
||||||
}
|
}
|
||||||
|
|
||||||
return true;
|
forAll(phaseNames_, i)
|
||||||
|
{
|
||||||
|
phases_.set
|
||||||
|
(
|
||||||
|
i,
|
||||||
|
mesh_.getObjectPtr<volScalarField>(phaseNames_[i])
|
||||||
|
);
|
||||||
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
return false;
|
return false;
|
||||||
@ -318,11 +318,11 @@ bool Foam::functionObjects::electricPotential::execute()
|
|||||||
Log << type() << " execute: " << name() << endl;
|
Log << type() << " execute: " << name() << endl;
|
||||||
|
|
||||||
tmp<volScalarField> tsigma = this->sigma();
|
tmp<volScalarField> tsigma = this->sigma();
|
||||||
const volScalarField& sigma = tsigma();
|
const auto& sigma = tsigma();
|
||||||
|
|
||||||
volScalarField& eV = getOrReadField(fieldName_);
|
volScalarField& eV = getOrReadField(Vname_);
|
||||||
|
|
||||||
for (label i = 1; i <= nCorr_; ++i)
|
for (int i = 1; i <= nCorr_; ++i)
|
||||||
{
|
{
|
||||||
fvScalarMatrix eVEqn
|
fvScalarMatrix eVEqn
|
||||||
(
|
(
|
||||||
@ -343,10 +343,10 @@ bool Foam::functionObjects::electricPotential::execute()
|
|||||||
bool Foam::functionObjects::electricPotential::write()
|
bool Foam::functionObjects::electricPotential::write()
|
||||||
{
|
{
|
||||||
Log << type() << " write: " << name() << nl
|
Log << type() << " write: " << name() << nl
|
||||||
<< tab << fieldName_
|
<< tab << Vname_
|
||||||
<< endl;
|
<< endl;
|
||||||
|
|
||||||
volScalarField& eV = getOrReadField(fieldName_);
|
volScalarField& eV = getOrReadField(Vname_);
|
||||||
|
|
||||||
if (writeDerivedFields_)
|
if (writeDerivedFields_)
|
||||||
{
|
{
|
||||||
|
|||||||
@ -118,9 +118,9 @@ Usage
|
|||||||
}
|
}
|
||||||
|
|
||||||
// Optional entries
|
// Optional entries
|
||||||
nCorr <label>;
|
nCorr <int>;
|
||||||
writeDerivedFields <bool>;
|
writeDerivedFields <bool>;
|
||||||
fieldName <word>;
|
V <word>;
|
||||||
|
|
||||||
// Inherited entries
|
// Inherited entries
|
||||||
...
|
...
|
||||||
@ -134,9 +134,9 @@ Usage
|
|||||||
libs | Library name: solverFunctionObjects | word | yes | -
|
libs | Library name: solverFunctionObjects | word | yes | -
|
||||||
sigma | Isotropic electrical conductivity of phase | scalar | yes | -
|
sigma | Isotropic electrical conductivity of phase | scalar | yes | -
|
||||||
epsilonr | Isotropic relative permittivity of phase | scalar | no | -
|
epsilonr | Isotropic relative permittivity of phase | scalar | no | -
|
||||||
nCorr | Number of corrector iterations | label | no | 1
|
nCorr | Number of corrector iterations | int | no | 1
|
||||||
writeDerivedFields | Flag to write extra fields | bool | no | false
|
writeDerivedFields | Flag to write extra fields | bool | no | false
|
||||||
fieldName | Name of operand field | word | no | electricPotential:V
|
V | Name of electric potential field | word | no | electricPotential:V
|
||||||
\endtable
|
\endtable
|
||||||
|
|
||||||
The inherited entries are elaborated in:
|
The inherited entries are elaborated in:
|
||||||
@ -199,11 +199,11 @@ class electricPotential
|
|||||||
//- Isotropic relative permittivity of a single phase
|
//- Isotropic relative permittivity of a single phase
|
||||||
dimensionedScalar epsilonr_;
|
dimensionedScalar epsilonr_;
|
||||||
|
|
||||||
//- Name of the operand field
|
//- Name of electric potential field
|
||||||
word fieldName_;
|
word Vname_;
|
||||||
|
|
||||||
//- Number of corrector iterations
|
//- Number of corrector iterations
|
||||||
label nCorr_;
|
int nCorr_;
|
||||||
|
|
||||||
//- Flag to write derived fields of
|
//- Flag to write derived fields of
|
||||||
//- electric field, current density and free-charge density
|
//- electric field, current density and free-charge density
|
||||||
|
|||||||
Reference in New Issue
Block a user