version 0.4.4
pathfollowingfunctions.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2025 The Ikarus Developers ikarus@ibb.uni-stuttgart.de
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
13#pragma once
14
15#include <cmath>
16#include <optional>
17#include <string>
18#include <utility>
19#include <vector>
20
21#include <dune/common/exceptions.hh>
22#include <dune/common/float_cmp.hh>
23
24#include <Eigen/Core>
25
30
31namespace Ikarus {
40{
41 double stepSize;
42 Eigen::VectorX<double> DD;
43 double Dlambda;
44 double f;
45 Eigen::VectorX<double> dfdDD;
46 double dfdDlambda;
48
49 void setZero(const Concepts::EigenType auto& firstParameter) {
50 stepSize = 0.0;
51 DD.resizeLike(firstParameter);
52 DD.setZero();
53 Dlambda = 0.0;
54 f = 0.0;
55 dfdDD.resizeLike(firstParameter);
56 dfdDD.setZero();
57 dfdDlambda = 0.0;
58 currentStep = 0;
59 }
60};
61
80{
89 void operator()(SubsidiaryArgs& args) const {
90 const auto root = sqrt(args.DD.squaredNorm() + psi * psi * args.Dlambda * args.Dlambda);
91 args.f = root - args.stepSize;
92 if (not computedInitialPredictor) {
93 args.dfdDD.setZero();
94 args.dfdDlambda = 0.0;
95 } else {
96 args.dfdDD = args.DD / root;
97 args.dfdDlambda = (psi * psi * args.Dlambda) / root;
98 }
99 }
100
113 template <typename F>
115 typename F::Domain::SolutionVectorType>)
116 void initialPrediction(typename F::Domain& req, F& residual, SubsidiaryArgs& args) {
117 SolverTypeTag solverTag;
118 using JacobianType = std::remove_cvref_t<typename F::Traits::template Range<1>>;
119 static_assert((traits::isSpecializationTypeAndNonTypes<Eigen::Matrix, JacobianType>::value) or
120 (traits::isSpecializationTypeNonTypeAndType<Eigen::SparseMatrix, JacobianType>::value),
121 "Linear solver not implemented for the chosen derivative type of the non-linear operator");
122
123 if constexpr (traits::isSpecializationTypeAndNonTypes<Eigen::Matrix, JacobianType>::value)
124 solverTag = SolverTypeTag::d_LDLT;
125 else
127
128 req.parameter() = 1.0;
129 decltype(auto) R = residual(req);
130 decltype(auto) K = derivative(residual)(req);
131
132 auto linearSolver = LinearSolver(solverTag); // for the linear predictor step
133 linearSolver.analyzePattern(K);
134 linearSolver.factorize(K);
135 linearSolver.solve(args.DD, -R);
136
137 const auto DD2 = args.DD.squaredNorm();
138
139 psi = sqrt(DD2);
140 auto s = sqrt(psi * psi + DD2);
141
142 args.DD = args.DD * args.stepSize / s;
143 args.Dlambda = args.stepSize / s;
144
145 req.globalSolution() = args.DD;
146 req.parameter() = args.Dlambda;
147 computedInitialPredictor = true;
148 }
149
160 template <typename F>
161 void intermediatePrediction(typename F::Domain& req, F& residual, SubsidiaryArgs& args) {
162 if (not computedInitialPredictor)
163 DUNE_THROW(Dune::InvalidStateException, "initialPrediction has to be called before intermediatePrediction.");
164 req.globalSolution() += args.DD;
165 req.parameter() += args.Dlambda;
166 }
167
169 constexpr std::string name() const { return "Arc length"; }
170
171private:
172 double psi{0.0};
173 bool computedInitialPredictor{false};
174};
175
189{
197 void operator()(SubsidiaryArgs& args) const {
198 args.f = args.Dlambda - args.stepSize;
199 args.dfdDD.setZero();
200 args.dfdDlambda = 1.0;
201 }
202
213 template <typename F>
214 void initialPrediction(typename F::Domain& req, F& residual, SubsidiaryArgs& args) {
215 args.Dlambda = args.stepSize;
216 req.parameter() = args.Dlambda;
217 computedInitialPredictor = true;
218 }
219
230 template <typename F>
231 void intermediatePrediction(typename F::Domain& req, F& residual, SubsidiaryArgs& args) {
232 if (not computedInitialPredictor)
233 DUNE_THROW(Dune::InvalidStateException, "initialPrediction has to be called before intermediatePrediction.");
234 req.parameter() += args.Dlambda;
235 }
236
238 constexpr std::string name() const { return "Load Control"; }
239
240private:
241 bool computedInitialPredictor{false};
242};
243
257{
263 explicit DisplacementControl(std::vector<int> p_controlledIndices)
264 : controlledIndices{std::move(p_controlledIndices)} {}
265
273 void operator()(SubsidiaryArgs& args) const {
274 const auto controlledDOFsNorm = args.DD(controlledIndices).norm();
275 args.f = controlledDOFsNorm - args.stepSize;
276 args.dfdDlambda = 0.0;
277 args.dfdDD.setZero();
278 args.dfdDD(controlledIndices) = args.DD(controlledIndices) / controlledDOFsNorm;
279 }
280
291 template <typename F>
292 void initialPrediction(typename F::Domain& req, F& residual, SubsidiaryArgs& args) {
293 args.DD(controlledIndices).array() = args.stepSize;
294 req.globalSolution() = args.DD;
295 computedInitialPredictor = true;
296 }
297
308 template <typename F>
309 void intermediatePrediction(typename F::Domain& req, F& residual, SubsidiaryArgs& args) {
310 if (not computedInitialPredictor)
311 DUNE_THROW(Dune::InvalidStateException, "initialPrediction has to be called before intermediatePrediction.");
312 req.globalSolution() += args.DD;
313 }
314
316 constexpr std::string name() const { return "Displacement Control"; }
317
318private:
319 std::vector<int> controlledIndices;
320 bool computedInitialPredictor{false};
321};
322} // namespace Ikarus
Contains stl-like type traits.
Collection of fallback default functions.
Type-erased linear solver with templated scalar type.
void initialPrediction(typename F::Domain &req, F &residual, SubsidiaryArgs &args)
Performs the initial prediction for the standard arc-length method.
Definition: pathfollowingfunctions.hh:116
Definition: assemblermanipulatorbuildingblocks.hh:22
LinearSolverTemplate< double > LinearSolver
Definition: linearsolver.hh:235
SolverTypeTag
Enumeration representing different solver types.
Definition: linearsolver.hh:31
Structure containing arguments for subsidiary functions.
Definition: pathfollowingfunctions.hh:40
double Dlambda
The increment in the load factor.
Definition: pathfollowingfunctions.hh:43
double dfdDlambda
The derivative of the subsidiary function with respect to Dlambda.
Definition: pathfollowingfunctions.hh:46
int currentStep
The current step index in the control routine.
Definition: pathfollowingfunctions.hh:47
void setZero(const Concepts::EigenType auto &firstParameter)
Definition: pathfollowingfunctions.hh:49
double f
The value of the subsidiary function.
Definition: pathfollowingfunctions.hh:44
double stepSize
The step size in the control routine.
Definition: pathfollowingfunctions.hh:41
Eigen::VectorX< double > dfdDD
The derivative of the subsidiary function with respect to DD.
Definition: pathfollowingfunctions.hh:45
Eigen::VectorX< double > DD
The vector representing the solution increment.
Definition: pathfollowingfunctions.hh:42
Structure representing the subsidiary function for the standard arc-length method.
Definition: pathfollowingfunctions.hh:80
void operator()(SubsidiaryArgs &args) const
Evaluates the subsidiary function for the standard arc-length method.
Definition: pathfollowingfunctions.hh:89
constexpr std::string name() const
The name of the PathFollowing method.
Definition: pathfollowingfunctions.hh:169
void intermediatePrediction(typename F::Domain &req, F &residual, SubsidiaryArgs &args)
Performs intermediate prediction for the standard arc-length method.
Definition: pathfollowingfunctions.hh:161
Structure representing the subsidiary function for the load control method.
Definition: pathfollowingfunctions.hh:189
void initialPrediction(typename F::Domain &req, F &residual, SubsidiaryArgs &args)
Performs initial prediction for the load control method.
Definition: pathfollowingfunctions.hh:214
void intermediatePrediction(typename F::Domain &req, F &residual, SubsidiaryArgs &args)
Performs intermediate prediction for the load control method.
Definition: pathfollowingfunctions.hh:231
constexpr std::string name() const
The name of the PathFollowing method.
Definition: pathfollowingfunctions.hh:238
void operator()(SubsidiaryArgs &args) const
Evaluates the subsidiary function for the load control method.
Definition: pathfollowingfunctions.hh:197
Structure representing the subsidiary function for the displacement control method.
Definition: pathfollowingfunctions.hh:257
void intermediatePrediction(typename F::Domain &req, F &residual, SubsidiaryArgs &args)
Performs intermediate prediction for the displacement control method.
Definition: pathfollowingfunctions.hh:309
void initialPrediction(typename F::Domain &req, F &residual, SubsidiaryArgs &args)
Performs initial prediction for the displacement control method.
Definition: pathfollowingfunctions.hh:292
constexpr std::string name() const
The name of the PathFollowing method.
Definition: pathfollowingfunctions.hh:316
void operator()(SubsidiaryArgs &args) const
Evaluates the subsidiary function for the displacement control method.
Definition: pathfollowingfunctions.hh:273
DisplacementControl(std::vector< int > p_controlledIndices)
Constructor for DisplacementControl.
Definition: pathfollowingfunctions.hh:263
Concept to check if a type is derived from Eigen::MatrixBase.
Definition: utils/concepts.hh:49
Concept to check if a linear solver implements all the needed functions for given vector and matrix t...
Definition: utils/concepts.hh:226
Several concepts.