Chi-Tech
lbsDO_sweepdata_00b_initFLUDS.cc
Go to the documentation of this file.
2
5
9
10#include "console/chi_console.h"
11
14
15#include "chi_runtime.h"
16#include "chi_log.h"
17#include "chi_mpi.h"
18
19#include <iomanip>
20
21#include "utils/chi_timer.h"
23
24// ###################################################################
25/**Initializes fluds_ data structures.*/
27{
28 namespace sweep_namespace = chi_mesh::sweep_management;
29 typedef sweep_namespace::AngleSetGroup TAngleSetGroup;
30 typedef sweep_namespace::AAH_AngleSet TAAH_AngleSet;
31
32 const auto& quadrature_sweep_info =
34
35 const auto& unique_so_groupings = quadrature_sweep_info.first;
36 const auto& dir_id_to_so_map = quadrature_sweep_info.second;
37
38 const size_t gs_num_grps = groupset.groups_.size();
39 const size_t gs_num_ss = groupset.grp_subset_infos_.size();
40
41 //=========================================== Passing the sweep boundaries
42 // to the angle aggregation
44 groupset.angle_agg_ = std::make_shared<AngleAgg>(
45 sweep_boundaries_, gs_num_grps, gs_num_ss, groupset.quadrature_, grid_ptr_);
46
47 TAngleSetGroup angle_set_group;
48 size_t angle_set_id = 0;
49 for (const auto& so_grouping : unique_so_groupings)
50 {
51 const size_t master_dir_id = so_grouping.front();
52 const size_t so_id = dir_id_to_so_map.at(master_dir_id);
53
54 const auto& sweep_ordering =
55 quadrature_spds_map_[groupset.quadrature_][so_id];
56 const auto& fluds_common_data =
58
59 // Compute direction subsets
60 const auto dir_subsets =
61 chi::MakeSubSets(so_grouping.size(), groupset.master_num_ang_subsets_);
62
63 for (size_t gs_ss = 0; gs_ss < gs_num_ss; gs_ss++)
64 {
65 const size_t gs_ss_size = groupset.grp_subset_infos_[gs_ss].ss_size;
66 for (const auto& dir_ss_info : dir_subsets)
67 {
68 const auto& dir_ss_begin = dir_ss_info.ss_begin;
69 const auto& dir_ss_end = dir_ss_info.ss_end;
70 const auto& dir_ss_size = dir_ss_info.ss_size;
71
72 std::vector<size_t> angle_indices(dir_ss_size, 0);
73 {
74 size_t k = 0;
75 for (size_t n = dir_ss_begin; n <= dir_ss_end; ++n)
76 angle_indices[k++] = so_grouping[n];
77 }
78
79 if (sweep_type_ == "AAH")
80 {
81 using namespace chi_mesh::sweep_management;
82 std::shared_ptr<FLUDS> fluds = std::make_shared<AAH_FLUDS>(
83 gs_ss_size,
84 angle_indices.size(),
85 dynamic_cast<const AAH_FLUDSCommonData&>(fluds_common_data));
86
87 auto angleSet =
88 std::make_shared<TAAH_AngleSet>(angle_set_id++,
89 gs_ss_size,
90 gs_ss,
91 *sweep_ordering,
92 fluds,
93 angle_indices,
97
98 angle_set_group.AngleSets().push_back(angleSet);
99 }
100 else if (sweep_type_ == "CBC")
101 {
103 "When using sweep_type \"CBC\" then "
104 "\"save_angular_flux\" must be true.");
105 using namespace chi_mesh::sweep_management;
106 std::shared_ptr<FLUDS> fluds = std::make_shared<CBC_FLUDS>(
107 gs_ss_size,
108 angle_indices.size(),
109 dynamic_cast<const CBC_FLUDSCommonData&>(fluds_common_data),
110 psi_new_local_[groupset.id_],
111 groupset.psi_uk_man_,
113
114 auto angleSet = std::make_shared<CBC_AngleSet>(angle_set_id++,
115 gs_ss_size,
116 *sweep_ordering,
117 fluds,
118 angle_indices,
120 gs_ss,
122
123 angle_set_group.AngleSets().push_back(angleSet);
124 }
125 else
126 ChiInvalidArgument("Unsupported sweeptype \"" + sweep_type_ + "\"");
127 } // for an_ss
128 } // for gs_ss
129 } // for so_grouping
130
131 groupset.angle_agg_->angle_set_groups.push_back(std::move(angle_set_group));
132
135 << " Initialized Angle Aggregation. "
136 << " Process memory = " << std::setprecision(3)
138
140}
#define ChiLogicalErrorIf(condition, message)
#define ChiInvalidArgument(message)
static chi::Timer program_timer
Definition: chi_runtime.h:79
static chi::ChiLog & log
Definition: chi_runtime.h:81
static chi::MPI_Info & mpi
Definition: chi_runtime.h:78
LogStream Log(LOG_LVL level=LOG_0)
Definition: chi_log.cc:35
static double GetMemoryUsageInMB()
Get current memory usage in Megabytes.
void Barrier() const
Definition: mpi_info.cc:38
std::string GetTimeString() const
Definition: chi_timer.cc:38
std::vector< std::shared_ptr< AngleSet > > & AngleSets()
Definition: anglesetgroup.h:11
void InitFluxDataStructures(LBSGroupset &groupset)
std::map< AngQuadPtr, FLUDSCommonDataPtrs > quadrature_fluds_commondata_map_
std::map< AngQuadPtr, SwpOrderGroupingInfo > quadrature_unq_so_grouping_map_
std::map< AngQuadPtr, SPDS_ptrs > quadrature_spds_map_
int master_num_ang_subsets_
Definition: lbs_groupset.h:48
AngleAggPtr angle_agg_
Definition: lbs_groupset.h:43
std::vector< chi::SubSetInfo > grp_subset_infos_
Definition: lbs_groupset.h:50
std::shared_ptr< chi_math::AngularQuadrature > quadrature_
Definition: lbs_groupset.h:42
chi_math::UnknownManager psi_uk_man_
Definition: lbs_groupset.h:83
std::vector< LBSGroup > groups_
Definition: lbs_groupset.h:41
std::vector< std::vector< double > > psi_new_local_
Definition: lbs_solver.h:96
chi_mesh::MeshContinuumPtr grid_ptr_
Definition: lbs_solver.h:75
std::shared_ptr< chi_math::SpatialDiscretization > discretization_
Definition: lbs_solver.h:74
std::map< uint64_t, std::shared_ptr< SweepBndry > > sweep_boundaries_
Definition: lbs_solver.h:86
lbs::Options options_
Definition: lbs_solver.h:61
MPILocalCommSetPtr grid_local_comm_set_
Definition: lbs_solver.h:78
std::vector< SubSetInfo > MakeSubSets(size_t num_items, size_t desired_num_subsets)
Definition: subsets.cc:10
int sweep_eager_limit
Definition: lbs_structs.h:127
bool verbose_inner_iterations
Definition: lbs_structs.h:143
bool save_angular_flux
Definition: lbs_structs.h:141
chi_mesh::sweep_management::AngleSetGroup TAngleSetGroup