examples.levelSet.electroChem package¶
Submodules¶
examples.levelSet.electroChem.adsorbingSurfactantEquation module¶
examples.levelSet.electroChem.adsorption module¶
This example tests 1D adsorption onto an interface and subsequent depletion from the bulk. The governing equations are given by,
and on the interface
There is a dimensionless number that governs whether the system is in
an interface limited (
) or diffusion limited (
)
regime. There are analytical solutions for both regimes. The dimensionless
number is given by:
The test solution provided here is for the case of interface limited kinetics. The analytical solutions are given by,
and
Make sure the dimensionless parameter is large enough
>>> (diffusion / cinf / L / L / rateConstant) > 100
True
Start time stepping:
>>> currentTime = 0.
>>> from builtins import range
>>> for i in range(totalTimeSteps):
... surfEqn.solve(surfactantVar, dt = dt)
... bulkEqn.solve(bulkVar, dt = dt)
... currentTime += dt
Compare the analytical and numerical results:
>>> theta = surfactantVar.interfaceVar[1]
>>> numerix.allclose(currentTimeFunc(theta), currentTime, rtol = 1e-4)()
1
>>> numerix.allclose(concentrationFunc(theta), bulkVar[1:], rtol = 1e-4)()
1
examples.levelSet.electroChem.gapFillDistanceVariable module¶
examples.levelSet.electroChem.gapFillMesh module¶
The gapFillMesh function glues 3 meshes together to form a composite mesh. The first mesh is a Grid2D object that is fine and deals with the area around the trench or via. The second mesh is a Gmsh2D object that forms a transition mesh from a fine to a course region. The third mesh is another Grid2D object that forms the boundary layer. This region consists of very large elements and is only used for the diffusion in the boundary layer.
examples.levelSet.electroChem.gold module¶
Model electrochemical superfill of gold using the CEAC mechanism.
This input file is a demonstration of the use of FiPy for modeling gold superfill. The material properties and experimental parameters used are roughly those that have been previously published [30].
To run this example from the base FiPy directory type:
$ python examples/levelSet/electroChem/gold.py
at the command line. The results of the simulation will be displayed and the word
finished
in the terminal at the end of the simulation. The simulation will
only run for 10 time steps. To run with a different number of time steps change
the numberOfSteps
argument as follows,
>>> runGold(numberOfSteps=10, displayViewers=False)
1
Change the displayViewers
argument to True
if you wish to see the
results displayed on the screen. This example has a more realistic
default boundary layer depth and thus requires gmsh to construct a
more complex mesh.
There are a few differences between the gold superfill model presented
in this example and in examples.levelSet.electroChem.simpleTrenchSystem
.
Most default
values have changed to account for a different metal ion (gold)
and catalyst (lead). In this system the catalyst is not present in
the electrolyte but instead has a non-zero initial coverage. Thus
quantities associated with bulk catalyst and catalyst accumulation
are not defined. The current density is given by,
The
more common representation of the current density includes an
exponential part. Here it is buried in and
. The
governing equation for catalyst evolution includes a term for
catalyst consumption on the interface and is given by
where is the consumption coefficient
(
consumptionRateConstant
). The trench geometry is also given a
slight taper, given by taperAngle
.
If the Mayavi plotting software is installed (see Installation) then a plot should appear that is updated every 10 time steps and will eventually resemble the image below.
examples.levelSet.electroChem.howToWriteAScript module¶
Tutorial for writing an electrochemical superfill script.
This input file demonstrates how to
create a new superfill script if the existing suite of scripts do
not meet the required needs. It provides the functionality of
examples.levelSet.electroChem.simpleTrenchSystem
.
To run this example from the base FiPy directory type:
$ python examples/levelSet/electroChem/howToWriteAScript.py --numberOfElements=10000 --numberOfSteps=800
at the command line. The results of the simulation will be displayed
and the word finished
in the terminal at the end of the
simulation. To obtain this example in a plain script file in order to
edit and run type:
$ python setup.py copy_script --From examples/levelSet/electroChem/howToWriteAScript.py --To myScript.py
in the base FiPy directory. The file myScript.py
will contain the
script.
The following is an explicit explanation of the input commands
required to set up and run the problem. At the top of the file all the
parameter values are set. Their use will be explained during the
instantiation of various objects and are the same as those explained in
examples.levelSet.electroChem.simpleTrenchSystem
.
The following parameters (all in S.I. units) represent,
physical constants,
>>> faradaysConstant = 9.6e4 >>> gasConstant = 8.314 >>> transferCoefficient = 0.5properties associated with the catalyst species,
>>> rateConstant0 = 1.76 >>> rateConstant3 = -245e-6 >>> catalystDiffusion = 1e-9 >>> siteDensity = 9.8e-6properties of the cupric ions,
>>> molarVolume = 7.1e-6 >>> charge = 2 >>> metalDiffusionCoefficient = 5.6e-10parameters dependent on experimental constraints,
>>> temperature = 298. >>> overpotential = -0.3 >>> bulkMetalConcentration = 250. >>> catalystConcentration = 5e-3 >>> catalystCoverage = 0.parameters obtained from experiments on flat copper electrodes,
>>> currentDensity0 = 0.26 >>> currentDensity1 = 45.general simulation control parameters,
>>> cflNumber = 0.2 >>> numberOfCellsInNarrowBand = 10 >>> cellsBelowTrench = 10 >>> cellSize = 0.1e-7parameters required for a trench geometry,
>>> trenchDepth = 0.5e-6 >>> aspectRatio = 2. >>> trenchSpacing = 0.6e-6 >>> boundaryLayerDepth = 0.3e-6
The hydrodynamic boundary layer depth (boundaryLayerDepth
) is
intentionally small in this example to keep the mesh at a reasonable
size.
Build the mesh:
>>> from fipy.tools.parser import parse
>>> numberOfElements = parse('--numberOfElements', action='store',
... type='int', default=-1)
>>> numberOfSteps = parse('--numberOfSteps', action='store',
... type='int', default=2)
>>> from fipy import *
>>> if numberOfElements != -1:
... pos = trenchSpacing * cellsBelowTrench / 4 / numberOfElements
... sqr = trenchSpacing * (trenchDepth + boundaryLayerDepth) \
... / (2 * numberOfElements)
... cellSize = pos + numerix.sqrt(pos**2 + sqr)
... else:
... cellSize = 0.1e-7
>>> yCells = cellsBelowTrench \
... + int((trenchDepth + boundaryLayerDepth) / cellSize)
>>> xCells = int(trenchSpacing / 2 / cellSize)
>>> from .metalIonDiffusionEquation import buildMetalIonDiffusionEquation
>>> from .adsorbingSurfactantEquation import AdsorbingSurfactantEquation
>>> from fipy import serialComm
>>> mesh = Grid2D(dx=cellSize,
... dy=cellSize,
... nx=xCells,
... ny=yCells,
... communicator=serialComm)
A distanceVariable
object,
, is required to store the position of the interface.
The distanceVariable
calculates its value so that it is a distance
function
(i.e. holds the distance at any point in the mesh from the electrolyte/metal
interface at ) and
.
First, create the variable, which is initially set to -1 everywhere.
Create an initial variable,
>>> narrowBandWidth = numberOfCellsInNarrowBand * cellSize
>>> distanceVar = DistanceVariable(
... name='distance variable',
... mesh= mesh,
... value=-1.,
... hasOld=1)
The electrolyte region will be the positive region of the domain while the metal region will be negative.
>>> bottomHeight = cellsBelowTrench * cellSize
>>> trenchHeight = bottomHeight + trenchDepth
>>> trenchWidth = trenchDepth / aspectRatio
>>> sideWidth = (trenchSpacing - trenchWidth) / 2
>>> x, y = mesh.cellCenters
>>> distanceVar.setValue(1., where=(y > trenchHeight)
... | ((y > bottomHeight)
... & (x < xCells * cellSize - sideWidth)))
>>> distanceVar.calcDistanceFunction(order=2)
The distanceVariable
has now been created to mark the interface. Some other
variables need to be created that govern the concentrations of various species.
Create the catalyst surfactant coverage, , variable.
This variable influences the deposition rate.
>>> catalystVar = SurfactantVariable(
... name="catalyst variable",
... value=catalystCoverage,
... distanceVar=distanceVar)
Create the bulk catalyst concentration, ,
in the electrolyte,
>>> bulkCatalystVar = CellVariable(
... name='bulk catalyst variable',
... mesh=mesh,
... value=catalystConcentration)
Create the bulk metal ion concentration,
, in the electrolyte.
>>> metalVar = CellVariable(
... name='metal variable',
... mesh=mesh,
... value=bulkMetalConcentration)
The following commands build the depositionRateVariable
,
. The
depositionRateVariable
is given by the following equation.
where is the metal molar volume,
is the metal ion
charge and
is Faraday’s constant. The current density is given
by
where is the metal ion concentration in the bulk at the
interface,
is the far-field bulk concentration of
metal ions,
is the transfer coefficient,
is the gas
constant,
is the temperature and
is the
overpotential. The exchange current density is an empirical
function of catalyst coverage,
The commands needed to build this equation are,
>>> expoConstant = -transferCoefficient * faradaysConstant \
... / (gasConstant * temperature)
>>> tmp = currentDensity1 \
... * catalystVar.interfaceVar
>>> exchangeCurrentDensity = currentDensity0 + tmp
>>> expo = numerix.exp(expoConstant * overpotential)
>>> currentDensity = expo * exchangeCurrentDensity * metalVar \
... / bulkMetalConcentration
>>> depositionRateVariable = currentDensity * molarVolume \
... / (charge * faradaysConstant)
Build the extension velocity variable . The extension
velocity uses the
extensionEquation
to spread the velocity at the interface to the
rest of the domain.
>>> extensionVelocityVariable = CellVariable(
... name='extension velocity',
... mesh=mesh,
... value=depositionRateVariable)
Using the variables created above the governing equations will be built. The governing equation for surfactant conservation is given by,
where is the coverage of catalyst at the interface,
is the curvature of the interface,
is the normal velocity
of the interface,
is the concentration of
catalyst in the bulk at the interface. The value
is given
by an empirical function of overpotential,
The above equation is represented by the
AdsorbingSurfactantEquation
in FiPy:
>>> surfactantEquation = AdsorbingSurfactantEquation(
... surfactantVar=catalystVar,
... distanceVar=distanceVar,
... bulkVar=bulkCatalystVar,
... rateConstant=rateConstant0 \
... + rateConstant3 * overpotential**3)
The variable is advected by the
advectionEquation
given by,
and is set up with the following commands:
>>> advectionEquation = TransientTerm() + AdvectionTerm(extensionVelocityVariable)
The diffusion of metal ions from the far field to the interface is governed by,
where,
The following boundary condition applies at ,
The metal ion diffusion equation is set up with the following commands.
>>> metalEquation = buildMetalIonDiffusionEquation(
... ionVar=metalVar,
... distanceVar=distanceVar,
... depositionRate=depositionRateVariable,
... diffusionCoeff=metalDiffusionCoefficient,
... metalIonMolarVolume=molarVolume,
... )
>>> metalVar.constrain(bulkMetalConcentration, mesh.facesTop)
The surfactant bulk diffusion equation solves the bulk diffusion of a species with a source term for the jump from the bulk to an interface. The governing equation is given by,
where,
The jump condition at the interface is defined by Langmuir
adsorption. Langmuir adsorption essentially states that the ability
for a species to jump from an electrolyte to an interface is
proportional to the concentration in the electrolyte, available site
density and a jump coefficient. The boundary condition
at is given by,
The surfactant bulk diffusion equation is set up with the following commands.
>>> from .surfactantBulkDiffusionEquation import buildSurfactantBulkDiffusionEquation
>>> bulkCatalystEquation = buildSurfactantBulkDiffusionEquation(
... bulkVar=bulkCatalystVar,
... distanceVar=distanceVar,
... surfactantVar=catalystVar,
... diffusionCoeff=catalystDiffusion,
... rateConstant=rateConstant0 * siteDensity
... )
>>> bulkCatalystVar.constrain(catalystConcentration, mesh.facesTop)
If running interactively, create viewers.
>>> if __name__ == '__main__':
... try:
... from .mayaviSurfactantViewer import MayaviSurfactantViewer
... viewer = MayaviSurfactantViewer(distanceVar,
... catalystVar.interfaceVar,
... zoomFactor=1e6,
... datamax=1.0,
... datamin=0.0,
... smooth=1)
... except:
... viewer = MultiViewer(viewers=(
... Viewer(distanceVar, datamin=-1e-9, datamax=1e-9),
... Viewer(catalystVar.interfaceVar)))
... from fipy.models.levelSet.surfactant.matplotlibSurfactantViewer import MatplotlibSurfactantViewer
... viewer = MatplotlibSurfactantViewer(catalystVar.interfaceVar)
... else:
... viewer = None
The levelSetUpdateFrequency
defines how often to call the
distanceEquation
to reinitialize the distanceVariable
to a
distance function.
>>> levelSetUpdateFrequency = int(0.8 * narrowBandWidth \
... / (cellSize * cflNumber * 2))
The following loop runs for numberOfSteps
time steps. The time step
is calculated with the CFL number and the maximum extension velocity.
to
throughout the whole domain using
.
>>> from builtins import range
>>> for step in range(numberOfSteps):
...
... if viewer is not None:
... viewer.plot()
...
... if step % levelSetUpdateFrequency == 0:
... distanceVar.calcDistanceFunction(order=2)
...
... extensionVelocityVariable.setValue(depositionRateVariable())
...
... distanceVar.updateOld()
... distanceVar.extendVariable(extensionVelocityVariable, order=2)
... dt = cflNumber * cellSize / extensionVelocityVariable.max()
... advectionEquation.solve(distanceVar, dt=dt)
... surfactantEquation.solve(catalystVar, dt=dt)
... metalEquation.solve(var=metalVar, dt=dt)
... bulkCatalystEquation.solve(var=bulkCatalystVar, dt=dt, solver=GeneralSolver())
The following is a short test case. It uses saved data from a simulation with 5 time steps. It is not a test for accuracy but a way to tell if something has changed or been broken.
>>> import os
>>> filepath = os.path.join(os.path.split(__file__)[0],
... "simpleTrenchSystem.gz")
>>> ##numerix.savetxt(filepath, numerix.array(catalystVar))
>>> print(catalystVar.allclose(numerix.loadtxt(filepath), rtol=1e-4))
1
>>> from fipy import input
>>> if __name__ == '__main__':
... input('finished')
examples.levelSet.electroChem.leveler module¶
Model electrochemical superfill of copper with leveler and accelerator additives.
This input file is a demonstration of the use of FiPy for modeling copper superfill with leveler and accelerator additives. The material properties and experimental parameters used are roughly those that have been previously published [31].
To run this example from the base FiPy directory type:
$ python examples/levelSet/electroChem/leveler.py
at the command line. The results of the simulation will be displayed and the word
finished
in the terminal at the end of the simulation. The simulation will
only run for 200 time steps. To run with a different number of time steps change
the numberOfSteps
argument as follows,
>>> runLeveler(numberOfSteps=10, displayViewers=False, cellSize=0.25e-7)
1
Change the displayViewers
argument to True
if you wish to see the
results displayed on the screen. This example requires gmsh to
construct the mesh.
This example models the case when suppressor, accelerator and leveler additives are present in the electrolyte. The suppressor is assumed to absorb quickly compared with the other additives. Any unoccupied surface sites are immediately covered with suppressor. The accelerator additive has more surface affinity than suppressor and is thus preferential adsorbed. The accelerator can also remove suppressor when the surface reaches full coverage. Similarly, the leveler additive has more surface affinity than both the suppressor and accelerator. This forms a simple set of assumptions for understanding the behavior of these additives.
The following is a complete description of the equations for the
model described here. Any equations that have been omitted are the
same as those given in
examples.levelSet.electroChem.simpleTrenchSystem
. The
current density is governed by
where represents
for suppressor,
for accelerator,
for leveler and
for vacant. This model
assumes a linear interpolation between the three cases of complete
coverage for each additive or vacant substrate. The governing
equations for the surfactants are given by,
It has been found experimentally that .
If the surface reaches full coverage, the equations do not naturally prevent the
coverage rising above full coverage due to the curvature terms. Thus, when
then the equation for accelerator becomes
and when
, the
equation for leveler becomes
The parameters ,
and
are both functions of
given by,
The following table shows the symbols used in the governing equations
and their corresponding arguments for the
runLeveler()
function.
The following images show accelerator and leveler contour plots that can be obtained by running this example.
examples.levelSet.electroChem.lines module¶
examples.levelSet.electroChem.matplotlibSurfactantViewer module¶
- class examples.levelSet.electroChem.matplotlibSurfactantViewer.MatplotlibSurfactantViewer(distanceVar, surfactantVar=None, levelSetValue=0.0, title=None, smooth=0, zoomFactor=1.0, animate=False, limits={}, **kwlimits)¶
Bases:
AbstractMatplotlibViewer
The MatplotlibSurfactantViewer creates a viewer with the Matplotlib python plotting package that displays a DistanceVariable.
Create a MatplotlibSurfactantViewer.
>>> from fipy import * >>> m = Grid2D(nx=100, ny=100) >>> x, y = m.cellCenters >>> v = CellVariable(mesh=m, value=x**2 + y**2 - 10**2) >>> s = CellVariable(mesh=m, value=sin(x / 10) * cos(y / 30)) >>> viewer = MatplotlibSurfactantViewer(distanceVar=v, surfactantVar=s) >>> from builtins import range >>> for r in range(1, 200): ... v.setValue(x**2 + y**2 - r**2) ... viewer.plot()
>>> from fipy import * >>> dx = 1. >>> dy = 1. >>> nx = 11 >>> ny = 11 >>> Lx = ny * dy >>> Ly = nx * dx >>> mesh = Grid2D(dx = dx, dy = dy, nx = nx, ny = ny) >>> # from fipy.models.levelSet.distanceFunction.distanceVariable import DistanceVariable >>> var = DistanceVariable(mesh = mesh, value = -1)
>>> x, y = mesh.cellCenters
>>> var.setValue(1, where=(x - Lx / 2.)**2 + (y - Ly / 2.)**2 < (Lx / 4.)**2) >>> var.calcDistanceFunction() >>> viewer = MatplotlibSurfactantViewer(var, smooth = 2) >>> viewer.plot() >>> viewer._promptForOpinion() >>> del viewer
>>> var = DistanceVariable(mesh = mesh, value = -1)
>>> var.setValue(1, where=(y > 2. * Ly / 3.) | ((x > Lx / 2.) & (y > Ly / 3.)) | ((y < Ly / 6.) & (x > Lx / 2))) >>> var.calcDistanceFunction() >>> viewer = MatplotlibSurfactantViewer(var) >>> viewer.plot() >>> viewer._promptForOpinion() >>> del viewer
>>> viewer = MatplotlibSurfactantViewer(var, smooth = 2) >>> viewer.plot() >>> viewer._promptForOpinion() >>> del viewer
- Parameters:
distanceVar (DistanceVariable) –
levelSetValue (float) – the value of the contour to be displayed
title (str) – displayed at the top of the Viewer window
animate (bool) – whether to show only the initial condition and the moving top boundary or to show all contours (Default)
limits (dict, optional) – a (deprecated) alternative to limit keyword arguments
xmin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
xmax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
ymin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
ymax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
zmin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
zmax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
datamin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
datamax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
- property axes¶
The Matplotlib
Axes
.
- property cmap¶
The Matplotlib
Colormap
.
- property colorbar¶
The Matplotlib
Colorbar
.
- property fig¶
The Matplotlib
Figure
.
- figaspect(figaspect, colorbar)¶
- property id¶
The Matplotlib
Figure
number.
- property limits¶
- plot(filename=None)¶
Update the display of the viewed variables.
- Parameters:
filename (str) – If not None, the name of a file to save the image into.
- plotMesh(filename=None)¶
Display a representation of the mesh
- Parameters:
filename (str) – If not None, the name of a file to save the image into.
- setLimits(limits={}, **kwlimits)¶
Update the limits.
- Parameters:
limits (dict, optional) – a (deprecated) alternative to limit keyword arguments
xmin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
xmax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
ymin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
ymax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
zmin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
zmax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
datamin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
datamax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
- property title¶
The text appearing at the top center.
(default: if
len(self.vars) == 1
, the name ofself.vars[0]
, otherwise""
.)
examples.levelSet.electroChem.mayaviSurfactantViewer module¶
- class examples.levelSet.electroChem.mayaviSurfactantViewer.MayaviSurfactantViewer(distanceVar, surfactantVar=None, levelSetValue=0.0, title=None, smooth=0, zoomFactor=1.0, animate=False, limits={}, **kwlimits)¶
Bases:
AbstractViewer
The MayaviSurfactantViewer creates a viewer with the Mayavi python plotting package that displays a DistanceVariable.
Create a MayaviSurfactantViewer.
>>> from fipy import * >>> dx = 1. >>> dy = 1. >>> nx = 11 >>> ny = 11 >>> Lx = ny * dy >>> Ly = nx * dx >>> mesh = Grid2D(dx = dx, dy = dy, nx = nx, ny = ny) >>> # from fipy.models.levelSet.distanceFunction.distanceVariable import DistanceVariable >>> var = DistanceVariable(mesh = mesh, value = -1)
>>> x, y = mesh.cellCenters
>>> var.setValue(1, where=(x - Lx / 2.)**2 + (y - Ly / 2.)**2 < (Lx / 4.)**2) >>> var.calcDistanceFunction() >>> viewer = MayaviSurfactantViewer(var, smooth = 2) >>> viewer.plot() >>> viewer._promptForOpinion() >>> del viewer
>>> var = DistanceVariable(mesh = mesh, value = -1)
>>> var.setValue(1, where=(y > 2. * Ly / 3.) | ((x > Lx / 2.) & (y > Ly / 3.)) | ((y < Ly / 6.) & (x > Lx / 2))) >>> var.calcDistanceFunction() >>> viewer = MayaviSurfactantViewer(var) >>> viewer.plot() >>> viewer._promptForOpinion() >>> del viewer
>>> viewer = MayaviSurfactantViewer(var, smooth = 2) >>> viewer.plot() >>> viewer._promptForOpinion() >>> del viewer
- Parameters:
distanceVar (DistanceVariable) –
levelSetValue (float) – the value of the contour to be displayed
title (str) – displayed at the top of the Viewer window
animate (bool) – whether to show only the initial condition and the moving top boundary or to show all contours (Default)
limits (dict, optional) – a (deprecated) alternative to limit keyword arguments
xmin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
xmax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
ymin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
ymax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
zmin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
zmax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
datamin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
datamax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
- property limits¶
- plotMesh(filename=None)¶
Display a representation of the mesh
- Parameters:
filename (str) – If not None, the name of a file to save the image into.
- setLimits(limits={}, **kwlimits)¶
Update the limits.
- Parameters:
limits (dict, optional) – a (deprecated) alternative to limit keyword arguments
xmin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
xmax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
ymin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
ymax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
zmin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
zmax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
datamin (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
datamax (float, optional) – displayed range of data. A 1D Viewer will only use xmin and xmax, a 2D viewer will also use ymin and ymax, and so on. All viewers will use datamin and datamax. Any limit set to a (default) value of None will autoscale.
- property title¶
The text appearing at the top center.
(default: if
len(self.vars) == 1
, the name ofself.vars[0]
, otherwise""
.)
examples.levelSet.electroChem.metalIonDiffusionEquation module¶
examples.levelSet.electroChem.simpleTrenchSystem module¶
Model electrochemical superfill using the CEAC mechanism.
This input file is a demonstration of the use of FiPy for modeling electrodeposition using the CEAC mechanism. The material properties and experimental parameters used are roughly those that have been previously published [32].
To run this example from the base FiPy directory type:
$ python examples/levelSet/electroChem/simpleTrenchSystem.py
at the command line. The results of the simulation will be displayed and the word
finished in the terminal at the end of the simulation. To run with a different
number of time steps change the numberOfSteps
argument as follows,
>>> runSimpleTrenchSystem(numberOfSteps=2, displayViewers=False)
1
Change the displayViewers
argument to True
if you wish to see the
results displayed on the
screen. Example examples.levelSet.electroChem.simpleTrenchSystem
gives explanation for
writing new scripts or modifying existing scripts that are
encapsulated by functions.
Any argument parameter can be changed. For example if the initial catalyst coverage is not 0, then it can be reset,
>>> runSimpleTrenchSystem(numberOfSteps=2, catalystCoverage=0.1, displayViewers=False)
0
The following image shows a schematic of a trench geometry along with
the governing equations for modeling electrodeposition with the CEAC
mechanism. All of the given equations are implemented in the
examples.levelSet.electroChem.simpleTrenchSystem.runSimpleTrenchSystem()
function. As stated above, all the parameters
in the equations can be changed with function arguments.
The following table shows the symbols used in the governing equations
and their corresponding arguments to the
runSimpleTrenchSystem()
function. The boundary layer depth is intentionally small in this
example in order not to complicate the mesh. Further examples will
simulate more realistic boundary layer depths but will also have more
complex meshes requiring the gmsh software.
If the Mayavi plotting software is installed (see Installation) then a plot should appear that is updated every 20 time steps and will eventually resemble the image below.