Chi-Tech
wgs_convergence_test.cc
Go to the documentation of this file.
2
3#include "wgs_context.h"
5
6#include "chi_runtime.h"
7#include "chi_log.h"
8
9#include "utils/chi_timer.h"
10
11#include <iomanip>
12
13namespace lbs
14{
15
16//###################################################################
17/**Customized convergence test.*/
18PetscErrorCode GSConvergenceTest(KSP ksp, PetscInt n, PetscReal rnorm,
19 KSPConvergedReason* convergedReason, void*)
20{
21 //======================================== Get data context
23 KSPGetApplicationContext(ksp,&context);
24
25 //======================================== Set rhs norm
26 double residual_scale = 1.0;
27 switch (context->residual_scale_type)
28 {
30 residual_scale = 1.0;
31 break;
33 if (context->rhs_norm > 1.0e-25)
34 residual_scale = 1.0/context->rhs_norm;
35 break;
37 if (context->rhs_preconditioned_norm > 1.0e-25)
38 residual_scale = 1.0/context->rhs_preconditioned_norm;
39 break;
41 if (context->custom_residual_scale > 1.0e-25)
42 residual_scale = 1.0/context->custom_residual_scale;
43 break;
44 }
45
46 //======================================== Compute test criterion
47 double tol;
48 int64_t maxIts;
49 KSPGetTolerances(ksp, nullptr,&tol, nullptr,&maxIts);
50
51 double scaled_residual = rnorm * residual_scale;
52
53 //======================================== Print iteration information
54 std::string offset;
55 if (context->groupset_.apply_wgdsa_ || context->groupset_.apply_tgdsa_)
56 offset = std::string(" ");
57
58 std::stringstream iter_info;
59 iter_info
61 << offset
62 << "WGS groups ["
63 << context->groupset_.groups_.front().id_
64 << "-"
65 << context->groupset_.groups_.back().id_
66 << "]"
67 << " Iteration " << std::setw(5) << n
68 << " Residual " << std::setw(9) << scaled_residual;
69
70 if (scaled_residual < tol)
71 {
72 *convergedReason = KSP_CONVERGED_RTOL;
73 iter_info << " CONVERGED\n";
74 }
75
76 if (context->log_info_) Chi::log.Log() << iter_info.str() << std::endl;
77
78 return KSP_CONVERGED_ITERATING;
79}
80
81}//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< LBSGroup > groups_
Definition: lbs_groupset.h:41
PetscErrorCode GSConvergenceTest(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *convergedReason, void *)
LBSGroupset & groupset_
Definition: wgs_context.h:28