Public Member Functions | List of all members
electrostaticDepositionFvPatchScalarField Class Reference

The electrostaticDeposition is a boundary condition to calculate electric potential (V) on a given boundary based on film thickness (h) and film resistance (R) fields which are updated based on a given patch-normal current density field (jn), Coulombic efficiency and film resistivity. More...

Inheritance diagram for electrostaticDepositionFvPatchScalarField:
Inheritance graph
[legend]
Collaboration diagram for electrostaticDepositionFvPatchScalarField:
Collaboration graph
[legend]

Public Member Functions

 TypeName ("electrostaticDeposition")
 Runtime type information. More...
 
 electrostaticDepositionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &)
 Construct from patch and internal field. More...
 
 electrostaticDepositionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
 Construct from patch, internal field and dictionary. More...
 
 electrostaticDepositionFvPatchScalarField (const electrostaticDepositionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &)
 Construct by mapping given electrostaticDepositionFvPatchScalarField onto a new patch. More...
 
 electrostaticDepositionFvPatchScalarField (const electrostaticDepositionFvPatchScalarField &)
 Construct as copy. More...
 
virtual tmp< fvPatchScalarFieldclone () const
 Construct and return a clone. More...
 
 electrostaticDepositionFvPatchScalarField (const electrostaticDepositionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &)
 Construct as copy setting internal field reference. More...
 
virtual tmp< fvPatchScalarFieldclone (const DimensionedField< scalar, volMesh > &iF) const
 Construct and return a clone setting internal field reference. More...
 
const scalarFieldh () const noexcept
 Return const access to film thickness patch field. More...
 
virtual void autoMap (const fvPatchFieldMapper &)
 Map (and resize as needed) from self given a mapping object. More...
 
virtual void rmap (const fvPatchScalarField &, const labelList &)
 Reverse map the given fvPatchField onto this fvPatchField. More...
 
tmp< scalarFieldsigma () const
 Return the isotropic electrical conductivity field of mixture. More...
 
virtual void updateCoeffs ()
 Update the coefficients associated with the patch field. More...
 
virtual void write (Ostream &) const
 Write. More...
 

Detailed Description

The electrostaticDeposition is a boundary condition to calculate electric potential (V) on a given boundary based on film thickness (h) and film resistance (R) fields which are updated based on a given patch-normal current density field (jn), Coulombic efficiency and film resistivity.

\[ j_n = - \sigma \nabla^\perp_p V = - \sigma (\vec{n}\cdot(\nabla V)_p) \]

\[ \frac{dh}{dt} = C_{eff} (j_n - j_{min}) \]

\[ \frac{dR}{dt} = \rho \frac{dh}{dt} = \rho C_{eff} (j_n - j_{min}) \]

\[ V_{film}^n = V_{film}^o + j_n R_\Delta \]

\[ V_{body} = j_n R_{body} \]

\[ V_p^n = V_i + V_{body} + V_{film}^n \]

where

$ j_n $ = Patch-normal current density [A/m^2]
$ V_p^n $ = Electric potential on film-fluid interface [volt = kg m^2/(A s^3)]
$ V_p^o $ = Previous time-step electric potential on the interface [volt]
$ V_{film} $ = Electric potential due to film resistance [volt]
$ V_{body} $ = Electric potential due to body resistance [volt]
$ V_i $ = Initial electric potential [volt]
$ R_\Delta$ = Film resistance (finite increment) [ohm m^2 = kg m^4/(A^2 s^3)]
$ R_{body} $ = Body resistance [ohm m^2 = kg m^4/(A^2 s^3)]
$ \rho $ = Isotropic film resistivity [ohm m = kg m^3/(A^2 s^3)]
$ h $ = Film thickness [m]
$ C_{eff} $ = Volumetric Coulombic efficiency [m^3/(A s)]
$ j_{min} $ = Minimum current density for deposition onset [A/m^2]
$ \sigma $ = Isotropic conductivity of mixture [S/m = A^2 s^3/(kg m^3)]
$ \vec{n} $ = Patch-normal unit vector [-]
Usage
Example of the boundary condition specification:
<patchName>
{
    // Mandatory entries
    type                    electrostaticDeposition;
    h                       <scalarField>;
    CoulombicEfficiency     <PatchFunction1>;
    resistivity             <PatchFunction1>;

    // Conditional mandatory entries

        // Option-1: single-phase
        sigma       <scalar>;

        // Option-2: multiphase
        phases
        {
            alpha.air
            {
                sigma       <scalar>;
            }
            alpha.water
            {
                sigma       <scalar>;
            }
            alpha.mercury
            {
                sigma       <scalar>;
            }
            ...
        }

    // Optional entries
    jMin                    <scalar>;
    qMin                    <scalar>;
    Rbody                   <scalar>;
    Vi                      <scalar>;
    Vanode                  <scalar>;
    qCumulative             <scalarField>;

    // Inherited entries
    ...
}

where the entries mean:

Property Description Type Reqd Deflt
type Type name: electrostaticDeposition word yes -
h Film thickness scalarField yes -
CoulombicEfficiency Coulombic efficiency PatchFunction1<scalar> yes -
resistivity Isotropic film resistivity PatchFunction1<scalar> yes -
sigma Isotropic electrical conductivity of phase scalar yes -
jMin Minimum current density for deposition onset scalar no 0
qMin Minimum accumulative specific charge for deposition onset scalar no 0
Rbody Resistance due to main body and/or pretreatment layers scalar no 0
Vi Initial electric potential scalar no 0
Vanode Anode electric potential scalar no GREAT
qCumulative Accumulative specific charge [A s/m^2] scalarField no 0

The inherited entries are elaborated in:

Note
  • Depletion or abrasion of material due to negative current is not allowed.
  • When accumulative specific charge (qCumulative) is less than minimum accumulative specific charge (qMin), no deposition occurs.
  • Boundary-condition updates are not allowed during outer corrections to prevent spurious accumulation of film thickness.
  • resistivity, jMin, qMin and Rbody are always non-negative.
Source files

Definition at line 315 of file electrostaticDepositionFvPatchScalarField.H.

Constructor & Destructor Documentation

◆ electrostaticDepositionFvPatchScalarField() [1/5]

electrostaticDepositionFvPatchScalarField ( const fvPatch p,
const DimensionedField< scalar, volMesh > &  iF 
)

Construct from patch and internal field.

Definition at line 134 of file electrostaticDepositionFvPatchScalarField.C.

Referenced by electrostaticDepositionFvPatchScalarField::clone().

Here is the caller graph for this function:

◆ electrostaticDepositionFvPatchScalarField() [2/5]

electrostaticDepositionFvPatchScalarField ( const fvPatch p,
const DimensionedField< scalar, volMesh > &  iF,
const dictionary dict 
)

◆ electrostaticDepositionFvPatchScalarField() [3/5]

Construct by mapping given electrostaticDepositionFvPatchScalarField onto a new patch.

Definition at line 261 of file electrostaticDepositionFvPatchScalarField.C.

◆ electrostaticDepositionFvPatchScalarField() [4/5]

Construct as copy.

Definition at line 291 of file electrostaticDepositionFvPatchScalarField.C.

◆ electrostaticDepositionFvPatchScalarField() [5/5]

Construct as copy setting internal field reference.

Definition at line 318 of file electrostaticDepositionFvPatchScalarField.C.

Member Function Documentation

◆ TypeName()

TypeName ( "electrostaticDeposition"  )

Runtime type information.

◆ clone() [1/2]

virtual tmp<fvPatchScalarField> clone ( ) const
inlinevirtual

Construct and return a clone.

Definition at line 495 of file electrostaticDepositionFvPatchScalarField.H.

References electrostaticDepositionFvPatchScalarField::electrostaticDepositionFvPatchScalarField().

Here is the call graph for this function:

◆ clone() [2/2]

virtual tmp<fvPatchScalarField> clone ( const DimensionedField< scalar, volMesh > &  iF) const
inlinevirtual

Construct and return a clone setting internal field reference.

Definition at line 516 of file electrostaticDepositionFvPatchScalarField.H.

References electrostaticDepositionFvPatchScalarField::electrostaticDepositionFvPatchScalarField().

Here is the call graph for this function:

◆ h()

const scalarField& h ( ) const
inlinenoexcept

Return const access to film thickness patch field.

Definition at line 534 of file electrostaticDepositionFvPatchScalarField.H.

◆ autoMap()

void autoMap ( const fvPatchFieldMapper m)
virtual

Map (and resize as needed) from self given a mapping object.

Definition at line 347 of file electrostaticDepositionFvPatchScalarField.C.

◆ rmap()

void rmap ( const fvPatchScalarField ptf,
const labelList addr 
)
virtual

Reverse map the given fvPatchField onto this fvPatchField.

Definition at line 370 of file electrostaticDepositionFvPatchScalarField.C.

◆ sigma()

Foam::tmp< Foam::scalarField > sigma ( ) const

Return the isotropic electrical conductivity field of mixture.

Definition at line 397 of file electrostaticDepositionFvPatchScalarField.C.

References tmp< T >::New(), Foam::foamVersion::patch, and tmp< T >::ref().

Here is the call graph for this function:

◆ updateCoeffs()

void updateCoeffs ( )
virtual

Update the coefficients associated with the patch field.

Definition at line 419 of file electrostaticDepositionFvPatchScalarField.C.

References Foam::endl(), forAll, Foam::gAverage(), Foam::gMax(), Foam::gMin(), Foam::Info, UPstream::master(), Foam::max(), Foam::min(), Foam::nl, Foam::operator==(), Foam::foamVersion::patch, tmp< T >::ref(), sigma(), and timeIndex.

Here is the call graph for this function:

◆ write()

void write ( Ostream os) const
virtual

The documentation for this class was generated from the following files: