Chi-Tech
lbs::LBSSolver Class Reference

#include <lbs_solver.h>

Inheritance diagram for lbs::LBSSolver:
chi_physics::Solver ChiObject lbs::DiffusionDFEMSolver lbs::DiscreteOrdinatesSolver lbs::DiscreteOrdinatesAdjointSolver lbs::DiscreteOrdinatesCurvilinearSolver

Public Types

typedef std::shared_ptr< AGSLinearSolver< Mat, Vec, KSP > > AGSLinSolverPtr
 
typedef std::shared_ptr< chi_math::LinearSolver< Mat, Vec, KSP > > LinSolvePtr
 

Public Member Functions

 LBSSolver (const std::string &text_name)
 
 LBSSolver (const chi::InputParameters &params)
 
 LBSSolver (const LBSSolver &)=delete
 
LBSSolveroperator= (const LBSSolver &)=delete
 
virtual ~LBSSolver ()=default
 
size_t GetSourceEventTag () const
 
double LastRestartWrite () const
 
double & LastRestartWrite ()
 
lbs::OptionsOptions ()
 
const lbs::OptionsOptions () const
 
void SetOptions (const chi::InputParameters &params)
 
void SetBoundaryOptions (const chi::InputParameters &params)
 
size_t NumMoments () const
 
size_t NumGroups () const
 
size_t NumPrecursors () const
 
size_t GetMaxPrecursorsPerMaterial () const
 
void AddGroup (int id)
 
const std::vector< LBSGroup > & Groups () const
 
void AddGroupset ()
 
std::vector< LBSGroupset > & Groupsets ()
 
const std::vector< LBSGroupset > & Groupsets () const
 
void AddPointSource (PointSource psrc)
 
void ClearPointSources ()
 
const std::vector< PointSource > & PointSources () const
 
const std::map< int, XSPtr > & GetMatID2XSMap () const
 
const std::map< int, IsotropicSrcPtr > & GetMatID2IsoSrcMap () const
 
const chi_math::SpatialDiscretizationSpatialDiscretization () const
 
const std::vector< UnitCellMatrices > & GetUnitCellMatrices () const
 
const chi_mesh::MeshContinuumGrid () const
 
const std::vector< lbs::CellLBSView > & GetCellTransportViews () const
 
std::map< uint64_t, BoundaryPreference > & BoundaryPreferences ()
 
const chi_math::UnknownManagerUnknownManager () const
 
size_t LocalNodeCount () const
 
size_t GlobalNodeCount () const
 
std::vector< double > & QMomentsLocal ()
 
const std::vector< double > & QMomentsLocal () const
 
std::vector< double > & ExtSrcMomentsLocal ()
 
const std::vector< double > & ExtSrcMomentsLocal () const
 
std::vector< double > & PhiOldLocal ()
 
const std::vector< double > & PhiOldLocal () const
 
std::vector< double > & PhiNewLocal ()
 
const std::vector< double > & PhiNewLocal () const
 
std::vector< double > & PrecursorsNewLocal ()
 
const std::vector< double > & PrecursorsNewLocal () const
 
std::vector< VecDbl > & PsiNewLocal ()
 
const std::vector< VecDbl > & PsiNewLocal () const
 
const std::map< uint64_t, std::shared_ptr< SweepBndry > > & SweepBoundaries () const
 
SetSourceFunction GetActiveSetSourceFunction () const
 
AGSLinSolverPtr GetPrimaryAGSSolver ()
 
std::vector< LinSolvePtr > & GetWGSSolvers ()
 
WGSContext< Mat, Vec, KSP > & GetWGSContext (int groupset_id)
 
virtual std::pair< size_t, size_t > GetNumPhiIterativeUnknowns ()
 
size_t MapPhiFieldFunction (size_t g, size_t m) const
 
size_t GetHandleToPowerGenFieldFunc () const
 
void Initialize () override
 
void InitMaterials ()
 
void InitializePointSources ()
 
void InitWGDSA (LBSGroupset &groupset, bool vaccum_bcs_are_dirichlet=true)
 
std::vector< double > WGSCopyOnlyPhi0 (const LBSGroupset &groupset, const std::vector< double > &phi_in)
 
void GSProjectBackPhi0 (const LBSGroupset &groupset, const std::vector< double > &input, std::vector< double > &output)
 
void AssembleWGDSADeltaPhiVector (const LBSGroupset &groupset, const std::vector< double > &phi_in, std::vector< double > &delta_phi_local)
 
void DisAssembleWGDSADeltaPhiVector (const LBSGroupset &groupset, const std::vector< double > &delta_phi_local, std::vector< double > &ref_phi_new)
 
void AssembleTGDSADeltaPhiVector (const LBSGroupset &groupset, const std::vector< double > &phi_in, std::vector< double > &delta_phi_local)
 
void DisAssembleTGDSADeltaPhiVector (const LBSGroupset &groupset, const std::vector< double > &delta_phi_local, std::vector< double > &ref_phi_new)
 
void WriteRestartData (const std::string &folder_name, const std::string &file_base)
 
void ReadRestartData (const std::string &folder_name, const std::string &file_base)
 
void WriteGroupsetAngularFluxes (const LBSGroupset &groupset, const std::string &file_base)
 
void ReadGroupsetAngularFluxes (LBSGroupset &groupset, const std::string &file_base)
 
std::vector< double > MakeSourceMomentsFromPhi ()
 
void WriteFluxMoments (const std::string &file_base, const std::vector< double > &flux_moments)
 
void ReadFluxMoments (const std::string &file_base, std::vector< double > &flux_moments, bool single_file=false)
 
void UpdateFieldFunctions ()
 
void SetPhiFromFieldFunctions (PhiSTLOption which_phi, const std::vector< size_t > &m_indices, const std::vector< size_t > &g_indices)
 
double ComputeFissionProduction (const std::vector< double > &phi)
 
double ComputeFissionRate (const std::vector< double > &phi)
 
void ComputePrecursors ()
 
virtual void SetPhiVectorScalarValues (std::vector< double > &phi_vector, double value)
 
virtual void ScalePhiVector (PhiSTLOption which_phi, double value)
 
virtual void SetGSPETScVecFromPrimarySTLvector (LBSGroupset &groupset, Vec x, PhiSTLOption which_phi)
 
virtual void SetPrimarySTLvectorFromGSPETScVec (LBSGroupset &groupset, Vec x_src, PhiSTLOption which_phi)
 
virtual void GSScopedCopyPrimarySTLvectors (LBSGroupset &groupset, const std::vector< double > &x_src, std::vector< double > &y)
 
virtual void GSScopedCopyPrimarySTLvectors (LBSGroupset &groupset, PhiSTLOption from_which_phi, PhiSTLOption to_which_phi)
 
virtual void SetGroupScopedPETScVecFromPrimarySTLvector (int first_group_id, int last_group_id, Vec x, const std::vector< double > &y)
 
virtual void SetPrimarySTLvectorFromGroupScopedPETScVec (int first_group_id, int last_group_id, Vec x_src, std::vector< double > &y)
 
virtual void SetMultiGSPETScVecFromPrimarySTLvector (const std::vector< int > &gs_ids, Vec x, PhiSTLOption which_phi)
 
virtual void SetPrimarySTLvectorFromMultiGSPETScVecFrom (const std::vector< int > &gs_ids, Vec x_src, PhiSTLOption which_phi)
 
- Public Member Functions inherited from chi_physics::Solver
 Solver (std::string in_text_name)
 
 Solver (std::string in_text_name, std::initializer_list< BasicOption > in_options)
 
 Solver (const chi::InputParameters &params)
 
virtual ~Solver ()=default
 
std::string TextName () const
 
BasicOptionsGetBasicOptions ()
 
const BasicOptionsGetBasicOptions () const
 
std::vector< std::shared_ptr< FieldFunctionGridBased > > & GetFieldFunctions ()
 
const std::vector< std::shared_ptr< FieldFunctionGridBased > > & GetFieldFunctions () const
 
TimeStepperGetTimeStepper ()
 
const TimeStepperGetTimeStepper () const
 
virtual void Initialize ()
 
virtual void Execute ()
 
virtual void Step ()
 
virtual void Advance ()
 
virtual chi::ParameterBlock GetInfo (const chi::ParameterBlock &params) const
 
virtual void SetProperties (const chi::ParameterBlock &params)
 
chi::ParameterBlock GetInfoWithPreCheck (const chi::ParameterBlock &params) const
 
- Public Member Functions inherited from ChiObject
 ChiObject ()
 
 ChiObject (const chi::InputParameters &params)
 
void SetStackID (size_t stack_id)
 
size_t StackID () const
 
virtual void PushOntoStack (std::shared_ptr< ChiObject > &new_object)
 
virtual ~ChiObject ()=default
 

Static Public Member Functions

static chi::InputParameters GetInputParameters ()
 
static chi::InputParameters OptionsBlock ()
 
static chi::InputParameters BoundaryOptionsBlock ()
 
- Static Public Member Functions inherited from chi_physics::Solver
static chi::InputParameters GetInputParameters ()
 
- Static Public Member Functions inherited from ChiObject
static chi::InputParameters GetInputParameters ()
 

Protected Types

typedef chi_mesh::sweep_management::CellFaceNodalMapping CellFaceNodalMapping
 

Protected Member Functions

virtual void PerformInputChecks ()
 
void PrintSimHeader ()
 
virtual void InitializeSpatialDiscretization ()
 
void ComputeUnitIntegrals ()
 
void InitializeGroupsets ()
 
void ComputeNumberOfMoments ()
 
virtual void InitializeParrays ()
 
void InitializeFieldFunctions ()
 
void InitializeBoundaries ()
 
virtual void InitializeSolverSchemes ()
 
virtual void InitializeWGSSolvers ()
 
void InitTGDSA (LBSGroupset &groupset)
 

Static Protected Member Functions

static void CleanUpWGDSA (LBSGroupset &groupset)
 
static void CleanUpTGDSA (LBSGroupset &groupset)
 

Protected Attributes

size_t source_event_tag_ = 0
 
double last_restart_write_ = 0.0
 
lbs::Options options_
 
size_t num_moments_ = 0
 
size_t num_groups_ = 0
 
size_t num_precursors_ = 0
 
size_t max_precursors_per_material_ = 0
 
std::vector< LBSGroupgroups_
 
std::vector< LBSGroupsetgroupsets_
 
std::vector< PointSourcepoint_sources_
 
std::map< int, XSPtrmatid_to_xs_map_
 
std::map< int, IsotropicSrcPtrmatid_to_src_map_
 
std::shared_ptr< chi_math::SpatialDiscretizationdiscretization_ = nullptr
 
chi_mesh::MeshContinuumPtr grid_ptr_
 
std::vector< CellFaceNodalMappinggrid_nodal_mappings_
 
MPILocalCommSetPtr grid_local_comm_set_ = nullptr
 
GridFaceHistogramPtr grid_face_histogram_ = nullptr
 
std::vector< UnitCellMatricesunit_cell_matrices_
 
std::map< uint64_t, UnitCellMatricesunit_ghost_cell_matrices_
 
std::vector< lbs::CellLBSViewcell_transport_views_
 
std::map< uint64_t, BoundaryPreferenceboundary_preferences_
 
std::map< uint64_t, std::shared_ptr< SweepBndry > > sweep_boundaries_
 
chi_math::UnknownManager flux_moments_uk_man_
 
size_t max_cell_dof_count_ = 0
 
uint64_t local_node_count_ = 0
 
uint64_t glob_node_count_ = 0
 
std::vector< double > q_moments_local_
 
std::vector< double > ext_src_moments_local_
 
std::vector< double > phi_new_local_
 
std::vector< double > phi_old_local_
 
std::vector< std::vector< double > > psi_new_local_
 
std::vector< double > precursor_new_local_
 
SetSourceFunction active_set_source_function_
 
std::vector< AGSLinSolverPtrags_solvers_
 
std::vector< LinSolvePtrwgs_solvers_
 
AGSLinSolverPtr primary_ags_solver_
 
std::map< std::pair< size_t, size_t >, size_t > phi_field_functions_local_map_
 
size_t power_gen_fieldfunc_local_handle_ = 0
 
std::shared_ptr< const chi_math::TimeIntegrationtime_integration_ = nullptr
 
- Protected Attributes inherited from chi_physics::Solver
BasicOptions basic_options_
 
std::vector< std::shared_ptr< FieldFunctionGridBased > > field_functions_
 
std::shared_ptr< TimeSteppertimestepper_ = nullptr
 

Detailed Description

Base class for all Linear Boltzmann Solvers.

Definition at line 49 of file lbs_solver.h.

Member Typedef Documentation

◆ AGSLinSolverPtr

typedef std::shared_ptr<AGSLinearSolver<Mat, Vec, KSP> > lbs::LBSSolver::AGSLinSolverPtr

Definition at line 52 of file lbs_solver.h.

◆ CellFaceNodalMapping

◆ LinSolvePtr

typedef std::shared_ptr<chi_math::LinearSolver<Mat, Vec, KSP> > lbs::LBSSolver::LinSolvePtr

Definition at line 53 of file lbs_solver.h.

Constructor & Destructor Documentation

◆ LBSSolver() [1/3]

lbs::LBSSolver::LBSSolver ( const std::string &  text_name)
explicit

Base class constructor.

Definition at line 13 of file lbs_00_constrdestr.cc.

◆ LBSSolver() [2/3]

lbs::LBSSolver::LBSSolver ( const chi::InputParameters params)
explicit

Input parameters based construction.

Definition at line 45 of file lbs_00_constrdestr.cc.

◆ LBSSolver() [3/3]

lbs::LBSSolver::LBSSolver ( const LBSSolver )
delete

◆ ~LBSSolver()

virtual lbs::LBSSolver::~LBSSolver ( )
virtualdefault

Member Function Documentation

◆ AddGroup()

void lbs::LBSSolver::AddGroup ( int  id)

Adds a group to the list of groups. If group id < 0, the id will be logically derived from the list size. If >= 0 the id will be set to the id specified.

Definition at line 117 of file lbs_00_constrdestr.cc.

◆ AddGroupset()

void lbs::LBSSolver::AddGroupset ( )

Adds a groupset to the list of groupsets. The groupset id will be logically derived from the list size.

Definition at line 129 of file lbs_00_constrdestr.cc.

◆ AddPointSource()

void lbs::LBSSolver::AddPointSource ( PointSource  psrc)

Adds a point source to the solver's point source list.

Definition at line 144 of file lbs_00_constrdestr.cc.

◆ AssembleTGDSADeltaPhiVector()

void lbs::LBSSolver::AssembleTGDSADeltaPhiVector ( const LBSGroupset groupset,
const std::vector< double > &  phi_in,
std::vector< double > &  delta_phi_local 
)

Assembles a delta-phi vector on the first moment.

Definition at line 87 of file lbs_03e_tgdsa.cc.

◆ AssembleWGDSADeltaPhiVector()

void lbs::LBSSolver::AssembleWGDSADeltaPhiVector ( const LBSGroupset groupset,
const std::vector< double > &  phi_in,
std::vector< double > &  delta_phi_local 
)

Assembles a delta-phi vector on the first moment.

Definition at line 139 of file lbs_03d_wgdsa.cc.

◆ BoundaryOptionsBlock()

chi::InputParameters lbs::LBSSolver::BoundaryOptionsBlock ( )
static

Definition at line 116 of file lbs_00a_setoptions.cc.

◆ BoundaryPreferences()

std::map< uint64_t, BoundaryPreference > & lbs::LBSSolver::BoundaryPreferences ( )

Read/Write access to the boundary preferences.

Definition at line 291 of file lbs_00_constrdestr.cc.

◆ CleanUpTGDSA()

void lbs::LBSSolver::CleanUpTGDSA ( LBSGroupset groupset)
staticprotected

Cleans up memory consuming items.

Definition at line 80 of file lbs_03e_tgdsa.cc.

◆ CleanUpWGDSA()

void lbs::LBSSolver::CleanUpWGDSA ( LBSGroupset groupset)
staticprotected

Cleans up memory consuming items.

Definition at line 60 of file lbs_03d_wgdsa.cc.

◆ ClearPointSources()

void lbs::LBSSolver::ClearPointSources ( )

Clears all the point sources from the solver's point source list.

Definition at line 150 of file lbs_00_constrdestr.cc.

◆ ComputeFissionProduction()

double LBSSolver::ComputeFissionProduction ( const std::vector< double > &  phi)

Compute the total fission production in the problem.

Author
Zachary Hardy.

Definition at line 11 of file lbs_06b_compute_fission_rates.cc.

◆ ComputeFissionRate()

double LBSSolver::ComputeFissionRate ( const std::vector< double > &  phi)

Computes the total fission rate in the problem.

Author
Zachary Hardy.

Definition at line 69 of file lbs_06b_compute_fission_rates.cc.

◆ ComputeNumberOfMoments()

void lbs::LBSSolver::ComputeNumberOfMoments ( )
protected

Computes the number of moments for the given mesher types

Definition at line 5 of file lbs_01f_compute_nummoments.cc.

◆ ComputePrecursors()

void lbs::LBSSolver::ComputePrecursors ( )

Compute the steady state delayed neutron precursor concentrations.

Definition at line 10 of file lbs_06c_compute_precursors.cc.

◆ ComputeUnitIntegrals()

void lbs::LBSSolver::ComputeUnitIntegrals ( )
protected

Definition at line 22 of file lbs_01d_init_spat_discr.cc.

◆ DisAssembleTGDSADeltaPhiVector()

void lbs::LBSSolver::DisAssembleTGDSADeltaPhiVector ( const LBSGroupset groupset,
const std::vector< double > &  delta_phi_local,
std::vector< double > &  ref_phi_new 
)

DAssembles a delta-phi vector on the first moment.

Definition at line 130 of file lbs_03e_tgdsa.cc.

◆ DisAssembleWGDSADeltaPhiVector()

void lbs::LBSSolver::DisAssembleWGDSADeltaPhiVector ( const LBSGroupset groupset,
const std::vector< double > &  delta_phi_local,
std::vector< double > &  ref_phi_new 
)

DAssembles a delta-phi vector on the first moment.

Definition at line 178 of file lbs_03d_wgdsa.cc.

◆ ExtSrcMomentsLocal() [1/2]

std::vector< double > & lbs::LBSSolver::ExtSrcMomentsLocal ( )

Read/write access to exterior src moments vector.

Definition at line 212 of file lbs_00_constrdestr.cc.

◆ ExtSrcMomentsLocal() [2/2]

const std::vector< double > & lbs::LBSSolver::ExtSrcMomentsLocal ( ) const

Read access to exterior src moments vector.

Definition at line 218 of file lbs_00_constrdestr.cc.

◆ GetActiveSetSourceFunction()

SetSourceFunction lbs::LBSSolver::GetActiveSetSourceFunction ( ) const

Definition at line 263 of file lbs_00_constrdestr.cc.

◆ GetCellTransportViews()

const std::vector< CellLBSView > & lbs::LBSSolver::GetCellTransportViews ( ) const

Returns a reference to the list of local cell transport views.

Definition at line 186 of file lbs_00_constrdestr.cc.

◆ GetHandleToPowerGenFieldFunc()

size_t lbs::LBSSolver::GetHandleToPowerGenFieldFunc ( ) const

Returns the local handle to the power generation field function, if enabled.

Definition at line 320 of file lbs_00_constrdestr.cc.

◆ GetInputParameters()

chi::InputParameters lbs::LBSSolver::GetInputParameters ( )
static

Returns the input parameters for this object.

Definition at line 19 of file lbs_00_constrdestr.cc.

◆ GetMatID2IsoSrcMap()

const std::map< int, IsotropicSrcPtr > & lbs::LBSSolver::GetMatID2IsoSrcMap ( ) const

Returns a reference to the map of material ids to Isotropic Srcs.

Definition at line 165 of file lbs_00_constrdestr.cc.

◆ GetMatID2XSMap()

const std::map< int, XSPtr > & lbs::LBSSolver::GetMatID2XSMap ( ) const

Returns a reference to the map of material ids to XSs.

Definition at line 159 of file lbs_00_constrdestr.cc.

◆ GetMaxPrecursorsPerMaterial()

size_t lbs::LBSSolver::GetMaxPrecursorsPerMaterial ( ) const

Returns the maximum number of precursors, for a material, as encountered accross all the materials. This will only be non-zero after initialization.

Definition at line 110 of file lbs_00_constrdestr.cc.

◆ GetNumPhiIterativeUnknowns()

std::pair< size_t, size_t > lbs::LBSSolver::GetNumPhiIterativeUnknowns ( )
virtual

Gets the local and global number of iterative unknowns. This normally is only the flux moments, however, the sweep based solvers might include delayed angular fluxes in this number.

Reimplemented in lbs::DiscreteOrdinatesSolver.

Definition at line 299 of file lbs_00_constrdestr.cc.

◆ GetPrimaryAGSSolver()

LBSSolver::AGSLinSolverPtr lbs::LBSSolver::GetPrimaryAGSSolver ( )

Definition at line 268 of file lbs_00_constrdestr.cc.

◆ GetSourceEventTag()

size_t lbs::LBSSolver::GetSourceEventTag ( ) const

Returns the source event tag used for logging the time it takes to set source moments.

Definition at line 81 of file lbs_00_constrdestr.cc.

◆ GetUnitCellMatrices()

const std::vector< UnitCellMatrices > & lbs::LBSSolver::GetUnitCellMatrices ( ) const

Returns read-only access to the unit cell matrices.

Definition at line 177 of file lbs_00_constrdestr.cc.

◆ GetWGSContext()

WGSContext< Mat, Vec, KSP > & lbs::LBSSolver::GetWGSContext ( int  groupset_id)

Definition at line 278 of file lbs_00_constrdestr.cc.

◆ GetWGSSolvers()

std::vector< LBSSolver::LinSolvePtr > & lbs::LBSSolver::GetWGSSolvers ( )

Definition at line 273 of file lbs_00_constrdestr.cc.

◆ GlobalNodeCount()

size_t lbs::LBSSolver::GlobalNodeCount ( ) const

Returns the global node count for the flux-moments data structures.

Definition at line 201 of file lbs_00_constrdestr.cc.

◆ Grid()

const chi_mesh::MeshContinuum & lbs::LBSSolver::Grid ( ) const

Obtains a reference to the grid.

Definition at line 183 of file lbs_00_constrdestr.cc.

◆ Groups()

const std::vector< LBSGroup > & lbs::LBSSolver::Groups ( ) const

Const accessor.

Definition at line 125 of file lbs_00_constrdestr.cc.

◆ Groupsets() [1/2]

std::vector< LBSGroupset > & lbs::LBSSolver::Groupsets ( )

Non-Const accessor.

Definition at line 135 of file lbs_00_constrdestr.cc.

◆ Groupsets() [2/2]

const std::vector< LBSGroupset > & lbs::LBSSolver::Groupsets ( ) const

Const accessor.

Definition at line 138 of file lbs_00_constrdestr.cc.

◆ GSProjectBackPhi0()

void lbs::LBSSolver::GSProjectBackPhi0 ( const LBSGroupset groupset,
const std::vector< double > &  input,
std::vector< double > &  output 
)

From the WGDSA DOFs, projects the scalar moments back into a primary STL vector.

Definition at line 107 of file lbs_03d_wgdsa.cc.

◆ GSScopedCopyPrimarySTLvectors() [1/2]

void lbs::LBSSolver::GSScopedCopyPrimarySTLvectors ( LBSGroupset groupset,
const std::vector< double > &  x_src,
std::vector< double > &  y 
)
virtual

Assembles a vector for a given groupset from a source vector.

Definition at line 143 of file lbs_07_vectorassembly.cc.

◆ GSScopedCopyPrimarySTLvectors() [2/2]

void lbs::LBSSolver::GSScopedCopyPrimarySTLvectors ( LBSGroupset groupset,
PhiSTLOption  from_which_phi,
PhiSTLOption  to_which_phi 
)
virtual

Assembles a vector for a given groupset from a source vector.

Reimplemented in lbs::DiscreteOrdinatesSolver.

Definition at line 171 of file lbs_07_vectorassembly.cc.

◆ Initialize()

void lbs::LBSSolver::Initialize ( )
overridevirtual

Initialize the solver.

Reimplemented from chi_physics::Solver.

Reimplemented in lbs::DiffusionDFEMSolver, lbs::DiscreteOrdinatesSolver, and lbs::DiscreteOrdinatesAdjointSolver.

Definition at line 8 of file lbs_01_main_initialize.cc.

◆ InitializeBoundaries()

void lbs::LBSSolver::InitializeBoundaries ( )
protected

Initializes transport related boundaries.

Definition at line 31 of file lbs_01h_init_boundaries.cc.

◆ InitializeFieldFunctions()

void lbs::LBSSolver::InitializeFieldFunctions ( )
protected

Definition at line 10 of file lbs_01ga_init_fieldfuncs.cc.

◆ InitializeGroupsets()

void lbs::LBSSolver::InitializeGroupsets ( )
protected

Initializes common groupset items.

Definition at line 5 of file lbs_01e_initgroupsets.cc.

◆ InitializeParrays()

void lbs::LBSSolver::InitializeParrays ( )
protectedvirtual

Initializes parallel arrays.

Definition at line 14 of file lbs_01g_initparrays.cc.

◆ InitializePointSources()

void lbs::LBSSolver::InitializePointSources ( )

Intializes all point sources.

ByRef

Definition at line 9 of file lbs_01i_init_point_sources.cc.

◆ InitializeSolverSchemes()

void lbs::LBSSolver::InitializeSolverSchemes ( )
protectedvirtual

Definition at line 9 of file lbs_01j_init_solver_schemes.cc.

◆ InitializeSpatialDiscretization()

void lbs::LBSSolver::InitializeSpatialDiscretization ( )
protectedvirtual

Reimplemented in lbs::DiscreteOrdinatesCurvilinearSolver.

Definition at line 13 of file lbs_01d_init_spat_discr.cc.

◆ InitializeWGSSolvers()

virtual void lbs::LBSSolver::InitializeWGSSolvers ( )
inlineprotectedvirtual

Reimplemented in lbs::DiffusionDFEMSolver, and lbs::DiscreteOrdinatesSolver.

Definition at line 236 of file lbs_solver.h.

◆ InitMaterials()

void lbs::LBSSolver::InitMaterials ( )

Initializes default materials and physics materials.

Definition at line 16 of file lbs_01c_initmaterials.cc.

◆ InitTGDSA()

void lbs::LBSSolver::InitTGDSA ( LBSGroupset groupset)
protected

Initializes the Within-Group DSA solver.

Definition at line 11 of file lbs_03e_tgdsa.cc.

◆ InitWGDSA()

void lbs::LBSSolver::InitWGDSA ( LBSGroupset groupset,
bool  vaccum_bcs_are_dirichlet = true 
)

Initializes the Within-Group DSA solver.

Definition at line 10 of file lbs_03d_wgdsa.cc.

◆ LastRestartWrite() [1/2]

double & lbs::LBSSolver::LastRestartWrite ( )

Returns a reference to the time at which the last restart was written.

Definition at line 87 of file lbs_00_constrdestr.cc.

◆ LastRestartWrite() [2/2]

double lbs::LBSSolver::LastRestartWrite ( ) const

Returns the time at which the last restart was written.

Definition at line 84 of file lbs_00_constrdestr.cc.

◆ LocalNodeCount()

size_t lbs::LBSSolver::LocalNodeCount ( ) const

Returns the local node count for the flux-moments data structures.

Definition at line 198 of file lbs_00_constrdestr.cc.

◆ MakeSourceMomentsFromPhi()

std::vector< double > lbs::LBSSolver::MakeSourceMomentsFromPhi ( )

Makes a source-moments vector from scattering and fission based on the latest phi-solution.

Definition at line 15 of file lbs_04c_io_flux_moments.cc.

◆ MapPhiFieldFunction()

size_t lbs::LBSSolver::MapPhiFieldFunction ( size_t  g,
size_t  m 
) const

Gets the local handle of a flux-moment based field function.

Definition at line 309 of file lbs_00_constrdestr.cc.

◆ NumGroups()

size_t lbs::LBSSolver::NumGroups ( ) const

Returns the number of groups for the solver. This will only be non-zero after initialization.

Definition at line 101 of file lbs_00_constrdestr.cc.

◆ NumMoments()

size_t lbs::LBSSolver::NumMoments ( ) const

Returns the number of moments for the solver. This will only be non-zero after initialization.

Definition at line 97 of file lbs_00_constrdestr.cc.

◆ NumPrecursors()

size_t lbs::LBSSolver::NumPrecursors ( ) const

Returns the number of precursors for the solver. This will only be non-zero after initialization.

Definition at line 105 of file lbs_00_constrdestr.cc.

◆ operator=()

LBSSolver & lbs::LBSSolver::operator= ( const LBSSolver )
delete

◆ Options() [1/2]

Options & lbs::LBSSolver::Options ( )

Returns a reference to the solver options.

Definition at line 90 of file lbs_00_constrdestr.cc.

◆ Options() [2/2]

const Options & lbs::LBSSolver::Options ( ) const

Returns a constant reference to the solver options.

Definition at line 93 of file lbs_00_constrdestr.cc.

◆ OptionsBlock()

chi::InputParameters lbs::LBSSolver::OptionsBlock ( )
static

Definition at line 13 of file lbs_00a_setoptions.cc.

◆ PerformInputChecks()

void lbs::LBSSolver::PerformInputChecks ( )
protectedvirtual

Performs general input checks before initialization continues.

Reimplemented in lbs::DiscreteOrdinatesCurvilinearSolver.

Definition at line 11 of file lbs_01a_inputchecks.cc.

◆ PhiNewLocal() [1/2]

std::vector< double > & lbs::LBSSolver::PhiNewLocal ( )

Read/write access to newest updated flux vector.

Definition at line 233 of file lbs_00_constrdestr.cc.

◆ PhiNewLocal() [2/2]

const std::vector< double > & lbs::LBSSolver::PhiNewLocal ( ) const

Read access to newest updated flux vector.

Definition at line 235 of file lbs_00_constrdestr.cc.

◆ PhiOldLocal() [1/2]

std::vector< double > & lbs::LBSSolver::PhiOldLocal ( )

Read/write access to last updated flux vector.

Definition at line 224 of file lbs_00_constrdestr.cc.

◆ PhiOldLocal() [2/2]

const std::vector< double > & lbs::LBSSolver::PhiOldLocal ( ) const

Read access to last updated flux vector.

Definition at line 227 of file lbs_00_constrdestr.cc.

◆ PointSources()

const std::vector< PointSource > & lbs::LBSSolver::PointSources ( ) const

Const accessor to the list of point sources.

Definition at line 153 of file lbs_00_constrdestr.cc.

◆ PrecursorsNewLocal() [1/2]

std::vector< double > & lbs::LBSSolver::PrecursorsNewLocal ( )

Read/write access to newest updated precursors vector.

Definition at line 241 of file lbs_00_constrdestr.cc.

◆ PrecursorsNewLocal() [2/2]

const std::vector< double > & lbs::LBSSolver::PrecursorsNewLocal ( ) const

Read access to newest updated precursors vector.

Definition at line 243 of file lbs_00_constrdestr.cc.

◆ PrintSimHeader()

void lbs::LBSSolver::PrintSimHeader ( )
protected

Prints header information of simulation.

Definition at line 9 of file lbs_01b_print_simheader.cc.

◆ PsiNewLocal() [1/2]

std::vector< VecDbl > & lbs::LBSSolver::PsiNewLocal ( )

Read/write access to newest updated angular flux vector.

Definition at line 249 of file lbs_00_constrdestr.cc.

◆ PsiNewLocal() [2/2]

const std::vector< VecDbl > & lbs::LBSSolver::PsiNewLocal ( ) const

Read access to newest updated angular flux vector.

Definition at line 251 of file lbs_00_constrdestr.cc.

◆ QMomentsLocal() [1/2]

std::vector< double > & lbs::LBSSolver::QMomentsLocal ( )

Read/write access to source moments vector.

Definition at line 204 of file lbs_00_constrdestr.cc.

◆ QMomentsLocal() [2/2]

const std::vector< double > & lbs::LBSSolver::QMomentsLocal ( ) const

Read access to source moments vector.

Definition at line 206 of file lbs_00_constrdestr.cc.

◆ ReadFluxMoments()

void lbs::LBSSolver::ReadFluxMoments ( const std::string &  file_base,
std::vector< double > &  flux_moments,
bool  single_file = false 
)

Reads a flux-moments vector from a file in the specified vector.

Definition at line 154 of file lbs_04c_io_flux_moments.cc.

◆ ReadGroupsetAngularFluxes()

void lbs::LBSSolver::ReadGroupsetAngularFluxes ( LBSGroupset groupset,
const std::string &  file_base 
)

Prints the groupset's angular fluxes to file.

Definition at line 107 of file lbs_04b_angular_fluxes.cc.

◆ ReadRestartData()

void lbs::LBSSolver::ReadRestartData ( const std::string &  folder_name,
const std::string &  file_base 
)

Read phi_old from restart file.

Definition at line 92 of file lbs_04a_restartdata.cc.

◆ ScalePhiVector()

void lbs::LBSSolver::ScalePhiVector ( PhiSTLOption  which_phi,
double  value 
)
virtual

Scales a flux moment vector. For sweep methods the delayed angular fluxes will also be scaled.

Reimplemented in lbs::DiscreteOrdinatesSolver.

Definition at line 38 of file lbs_07_vectorassembly.cc.

◆ SetBoundaryOptions()

void lbs::LBSSolver::SetBoundaryOptions ( const chi::InputParameters params)

Definition at line 250 of file lbs_00a_setoptions.cc.

◆ SetGroupScopedPETScVecFromPrimarySTLvector()

void lbs::LBSSolver::SetGroupScopedPETScVecFromPrimarySTLvector ( int  first_group_id,
int  last_group_id,
Vec  x,
const std::vector< double > &  y 
)
virtual

Assembles a vector for a given group span from a source vector.

Definition at line 217 of file lbs_07_vectorassembly.cc.

◆ SetGSPETScVecFromPrimarySTLvector()

void lbs::LBSSolver::SetGSPETScVecFromPrimarySTLvector ( LBSGroupset groupset,
Vec  x,
PhiSTLOption  which_phi 
)
virtual

Assembles a vector for a given groupset from a source vector.

Reimplemented in lbs::DiscreteOrdinatesSolver.

Definition at line 55 of file lbs_07_vectorassembly.cc.

◆ SetMultiGSPETScVecFromPrimarySTLvector()

void lbs::LBSSolver::SetMultiGSPETScVecFromPrimarySTLvector ( const std::vector< int > &  gs_ids,
Vec  x,
PhiSTLOption  which_phi 
)
virtual

Assembles a PETSc vector from multiple groupsets.

Reimplemented in lbs::DiscreteOrdinatesSolver.

Definition at line 290 of file lbs_07_vectorassembly.cc.

◆ SetOptions()

void lbs::LBSSolver::SetOptions ( const chi::InputParameters params)

Definition at line 151 of file lbs_00a_setoptions.cc.

◆ SetPhiFromFieldFunctions()

void lbs::LBSSolver::SetPhiFromFieldFunctions ( PhiSTLOption  which_phi,
const std::vector< size_t > &  m_indices,
const std::vector< size_t > &  g_indices 
)

Sets the internal phi vector to the value in the associated field function.

Definition at line 136 of file lbs_05a_update_fieldfuncs.cc.

◆ SetPhiVectorScalarValues()

void lbs::LBSSolver::SetPhiVectorScalarValues ( std::vector< double > &  phi_vector,
double  value 
)
virtual

Sets a value to the zeroth (scalar) moment of the vector.

Definition at line 8 of file lbs_07_vectorassembly.cc.

◆ SetPrimarySTLvectorFromGroupScopedPETScVec()

void lbs::LBSSolver::SetPrimarySTLvectorFromGroupScopedPETScVec ( int  first_group_id,
int  last_group_id,
Vec  x_src,
std::vector< double > &  y 
)
virtual

Assembles a vector for a given groupset from a source vector.

Definition at line 253 of file lbs_07_vectorassembly.cc.

◆ SetPrimarySTLvectorFromGSPETScVec()

void lbs::LBSSolver::SetPrimarySTLvectorFromGSPETScVec ( LBSGroupset groupset,
Vec  x_src,
PhiSTLOption  which_phi 
)
virtual

Assembles a vector for a given groupset from a source vector.

Reimplemented in lbs::DiscreteOrdinatesSolver.

Definition at line 99 of file lbs_07_vectorassembly.cc.

◆ SetPrimarySTLvectorFromMultiGSPETScVecFrom()

void lbs::LBSSolver::SetPrimarySTLvectorFromMultiGSPETScVecFrom ( const std::vector< int > &  gs_ids,
Vec  x_src,
PhiSTLOption  which_phi 
)
virtual

Disassembles a multiple Groupset PETSc vector STL vectors.

Reimplemented in lbs::DiscreteOrdinatesSolver.

Definition at line 339 of file lbs_07_vectorassembly.cc.

◆ SpatialDiscretization()

const chi_math::SpatialDiscretization & lbs::LBSSolver::SpatialDiscretization ( ) const

Obtains a reference to the spatial discretization.

Definition at line 171 of file lbs_00_constrdestr.cc.

◆ SweepBoundaries()

const std::map< uint64_t, std::shared_ptr< SweepBndry > > & lbs::LBSSolver::SweepBoundaries ( ) const

Returns the sweep boundaries as a read only reference

Definition at line 258 of file lbs_00_constrdestr.cc.

◆ UnknownManager()

const chi_math::UnknownManager & lbs::LBSSolver::UnknownManager ( ) const

Obtains a reference to the unknown manager for flux-moments.

Definition at line 192 of file lbs_00_constrdestr.cc.

◆ UpdateFieldFunctions()

void lbs::LBSSolver::UpdateFieldFunctions ( )

Copy relevant section of phi_old to the field functions.

Definition at line 13 of file lbs_05a_update_fieldfuncs.cc.

◆ WGSCopyOnlyPhi0()

std::vector< double > lbs::LBSSolver::WGSCopyOnlyPhi0 ( const LBSGroupset groupset,
const std::vector< double > &  phi_in 
)

Creates a vector from a lbs primary stl vector where only the scalar moments are mapped to the DOFs needed by WGDSA.

Definition at line 69 of file lbs_03d_wgdsa.cc.

◆ WriteFluxMoments()

void lbs::LBSSolver::WriteFluxMoments ( const std::string &  file_base,
const std::vector< double > &  flux_moments 
)

Writes a given flux-moments vector to file.

Definition at line 37 of file lbs_04c_io_flux_moments.cc.

◆ WriteGroupsetAngularFluxes()

void lbs::LBSSolver::WriteGroupsetAngularFluxes ( const LBSGroupset groupset,
const std::string &  file_base 
)

Writes the groupset's angular fluxes to file.

Definition at line 15 of file lbs_04b_angular_fluxes.cc.

◆ WriteRestartData()

void lbs::LBSSolver::WriteRestartData ( const std::string &  folder_name,
const std::string &  file_base 
)

Writes phi_old to restart file.

Definition at line 13 of file lbs_04a_restartdata.cc.

Field Documentation

◆ active_set_source_function_

SetSourceFunction lbs::LBSSolver::active_set_source_function_
protected

Definition at line 99 of file lbs_solver.h.

◆ ags_solvers_

std::vector<AGSLinSolverPtr> lbs::LBSSolver::ags_solvers_
protected

Definition at line 101 of file lbs_solver.h.

◆ boundary_preferences_

std::map<uint64_t, BoundaryPreference> lbs::LBSSolver::boundary_preferences_
protected

Definition at line 85 of file lbs_solver.h.

◆ cell_transport_views_

std::vector<lbs::CellLBSView> lbs::LBSSolver::cell_transport_views_
protected

Definition at line 83 of file lbs_solver.h.

◆ discretization_

std::shared_ptr<chi_math::SpatialDiscretization> lbs::LBSSolver::discretization_ = nullptr
protected

Definition at line 74 of file lbs_solver.h.

◆ ext_src_moments_local_

std::vector<double> lbs::LBSSolver::ext_src_moments_local_
protected

Definition at line 94 of file lbs_solver.h.

◆ flux_moments_uk_man_

chi_math::UnknownManager lbs::LBSSolver::flux_moments_uk_man_
protected

Definition at line 88 of file lbs_solver.h.

◆ glob_node_count_

uint64_t lbs::LBSSolver::glob_node_count_ = 0
protected

Definition at line 92 of file lbs_solver.h.

◆ grid_face_histogram_

GridFaceHistogramPtr lbs::LBSSolver::grid_face_histogram_ = nullptr
protected

Definition at line 79 of file lbs_solver.h.

◆ grid_local_comm_set_

MPILocalCommSetPtr lbs::LBSSolver::grid_local_comm_set_ = nullptr
protected

Definition at line 78 of file lbs_solver.h.

◆ grid_nodal_mappings_

std::vector<CellFaceNodalMapping> lbs::LBSSolver::grid_nodal_mappings_
protected

Definition at line 77 of file lbs_solver.h.

◆ grid_ptr_

chi_mesh::MeshContinuumPtr lbs::LBSSolver::grid_ptr_
protected

Definition at line 75 of file lbs_solver.h.

◆ groups_

std::vector<LBSGroup> lbs::LBSSolver::groups_
protected

Definition at line 67 of file lbs_solver.h.

◆ groupsets_

std::vector<LBSGroupset> lbs::LBSSolver::groupsets_
protected

Definition at line 68 of file lbs_solver.h.

◆ last_restart_write_

double lbs::LBSSolver::last_restart_write_ = 0.0
protected

Definition at line 59 of file lbs_solver.h.

◆ local_node_count_

uint64_t lbs::LBSSolver::local_node_count_ = 0
protected

Definition at line 91 of file lbs_solver.h.

◆ matid_to_src_map_

std::map<int, IsotropicSrcPtr> lbs::LBSSolver::matid_to_src_map_
protected

Definition at line 72 of file lbs_solver.h.

◆ matid_to_xs_map_

std::map<int, XSPtr> lbs::LBSSolver::matid_to_xs_map_
protected

Definition at line 71 of file lbs_solver.h.

◆ max_cell_dof_count_

size_t lbs::LBSSolver::max_cell_dof_count_ = 0
protected

Definition at line 90 of file lbs_solver.h.

◆ max_precursors_per_material_

size_t lbs::LBSSolver::max_precursors_per_material_ = 0
protected

Definition at line 65 of file lbs_solver.h.

◆ num_groups_

size_t lbs::LBSSolver::num_groups_ = 0
protected

Definition at line 63 of file lbs_solver.h.

◆ num_moments_

size_t lbs::LBSSolver::num_moments_ = 0
protected

Definition at line 62 of file lbs_solver.h.

◆ num_precursors_

size_t lbs::LBSSolver::num_precursors_ = 0
protected

Definition at line 64 of file lbs_solver.h.

◆ options_

lbs::Options lbs::LBSSolver::options_
protected

Definition at line 61 of file lbs_solver.h.

◆ phi_field_functions_local_map_

std::map<std::pair<size_t, size_t>, size_t> lbs::LBSSolver::phi_field_functions_local_map_
protected

Definition at line 105 of file lbs_solver.h.

◆ phi_new_local_

std::vector<double> lbs::LBSSolver::phi_new_local_
protected

Definition at line 95 of file lbs_solver.h.

◆ phi_old_local_

std::vector<double> lbs::LBSSolver::phi_old_local_
protected

Definition at line 95 of file lbs_solver.h.

◆ point_sources_

std::vector<PointSource> lbs::LBSSolver::point_sources_
protected

Definition at line 69 of file lbs_solver.h.

◆ power_gen_fieldfunc_local_handle_

size_t lbs::LBSSolver::power_gen_fieldfunc_local_handle_ = 0
protected

Definition at line 106 of file lbs_solver.h.

◆ precursor_new_local_

std::vector<double> lbs::LBSSolver::precursor_new_local_
protected

Definition at line 97 of file lbs_solver.h.

◆ primary_ags_solver_

AGSLinSolverPtr lbs::LBSSolver::primary_ags_solver_
protected

Definition at line 103 of file lbs_solver.h.

◆ psi_new_local_

std::vector<std::vector<double> > lbs::LBSSolver::psi_new_local_
protected

Definition at line 96 of file lbs_solver.h.

◆ q_moments_local_

std::vector<double> lbs::LBSSolver::q_moments_local_
protected

Definition at line 94 of file lbs_solver.h.

◆ source_event_tag_

size_t lbs::LBSSolver::source_event_tag_ = 0
protected

Definition at line 58 of file lbs_solver.h.

◆ sweep_boundaries_

std::map<uint64_t, std::shared_ptr<SweepBndry> > lbs::LBSSolver::sweep_boundaries_
protected

Definition at line 86 of file lbs_solver.h.

◆ time_integration_

std::shared_ptr<const chi_math::TimeIntegration> lbs::LBSSolver::time_integration_ = nullptr
protected

Time integration parameter meant to be set by an executor

Definition at line 109 of file lbs_solver.h.

◆ unit_cell_matrices_

std::vector<UnitCellMatrices> lbs::LBSSolver::unit_cell_matrices_
protected

Definition at line 81 of file lbs_solver.h.

◆ unit_ghost_cell_matrices_

std::map<uint64_t, UnitCellMatrices> lbs::LBSSolver::unit_ghost_cell_matrices_
protected

Definition at line 82 of file lbs_solver.h.

◆ wgs_solvers_

std::vector<LinSolvePtr> lbs::LBSSolver::wgs_solvers_
protected

Definition at line 102 of file lbs_solver.h.


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