43 #ifndef EPOT_GSSOLVER_HPP
44 #define EPOT_GSSOLVER_HPP 1
63 double _local_Ulim_fac;
68 double solve_pexp_potential(
double epf,
double cof,
double rhs,
double p )
const;
69 double solve_nsimp_potential(
double epf,
double cof,
double rhs,
double p )
const;
70 void prepare_local_gnewton_settings(
void );
72 double gs_loop_3d(
void )
const;
73 double gs_process_near_solid_3d(
const uint8_t *nearsolid_ptr,
74 uint32_t a, uint32_t dj, uint32_t dk, uint8_t bindex )
const;
75 double gs_process_pure_vacuum_3d( uint32_t a, uint32_t dj, uint32_t dk )
const;
76 double gs_process_neumann_3d( uint32_t a, uint32_t dj, uint32_t dk, uint8_t bindex )
const;
79 double gs_loop_cyl(
void )
const;
80 double gs_process_near_solid_cyl(
const uint8_t *nearsolid_ptr,
81 uint32_t i, uint32_t j, uint8_t bindex )
const;
82 double gs_process_pure_vacuum_cyl( uint32_t i, uint32_t j )
const;
83 double gs_process_neumann_cyl( uint32_t i, uint32_t j, uint8_t bindex )
const;
86 double gs_loop_2d(
void )
const;
87 double gs_process_near_solid_2d(
const uint8_t *nearsolid_ptr,
88 uint32_t a, uint32_t dj,
89 uint8_t bindex )
const;
90 double gs_process_pure_vacuum_2d( uint32_t a, uint32_t dj )
const;
91 double gs_process_neumann_2d( uint32_t a, uint32_t dj, uint8_t bindex )
const;
94 double gs_loop_1d(
void )
const;
95 double gs_process_near_solid_1d(
const uint8_t *nearsolid_ptr,
96 uint32_t i, uint8_t bindex )
const;
97 double gs_process_pure_vacuum_1d( uint32_t i )
const;
98 double gs_process_neumann_1d( uint32_t i, uint8_t bindex )
const;
101 double near_solid_neumann_rhs_contribution( uint32_t i, uint32_t j, uint32_t k,
102 uint8_t bindex,
const Vec3D &x )
const;
104 void postprocess(
void );
111 double error_scale(
void );
115 virtual void reset_problem(
void );
151 void set_w(
double w );
183 virtual void debug_print( std::ostream &os )
const;
187 virtual void save( std::ostream &s )
const;
double get_error_estimate(void) const
Get estimate of relative solution error.
Definition: epot_gssolver.cpp:89
virtual void debug_print(std::ostream &os) const
Print debugging information to os.
Definition: epot_gssolver.cpp:1150
Scalar field class.
Definition: meshscalarfield.hpp:70
Geometry defining class.
Definition: geometry.hpp:179
Class for constructing the linear/nonlinear problem for the solver.
Definition: epot_solver.hpp:207
double get_potential_change_norm(void) const
Get potential change norm.
Definition: epot_gssolver.cpp:83
void set_w(double w)
Sets relaxation parameter.
Definition: epot_gssolver.cpp:107
virtual ~EpotGSSolver()
Destructor.
Definition: epot_gssolver.hpp:133
EpotGSSolver(Geometry &geom)
Constructor.
Definition: epot_gssolver.cpp:55
void set_imax(uint32_t imax)
Sets maximum iteration count.
Definition: epot_gssolver.cpp:101
Poisson equation problem for solving electric potential.
void set_plasma_solver_parameters(double Ulim_fac, uint32_t imax, double eps)
Set node wise plasma solver parameters.
Definition: epot_gssolver.cpp:113
uint32_t get_iter(void) const
Get number of iteration rounds done with last solve().
Definition: epot_gssolver.cpp:95
void set_eps(double eps)
Sets the accuracy request.
Definition: epot_gssolver.cpp:77
virtual void save(std::ostream &s) const
Saves problem data to stream.
Definition: epot_gssolver.cpp:1144
Three dimensional vector.
Definition: vec3d.hpp:58
Gauss-Seidel solver for Electric potential problem.
Definition: epot_gssolver.hpp:52