PRISMS-PF Manual
Loading...
Searching...
No Matches
linear_helper.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: © 2025 PRISMS Center at the University of Michigan
2// SPDX-License-Identifier: GNU Lesser General Public Version 2.1
3
4#pragma once
5
6#include <deal.II/lac/precondition.h>
7#include <deal.II/lac/solver_cg.h>
8#include <deal.II/lac/solver_control.h>
9
13#include <prismspf/core/timer.h>
15#include <prismspf/core/types.h>
16
19
20#include <prismspf/config.h>
21
23
32template <typename number>
33class CGSolver : public dealii::SolverCG<BlockVector<number>>
34{
35public:
36 using AddData = dealii::SolverCG<BlockVector<number>>::AdditionalData;
37 using BaseClass = dealii::SolverCG<BlockVector<number>>;
38 using BaseClass::all_condition_numbers_signal;
39 using BaseClass::all_eigenvalues_signal;
40 using BaseClass::condition_number_signal;
41 using BaseClass::determine_beta_by_flexible_formula;
42 using BaseClass::eigenvalues_signal;
43
47 CGSolver(dealii::SolverControl &cn,
48 dealii::VectorMemory<BlockVector<number>> &mem,
49 const AddData &data = AddData())
50 : dealii::SolverCG<BlockVector<number>>(cn, mem, data)
51 {}
52
57 explicit CGSolver(dealii::SolverControl &cn, const AddData &data = AddData())
58 : dealii::SolverCG<BlockVector<number>>(cn, data)
59 {}
60
65 template <typename MatrixType, typename PreconditionerType>
66 void
70 const PreconditionerType &preconditioner)
71 {
72 using namespace dealii;
73 using VectorType = BlockVector<number>;
74
75 SolverControl::State solver_state = SolverControl::iterate;
76
77 LogStream::Prefix prefix("cg");
78
79 // Should we build the matrix for eigenvalue computations?
80 const bool do_eigenvalues =
82 !eigenvalues_signal.empty() || !all_eigenvalues_signal.empty();
83
84 // vectors used for eigenvalue computations
85 std::vector<typename VectorType::value_type> diagonal;
86 std::vector<typename VectorType::value_type> offdiagonal;
87
88 typename VectorType::value_type eigen_beta_alpha = 0;
89
90 int it = 0;
91
92 internal::SolverCG::IterationWorker<VectorType, MatrixType, PreconditionerType>
93 worker(A, preconditioner, determine_beta_by_flexible_formula, this->memory, x);
94
95 worker.startup(b);
96
97 solver_state = this->iteration_status(0, worker.residual_norm, x);
98 if (solver_state != SolverControl::iterate)
99 {
100 return;
101 }
102
103 while (solver_state == SolverControl::iterate)
104 {
105 ++it;
106
107 worker.do_iteration(it);
108
110
111 if (it > 1)
112 {
113 this->coefficients_signal(worker.previous_alpha, worker.beta);
114 // set up the vectors containing the diagonal and the off diagonal
115 // of the projected matrix.
116 if (do_eigenvalues)
117 {
118 diagonal.push_back(number(1.) / worker.previous_alpha + eigen_beta_alpha);
119 eigen_beta_alpha = worker.beta / worker.previous_alpha;
120 offdiagonal.push_back(std::sqrt(worker.beta) / worker.previous_alpha);
121 }
126 }
127
128 solver_state = this->iteration_status(it, worker.residual_norm, x);
129 }
130
131 worker.finalize_after_convergence(it);
132
137
138 AssertThrow(solver_state == SolverControl::success,
139 SolverControl::NoConvergence(it, worker.residual_norm));
140 }
141};
142
143// in linearsolver
144/* dealii::SolverControl::State
145rmse(const BlockVector<number> &residual,
146 number tolerance,
147 const InvMManager<dim, degree, number> &invm,
148 unsigned int relative_level)
149{
150 const auto &field_to_block = solutions.get_global_to_block_index();
151 number sum = 0.0;
152 for (unsigned int field_index : solve_group.field_indices)
153 {
154 const unsigned int block_index = field_to_block[field_index];
155 TensorRank rank = solve_context->get_field_attributes()[field_index].field_type;
156
157 SolutionVector<number> vec = residual.block(block_index); // residual
158 vec.scale(vec); // l2
159
160 vec.scale(invm.get_jxw(rank, relative_level));
161 sum += vec.mean_value();
162 }
163 using std::sqrt;
164 if (sqrt(sum) < tolerance)
165 {
166 return dealii::SolverControl::success;
167 }
168 return dealii::SolverControl::iterate;
169} */
170
171// before solving
172/* cg_solver.connect(
173 [&](const unsigned int, const double, const BlockVector<number> &)
174 {
175 return rmse(*(cg_solver.r_ptr),
176 params.tolerance,
177 solve_context->get_invm_manager(),
178 relative_level);
179 }); */
This class exists as a hack to get access to the residual vector for the custom convergence criterion...
Definition linear_helper.h:34
void solve_better(const MatrixType &A, BlockVector< number > &x, const BlockVector< number > &b, const PreconditionerType &preconditioner)
Solve method we are replacing in prisms-pf because we need access to the residual.
Definition linear_helper.h:67
dealii::SolverCG< BlockVector< number > > BaseClass
Definition linear_helper.h:37
dealii::SolverCG< BlockVector< number > >::AdditionalData AddData
Definition linear_helper.h:36
CGSolver(dealii::SolverControl &cn, const AddData &data=AddData())
Constructor. Use an object of type GrowingVectorMemory as a default to allocate memory.
Definition linear_helper.h:57
CGSolver(dealii::SolverControl &cn, dealii::VectorMemory< BlockVector< number > > &mem, const AddData &data=AddData())
Constructor.
Definition linear_helper.h:47
@ Value
Use value of the variable as a criterion for refinement.
Definition grid_refiner_criterion.h:31
dealii::LinearAlgebra::distributed::BlockVector< number > BlockVector
Typedef for solution block vector.
Definition group_solution_handler.h:29
Definition conditional_ostreams.cc:20
PreconditionerType
Preconditioner type.
Definition type_enums.h:100