chombo-discharge
Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes | List of all members
FieldSolver Class Referenceabstract

Base class for electrostatic solvers. More...

#include <CD_FieldSolver.H>

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

Public Member Functions

 FieldSolver ()
 Constructor.
 
 FieldSolver (const FieldSolver &a_other)=delete
 Disallowed copy constructor. More...
 
 FieldSolver (const FieldSolver &&a_other)=delete
 Disallowed move constructor. More...
 
FieldSolveroperator= (const FieldSolver &a_other)=delete
 Disallowed copy assignment operator. More...
 
FieldSolveroperator= (const FieldSolver &&a_other)=delete
 Disallowed move assignment operator. More...
 
virtual ~FieldSolver ()
 Constructor.
 
virtual void setupSolver ()=0
 Set up solver routines.
 
virtual void setSolverPermittivities (const MFAMRCellData &a_permittivityCell, const MFAMRFluxData &a_permittivityFace, const MFAMRIVData &a_permittivityEB)
 A special routine for when solver permittivities need to change but solver does not. More...
 
virtual bool solve (const bool a_zeroPhi=false)
 Solve Poisson equation using m_potential, m_rho, and m_sigma. More...
 
virtual bool solve (MFAMRCellData &a_potential, const bool a_zerophi=false)
 Solve Poisson equation onto a_phi using m_rho, and m_sigma as right-hand sides. More...
 
virtual bool solve (MFAMRCellData &a_phi, const MFAMRCellData &a_rho, const EBAMRIVData &a_sigma, const bool a_zerophi=false)=0
 Solves Poisson equation onto a_phi using a_rho and a_sigma as right-hand sides. More...
 
virtual void computeElectricField ()
 Compute the cell-centered electric field. More...
 
virtual void computeElectricField (MFAMRCellData &a_E, const MFAMRCellData &a_potential) const =0
 Compute the cell-centered electric field. More...
 
virtual void computeElectricField (MFAMRFluxData &a_E, const MFAMRCellData &a_potential) const =0
 Compute the face-centered electric field. More...
 
virtual void computeElectricField (EBAMRCellData &a_E, const phase::which_phase a_phase, const MFAMRCellData &a_potential) const =0
 Compute the cell-centered electric field on a specific phase. More...
 
virtual void computeElectricField (EBAMRFluxData &a_E, const phase::which_phase a_phase, const MFAMRCellData &a_potential) const =0
 Compute the face-centered electric field on a specific phase. More...
 
virtual void parseOptions ()=0
 Parse options (for derived class)
 
virtual void parseRuntimeOptions ()=0
 Parse runtime options (for derived class)
 
virtual void allocate ()
 Allocates internal storage for FieldSolver. Derived classes may want to overwrite.
 
virtual void preRegrid (const int a_lbase, const int a_oldFinestLevel)
 Cache state before regridding. More...
 
virtual void computeDisplacementField (MFAMRCellData &a_displacementField, const MFAMRCellData &a_electricField)
 Compute displacement field from the electric field. More...
 
virtual void deallocate ()
 Deallocate internal storage.
 
virtual void regrid (const int a_lmin, const int a_oldFinestLevel, const int a_newFinestLevel)
 Regrid method. More...
 
virtual void registerOperators ()=0
 Register operators for AMR. Derived classes have to implement these.
 
void setRho (const Real a_rho)
 Set space charge to constant value everywhere. More...
 
void setRho (const std::function< Real(const RealVect)> &a_rho)
 Set space charge to spatially dependent function. More...
 
void setSigma (const Real a_sigma)
 Set surface charge to specified value.
 
void setSigma (const std::function< Real(const RealVect)> &a_sigma)
 Set surface charge density to spatially dependent function. More...
 
void setComputationalGeometry (const RefCountedPtr< ComputationalGeometry > &a_computationalGeometry)
 Set the computational geometry. More...
 
void setAmr (const RefCountedPtr< AmrMesh > &a_amr)
 Set the amr object. More...
 
virtual void setPermittivities ()
 Set the permittivities. More...
 
virtual void writePlotFile ()
 Write plot file. More...
 
virtual void postCheckpoint ()
 Write checkpoint data for a level @paramo[out] a_handle HDF5 handle. More...
 
virtual void writePlotData (LevelData< EBCellFAB > &a_output, int &a_comp, const std::string a_outputRealm, const int a_level, const bool a_forceNoInterp=false) const noexcept
 Write output data to a_output. More...
 
virtual void writeMultifluidData (LevelData< EBCellFAB > &a_output, int &a_comp, const MFAMRCellData &a_data, const phase::which_phase a_phase, const std::string a_outputRealm, const int a_level, const bool a_interp) const noexcept
 Write multifluid data to single-fluid data holders. More...
 
virtual void writeSurfaceData (LevelData< EBCellFAB > &a_output, int &a_comp, const LevelData< BaseIVFAB< Real >> &a_data, const std::string a_outputRealm, const int a_level) const noexcept
 Write surface data to volume data holder. More...
 
void setRealm (const std::string a_realm)
 Set Realm. More...
 
void setTime (const int a_timeStep, const Real a_time, const Real a_dt)
 Set time for this solver. More...
 
void setVerbosity (const int a_verbosity)
 Set verbosity. More...
 
virtual void setVoltage (std::function< Real(const Real a_time)> a_voltage)
 Set potential dependence in time. More...
 
virtual void setDomainSideBcFunction (const int a_dir, const Side::LoHiSide a_side, const ElectrostaticDomainBc::BcFunction &a_function)
 Boundary condition function on a domain side. More...
 
virtual void setElectrodeDirichletFunction (const int a_electrode, const ElectrostaticEbBc::BcFunction &a_function)
 Set embedded boundary Dirichlet function on a specific electrode. More...
 
virtual int getNumberOfPlotVariables () const
 Get number of output fields. More...
 
const std::function< Real(const Real a_time)> & getVoltage () const
 Get voltage function. More...
 
Real getCurrentVoltage () const
 Get current voltage. More...
 
Real getTime () const
 Get time. More...
 
Real computeCapacitance ()
 Compute the capacitance. More...
 
Real computeEnergy (const MFAMRCellData &a_electricField)
 Compute energy density U = 0.5*int(E.dot.D dV) More...
 
virtual Vector< std::string > getPlotVariableNames () const
 Get output plot names. More...
 
virtual Vector< long long > computeLoads (const DisjointBoxLayout &a_dbl, const int a_level)
 Get computational loads for a level. More...
 
void setDataLocation (const Location::Cell a_dataLocation)
 Set the data location for the solver. More...
 
std::string getRealm () const
 Get the Realm where this solver is registered. More...
 
MFAMRCellDatagetPotential ()
 Get potential on both phases.
 
MFAMRCellDatagetElectricField ()
 Get electric field on both phases.
 
MFAMRCellDatagetRho ()
 Get storage for the space charge density. More...
 
MFAMRCellDatagetResidue ()
 Get the residue.
 
MFAMRCellDatagetPermittivityCell ()
 Get cell-centered permittivity.
 
MFAMRFluxDatagetPermittivityFace ()
 Get face-centered permittivity.
 
MFAMRIVDatagetPermittivityEB ()
 Get irregular b coefficient.
 
EBAMRIVDatagetSigma ()
 Get m_sigma.
 
Location::Cell getDataLocation () const
 Get data location.
 

Protected Member Functions

virtual void parseVerbosity ()
 Parse solver class verbosity.
 
virtual void parsePlotVariables ()
 Function which parses which plot variables to write to plot files.
 
virtual void parseDomainBc ()
 Parse domain boundary conditions.
 
virtual void parseRegridSlopes ()
 Parse slope regrid.
 
virtual void setDefaultDomainBcFunctions ()
 Set default BC functions. This sets all the m_domainBcFunction objects to s_defaultDomainBcFunction, which return 1 everywhere.
 
virtual void setDefaultEbBcFunctions ()
 Set default Dirichlet boundary conditions on the embedded boundaries. More...
 
Real getDielectricPermittivity (const RealVect &a_pos) const
 Get relative permittivity at some point in space. More...
 
virtual void setCellPermittivities (EBCellFAB &a_perm, const Box &a_cellBox, const EBISBox &a_ebisbox, const RealVect &a_probLo, const Real &a_dx)
 Set cell-centered permittivities. More...
 
virtual void setFacePermittivities (EBFluxFAB &a_perm, const Box &a_cellBox, const EBISBox &a_ebisbox, const RealVect &a_probLo, const Real &a_dx)
 Set face-centered permittivities. More...
 
virtual void setEbPermittivities (BaseIVFAB< Real > &a_perm, const Box &a_cellBox, const EBISBox &a_ebisbox, const RealVect &a_origin, const Real &a_dx)
 Set EB-centered permittivities. More...
 
virtual ElectrostaticDomainBc::BcType parseBcString (const std::string a_str) const
 Returns BC type based on string. More...
 
virtual std::string makeBcString (const int a_dir, const Side::LoHiSide a_side) const
 Shortcut for making a boundary condition string. More...
 

Protected Attributes

Location::Cell m_dataLocation
 Flag which specifies that data location.
 
Location::Face m_faceLocation
 Flag which specifies where the permittivities are stored.
 
std::string m_realm
 Realm where this solver is registered.
 
std::string m_className
 Class name (i.e., "FieldSolver" for the base class)
 
RefCountedPtr< MultiFluidIndexSpacem_multifluidIndexSpace
 Multifluid index space.
 
RefCountedPtr< ComputationalGeometrym_computationalGeometry
 Computational geometry.
 
RefCountedPtr< AmrMeshm_amr
 AMR - needed for pretty much everything.
 
MFAMRCellData m_cache
 Cached state used for regridding.
 
MFAMRCellData m_potential
 State data, i.e. the potential. The centering of this is the same as m_dataLocation.
 
MFAMRCellData m_electricField
 Electric field. The centering of this is the same as m_dataLocation.
 
MFAMRCellData m_rho
 Storage for space charge density.
 
EBAMRIVData m_sigma
 Storage for surface charge density.
 
MFAMRCellData m_residue
 Residue, e.g. used after solving Poisson equation.
 
MFAMRCellData m_permittivityCell
 Cell permittivity.
 
MFAMRFluxData m_permittivityFace
 Face permittivity.
 
MFAMRIVData m_permittivityEB
 EB permittivity.
 
bool m_isVoltageSet
 Flag for checking if voltage has been set.
 
bool m_plotPotential
 If true, potential will be added to plot files.
 
bool m_plotRho
 If true, space charge will be added to plot files.
 
bool m_plotElectricField
 If true, the electric field will be added to plot files.
 
bool m_plotElectricFieldSolid
 If true, the electric field on the inside of dielectrics will be added to plot files.
 
bool m_plotResidue
 If true, the residue will be added to plot files.
 
bool m_plotPermittivity
 If true, the permittivity will be added to plot files.
 
bool m_plotSigma
 If true, m_sigma will be added to plot files.
 
bool m_regridSlopes
 Use slopes when regridding or ont.
 
int m_verbosity
 Verbosity for this calss.
 
int m_timeStep
 Current time step.
 
Real m_dt
 Last time step increment.
 
Real m_time
 Current time.
 
ElectrostaticDomainBc m_domainBc
 Domain boundary conditions for FieldSolver.
 
ElectrostaticEbBc m_ebBc
 Dirichlet boundary conditions on electrodes.
 
std::map< ElectrostaticDomainBc::DomainSide, ElectrostaticDomainBc::BcFunctionm_domainBcFunctions
 Domain BC functions. This is used to map space/time to a voltage/field at the domain faces.
 
std::vector< std::pair< Electrode, ElectrostaticEbBc::BcFunction > > m_electrodeBcFunctions
 BC functions (Dirichlet) on the electrodes. Used to map space/time to a voltage on the electrodes.
 
std::function< Real(const Real a_time)> m_voltage
 Voltage function.
 

Static Protected Attributes

constexpr static int m_comp = 0
 Component number where data is stored.
 
constexpr static int m_nComp = 1
 Number of components in data holders.
 

Detailed Description

Base class for electrostatic solvers.

This class contains an interface to a solve routine, data holders, and interface to boundary conditions. Solver routines are supposed to be implemented by derived classes.

Constructor & Destructor Documentation

◆ FieldSolver() [1/2]

FieldSolver::FieldSolver ( const FieldSolver a_other)
delete

Disallowed copy constructor.

Parameters
[in]a_otherOther solver

◆ FieldSolver() [2/2]

FieldSolver::FieldSolver ( const FieldSolver &&  a_other)
delete

Disallowed move constructor.

Parameters
[in]a_otherOther solver

Member Function Documentation

◆ computeCapacitance()

Real FieldSolver::computeCapacitance ( )

Compute the capacitance.

This will first obtain a solution without any sources, and then compute the energy density. The capacitance is then C = 2*EnergyDensity/(V*V)

Returns
Capacitance for system.

◆ computeDisplacementField()

void FieldSolver::computeDisplacementField ( MFAMRCellData a_displacementField,
const MFAMRCellData a_electricField 
)
virtual

Compute displacement field from the electric field.

Parameters
[out]a_displacementFieldDisplacement field
[in]a_electricFieldElectric field

◆ computeElectricField() [1/5]

void FieldSolver::computeElectricField ( )
virtual

Compute the cell-centered electric field.

This uses m_potential for computing the electric field and puts the result into m_electricField.

◆ computeElectricField() [2/5]

virtual void FieldSolver::computeElectricField ( EBAMRCellData a_E,
const phase::which_phase  a_phase,
const MFAMRCellData a_potential 
) const
pure virtual

Compute the cell-centered electric field on a specific phase.

Parameters
[out]a_electricFieldFace-centered electric field
[in]a_phasePhase
[in]a_potentialCell-centered potential

Implemented in FieldSolverMultigrid.

◆ computeElectricField() [3/5]

virtual void FieldSolver::computeElectricField ( EBAMRFluxData a_E,
const phase::which_phase  a_phase,
const MFAMRCellData a_potential 
) const
pure virtual

Compute the face-centered electric field on a specific phase.

Parameters
[out]a_electricFieldFace-centered electric field
[in]a_phasePhase
[in]a_potentialCell-centered potential

Implemented in FieldSolverMultigrid.

◆ computeElectricField() [4/5]

virtual void FieldSolver::computeElectricField ( MFAMRCellData a_E,
const MFAMRCellData a_potential 
) const
pure virtual

Compute the cell-centered electric field.

Parameters
[out]a_electricFieldCell-centered electric field
[in]a_potentialCell-centered potential

Implemented in FieldSolverMultigrid.

◆ computeElectricField() [5/5]

virtual void FieldSolver::computeElectricField ( MFAMRFluxData a_E,
const MFAMRCellData a_potential 
) const
pure virtual

Compute the face-centered electric field.

Parameters
[out]a_electricFieldFace-centered electric field
[in]a_potentialCell-centered potential

Implemented in FieldSolverMultigrid.

◆ computeEnergy()

Real FieldSolver::computeEnergy ( const MFAMRCellData a_electricField)

Compute energy density U = 0.5*int(E.dot.D dV)

Parameters
[in]a_electricFieldThe electric field.
Returns
Energy density

◆ computeLoads()

Vector< long long > FieldSolver::computeLoads ( const DisjointBoxLayout &  a_dbl,
const int  a_level 
)
virtual

Get computational loads for a level.

Returns
Loads for each box on a grid level.
Note
The return vector should have the same order as the boxes in a_dbl. E.g. ret[0] must be the load for a_dbl.boxArray()[0];

The default implementation returns the number of cells in the grid patch as a proxy for the load.

Reimplemented in FieldSolverMultigrid.

◆ getCurrentVoltage()

Real FieldSolver::getCurrentVoltage ( ) const

Get current voltage.

Returns
Evaluates m_voltage(m_time) and returns the result.

◆ getDielectricPermittivity()

Real FieldSolver::getDielectricPermittivity ( const RealVect &  a_pos) const
inlineprotected

Get relative permittivity at some point in space.

Parameters
[in]a_positionPhysical position
Note
This routine is used for computing the permittivity inside dielectrics.

◆ getNumberOfPlotVariables()

int FieldSolver::getNumberOfPlotVariables ( ) const
virtual

Get number of output fields.

Returns
Number of plot variables

◆ getPlotVariableNames()

Vector< std::string > FieldSolver::getPlotVariableNames ( ) const
virtual

Get output plot names.

Returns
All the plot variable names.

◆ getRealm()

std::string FieldSolver::getRealm ( ) const

Get the Realm where this solver is registered.

Returns
Returns realm (as string)

◆ getRho()

MFAMRCellData & FieldSolver::getRho ( )

Get storage for the space charge density.

This lives in the FieldSolver so users don't need to allocate their own storage for the space charge.

◆ getTime()

Real FieldSolver::getTime ( ) const

Get time.

Returns
m_time

◆ getVoltage()

const std::function< Real(const Real a_time)> & FieldSolver::getVoltage ( ) const

Get voltage function.

Returns
Returns m_voltage

◆ makeBcString()

std::string FieldSolver::makeBcString ( const int  a_dir,
const Side::LoHiSide  a_side 
) const
protectedvirtual

Shortcut for making a boundary condition string.

Parameters
[in]a_dirDirection.
[in]a_sideCoordinate side.
Returns
Returns string of type m_className.bc.direction.side.

◆ operator=() [1/2]

FieldSolver& FieldSolver::operator= ( const FieldSolver &&  a_other)
delete

Disallowed move assignment operator.

Parameters
[in]a_otherOther solver

◆ operator=() [2/2]

FieldSolver& FieldSolver::operator= ( const FieldSolver a_other)
delete

Disallowed copy assignment operator.

Parameters
[in]a_otherOther solver

◆ parseBcString()

ElectrostaticDomainBc::BcType FieldSolver::parseBcString ( const std::string  a_str) const
protectedvirtual

Returns BC type based on string.

Parameters
[in]a_strBoundary condition string. Must be "dirichlet number", "neumann number", "dirichlet_custom" or "neumann_custom"
Returns
This returns the boundary condition type, either ElectrostaticDomainBc::Dirichlet or ElectrostaticDomainBc::Neumann

◆ postCheckpoint()

void FieldSolver::postCheckpoint ( )
virtual

Write checkpoint data for a level @paramo[out] a_handle HDF5 handle.

Parameters
[in]a_levelGrid level

This writes m_potential[a_level] to the checkpoint file.

Read checkpoint data onto a level

Parameters
[in]a_handleHDF5 handle.
[in]constint a_level Grid level

This fills m_potential[a_level] with data from a_handle.

Post checkpoint method.

Default implementation does not do anything, but derived classes will usually need it.

◆ preRegrid()

void FieldSolver::preRegrid ( const int  a_lbase,
const int  a_oldFinestLevel 
)
virtual

Cache state before regridding.

Parameters
[in]a_lbaseCoarsest level which changes during regrid.
[in]a_oldFinestLevelFinest AMR level before regrid.

This allocates m_cache which gets a copy of m_potential

Reimplemented in FieldSolverMultigrid.

◆ regrid()

void FieldSolver::regrid ( const int  a_lmin,
const int  a_oldFinestLevel,
const int  a_newFinestLevel 
)
virtual

Regrid method.

Parameters
[in]a_lminCoarsest level allowed to change.
[in]a_oldFinestLevelFinest level before the regrid.
[in]a_newFinestLevelFinest level after the regrid.

This linearly interpolates (with limiters) m_potential to the new grids and recomputes the electric field (from the interpolated potential).

Reimplemented in FieldSolverMultigrid.

◆ setAmr()

void FieldSolver::setAmr ( const RefCountedPtr< AmrMesh > &  a_amr)

Set the amr object.

Parameters
[in]a_amrAmrMesh object.

◆ setCellPermittivities()

void FieldSolver::setCellPermittivities ( EBCellFAB &  a_perm,
const Box &  a_cellBox,
const EBISBox &  a_ebisbox,
const RealVect &  a_probLo,
const Real &  a_dx 
)
protectedvirtual

Set cell-centered permittivities.

Parameters
[out]a_permPermittivity (on either cell center or centroid)
[in]a_cellBoxComputation box
[in]a_ebisboxEBIS box
[in]a_probLoLower-left corner off comptuational domain
[in]a_dxResolution

◆ setComputationalGeometry()

void FieldSolver::setComputationalGeometry ( const RefCountedPtr< ComputationalGeometry > &  a_computationalGeometry)

Set the computational geometry.

Parameters
[in]a_computationalGeometryComputational geometry.

◆ setDataLocation()

void FieldSolver::setDataLocation ( const Location::Cell  a_dataLocation)

Set the data location for the solver.

Parameters
[in]a_dataLocationDesired data location.

Only Location::Cell::Center and Location::Cell::Centroid are accepted arguments.

◆ setDefaultEbBcFunctions()

void FieldSolver::setDefaultEbBcFunctions ( )
protectedvirtual

Set default Dirichlet boundary conditions on the embedded boundaries.

For each electrode the default boundary condition is set to m_voltage*electrode.getFraction().

◆ setDomainSideBcFunction()

void FieldSolver::setDomainSideBcFunction ( const int  a_dir,
const Side::LoHiSide  a_side,
const ElectrostaticDomainBc::BcFunction a_function 
)
virtual

Boundary condition function on a domain side.

Parameters
[in]a_dirCoordinate direction.
[in]a_sideSide in the coordinate direction.
[in]a_functionBoundary condition function.

This sets a boundary condition for a particular side on a wall. The user must also specify how to use this BC in the input script.

◆ setEbPermittivities()

void FieldSolver::setEbPermittivities ( BaseIVFAB< Real > &  a_perm,
const Box &  a_cellBox,
const EBISBox &  a_ebisbox,
const RealVect &  a_origin,
const Real &  a_dx 
)
protectedvirtual

Set EB-centered permittivities.

Parameters
[out]a_permPermittivity (on face centroid or center)
[in]a_cellBoxComputation box
[in]a_ebisboxEBIS box
[in]a_probLoLower-left corner off comptuational domain
[in]a_dxResolution

◆ setElectrodeDirichletFunction()

void FieldSolver::setElectrodeDirichletFunction ( const int  a_electrode,
const ElectrostaticEbBc::BcFunction a_function 
)
virtual

Set embedded boundary Dirichlet function on a specific electrode.

Parameters
[in]a_electrodeelectrode index. Follows the same order as ComputationalGeometry.
[in]a_functionVoltage on the electrode.

◆ setFacePermittivities()

void FieldSolver::setFacePermittivities ( EBFluxFAB &  a_perm,
const Box &  a_cellBox,
const EBISBox &  a_ebisbox,
const RealVect &  a_probLo,
const Real &  a_dx 
)
protectedvirtual

Set face-centered permittivities.

Parameters
[out]a_permPermittivity (on face centroid or center)
[in]a_cellBoxComputation box
[in]a_ebisboxEBIS box
[in]a_probLoLower-left corner off comptuational domain
[in]a_dxResolution

◆ setPermittivities()

void FieldSolver::setPermittivities ( )
virtual

Set the permittivities.

This sets m_permittivityCell, m_permittivityFace, and m_permittivityEB.

Reimplemented in FieldSolverMultigrid.

◆ setRealm()

void FieldSolver::setRealm ( const std::string  a_realm)

Set Realm.

Parameters
[in]a_realmRealm identifier.

◆ setRho() [1/2]

void FieldSolver::setRho ( const Real  a_rho)

Set space charge to constant value everywhere.

Parameters
[in]a_rhoValue of space charge

◆ setRho() [2/2]

void FieldSolver::setRho ( const std::function< Real(const RealVect)> &  a_rho)

Set space charge to spatially dependent function.

Parameters
[in]a_rhoSpatially dependent value of space charge

◆ setSigma()

void FieldSolver::setSigma ( const std::function< Real(const RealVect)> &  a_sigma)

Set surface charge density to spatially dependent function.

Parameters
[in]a_sigmaSpatially dependent value of the surface charge density.

◆ setSolverPermittivities()

void FieldSolver::setSolverPermittivities ( const MFAMRCellData a_permittivityCell,
const MFAMRFluxData a_permittivityFace,
const MFAMRIVData a_permittivityEB 
)
virtual

A special routine for when solver permittivities need to change but solver does not.

Parameters
[in]a_permittivityCellPermittivity on cell center
[in]a_permittivityFacePermittivity on faces
[in]a_permittivityEBPermittivity on EB
Note
Used by FieldSolverMultigrid to change permittivities/conductivities under the hood of the multigrid solver.

Reimplemented in FieldSolverMultigrid.

◆ setTime()

void FieldSolver::setTime ( const int  a_timeStep,
const Real  a_time,
const Real  a_dt 
)

Set time for this solver.

Parameters
[in]a_timeStepTime step
[in]a_timeTime (in seconds).
[in]a_dtTime step size (in seconds).

This sets m_timeStep, m_time, and m_dt.

◆ setVerbosity()

void FieldSolver::setVerbosity ( const int  a_verbosity)

Set verbosity.

Parameters
[in]a_verbosityVerbosity factor (lower yields less printed output).

◆ setVoltage()

void FieldSolver::setVoltage ( std::function< Real(const Real a_time)>  a_voltage)
virtual

Set potential dependence in time.

Parameters
[in]a_voltageFunction pointer which sets the voltage travel curve.

If you want something more complex, the voltage can be set individually for each electrode using setElectrodeDirichletFunction.

◆ solve() [1/3]

bool FieldSolver::solve ( const bool  a_zeroPhi = false)
virtual

Solve Poisson equation using m_potential, m_rho, and m_sigma.

Parameters
[in]a_zeroPhiSet m_potential to zero before calling other function.
Returns
True if we found a solution and false otherwise.

This calls the other version.

◆ solve() [2/3]

virtual bool FieldSolver::solve ( MFAMRCellData a_phi,
const MFAMRCellData a_rho,
const EBAMRIVData a_sigma,
const bool  a_zerophi = false 
)
pure virtual

Solves Poisson equation onto a_phi using a_rho and a_sigma as right-hand sides.

Parameters
[in,out]a_potentialPotential
[in]a_rhoSpace charge density
[in]a_sigmaSurface charge density.
[in]a_zeroPhiSet a_potential to zero first.
Returns
True if we found a solution and false otherwise.
Note
a_sigma must be defined on the gas phase.

Implemented in FieldSolverMultigrid.

◆ solve() [3/3]

bool FieldSolver::solve ( MFAMRCellData a_potential,
const bool  a_zerophi = false 
)
virtual

Solve Poisson equation onto a_phi using m_rho, and m_sigma as right-hand sides.

Parameters
[in,out]a_potentialPotential
[in]a_zeroPhiSet a_potential to zero before calling other function.
Returns
True if we found a solution and false otherwise.
Note
This calls the other (pure) version.

◆ writeMultifluidData()

void FieldSolver::writeMultifluidData ( LevelData< EBCellFAB > &  a_output,
int &  a_comp,
const MFAMRCellData a_data,
const phase::which_phase  a_phase,
const std::string  a_outputRealm,
const int  a_level,
const bool  a_interp 
) const
virtualnoexcept

Write multifluid data to single-fluid data holders.

This takes the valid data from each phase and writes it to a_output. On cut-cells, the a_phase flag determines which phase will be plotted.

Parameters
[in,out]a_outputOutput data holder
[in,out]a_compCurrent component in a_output.
[in]a_dataMultifluid data holder
[in]a_phaseMain phase in cut-cells.
[in]a_outputRealmRealm where a_output belongs
[in]a_levelGrid level
[in]a_interpSpecial flag which, if true, tells AmrMesh to interpolate the data to centroids (in irregular cells).

◆ writePlotData()

void FieldSolver::writePlotData ( LevelData< EBCellFAB > &  a_output,
int &  a_comp,
const std::string  a_outputRealm,
const int  a_level,
const bool  a_forceNoInterp = false 
) const
virtualnoexcept

Write output data to a_output.

Parameters
[in,out]a_outputOutput data holder
[in,out]a_compCurrent component in a_output.
[in]a_outputRealmRealm where a_output belongs
[in]a_levelGrid level
[in]a_forceNoInterpForces data to be written using native centering
Returns
On output, the solver will have written its plot variables to a_output and increment a_comp by the number of plotted variables.

◆ writePlotFile()

void FieldSolver::writePlotFile ( )
virtual

Write plot file.

This writes a plot file named FieldSolver.stepXXXXXXX.2d.hdf5 using the specified plot variables.

◆ writeSurfaceData()

void FieldSolver::writeSurfaceData ( LevelData< EBCellFAB > &  a_output,
int &  a_comp,
const LevelData< BaseIVFAB< Real >> &  a_data,
const std::string  a_outputRealm,
const int  a_level 
) const
virtualnoexcept

Write surface data to volume data holder.

This takes the valid data from a_data and puts in onto a_output.

Parameters
[in,out]a_outputOutput data holder
[in,out]a_compCurrent component in a_output.
[in]a_dataSurface data.
[in]a_outputRealmRealm where a_output belongs
[in]a_levelGrid level

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