version 0.4.1
pathfollowingfunctions.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2025 The Ikarus Developers mueller@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
23#include <Eigen/Core>
24
29
30namespace Ikarus {
39{
40 double stepSize;
41 Eigen::VectorX<double> DD;
42 double Dlambda;
43 double f;
44 Eigen::VectorX<double> dfdDD;
45 double dfdDlambda;
47
48 void setZero(const Concepts::EigenType auto& firstParameter) {
49 stepSize = 0.0;
50 DD.resizeLike(firstParameter);
51 DD.setZero();
52 Dlambda = 0.0;
53 f = 0.0;
54 dfdDD.resizeLike(firstParameter);
55 dfdDD.setZero();
56 dfdDlambda = 0.0;
57 currentStep = 0;
58 }
59};
60
79{
88 void operator()(SubsidiaryArgs& args) const {
89 if (psi) {
90 const auto root = sqrt(args.DD.squaredNorm() + psi.value() * psi.value() * args.Dlambda * args.Dlambda);
91 args.f = root - args.stepSize;
92 args.dfdDD = args.DD / root;
93 args.dfdDlambda = (psi.value() * psi.value() * args.Dlambda) / root;
94 } else {
95 DUNE_THROW(Dune::InvalidStateException, "You have to call initialPrediction first. Otherwise psi is not defined");
96 }
97 }
98
111 template <typename F>
113 typename F::Domain::SolutionVectorType>)
114 void initialPrediction(typename F::Domain& req, F& residual, SubsidiaryArgs& args) {
115 SolverTypeTag solverTag;
116 using JacobianType = std::remove_cvref_t<typename F::Traits::template Range<1>>;
117 static_assert((traits::isSpecializationTypeAndNonTypes<Eigen::Matrix, JacobianType>::value) or
118 (traits::isSpecializationTypeNonTypeAndType<Eigen::SparseMatrix, JacobianType>::value),
119 "Linear solver not implemented for the chosen derivative type of the non-linear operator");
120
121 if constexpr (traits::isSpecializationTypeAndNonTypes<Eigen::Matrix, JacobianType>::value)
122 solverTag = SolverTypeTag::d_LDLT;
123 else
125
126 req.parameter() = 1.0;
127 decltype(auto) R = residual(req);
128 decltype(auto) K = derivative(residual)(req);
129
130 auto linearSolver = LinearSolver(solverTag); // for the linear predictor step
131 linearSolver.analyzePattern(K);
132 linearSolver.factorize(K);
133 linearSolver.solve(args.DD, -R);
134
135 const auto DD2 = args.DD.squaredNorm();
136
137 psi = sqrt(DD2);
138 auto s = sqrt(psi.value() * psi.value() + DD2);
139
140 args.DD = args.DD * args.stepSize / s;
141 args.Dlambda = args.stepSize / s;
142
143 req.globalSolution() = args.DD;
144 req.parameter() = args.Dlambda;
145 }
146
157 template <typename F>
158 void intermediatePrediction(typename F::Domain& req, F& residual, SubsidiaryArgs& args) {
159 req.globalSolution() += args.DD;
160 req.parameter() += args.Dlambda;
161 }
162
164 constexpr auto name() const { return std::string("Arc length"); }
165
166private:
167 std::optional<double> psi;
168};
169
183{
191 void operator()(SubsidiaryArgs& args) const {
192 args.f = args.Dlambda - args.stepSize;
193 args.dfdDD.setZero();
194 args.dfdDlambda = 1.0;
195 }
196
207 template <typename F>
208 void initialPrediction(typename F::Domain& req, F& residual, SubsidiaryArgs& args) {
209 args.Dlambda = args.stepSize;
210 req.parameter() = args.Dlambda;
211 }
212
223 template <typename F>
224 void intermediatePrediction(typename F::Domain& req, F& residual, SubsidiaryArgs& args) {
225 req.parameter() += args.Dlambda;
226 }
227
229 constexpr auto name() const { return std::string("Load Control"); }
230};
231
245{
251 explicit DisplacementControl(std::vector<int> p_controlledIndices)
252 : controlledIndices{std::move(p_controlledIndices)} {}
253
261 void operator()(SubsidiaryArgs& args) const {
262 const auto controlledDOFsNorm = args.DD(controlledIndices).norm();
263 args.f = controlledDOFsNorm - args.stepSize;
264 args.dfdDlambda = 0.0;
265 args.dfdDD.setZero();
266 args.dfdDD(controlledIndices) = args.DD(controlledIndices) / controlledDOFsNorm;
267 }
268
279 template <typename F>
280 void initialPrediction(typename F::Domain& req, F& residual, SubsidiaryArgs& args) {
281 args.DD(controlledIndices).array() = args.stepSize;
282 req.globalSolution() = args.DD;
283 }
284
295 template <typename F>
296 void intermediatePrediction(typename F::Domain& req, F& residual, SubsidiaryArgs& args) {
297 req.globalSolution() += args.DD;
298 }
299
301 constexpr auto name() const { return std::string("Displacement Control"); }
302
303private:
304 std::vector<int> controlledIndices;
305};
306} // 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:114
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:39
double Dlambda
The increment in the load factor.
Definition: pathfollowingfunctions.hh:42
double dfdDlambda
The derivative of the subsidiary function with respect to Dlambda.
Definition: pathfollowingfunctions.hh:45
int currentStep
The current step index in the control routine.
Definition: pathfollowingfunctions.hh:46
void setZero(const Concepts::EigenType auto &firstParameter)
Definition: pathfollowingfunctions.hh:48
double f
The value of the subsidiary function.
Definition: pathfollowingfunctions.hh:43
double stepSize
The step size in the control routine.
Definition: pathfollowingfunctions.hh:40
Eigen::VectorX< double > dfdDD
The derivative of the subsidiary function with respect to DD.
Definition: pathfollowingfunctions.hh:44
Eigen::VectorX< double > DD
The vector representing the solution increment.
Definition: pathfollowingfunctions.hh:41
Structure representing the subsidiary function for the standard arc-length method.
Definition: pathfollowingfunctions.hh:79
void operator()(SubsidiaryArgs &args) const
Evaluates the subsidiary function for the standard arc-length method.
Definition: pathfollowingfunctions.hh:88
void intermediatePrediction(typename F::Domain &req, F &residual, SubsidiaryArgs &args)
Performs intermediate prediction for the standard arc-length method.
Definition: pathfollowingfunctions.hh:158
constexpr auto name() const
The name of the PathFollowing method.
Definition: pathfollowingfunctions.hh:164
Structure representing the subsidiary function for the load control method.
Definition: pathfollowingfunctions.hh:183
constexpr auto name() const
The name of the PathFollowing method.
Definition: pathfollowingfunctions.hh:229
void initialPrediction(typename F::Domain &req, F &residual, SubsidiaryArgs &args)
Performs initial prediction for the load control method.
Definition: pathfollowingfunctions.hh:208
void intermediatePrediction(typename F::Domain &req, F &residual, SubsidiaryArgs &args)
Performs intermediate prediction for the load control method.
Definition: pathfollowingfunctions.hh:224
void operator()(SubsidiaryArgs &args) const
Evaluates the subsidiary function for the load control method.
Definition: pathfollowingfunctions.hh:191
Structure representing the subsidiary function for the displacement control method.
Definition: pathfollowingfunctions.hh:245
constexpr auto name() const
The name of the PathFollowing method.
Definition: pathfollowingfunctions.hh:301
void intermediatePrediction(typename F::Domain &req, F &residual, SubsidiaryArgs &args)
Performs intermediate prediction for the displacement control method.
Definition: pathfollowingfunctions.hh:296
void initialPrediction(typename F::Domain &req, F &residual, SubsidiaryArgs &args)
Performs initial prediction for the displacement control method.
Definition: pathfollowingfunctions.hh:280
void operator()(SubsidiaryArgs &args) const
Evaluates the subsidiary function for the displacement control method.
Definition: pathfollowingfunctions.hh:261
DisplacementControl(std::vector< int > p_controlledIndices)
Constructor for DisplacementControl.
Definition: pathfollowingfunctions.hh:251
Concept to check if a type is derived from Eigen::MatrixBase.
Definition: utils/concepts.hh:46
Concept to check if a linear solver implements all the needed functions for given vector and matrix t...
Definition: utils/concepts.hh:223
Several concepts.