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 <string>
17#include <utility>
18#include <vector>
19
20#include <dune/common/exceptions.hh>
21#include <dune/common/float_cmp.hh>
22
23#include <Eigen/Core>
24
28
29namespace Ikarus {
30
31namespace Impl {
32
33 template <typename NLS>
34 void checkHasNoValidIDBCForceFunction() {
35 static_assert(not Concepts::HasValidIDBCForceFunction<NLS>,
36 "The given path following function is not implemented to handle inhomogeneous Dirichlet BCs.");
37 }
38
39} // namespace Impl
40
49{
50 double stepSize;
51 Eigen::VectorX<double> DD;
52 double Dlambda;
53 double f;
54 Eigen::VectorX<double> dfdDD;
55 double dfdDlambda;
57
58 void setZero(const Concepts::EigenType auto& firstParameter) {
59 stepSize = 0.0;
60 DD.resizeLike(firstParameter);
61 DD.setZero();
62 Dlambda = 0.0;
63 f = 0.0;
64 dfdDD.resizeLike(firstParameter);
65 dfdDD.setZero();
66 dfdDlambda = 0.0;
67 currentStep = 0;
68 }
69};
70
89{
102 template <typename NLS>
103 void operator()(typename NLS::Domain& req, NLS& nonlinearSolver, SubsidiaryArgs& args) const {
104 const auto idbcDelta = idbcIncrement(req, nonlinearSolver, args.Dlambda);
105 const auto root = sqrt(args.DD.squaredNorm() + idbcDelta.squaredNorm() + psi * psi * args.Dlambda * args.Dlambda);
106 args.f = root - args.stepSize;
107 if (not computedInitialPredictor) {
108 args.dfdDD.setZero();
109 args.dfdDlambda = 0.0;
110 } else {
111 args.dfdDD = args.DD / root;
112 const auto Dhat0 = idbcDelta / args.Dlambda; // extracting Dhat0, assuming IDBC: Dhat = lambda * Dhat0
113 args.dfdDlambda = ((Dhat0.dot(Dhat0) + psi * psi) * args.Dlambda) / root;
114 }
115 }
116
132 template <typename NLS>
133 void initialPrediction(typename NLS::Domain& req, NLS& nonlinearSolver, SubsidiaryArgs& args) {
134 auto& residual = nonlinearSolver.residual();
135 double dlambda = 1.0;
136 auto idbcDelta = idbcIncrement(req, nonlinearSolver, dlambda);
137 const auto Dhat0 = idbcDelta / dlambda;
138 const auto idbcForceFunction = nonlinearSolver.idbcForceFunction();
139 const auto updateFunction = nonlinearSolver.updateFunction();
140 req.parameter() = dlambda;
141 decltype(auto) R = residual(req);
142 decltype(auto) K = derivative(residual)(req);
143 auto linearSolver = createSPDLinearSolverFromNonLinearSolver(nonlinearSolver);
144
146 R += idbcForceFunction(req);
147
148 linearSolver.analyzePattern(K);
149 linearSolver.factorize(K);
150 linearSolver.solve(args.DD, -R);
151
152 const auto DD2 = args.DD.squaredNorm();
153 const auto Dhat02 = Dhat0.squaredNorm();
154
155 psi = sqrt(DD2 + Dhat02);
156 auto s = sqrt(DD2 + Dhat02 + psi * psi * dlambda * dlambda);
157
158 args.DD = args.DD * args.stepSize / s;
159 args.Dlambda = args.stepSize / s;
160
161 idbcDelta *= args.Dlambda / dlambda;
162
163 updateFunction(req.globalSolution(), args.DD);
164 req.parameter() = args.Dlambda;
165
166 // modify the globalSolution() considering inhomogeneous Dirichlet BCs
168 req.syncParameterAndGlobalSolution(updateFunction);
169
170 // rescaling of psi
171 psi = sqrt((idbcDelta.squaredNorm() + args.DD.squaredNorm()) / (args.Dlambda * args.Dlambda));
172
173 computedInitialPredictor = true;
174 }
175
187 template <typename NLS>
188 void intermediatePrediction(typename NLS::Domain& req, NLS& nonlinearSolver, SubsidiaryArgs& args) {
189 if (not computedInitialPredictor)
190 DUNE_THROW(Dune::InvalidStateException, "initialPrediction has to be called before intermediatePrediction.");
191 nonlinearSolver.updateFunction()(req.globalSolution(), args.DD);
192 req.parameter() += args.Dlambda;
194 req.syncParameterAndGlobalSolution(nonlinearSolver.updateFunction());
195 }
196
198 constexpr std::string name() const { return "Arc length"; }
199
200private:
201 double psi{0.0};
202 bool computedInitialPredictor{false};
203};
204
218{
230 template <typename NLS>
231 void operator()(typename NLS::Domain& req, NLS& nonlinearSolver, SubsidiaryArgs& args) const {
232 Impl::checkHasNoValidIDBCForceFunction<NLS>();
233 args.f = args.Dlambda - args.stepSize;
234 args.dfdDD.setZero();
235 args.dfdDlambda = 1.0;
236 }
237
249 template <typename NLS>
250 void initialPrediction(typename NLS::Domain& req, NLS& nonlinearSolver, SubsidiaryArgs& args) {
251 Impl::checkHasNoValidIDBCForceFunction<NLS>();
252 args.Dlambda = args.stepSize;
253 req.parameter() = args.Dlambda;
254 computedInitialPredictor = true;
255 }
256
268 template <typename NLS>
269 void intermediatePrediction(typename NLS::Domain& req, NLS& nonlinearSolver, SubsidiaryArgs& args) {
270 Impl::checkHasNoValidIDBCForceFunction<NLS>();
271 if (not computedInitialPredictor)
272 DUNE_THROW(Dune::InvalidStateException, "initialPrediction has to be called before intermediatePrediction.");
273 req.parameter() += args.Dlambda;
274 }
275
277 constexpr std::string name() const { return "Load Control"; }
278
279private:
280 bool computedInitialPredictor{false};
281};
282
296{
302 explicit DisplacementControl(std::vector<int> p_controlledIndices)
303 : controlledIndices{std::move(p_controlledIndices)} {}
304
316 template <typename NLS>
317 void operator()(typename NLS::Domain& req, NLS& nonlinearSolver, SubsidiaryArgs& args) const {
318 Impl::checkHasNoValidIDBCForceFunction<NLS>();
319 const auto controlledDOFsNorm = args.DD(controlledIndices).norm();
320 args.f = controlledDOFsNorm - args.stepSize;
321 args.dfdDlambda = 0.0;
322 args.dfdDD.setZero();
323 args.dfdDD(controlledIndices) = args.DD(controlledIndices) / controlledDOFsNorm;
324 }
325
337 template <typename NLS>
338 void initialPrediction(typename NLS::Domain& req, NLS& nonlinearSolver, SubsidiaryArgs& args) {
339 Impl::checkHasNoValidIDBCForceFunction<NLS>();
340 args.DD(controlledIndices).array() = args.stepSize;
341 req.globalSolution() = args.DD;
342 computedInitialPredictor = true;
343 }
344
356 template <typename NLS>
357 void intermediatePrediction(typename NLS::Domain& req, NLS& nonlinearSolver, SubsidiaryArgs& args) {
358 Impl::checkHasNoValidIDBCForceFunction<NLS>();
359 if (not computedInitialPredictor)
360 DUNE_THROW(Dune::InvalidStateException, "initialPrediction has to be called before intermediatePrediction.");
361 req.globalSolution() += args.DD;
362 }
363
365 constexpr std::string name() const { return "Displacement Control"; }
366
367private:
368 std::vector<int> controlledIndices;
369 bool computedInitialPredictor{false};
370};
371} // namespace Ikarus
Collection of fallback default functions.
Common things from the control routines.
void initialPrediction(typename NLS::Domain &req, NLS &nonlinearSolver, SubsidiaryArgs &args)
Performs the initial prediction for the standard arc-length method.
Definition: pathfollowingfunctions.hh:133
Definition: assemblermanipulatorbuildingblocks.hh:22
auto createSPDLinearSolverFromNonLinearSolver(const NLS &nls)
A helper function that returns a linear solver suitable for symmetric, positive-definite sparse or de...
Definition: common.hh:40
NLS::Domain::SolutionVectorType idbcIncrement(typename NLS::Domain &x, const NLS &nls, double Dlambda)
A helper function to calculate the increment in the solution vector based on inhomogeneous Dirichlet ...
Definition: common.hh:64
Structure containing arguments for subsidiary functions.
Definition: pathfollowingfunctions.hh:49
double Dlambda
The increment in the load factor.
Definition: pathfollowingfunctions.hh:52
double dfdDlambda
The derivative of the subsidiary function with respect to Dlambda.
Definition: pathfollowingfunctions.hh:55
int currentStep
The current step index in the control routine.
Definition: pathfollowingfunctions.hh:56
void setZero(const Concepts::EigenType auto &firstParameter)
Definition: pathfollowingfunctions.hh:58
double f
The value of the subsidiary function.
Definition: pathfollowingfunctions.hh:53
double stepSize
The step size in the control routine.
Definition: pathfollowingfunctions.hh:50
Eigen::VectorX< double > dfdDD
The derivative of the subsidiary function with respect to DD.
Definition: pathfollowingfunctions.hh:54
Eigen::VectorX< double > DD
The vector representing the solution increment.
Definition: pathfollowingfunctions.hh:51
Structure representing the subsidiary function for the standard arc-length method.
Definition: pathfollowingfunctions.hh:89
void operator()(typename NLS::Domain &req, NLS &nonlinearSolver, SubsidiaryArgs &args) const
Evaluates the subsidiary function for the standard arc-length method.
Definition: pathfollowingfunctions.hh:103
constexpr std::string name() const
The name of the PathFollowing method.
Definition: pathfollowingfunctions.hh:198
void intermediatePrediction(typename NLS::Domain &req, NLS &nonlinearSolver, SubsidiaryArgs &args)
Performs intermediate prediction for the standard arc-length method.
Definition: pathfollowingfunctions.hh:188
Structure representing the subsidiary function for the load control method.
Definition: pathfollowingfunctions.hh:218
void operator()(typename NLS::Domain &req, NLS &nonlinearSolver, SubsidiaryArgs &args) const
Evaluates the subsidiary function for the load control method.
Definition: pathfollowingfunctions.hh:231
void initialPrediction(typename NLS::Domain &req, NLS &nonlinearSolver, SubsidiaryArgs &args)
Performs initial prediction for the load control method.
Definition: pathfollowingfunctions.hh:250
void intermediatePrediction(typename NLS::Domain &req, NLS &nonlinearSolver, SubsidiaryArgs &args)
Performs intermediate prediction for the load control method.
Definition: pathfollowingfunctions.hh:269
constexpr std::string name() const
The name of the PathFollowing method.
Definition: pathfollowingfunctions.hh:277
Structure representing the subsidiary function for the displacement control method.
Definition: pathfollowingfunctions.hh:296
void operator()(typename NLS::Domain &req, NLS &nonlinearSolver, SubsidiaryArgs &args) const
Evaluates the subsidiary function for the displacement control method.
Definition: pathfollowingfunctions.hh:317
constexpr std::string name() const
The name of the PathFollowing method.
Definition: pathfollowingfunctions.hh:365
void intermediatePrediction(typename NLS::Domain &req, NLS &nonlinearSolver, SubsidiaryArgs &args)
Performs intermediate prediction for the displacement control method.
Definition: pathfollowingfunctions.hh:357
DisplacementControl(std::vector< int > p_controlledIndices)
Constructor for DisplacementControl.
Definition: pathfollowingfunctions.hh:302
void initialPrediction(typename NLS::Domain &req, NLS &nonlinearSolver, SubsidiaryArgs &args)
Performs initial prediction for the displacement control method.
Definition: pathfollowingfunctions.hh:338
A concept to check if the underlying solver has a valid function to handle inhomogeneous Dirichlet BC...
Definition: common.hh:26
Concept to check if a type is derived from Eigen::MatrixBase.
Definition: utils/concepts.hh:49
Several concepts.