Chi-Tech
pi_keigen_02_exec.cc
Go to the documentation of this file.
1#include "pi_keigen.h"
2
3#include "chi_runtime.h"
4#include "chi_log.h"
5#include "utils/chi_timer.h"
6
8
9#include <iomanip>
10
11namespace lbs
12{
13
14// ##################################################################
15/**Executes the solver.*/
17{
18 using namespace chi_math;
19
20 double F_prev = 1.0;
21 k_eff_ = 1.0;
22 double k_eff_prev = 1.0;
23 double k_eff_change = 1.0;
24
25 //================================================== Start power iterations
26 int nit = 0;
27 bool converged = false;
28 while (nit < max_iters_)
29 {
30 //================================= Set the fission source
31 SetLBSFissionSource(phi_old_local_, /*additive=*/false);
33
34 //================================= This solves the inners for transport
35 primary_ags_solver_->Setup();
36 primary_ags_solver_->Solve();
37
38 //================================= Recompute k-eigenvalue
40 k_eff_ = F_new / F_prev * k_eff_;
41 double reactivity = (k_eff_ - 1.0) / k_eff_;
42
43 //================================= Check convergence, bookkeeping
44 k_eff_change = fabs(k_eff_ - k_eff_prev) / k_eff_;
45 k_eff_prev = k_eff_;
46 F_prev = F_new;
47 nit += 1;
48
49 if (k_eff_change < std::max(k_tolerance_, 1.0e-12)) converged = true;
50
51 //================================= Print iteration summary
53 {
54 std::stringstream k_iter_info;
55 k_iter_info << Chi::program_timer.GetTimeString() << " "
56 << " Iteration " << std::setw(5) << nit << " k_eff "
57 << std::setw(11) << std::setprecision(7) << k_eff_
58 << " k_eff change " << std::setw(12) << k_eff_change
59 << " reactivity " << std::setw(10) << reactivity * 1e5;
60 if (converged) k_iter_info << " CONVERGED\n";
61
62 Chi::log.Log() << k_iter_info.str();
63 }
64
65 if (converged) break;
66 } // for k iterations
67
68 //================================================== Print summary
69 Chi::log.Log() << "\n";
70 Chi::log.Log() << " Final k-eigenvalue : "
71 << std::setprecision(7) << k_eff_;
72 Chi::log.Log() << " Final change : "
73 << std::setprecision(6) << k_eff_change << " (num_TrOps:"
74 << front_wgs_context_->counter_applications_of_inv_op_ << ")"
75 << "\n";
76 Chi::log.Log() << "\n";
77
79 {
82 }
83
85
87 << "LinearBoltzmann::KEigenvalueSolver execution completed\n\n";
88}
89
90} // namespace lbs
static chi::Timer program_timer
Definition: chi_runtime.h:79
static chi::ChiLog & log
Definition: chi_runtime.h:81
LogStream Log(LOG_LVL level=LOG_0)
Definition: chi_log.cc:35
std::string GetTimeString() const
Definition: chi_timer.cc:38
std::vector< double > & PrecursorsNewLocal()
double ComputeFissionProduction(const std::vector< double > &phi)
lbs::Options & Options()
std::shared_ptr< lbs::WGSContext< Mat, Vec, KSP > > front_wgs_context_
Definition: pi_keigen.h:26
std::shared_ptr< AGSLinearSolver< Mat, Vec, KSP > > primary_ags_solver_
Definition: pi_keigen.h:22
void SetLBSFissionSource(const VecDbl &input, bool additive)
void Scale(VecDbl &x, const double &val)
bool verbose_outer_iterations
Definition: lbs_structs.h:145
bool use_precursors
Definition: lbs_structs.h:138