Class DiffusionSolver
Defined in File diffusion.h
Nested Relationships
Nested Types
Inheritance Relationships
Derived Types
public opensn::DiffusionMIPSolver
(Class DiffusionMIPSolver)public opensn::DiffusionPWLCSolver
(Class DiffusionPWLCSolver)
Class Documentation
-
class DiffusionSolver
Generic diffusion solver for acceleration.
Subclassed by opensn::DiffusionMIPSolver, opensn::DiffusionPWLCSolver
Public Functions
-
DiffusionSolver(std::string name, const SpatialDiscretization &sdm, const UnknownManager &uk_man, std::map<uint64_t, BoundaryCondition> bcs, MatID2XSMap map_mat_id_2_xs, const std::vector<UnitCellMatrices> &unit_cell_matrices, bool suppress_bcs, bool requires_ghosts, bool verbose)
-
std::string GetName() const
Returns the assigned name.
-
const Vec &GetRHS() const
Returns the right-hand side petsc vector.
-
inline const std::map<uint64_t, BoundaryCondition> &GetBCS() const
Returns the assigned unknown structure.
-
const UnknownManager &GetUnknownStructure() const
-
const class SpatialDiscretization &GetSpatialDiscretization() const
Returns the associated spatial discretization.
-
std::pair<size_t, size_t> GetNumPhiIterativeUnknowns()
-
virtual ~DiffusionSolver()
-
void Initialize()
Initializes the diffusion solver. This involves creating the sparse matrix with the appropriate sparsity pattern. Creating the RHS vector. Creating the KSP solver. Setting the very specialized parameters for Hypre’s BooomerAMG. Note:
PCSetFromOptions
andKSPSetFromOptions
are called at the end. Therefore, any number of additional PETSc options can be passed via the commandline.
-
virtual void AssembleAand_b(const std::vector<double> &q_vector) = 0
-
virtual void Assemble_b(const std::vector<double> &q_vector) = 0
-
virtual void Assemble_b(Vec petsc_q_vector) = 0
-
void AddToRHS(const std::vector<double> &values)
Adds to the right-hand side without applying spatial discretization.
-
void AddToMatrix(const std::vector<int64_t> &rows, const std::vector<int64_t> &cols, const std::vector<double> &vals)
Adds to the entries into the matrix without applying spatial discretization.
-
void Solve(std::vector<double> &solution, bool use_initial_guess = false)
Solves the system and stores the local solution in the vector provide.
- Parameters:
solution – Vector in to which the solution will be parsed.
use_initial_guess – bool [Default:False] Flag, when set, will use the values of the output solution as initial guess.
-
void Solve(Vec petsc_solution, bool use_initial_guess = false)
Solves the system and stores the local solution in the vector provide.
- Parameters:
petsc_solution – Vector in to which the solution will be parsed.
use_initial_guess – bool [Default:False] Flag, when set, will use the values of the output solution as initial guess.
Public Members
-
struct opensn::DiffusionSolver::Options options
Protected Types
-
using MatID2XSMap = std::map<int, Multigroup_D_and_sigR>
Protected Attributes
-
const std::string name_
-
const std::shared_ptr<MeshContinuum> grid_
-
const class SpatialDiscretization &sdm_
-
const UnknownManager uk_man_
-
const std::map<uint64_t, BoundaryCondition> bcs_
-
const MatID2XSMap mat_id_2_xs_map_
-
const std::vector<UnitCellMatrices> &unit_cell_matrices_
-
const int64_t num_local_dofs_
-
const int64_t num_global_dofs_
-
Mat A_ = nullptr
-
Vec rhs_ = nullptr
-
KSP ksp_ = nullptr
-
const bool requires_ghosts_
-
const bool suppress_bcs_
-
struct Options
Public Members
-
double residual_tolerance = 1.0e-4
Residual tol. relative to rhs.
-
int max_iters = 100
Maximum iterations.
-
bool verbose = false
Verbosity flag.
-
bool perform_symmetry_check = false
For debugging only (very expensive)
-
std::string additional_options_string
-
double penalty_factor = 4.0
-
double residual_tolerance = 1.0e-4
-
DiffusionSolver(std::string name, const SpatialDiscretization &sdm, const UnknownManager &uk_man, std::map<uint64_t, BoundaryCondition> bcs, MatID2XSMap map_mat_id_2_xs, const std::vector<UnitCellMatrices> &unit_cell_matrices, bool suppress_bcs, bool requires_ghosts, bool verbose)