20#include <dune/common/exceptions.hh>
21#include <dune/common/float_cmp.hh>
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.");
51 Eigen::VectorX<double>
DD;
60 DD.resizeLike(firstParameter);
64 dfdDD.resizeLike(firstParameter);
102 template <
typename NLS>
105 const auto root = sqrt(args.
DD.squaredNorm() + idbcDelta.squaredNorm() + psi * psi * args.
Dlambda * args.
Dlambda);
107 if (not computedInitialPredictor) {
108 args.
dfdDD.setZero();
112 const auto Dhat0 = idbcDelta / args.
Dlambda;
132 template <
typename NLS>
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);
146 R += idbcForceFunction(req);
148 linearSolver.analyzePattern(K);
149 linearSolver.factorize(K);
150 linearSolver.solve(args.
DD, -R);
152 const auto DD2 = args.
DD.squaredNorm();
153 const auto Dhat02 = Dhat0.squaredNorm();
155 psi = sqrt(DD2 + Dhat02);
156 auto s = sqrt(DD2 + Dhat02 + psi * psi * dlambda * dlambda);
161 idbcDelta *= args.
Dlambda / dlambda;
163 updateFunction(req.globalSolution(), args.
DD);
164 req.parameter() = args.
Dlambda;
168 req.syncParameterAndGlobalSolution(updateFunction);
171 psi = sqrt((idbcDelta.squaredNorm() + args.
DD.squaredNorm()) / (args.
Dlambda * args.
Dlambda));
173 computedInitialPredictor =
true;
187 template <
typename NLS>
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());
198 constexpr std::string
name()
const {
return "Arc length"; }
202 bool computedInitialPredictor{
false};
230 template <
typename NLS>
232 Impl::checkHasNoValidIDBCForceFunction<NLS>();
234 args.
dfdDD.setZero();
249 template <
typename NLS>
251 Impl::checkHasNoValidIDBCForceFunction<NLS>();
253 req.parameter() = args.
Dlambda;
254 computedInitialPredictor =
true;
268 template <
typename NLS>
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;
277 constexpr std::string
name()
const {
return "Load Control"; }
280 bool computedInitialPredictor{
false};
303 : controlledIndices{std::move(p_controlledIndices)} {}
316 template <
typename NLS>
318 Impl::checkHasNoValidIDBCForceFunction<NLS>();
319 const auto controlledDOFsNorm = args.
DD(controlledIndices).norm();
320 args.
f = controlledDOFsNorm - args.
stepSize;
322 args.
dfdDD.setZero();
323 args.
dfdDD(controlledIndices) = args.
DD(controlledIndices) / controlledDOFsNorm;
337 template <
typename NLS>
339 Impl::checkHasNoValidIDBCForceFunction<NLS>();
340 args.
DD(controlledIndices).array() = args.
stepSize;
341 req.globalSolution() = args.
DD;
342 computedInitialPredictor =
true;
356 template <
typename NLS>
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;
365 constexpr std::string
name()
const {
return "Displacement Control"; }
368 std::vector<int> controlledIndices;
369 bool computedInitialPredictor{
false};
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