Chi-Tech
lbs_00_constrdestr.cc
Go to the documentation of this file.
1#include "lbs_solver.h"
2
3#include "chi_log.h"
4
7
8namespace lbs
9{
10//RegisterChiObject(lbs, LBSSolver); Should not be constructible
11
12/**Base class constructor.*/
13LBSSolver::LBSSolver(const std::string& text_name)
14 : chi_physics::Solver(text_name)
15{
16}
17
18/**Returns the input parameters for this object.*/
20{
23
24 // clang-format off
25 params.ChangeExistingParamToOptional("name", "LBSDatablock");
26
27 params.AddRequiredParameter<size_t>(
28 "num_groups", "The total number of groups within the solver");
29
31 "groupsets",
32 "An array of blocks each specifying the input parameters for a "
33 "<TT>lbs::LBSGroupset</TT>.");
34 params.LinkParameterToBlock("groupsets", "lbs::LBSGroupset");
35
37 "Block of options. See <TT>lbs::OptionsBlock</TT>.");
38 params.LinkParameterToBlock("options", "lbs::OptionsBlock");
39 // clang-format on
40
41 return params;
42}
43
44/**Input parameters based construction.*/
46 : chi_physics::Solver(params)
47{
48 //=================================== Make groups
49 const size_t num_groups = params.GetParamValue<size_t>("num_groups");
50 for (size_t g = 0; g < num_groups; ++g)
51 groups_.push_back(LBSGroup(static_cast<int>(g)));
52
53 //=================================== Make groupsets
54 const auto& groupsets_array = params.GetParam("groupsets");
55
56 const size_t num_gs = groupsets_array.NumParameters();
57 for (size_t gs = 0; gs < num_gs; ++gs)
58 {
59 const auto& groupset_params = groupsets_array.GetParam(gs);
60
61 chi::InputParameters gs_input_params =
63 gs_input_params.SetObjectType("LBSSolver:LBSGroupset");
64 gs_input_params.AssignParameters(groupset_params);
65
66 groupsets_.emplace_back(gs_input_params, gs, *this);
67 } // for gs
68
69 //=================================== Options
70 if (params.ParametersAtAssignment().Has("options"))
71 {
72 auto options_params = LBSSolver::OptionsBlock();
73 options_params.AssignParameters(params.GetParam("options"));
74
75 this->SetOptions(options_params);
76 }
77}
78
79/**Returns the source event tag used for logging the time it
80 * takes to set source moments.*/
82
83/**Returns the time at which the last restart was written.*/
85
86/**Returns a reference to the time at which the last restart was written.*/
88
89/**Returns a reference to the solver options.*/
91
92/**Returns a constant reference to the solver options.*/
93const Options& LBSSolver::Options() const { return options_; }
94
95/**Returns the number of moments for the solver. This will only be non-zero
96 * after initialization.*/
97size_t LBSSolver::NumMoments() const { return num_moments_; }
98
99/**Returns the number of groups for the solver. This will only be non-zero
100 * after initialization.*/
101size_t LBSSolver::NumGroups() const { return num_groups_; }
102
103/**Returns the number of precursors for the solver. This will only be non-zero
104 * after initialization.*/
106
107/**Returns the maximum number of precursors, for a material, as encountered
108 * accross all the materials. This will only be non-zero
109 * after initialization.*/
111{
113}
114
115/**Adds a group to the list of groups. If group id < 0, the id will be logically
116 * derived from the list size. If >= 0 the id will be set to the id specified.*/
118{
119 if (id < 0) groups_.emplace_back(static_cast<int>(groups_.size()));
120 else
121 groups_.emplace_back(id);
122}
123
124/**Const accessor.*/
125const std::vector<LBSGroup>& LBSSolver::Groups() const { return groups_; }
126
127/**Adds a groupset to the list of groupsets. The groupset id will be logically
128 * derived from the list size.*/
130{
131 groupsets_.emplace_back(static_cast<int>(groupsets_.size()));
132}
133
134/**Non-Const accessor.*/
135std::vector<LBSGroupset>& LBSSolver::Groupsets() { return groupsets_; }
136
137/**Const accessor.*/
138const std::vector<LBSGroupset>& LBSSolver::Groupsets() const
139{
140 return groupsets_;
141}
142
143/**Adds a point source to the solver's point source list.*/
145{
146 point_sources_.push_back(std::move(psrc));
147}
148
149/**Clears all the point sources from the solver's point source list.*/
151
152/**Const accessor to the list of point sources.*/
153const std::vector<PointSource>& LBSSolver::PointSources() const
154{
155 return point_sources_;
156}
157
158/**Returns a reference to the map of material ids to XSs.*/
159const std::map<int, XSPtr>& LBSSolver::GetMatID2XSMap() const
160{
161 return matid_to_xs_map_;
162}
163
164/**Returns a reference to the map of material ids to Isotropic Srcs.*/
165const std::map<int, IsotropicSrcPtr>& LBSSolver::GetMatID2IsoSrcMap() const
166{
167 return matid_to_src_map_;
168}
169
170/**Obtains a reference to the spatial discretization.*/
172{
173 return *discretization_;
174}
175
176/**Returns read-only access to the unit cell matrices.*/
177const std::vector<UnitCellMatrices>& LBSSolver::GetUnitCellMatrices() const
178{
179 return unit_cell_matrices_;
180}
181
182/**Obtains a reference to the grid.*/
184
185/**Returns a reference to the list of local cell transport views.*/
186const std::vector<CellLBSView>& LBSSolver::GetCellTransportViews() const
187{
189}
190
191/**Obtains a reference to the unknown manager for flux-moments.*/
193{
195}
196
197/**Returns the local node count for the flux-moments data structures.*/
199
200/**Returns the global node count for the flux-moments data structures.*/
202
203/**Read/write access to source moments vector.*/
204std::vector<double>& LBSSolver::QMomentsLocal() { return q_moments_local_; }
205/**Read access to source moments vector.*/
206const std::vector<double>& LBSSolver::QMomentsLocal() const
207{
208 return q_moments_local_;
209}
210
211/**Read/write access to exterior src moments vector.*/
212std::vector<double>& LBSSolver::ExtSrcMomentsLocal()
213{
215}
216
217/**Read access to exterior src moments vector.*/
218const std::vector<double>& LBSSolver::ExtSrcMomentsLocal() const
219{
221}
222
223/**Read/write access to last updated flux vector.*/
224std::vector<double>& LBSSolver::PhiOldLocal() { return phi_old_local_; }
225
226/**Read access to last updated flux vector.*/
227const std::vector<double>& LBSSolver::PhiOldLocal() const
228{
229 return phi_old_local_;
230}
231
232/**Read/write access to newest updated flux vector.*/
233std::vector<double>& LBSSolver::PhiNewLocal() { return phi_new_local_; }
234/**Read access to newest updated flux vector.*/
235const std::vector<double>& LBSSolver::PhiNewLocal() const
236{
237 return phi_new_local_;
238}
239
240/**Read/write access to newest updated precursors vector.*/
241std::vector<double>& LBSSolver::PrecursorsNewLocal() { return phi_new_local_; }
242/**Read access to newest updated precursors vector.*/
243const std::vector<double>& LBSSolver::PrecursorsNewLocal() const
244{
245 return phi_new_local_;
246}
247
248/**Read/write access to newest updated angular flux vector.*/
249std::vector<VecDbl>& LBSSolver::PsiNewLocal() { return psi_new_local_; }
250/**Read access to newest updated angular flux vector.*/
251const std::vector<VecDbl>& LBSSolver::PsiNewLocal() const
252{
253 return psi_new_local_;
254}
255
256/**Returns the sweep boundaries as a read only reference*/
257const std::map<uint64_t, std::shared_ptr<SweepBndry>>&
259{
260 return sweep_boundaries_;
261}
262
264{
266}
267
269{
270 return primary_ags_solver_;
271}
272
273std::vector<LBSSolver::LinSolvePtr>& LBSSolver::GetWGSSolvers()
274{
275 return wgs_solvers_;
276}
277
279{
280 auto& wgs_solver = wgs_solvers_[groupset_id];
281 auto& raw_context = wgs_solver->GetContext();
282
283 typedef WGSContext<Mat, Vec, KSP> LBSWGSContext;
284 auto wgs_context_ptr = std::dynamic_pointer_cast<LBSWGSContext>(raw_context);
285
286 ChiLogicalErrorIf(not wgs_context_ptr, "Failed to cast WGSContext");
287 return *wgs_context_ptr;
288}
289
290/**Read/Write access to the boundary preferences.*/
291std::map<uint64_t, BoundaryPreference>& LBSSolver::BoundaryPreferences()
292{
294}
295
296/**Gets the local and global number of iterative unknowns. This normally is
297 * only the flux moments, however, the sweep based solvers might include
298 * delayed angular fluxes in this number.*/
299std::pair<size_t, size_t> LBSSolver::GetNumPhiIterativeUnknowns()
300{
301 const auto& sdm = *discretization_;
302 const size_t num_local_phi_dofs = sdm.GetNumLocalDOFs(flux_moments_uk_man_);
303 const size_t num_globl_phi_dofs = sdm.GetNumGlobalDOFs(flux_moments_uk_man_);
304
305 return {num_local_phi_dofs, num_globl_phi_dofs};
306}
307
308/**Gets the local handle of a flux-moment based field function.*/
309size_t LBSSolver::MapPhiFieldFunction(size_t g, size_t m) const
310{
312 std::string("Failure to map phi field function g") +
313 std::to_string(g) + " m" + std::to_string(m));
314
315 return phi_field_functions_local_map_.at({g, m});
316}
317
318/**Returns the local handle to the power generation field function, if
319 * enabled.*/
321{
323 "Called when options_.power_field_function_on == false");
324
326}
327
328} // namespace lbs
#define ChiLogicalErrorIf(condition, message)
void AssignParameters(const ParameterBlock &params)
Assigns parameters with thorough type checks, deprecation checks, unused parameter checks.
void AddRequiredParameter(const std::string &name, const std::string &doc_string)
void AddRequiredParameterArray(const std::string &name, const std::string &doc_string)
const ParameterBlock & ParametersAtAssignment() const
void SetObjectType(const std::string &obj_type)
void AddOptionalParameterBlock(const std::string &name, const ParameterBlock &block, const std::string &doc_string)
void LinkParameterToBlock(const std::string &param_name, const std::string &block_name)
void ChangeExistingParamToOptional(const std::string &name, T value, const std::string &doc_string="")
bool Has(const std::string &param_name) const
size_t NumParameters() const
T GetParamValue(const std::string &param_name) const
ParameterBlock & GetParam(const std::string &param_name)
static chi::InputParameters GetInputParameters()
Definition: chi_solver.cc:16
static chi::InputParameters GetInputParameters()
Definition: lbs_groupset.cc:20
AGSLinSolverPtr GetPrimaryAGSSolver()
std::vector< double > & PhiNewLocal()
static chi::InputParameters GetInputParameters()
std::map< uint64_t, BoundaryPreference > & BoundaryPreferences()
std::vector< std::vector< double > > psi_new_local_
Definition: lbs_solver.h:96
std::vector< PointSource > point_sources_
Definition: lbs_solver.h:69
chi_mesh::MeshContinuumPtr grid_ptr_
Definition: lbs_solver.h:75
std::shared_ptr< chi_math::SpatialDiscretization > discretization_
Definition: lbs_solver.h:74
std::vector< lbs::CellLBSView > cell_transport_views_
Definition: lbs_solver.h:83
std::vector< UnitCellMatrices > unit_cell_matrices_
Definition: lbs_solver.h:81
virtual std::pair< size_t, size_t > GetNumPhiIterativeUnknowns()
size_t MapPhiFieldFunction(size_t g, size_t m) const
std::vector< double > & QMomentsLocal()
std::vector< LBSGroupset > & Groupsets()
std::vector< double > & PrecursorsNewLocal()
size_t GetHandleToPowerGenFieldFunc() const
const chi_math::SpatialDiscretization & SpatialDiscretization() const
std::vector< double > phi_new_local_
Definition: lbs_solver.h:95
std::vector< LinSolvePtr > & GetWGSSolvers()
size_t NumGroups() const
std::vector< LBSGroup > groups_
Definition: lbs_solver.h:67
const std::vector< PointSource > & PointSources() const
std::map< uint64_t, std::shared_ptr< SweepBndry > > sweep_boundaries_
Definition: lbs_solver.h:86
SetSourceFunction GetActiveSetSourceFunction() const
double last_restart_write_
Definition: lbs_solver.h:59
size_t LocalNodeCount() const
size_t num_precursors_
Definition: lbs_solver.h:64
chi_math::UnknownManager flux_moments_uk_man_
Definition: lbs_solver.h:88
const std::vector< lbs::CellLBSView > & GetCellTransportViews() const
std::map< uint64_t, BoundaryPreference > boundary_preferences_
Definition: lbs_solver.h:85
const std::vector< LBSGroup > & Groups() const
size_t power_gen_fieldfunc_local_handle_
Definition: lbs_solver.h:106
std::vector< VecDbl > & PsiNewLocal()
std::vector< double > & ExtSrcMomentsLocal()
std::vector< double > q_moments_local_
Definition: lbs_solver.h:94
const std::map< uint64_t, std::shared_ptr< SweepBndry > > & SweepBoundaries() const
size_t num_groups_
Definition: lbs_solver.h:63
uint64_t glob_node_count_
Definition: lbs_solver.h:92
const chi_math::UnknownManager & UnknownManager() const
lbs::Options options_
Definition: lbs_solver.h:61
size_t GetMaxPrecursorsPerMaterial() const
double LastRestartWrite() const
std::vector< double > ext_src_moments_local_
Definition: lbs_solver.h:94
std::vector< double > phi_old_local_
Definition: lbs_solver.h:95
std::map< std::pair< size_t, size_t >, size_t > phi_field_functions_local_map_
Definition: lbs_solver.h:105
const chi_mesh::MeshContinuum & Grid() const
std::map< int, IsotropicSrcPtr > matid_to_src_map_
Definition: lbs_solver.h:72
size_t source_event_tag_
Definition: lbs_solver.h:58
size_t NumMoments() const
const std::map< int, XSPtr > & GetMatID2XSMap() const
const std::map< int, IsotropicSrcPtr > & GetMatID2IsoSrcMap() const
WGSContext< Mat, Vec, KSP > & GetWGSContext(int groupset_id)
void AddGroup(int id)
std::vector< LBSGroupset > groupsets_
Definition: lbs_solver.h:68
size_t GlobalNodeCount() const
const std::vector< UnitCellMatrices > & GetUnitCellMatrices() const
std::vector< LinSolvePtr > wgs_solvers_
Definition: lbs_solver.h:102
std::shared_ptr< AGSLinearSolver< Mat, Vec, KSP > > AGSLinSolverPtr
Definition: lbs_solver.h:52
AGSLinSolverPtr primary_ags_solver_
Definition: lbs_solver.h:103
lbs::Options & Options()
uint64_t local_node_count_
Definition: lbs_solver.h:91
size_t NumPrecursors() const
std::vector< double > & PhiOldLocal()
std::map< int, XSPtr > matid_to_xs_map_
Definition: lbs_solver.h:71
size_t max_precursors_per_material_
Definition: lbs_solver.h:65
SetSourceFunction active_set_source_function_
Definition: lbs_solver.h:99
void AddPointSource(PointSource psrc)
static chi::InputParameters OptionsBlock()
void SetOptions(const chi::InputParameters &params)
size_t num_moments_
Definition: lbs_solver.h:62
LBSSolver(const std::string &text_name)
size_t GetSourceEventTag() const
std::function< void(LBSGroupset &groupset, std::vector< double > &destination_q, const std::vector< double > &phi, SourceFlags source_flags)> SetSourceFunction
Definition: lbs_structs.h:110
bool power_field_function_on
Definition: lbs_structs.h:147