Navigation

Main Page
Download
Support
Installation
Tutorial
Examples
Reference Manual
   Version 1.0.6dev
      Class Index
      File List
   Version 1.0.6
   Version 1.0.5new_solver
   Version 1.0.5dev
   Version 1.0.5b
   Version 1.0.4dev
   Version 1.0.4
Publications


Hosted by Get Ion Beam Simulator at SourceForge.net. Fast, secure and Free Open Source software downloads

Class for constructing the linear/nonlinear problem for the solver. More...

#include <epot_solver.hpp>

Inheritance diagram for EpotSolver:
EpotGSSolver EpotMGSolver EpotMGSubSolver EpotMatrixSolver EpotBiCGSTABSolver EpotUMFPACKSolver

Public Member Functions

 EpotSolver (Geometry &geom)
 Constructor for solver from geom. More...
 
 EpotSolver (const EpotSolver &epsolver, Geometry &geom)
 Constructor for solver from geom. Parameters from epsolver are copied to new solver. More...
 
 EpotSolver (Geometry &geom, std::istream &s)
 Construct from file. More...
 
virtual ~EpotSolver ()
 Destructor. More...
 
void set_parameters (const EpotSolver &epsolver)
 Copy parameters from solver epsolver. More...
 
void set_forced_potential_volume (double force_pot, CallbackFunctorB_V *force_pot_func)
 Define forced potential volume. More...
 
void set_forced_potential_volume (CallbackFunctorD_V *force_pot_func)
 Define forced potential volume. More...
 
void set_initial_plasma (double Up, CallbackFunctorB_V *init_plasma_func)
 Define initial plasma to the problem. More...
 
void set_plasma_calc_region (CallbackFunctorB_V *plasma_calc_func)
 Define plasma calculation region. More...
 
void set_pexp_plasma (double rhoe, double Te, double Up)
 Enable plasma model for positive ion extraction problem. More...
 
void set_nsimp_initial_plasma (CallbackFunctorB_V *init_plasma_func)
 Define initial plasma boundary location to negative ion extraction problem. More...
 
void set_nsimp_plasma (double rhop, double Ep, std::vector< double > rhoi, std::vector< double > Ei)
 Enable plasma model for negative ion extraction problem. More...
 

Protected Member Functions

void shield_newton (double &rhs, double &drhs, double epot) const
 Return non-linear right-hand-side and it's derivative for vacuum node in shielding model plasma. More...
 
void pexp_newton (double &rhs, double &drhs, double epot) const
 Return non-linear right-hand-side and it's derivative for vacuum node in positive ion plasma. More...
 
void nsimp_newton (double &rhs, double &drhs, double epot) const
 Return non-linear right-hand-side and it's derivative for vacuum node in negative ion plasma. More...
 
uint8_t boundary_index (uint32_t i) const
 Return bitmask indicating to which boundaries the node belongs to. More...
 
uint8_t boundary_index (uint32_t i, uint32_t j) const
 Return bitmask indicating to which boundaries the node belongs to. More...
 
uint8_t boundary_index (uint32_t i, uint32_t j, uint32_t k) const
 Return bitmask indicating to which boundaries the node belongs to. More...
 
uint8_t boundary_index_general (uint32_t i, uint32_t j, uint32_t k) const
 Return bitmask indicating to which boundaries the node belongs to. More...
 
void preprocess (MeshScalarField &epot)
 Do preprocessing action before solving. More...
 
void postprocess (void)
 Do postprocessing action after solving. More...
 
virtual void reset_problem (void)=0
 Reset solver/problem settings. More...
 
MeshScalarFieldevaluate_scharge (const ScalarField &__scharge) const
 
virtual void subsolve (MeshScalarField &epot, const MeshScalarField &scharge)=0
 Solve problem with given mesh based space charge. More...
 

Protected Attributes

Geometry_geom
 Geometry reference. More...
 
plasma_mode_e _plasma
 Plasma simulation mode. More...
 
double _rhoe
 Electron charge density (C/m3), < 0. More...
 
double _Te
 Electron thermal energy, > 0. More...
 
double _Up
 Plasma potential, > 0. More...
 
std::vector< double > _rhoi
 Charge density for positive ions, first fast protons, then thermal ions. More...
 
std::vector< double > _Ei
 Energy for positive ions, first fast protons, then thermal ions. More...
 
double _force_pot
 Potential to be forced. More...
 
CallbackFunctorB_V_force_pot_func
 Force region potential function. More...
 
CallbackFunctorD_V_force_pot_func2
 Force region potential function. More...
 
CallbackFunctorB_V_init_plasma_func
 Initial plasma region function. More...
 
CallbackFunctorB_V_plasma_calc_func
 Definition of plasma calculation region. More...
 
double _plA
 Plasma parameter. For positive ion extraction: rho_th * h^2 / epsilon_0, for negative ion extraction: rho_f * h^2 / epsilon_0 for shield1: 1/Tm. More...
 
double _plB
 Plasma parameter. For positive ion extraction: 1/Te, for negative ion extraction: E_f,i for shield1: phi_m/Tm. More...
 
double _plC
 Plasma parameter for positive ion extraction. Up/Te. More...
 
std::vector< double > _plD
 Plasma parameter for negative ion extraction. rho_th,i * h^2 / epsilon_0. More...
 
std::vector< double > _plE
 Plasma parameter for negative ion extraction. 1/Ti. More...
 

Detailed Description

Class for constructing the linear/nonlinear problem for the solver.

EpotProblem class constructs the Poisson (equation) problem in finite difference form from Geometry (mesh) and various parameters and it presents the problem to the Solver via matrix/vector representation. In case of linear problem (no plasma model), the Poisson equation

\[ \nabla^2 \phi = -\frac{\rho}{\epsilon_0} \]

is

\[ \frac{\partial^2 \phi}{\partial x^2} = -\frac{\rho}{\epsilon_0} \]

in 1D, which is discretized into

\[ \phi_{i-1} + \phi_{i+1} - 2\phi_{i} = -h^2 \frac{\rho_{i}}{\epsilon_0} \]

using finite differences. In 2D coordinates the discretized form is

\[ \phi_{i-1,j} + \phi_{i+1,j} + \phi_{i,j-1} + \phi_{i,j+1} - 4\phi_{i,j} = -h^2 \frac{\rho_{i,j}}{\epsilon_0} \]

and in 3D it is

\[ \phi_{i-1,j,k} + \phi_{i+1,j,k} + \phi_{i,j-1,k} + \phi_{i,j+1,k} + \phi_{i,j,k-1} + \phi_{i,j,k+1} - 6\phi_{i,j} = -h^2 \frac{\rho_{i,j,k}}{\epsilon_0}. \]

In cylindrical coordinates the Poisson equation is

\[ \frac{\partial^2 \phi}{\partial r^2} + \frac{1}{r} \frac{\partial \phi}{\partial r} + \frac{1}{r^2} \frac{\partial^2 \phi}{\partial \theta^2} + \frac{\partial^2 \phi}{\partial z^2} = -\frac{\rho}{\epsilon_0}, \]

where $ \frac{\partial^2 \phi}{\partial \theta^2} = 0 $ because of cylindrical symmetry of the simulations. Therefore the discretized form becomes

\[ \phi_{i-1,j} + \phi_{i+1,j} + \left( 1 - \frac{h}{2r_j} \right) \phi_{i,j-1} + \left( 1 + \frac{h}{2r_j} \right) \phi_{i,j+1} - 4\phi_{i,j} = -h^2 \frac{\rho_{i,j}}{\epsilon_0}, \]

where $ \frac{h}{2r_j} = \frac{h}{2jh} = \frac{1}{2j} $ because the radius $ r_j = hj $. At the symmetry axis there is an exception because both $ r $ and $ \frac{\partial \phi}{\partial r} $ approach zero. By using Bernoulli-L'Hopital rule we can evaluate

\[ \lim_{r \rightarrow 0} \frac{1}{r} \frac{\partial \phi}{\partial r} = \lim_{r \rightarrow 0} \frac{\partial^2 \phi}{\partial r^2}, \]

which turns the Poisson equation to

\[ 2 \frac{\partial^2 \phi}{\partial r^2} + \frac{\partial^2 \phi}{\partial z^2} = -\frac{\rho}{\epsilon_0} \]

on axis. Discretation of the equation gives us

\[ \phi_{i-1,j} + \phi_{i+1,j} + 2\phi_{i,j+1} + 2\phi_{i,j-1} - 6\phi_{i,j} = -h^2 \frac{\rho_{i,j}}{\epsilon_0}. \]

Here $ \phi_{i,j-1} = \phi_{i,j+1} $ so that the final form is

\[ \phi_{i-1,j} + \phi_{i+1,j} + 4\phi_{i,j+1} - 6\phi_{i,j} = -h^2 \frac{\rho_{i,j}}{\epsilon_0}. \]

In addition to the Poisson equation the problem matrix and vector also contain finite difference representations of the boundary conditions. The Dirichlet boundary condition is defined by constant potential at boundary, i.e. $ \phi_{i,j} = \phi_{\mathrm{const}} $. The Neumann boundary condition can be defined as first order discretation

\[ \frac{\phi_{i+1}-\phi_{i}}{h} = N_{\mathrm{const}} \]

or second order discretation

\[ \frac{-\phi_{i+2}+4\phi_{i+1}-3\phi_{i}}{2h} = N_{\mathrm{const}} \]

selected by the user.

The plasma problems are described by using nonlinear models for screening charges in the plasma. For positive ion extraction for example, the screening charge is an electron population at the plasma potential $ \phi_p $ with a thermal energy distribution with temperature $ T_e $. The screening charge is therefore

\[ \rho_e = \rho_{e0} \exp \left( \frac{\phi-\phi_p}{kT_e/e} \right), \]

where electron charge density at plasma potential $ (\rho_{e0}) $ is the same as the total positive beam space charge density for enabling plasma neutrality.

Constructor & Destructor Documentation

◆ EpotSolver() [1/3]

EpotSolver::EpotSolver ( Geometry geom)

Constructor for solver from geom.

◆ EpotSolver() [2/3]

EpotSolver::EpotSolver ( const EpotSolver epsolver,
Geometry geom 
)

Constructor for solver from geom. Parameters from epsolver are copied to new solver.

◆ EpotSolver() [3/3]

EpotSolver::EpotSolver ( Geometry geom,
std::istream &  s 
)

Construct from file.

◆ ~EpotSolver()

virtual EpotSolver::~EpotSolver ( )
inlinevirtual

Destructor.

Member Function Documentation

◆ boundary_index() [1/3]

uint8_t EpotSolver::boundary_index ( uint32_t  i) const
protected

Return bitmask indicating to which boundaries the node belongs to.

1D version. The node i is checked. The return value has bits set accoding to, which boundaries the node belongs to. The lowest bit is xmin, second of xmax, third is ymin, fourth is ymax, fifth is zmin, sixth is zmax. If the node isn't a boundary node, a 0 is returned.

◆ boundary_index() [2/3]

uint8_t EpotSolver::boundary_index ( uint32_t  i,
uint32_t  j 
) const
protected

Return bitmask indicating to which boundaries the node belongs to.

2D version. The node (i,j) is checked. The return value has bits set accoding to, which boundaries the node belongs to. The lowest bit is xmin, second of xmax, third is ymin, fourth is ymax, fifth is zmin, sixth is zmax. If the node isn't a boundary node, a 0 is returned.

◆ boundary_index() [3/3]

uint8_t EpotSolver::boundary_index ( uint32_t  i,
uint32_t  j,
uint32_t  k 
) const
protected

Return bitmask indicating to which boundaries the node belongs to.

3D version. The node (i,j,k) is checked. The return value has bits set accoding to, which boundaries the node belongs to. The lowest bit is xmin, second of xmax, third is ymin, fourth is ymax, fifth is zmin, sixth is zmax. If the node isn't a boundary node, a 0 is returned.

◆ boundary_index_general()

uint8_t EpotSolver::boundary_index_general ( uint32_t  i,
uint32_t  j,
uint32_t  k 
) const
protected

Return bitmask indicating to which boundaries the node belongs to.

Version applicable to all dimensions. The node (i,j,k) is checked. The return value has bits set accoding to, which boundaries the node belongs to. The lowest bit is xmin, second of xmax, third is ymin, fourth is ymax, fifth is zmin, sixth is zmax. If the node isn't a boundary node, a 0 is returned.

◆ nsimp_newton()

void EpotSolver::nsimp_newton ( double &  rhs,
double &  drhs,
double  epot 
) const
protected

Return non-linear right-hand-side and it's derivative for vacuum node in negative ion plasma.

Returns the non-linear part of the right-hand-side and it's derivative resulting from the compensating space charge in negative ion extraction case.

◆ pexp_newton()

void EpotSolver::pexp_newton ( double &  rhs,
double &  drhs,
double  epot 
) const
protected

Return non-linear right-hand-side and it's derivative for vacuum node in positive ion plasma.

Returns the non-linear part of the right-hand-side and it's derivative resulting from the compensating space charge in positive ion extraction case.

◆ postprocess()

void EpotSolver::postprocess ( void  )
protected

Do postprocessing action after solving.

Return near solid neumann points as near solid. Restore distance data. Remove fixed node tags.

◆ preprocess()

void EpotSolver::preprocess ( MeshScalarField epot)
protected

Do preprocessing action before solving.

Does precalculation of plasma paramters. Sets near solid neumann points as neumann and saves distance to a vector. Mark fixed and initial plasma nodes and set potentials.

◆ reset_problem()

virtual void EpotSolver::reset_problem ( void  )
protectedpure virtual

Reset solver/problem settings.

◆ set_forced_potential_volume() [1/2]

void EpotSolver::set_forced_potential_volume ( CallbackFunctorD_V force_pot_func)

Define forced potential volume.

The force_pot_func functor will be called for every vacuum grid node at preprocessing stage of solve. If it returns a finite value, the node will be forced to the potential defined by the returned value. Otherwise the node will be treated as a regular vacuum node.

◆ set_forced_potential_volume() [2/2]

void EpotSolver::set_forced_potential_volume ( double  force_pot,
CallbackFunctorB_V force_pot_func 
)

Define forced potential volume.

Vacuum volume inside the volume defined by function force_pot_func will be forced to potential force_pot_func. This function is designed to be used with negative ion plasma extraction to stabilize plasma close non-physical boundaries.

This function may be removed in future.

◆ set_initial_plasma()

void EpotSolver::set_initial_plasma ( double  Up,
CallbackFunctorB_V init_plasma_func 
)

Define initial plasma to the problem.

Initial plasma volume is defined in the area given by callback functor init_plasma_func.

◆ set_nsimp_initial_plasma()

void EpotSolver::set_nsimp_initial_plasma ( CallbackFunctorB_V init_plasma_func)

Define initial plasma boundary location to negative ion extraction problem.

Initial plasma volume is defined in the area given by plasma_func.

◆ set_nsimp_plasma()

void EpotSolver::set_nsimp_plasma ( double  rhop,
double  Ep,
std::vector< double >  rhoi,
std::vector< double >  Ei 
)

Enable plasma model for negative ion extraction problem.

The positive (analytic) space charges for the negative ion plasma extraction are set using this function. The positive ions consist of fast (directed) protons and any number of thermal positive ions trapped at the plasma boundary in the zero potential.

The parameters set are rhop, the space charge density of protons and Ep, the energy of protons at zero potential. Vectors rhoi and Ei are used to set the space charge densities and thermal energies of the trapped ions.

For a typical negative ion extraction case the compensating charges are positive (rhoi>0). This function allows also negative compensating species to be simulated with rhoi<0 and Ei<0.

The boundary condition for the plasma should be BOUND_DIRICHLET with zero volts.

◆ set_parameters()

void EpotSolver::set_parameters ( const EpotSolver epsolver)

Copy parameters from solver epsolver.

◆ set_pexp_plasma()

void EpotSolver::set_pexp_plasma ( double  rhoe,
double  Te,
double  Up 
)

Enable plasma model for positive ion extraction problem.

Enable plasma model with background electron charge density of rhoe and electron temperature Te. The plasma potential is set to Up.

The boundary condition for the plasma should be BOUND_NEUMANN with gradient of zero V/m.

◆ set_plasma_calc_region()

void EpotSolver::set_plasma_calc_region ( CallbackFunctorB_V plasma_calc_func)

Define plasma calculation region.

By default the plasma compensation function takes place in the whole geometry. With this function one can define a functor plasma_calc_func, which is called for all calculation nodes during Poisson solve. The function should return true if plasma calculation is enabled at that point and false if not.

Implemented only in matrix solvers.

◆ shield_newton()

void EpotSolver::shield_newton ( double &  rhs,
double &  drhs,
double  epot 
) const
protected

Return non-linear right-hand-side and it's derivative for vacuum node in shielding model plasma.

Returns the non-linear part of the right-hand-side and it's derivative.

◆ subsolve()

virtual void EpotSolver::subsolve ( MeshScalarField epot,
const MeshScalarField scharge 
)
protectedpure virtual

Solve problem with given mesh based space charge.

Member Data Documentation

◆ _Ei

std::vector<double> EpotSolver::_Ei
protected

Energy for positive ions, first fast protons, then thermal ions.

◆ _force_pot

double EpotSolver::_force_pot
protected

Potential to be forced.

◆ _force_pot_func

CallbackFunctorB_V* EpotSolver::_force_pot_func
protected

Force region potential function.

◆ _force_pot_func2

CallbackFunctorD_V* EpotSolver::_force_pot_func2
protected

Force region potential function.

◆ _geom

Geometry& EpotSolver::_geom
protected

Geometry reference.

◆ _init_plasma_func

CallbackFunctorB_V* EpotSolver::_init_plasma_func
protected

Initial plasma region function.

◆ _plA

double EpotSolver::_plA
protected

Plasma parameter. For positive ion extraction: rho_th * h^2 / epsilon_0, for negative ion extraction: rho_f * h^2 / epsilon_0 for shield1: 1/Tm.

◆ _plasma

plasma_mode_e EpotSolver::_plasma
protected

Plasma simulation mode.

◆ _plasma_calc_func

CallbackFunctorB_V* EpotSolver::_plasma_calc_func
protected

Definition of plasma calculation region.

◆ _plB

double EpotSolver::_plB
protected

Plasma parameter. For positive ion extraction: 1/Te, for negative ion extraction: E_f,i for shield1: phi_m/Tm.

◆ _plC

double EpotSolver::_plC
protected

Plasma parameter for positive ion extraction. Up/Te.

◆ _plD

std::vector<double> EpotSolver::_plD
protected

Plasma parameter for negative ion extraction. rho_th,i * h^2 / epsilon_0.

◆ _plE

std::vector<double> EpotSolver::_plE
protected

Plasma parameter for negative ion extraction. 1/Ti.

◆ _rhoe

double EpotSolver::_rhoe
protected

Electron charge density (C/m3), < 0.

◆ _rhoi

std::vector<double> EpotSolver::_rhoi
protected

Charge density for positive ions, first fast protons, then thermal ions.

◆ _Te

double EpotSolver::_Te
protected

Electron thermal energy, > 0.

◆ _Up

double EpotSolver::_Up
protected

Plasma potential, > 0.


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


Reference manual for Ion Beam Simulator 1.0.6dev
Generated by Doxygen 1.9.1 on Thu Sep 11 2025 09:37:24.