G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsIterative.hpp
Go to the documentation of this file.
1
16#pragma once
17
19
21
22#include <sstream>
23
24namespace gismo
25{
26
27template <class T>
29 : assembler(assembler_),
30 m_options(defaultOptions())
31{
32 solVector.setZero(assembler.numDofs(),1);
33 fixedDoFs = assembler.allFixedDofs();
34 for (index_t d = 0; d < (index_t)(fixedDoFs.size()); ++d)
35 fixedDoFs[d].setZero();
36 reset();
37}
38
39template <class T>
41 const gsMatrix<T> & initFreeDoFs)
42 : assembler(assembler_),
43 solVector(initFreeDoFs),
44 m_options(defaultOptions())
45{
46 fixedDoFs = assembler.allFixedDofs();
47 assembler.homogenizeFixedDofs(-1);
48 reset();
49}
50
51template <class T>
53 const gsMatrix<T> & initFreeDoFs,
54 const std::vector<gsMatrix<T> > & initFixedDoFs)
55 : assembler(assembler_),
56 solVector(initFreeDoFs),
57 fixedDoFs(initFixedDoFs),
58 m_options(defaultOptions())
59{
60 reset();
61}
62
63template <class T>
65{
66 m_status = working;
67 numIterations = 0;
68 residualNorm = 0.;
69 initResidualNorm = 1.;
70 updateNorm = 0.;
71 initUpdateNorm = 1.;
72}
73
74template <class T>
76{
77 gsOptionList opt;
79 opt.addInt("Solver","Linear solver to use",linear_solver::LU);
81 opt.addInt("MaxIters","Maximum number of iterations per loop",50);
82 opt.addReal("AbsTol","Absolute tolerance for the convergence cretiria",1e-12);
83 opt.addReal("RelTol","Relative tolerance for the stopping criteria",1e-9);
85 opt.addInt("Verbosity","Amount of information printed to the terminal: none, some, all",solver_verbosity::none);
86 opt.addInt("IterType","Type of iteration: update or next/full",iteration_type::update);
87 return opt;
88}
89
90template <class T>
92{
93 while (m_status == working)
94 {
95 if (!compute())
96 {
97 m_status = bad_solution;
98 goto abort;
99 }
100 if (m_options.getInt("Verbosity") == solver_verbosity::all)
101 gsInfo << status() << std::endl;
102 if (residualNorm < m_options.getReal("AbsTol") ||
103 updateNorm < m_options.getReal("AbsTol") ||
104 residualNorm/initResidualNorm < m_options.getReal("RelTol") ||
105 updateNorm/initUpdateNorm < m_options.getReal("RelTol"))
106 m_status = converged;
107 else if (numIterations == m_options.getInt("MaxIters"))
108 m_status = interrupted;
109 }
110
111 abort:;
112 if (m_options.getInt("Verbosity") != solver_verbosity::none)
113 gsInfo << status() << std::endl;
114}
115
116template <class T>
118{
119 // update mode: set Dirichlet BC to zero after the first iteration
120 if (numIterations == 1 && m_options.getInt("IterType") == iteration_type::update)
121 assembler.homogenizeFixedDofs(-1);
122
123 if (!assembler.assemble(solVector,fixedDoFs))
124 return false;
125
126 gsVector<T> solutionVector;
127 if (m_options.getInt("Solver") == linear_solver::LU)
128 {
129#ifdef GISMO_WITH_PARDISO
130 gsSparseSolver<>::PardisoLU solver(assembler.matrix());
131 solutionVector = solver.solve(assembler.rhs());
132#else
133 gsSparseSolver<>::LU solver(assembler.matrix());
134 solutionVector = solver.solve(assembler.rhs());
135#endif
136 }
137 if (m_options.getInt("Solver") == linear_solver::LDLT)
138 {
139#ifdef GISMO_WITH_PARDISO
140 gsSparseSolver<>::PardisoLDLT solver(assembler.matrix());
141 solutionVector = solver.solve(assembler.rhs());
142#else
143 gsSparseSolver<>::SimplicialLDLT solver(assembler.matrix());
144 solutionVector = solver.solve(assembler.rhs());
145#endif
146 }
147 if (m_options.getInt("Solver") == linear_solver::BiCGSTABDiagonal)
148 {
149 gsSparseSolver<>::BiCGSTABDiagonal solver(assembler.matrix());
150 solutionVector = solver.solve(assembler.rhs());
151 }
152 if (m_options.getInt("Solver") == linear_solver::CGDiagonal)
153 {
154 gsSparseSolver<>::CGDiagonal solver(assembler.matrix());
155 solutionVector = solver.solve(assembler.rhs());
156 }
157
158 if (m_options.getInt("IterType") == iteration_type::update)
159 {
160 updateNorm = solutionVector.norm();
161 residualNorm = assembler.rhs().norm();
162 solVector += solutionVector;
163 // update fixed degrees fo freedom at the first iteration only (they are zero afterwards)
164 if (numIterations == 0)
165 for (index_t d = 0; d < (index_t)(fixedDoFs.size()); ++d)
166 fixedDoFs[d] += assembler.fixedDofs(d);
167 }
168 else if (m_options.getInt("IterType") == iteration_type::next)
169 {
170 updateNorm = (solutionVector-solVector).norm();
171 residualNorm = 1.; // residual is not defined
172 solVector = solutionVector;
173 // copy the fixed degrees of freedom
174 if (numIterations == 0)
175 fixedDoFs = assembler.allFixedDofs();
176 }
177
178 if (numIterations == 0)
179 {
180 initUpdateNorm = updateNorm;
181 initResidualNorm = residualNorm;
182 }
183 numIterations++;
184
185 return true;
186}
187
188template <class T>
190{
191 std::string statusString;
192 if (m_status == converged)
193 statusString = "Iterative solver converged after " +
194 util::to_string(numIterations) + " iteration(s).";
195 else if (m_status == interrupted)
196 statusString = "Iterative solver was interrupted after " +
197 util::to_string(numIterations) + " iteration(s).";
198 else if (m_status == bad_solution)
199 statusString = "Iterative solver was interrupted after " +
200 util::to_string(numIterations) + " iteration(s) due to an invalid solution";
201 else if (m_status == working)
202 statusString = "It: " + util::to_string(numIterations) +
203 ", updAbs: " + util::to_string(updateNorm) +
204 ", updRel: " + util::to_string(updateNorm/initUpdateNorm) +
205 ", resAbs: " + util::to_string(residualNorm) +
206 ", resRel: " + util::to_string(residualNorm/initResidualNorm);
207 return statusString;
208}
209
210template <class T>
211void gsIterative<T>::setFixedDofs(const std::vector<gsMatrix<T> > & ddofs)
212{
213 GISMO_ENSURE(ddofs.size() >= fixedDoFs.size(), "Wrong size of the container with fixed DoFs: " + util::to_string(ddofs.size()) +
214 ". Must be at least: " + util::to_string(fixedDoFs.size()));
215
216 for (short_t d = 0; d < (short_t)(fixedDoFs.size()); ++d)
217 {
218 GISMO_ENSURE(fixedDoFs[d].rows() == ddofs[d].rows(),"Wrong number of fixed DoFs for " + util::to_string(d) + "component: " +
219 util::to_string(ddofs[d].rows()) + ". Must be: " + util::to_string(fixedDoFs[d].rows()));
220 fixedDoFs[d] = ddofs[d];
221 }
222}
223
224template <class T>
226{
227 solVecSaved = solVector;
228 ddofsSaved = fixedDoFs;
229}
230
231template <class T>
233{
234 GISMO_ASSERT(solVecSaved.rows() == solVector.rows(),"No state saved!");
235 solVector = solVecSaved;
236 fixedDoFs = ddofsSaved;
237}
238
239} // namespace ends
Extends the gsAssembler class by adding functionality necessary for a general nonlinear solver....
Definition gsBaseAssembler.h:27
bool compute()
computes update or the next solution
Definition gsIterative.hpp:117
virtual void setFixedDofs(const std::vector< gsMatrix< T > > &ddofs)
set all fixed degrees of freedom
Definition gsIterative.hpp:211
gsBaseAssembler< T > & assembler
assembler object that generates the linear system
Definition gsIterative.h:106
std::string status()
return solver status as a string
Definition gsIterative.hpp:189
gsIterative(gsBaseAssembler< T > &assembler_)
constructor without an initial guess. Assumes a zero initial guess.
Definition gsIterative.hpp:28
void recoverState()
recover solver state from saved state
Definition gsIterative.hpp:232
gsMatrix< T > solVector
solution vector
Definition gsIterative.h:108
void saveState()
save solver state
Definition gsIterative.hpp:225
void reset()
reset the solver state
Definition gsIterative.hpp:64
static gsOptionList defaultOptions()
default option list. used for initialization
Definition gsIterative.hpp:75
std::vector< gsMatrix< T > > fixedDoFs
current Dirichlet DoFs that the solution satisfies
Definition gsIterative.h:110
void solve()
solution procedure
Definition gsIterative.hpp:91
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Class which holds a list of parameters/options, and provides easy access to them.
Definition gsOptionList.h:33
void addInt(const std::string &label, const std::string &desc, const index_t &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:201
void addReal(const std::string &label, const std::string &desc, const Real &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:211
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
std::string to_string(const C &value)
Converts value to string, assuming "operator<<" defined on C.
Definition gsUtils.h:56
Base class for assemblers of gsElasticity.
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
#define gsInfo
Definition gsDebug.h:43
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
A class providing an iterative solver for nonlinear problems.
The G+Smo namespace, containing all definitions for the library.
@ interrupted
method successfully converged
Definition gsBaseUtils.h:100
@ working
solver was interrupted after exceeding the limit of iterations
Definition gsBaseUtils.h:101
@ bad_solution
solver working
Definition gsBaseUtils.h:102
@ next
each iteration yields an update
Definition gsBaseUtils.h:94
@ BiCGSTABDiagonal
Conjugate gradient solver with diagonal (a.k.a. Jacobi) preconditioning: iterative(!...
Definition gsBaseUtils.h:72
@ LDLT
LU decomposition: direct, no matrix requirements, robust but a bit slow, Eigen and Pardiso available.
Definition gsBaseUtils.h:70
@ CGDiagonal
Cholesky decomposition pivoting: direct, simmetric positive or negative semidefinite,...
Definition gsBaseUtils.h:71
@ all
only essential output
Definition gsBaseUtils.h:84