version 0.4.1
pathfollowingfunctions.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2024 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
67{
76 void operator()(SubsidiaryArgs& args) const {
77 if (psi) {
78 const auto root = sqrt(args.DD.squaredNorm() + psi.value() * psi.value() * args.Dlambda * args.Dlambda);
79 args.f = root - args.stepSize;
80 args.dfdDD = args.DD / root;
81 args.dfdDlambda = (psi.value() * psi.value() * args.Dlambda) / root;
82 } else {
83 DUNE_THROW(Dune::InvalidStateException, "You have to call initialPrediction first. Otherwise psi is not defined");
84 }
85 }
86
98 template <typename NLO>
99 void initialPrediction(NLO& nonLinearOperator, SubsidiaryArgs& args) {
100 SolverTypeTag solverTag;
101 using JacobianType = std::remove_cvref_t<typename NLO::DerivativeType>;
102 static_assert((traits::isSpecializationTypeAndNonTypes<Eigen::Matrix, JacobianType>::value) or
103 (traits::isSpecializationTypeNonTypeAndType<Eigen::SparseMatrix, JacobianType>::value),
104 "Linear solver not implemented for the chosen derivative type of the non-linear operator");
105
106 if constexpr (traits::isSpecializationTypeAndNonTypes<Eigen::Matrix, JacobianType>::value)
107 solverTag = SolverTypeTag::d_LDLT;
108 else
110
111 nonLinearOperator.lastParameter() = 1.0; // lambda =1.0
112 nonLinearOperator.template update<0>();
113 const auto& R = nonLinearOperator.value();
114 const auto& K = nonLinearOperator.derivative();
115
116 static constexpr bool isLinearSolver =
117 Ikarus::Concepts::LinearSolverCheck<decltype(LinearSolver(solverTag)), typename NLO::DerivativeType,
118 typename NLO::ValueType>;
119 static_assert(isLinearSolver,
120 "Initial predictor step in the standard arc-length method doesn't have a linear solver");
121
122 auto linearSolver = LinearSolver(solverTag); // for the linear predictor step
123 linearSolver.analyzePattern(K);
124 linearSolver.factorize(K);
125 linearSolver.solve(args.DD, -R);
126
127 const auto DD2 = args.DD.squaredNorm();
128
129 psi = sqrt(DD2);
130 auto s = sqrt(psi.value() * psi.value() + DD2);
131
132 args.DD = args.DD * args.stepSize / s;
133 args.Dlambda = args.stepSize / s;
134
135 nonLinearOperator.firstParameter() = args.DD;
136 nonLinearOperator.lastParameter() = args.Dlambda;
137 }
138
148 template <typename NLO>
149 void intermediatePrediction(NLO& nonLinearOperator, SubsidiaryArgs& args) {
150 nonLinearOperator.firstParameter() += args.DD;
151 nonLinearOperator.lastParameter() += args.Dlambda;
152 }
153
155 constexpr auto name() const { return std::string("Arc length"); }
156
157private:
158 std::optional<double> psi;
159};
160
174{
182 void operator()(SubsidiaryArgs& args) const {
183 args.f = args.Dlambda - args.stepSize;
184 args.dfdDD.setZero();
185 args.dfdDlambda = 1.0;
186 }
187
197 template <typename NLO>
198 void initialPrediction(NLO& nonLinearOperator, SubsidiaryArgs& args) {
199 args.Dlambda = args.stepSize;
200 nonLinearOperator.lastParameter() = args.Dlambda;
201 }
202
212 template <typename NLO>
213 void intermediatePrediction(NLO& nonLinearOperator, SubsidiaryArgs& args) {
214 nonLinearOperator.lastParameter() += args.Dlambda;
215 }
216
218 constexpr auto name() const { return std::string("Load Control"); }
219};
220
234{
240 explicit DisplacementControl(std::vector<int> p_controlledIndices)
241 : controlledIndices{std::move(p_controlledIndices)} {}
242
250 void operator()(SubsidiaryArgs& args) const {
251 const auto controlledDOFsNorm = args.DD(controlledIndices).norm();
252 args.f = controlledDOFsNorm - args.stepSize;
253 args.dfdDlambda = 0.0;
254 args.dfdDD.setZero();
255 args.dfdDD(controlledIndices) = args.DD(controlledIndices) / controlledDOFsNorm;
256 }
257
267 template <typename NLO>
268 void initialPrediction(NLO& nonLinearOperator, SubsidiaryArgs& args) {
269 args.DD(controlledIndices).array() = args.stepSize;
270 nonLinearOperator.firstParameter() = args.DD;
271 }
272
282 template <typename NLO>
283 void intermediatePrediction(NLO& nonLinearOperator, SubsidiaryArgs& args) {
284 nonLinearOperator.firstParameter() += args.DD;
285 }
286
288 constexpr auto name() const { return std::string("Displacement Control"); }
289
290private:
291 std::vector<int> controlledIndices;
292};
293} // namespace Ikarus
Several concepts.
Contains stl-like type traits.
Collection of fallback default functions.
Type-erased linear solver with templated scalar type.
void initialPrediction(NLO &nonLinearOperator, SubsidiaryArgs &args)
Performs the initial prediction for the standard arc-length method.
Definition: pathfollowingfunctions.hh:99
Definition: simpleassemblers.hh:22
LinearSolverTemplate< double > LinearSolver
Definition: linearsolver.hh:245
SolverTypeTag
Enumeration representing different solver types.
Definition: linearsolver.hh:28
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
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:67
void operator()(SubsidiaryArgs &args) const
Evaluates the subsidiary function for the standard arc-length method.
Definition: pathfollowingfunctions.hh:76
void intermediatePrediction(NLO &nonLinearOperator, SubsidiaryArgs &args)
Performs intermediate prediction for the standard arc-length method.
Definition: pathfollowingfunctions.hh:149
constexpr auto name() const
The name of the PathFollowing method.
Definition: pathfollowingfunctions.hh:155
Structure representing the subsidiary function for the load control method.
Definition: pathfollowingfunctions.hh:174
void initialPrediction(NLO &nonLinearOperator, SubsidiaryArgs &args)
Performs initial prediction for the load control method.
Definition: pathfollowingfunctions.hh:198
constexpr auto name() const
The name of the PathFollowing method.
Definition: pathfollowingfunctions.hh:218
void intermediatePrediction(NLO &nonLinearOperator, SubsidiaryArgs &args)
Performs intermediate prediction for the load control method.
Definition: pathfollowingfunctions.hh:213
void operator()(SubsidiaryArgs &args) const
Evaluates the subsidiary function for the load control method.
Definition: pathfollowingfunctions.hh:182
Structure representing the subsidiary function for the displacement control method.
Definition: pathfollowingfunctions.hh:234
constexpr auto name() const
The name of the PathFollowing method.
Definition: pathfollowingfunctions.hh:288
void operator()(SubsidiaryArgs &args) const
Evaluates the subsidiary function for the displacement control method.
Definition: pathfollowingfunctions.hh:250
void intermediatePrediction(NLO &nonLinearOperator, SubsidiaryArgs &args)
Performs intermediate prediction for the displacement control method.
Definition: pathfollowingfunctions.hh:283
DisplacementControl(std::vector< int > p_controlledIndices)
Constructor for DisplacementControl.
Definition: pathfollowingfunctions.hh:240
void initialPrediction(NLO &nonLinearOperator, SubsidiaryArgs &args)
Performs initial prediction for the displacement control method.
Definition: pathfollowingfunctions.hh:268
Concept to check if a linear solver implements all the needed functions for given vector and matrix t...
Definition: concepts.hh:210