rsvs3D
0.0.0
Codes for the c++ implementation of the 3D RSVS
|
Provides the infrastructure for calculation of the RSVS equations. More...
#include <Eigen>
#include <fstream>
#include <iostream>
#include <vector>
#include "RSVScalc.hpp"
#include "arraystructures.hpp"
#include "mesh.hpp"
#include "snake.hpp"
#include "triangulate.hpp"
#include "warning.hpp"
Go to the source code of this file.
Functions | |
void | ResizeLagrangianMultiplier (const RSVScalc &calcobj, Eigen::VectorXd &lagMultAct, bool &isLarge, bool &isNan) |
Resizes the lagrangian multiplier LagMultAct based on whether any of its values are nan or too large. More... | |
template<class T > | |
bool | SQPstep (const RSVScalc &calcobj, const Eigen::MatrixXd &dConstrAct, const Eigen::RowVectorXd &dObjAct, const Eigen::VectorXd &constrAct, Eigen::VectorXd &lagMultAct, Eigen::VectorXd &deltaDVAct, bool &isNan, bool &isLarge, bool attemptConstrOnly=true) |
Template for calculation of an SQP step. More... | |
template<class T > | |
bool | SQPsens_sparse (const Eigen::MatrixXd &sensMult, const MatrixXd_sparse &sensInv, Eigen::MatrixXd &sensRes) |
Calculation for the sparse SQP sensitivity. More... | |
template<class T > | |
bool | SQPsens (const Eigen::MatrixXd &sensMult, const Eigen::MatrixXd &sensInv, Eigen::MatrixXd &sensRes) |
template<class T > | |
bool | SQPstep_sparse (const RSVScalc &calcobj, const MatrixXd_sparse &dConstrAct, const Eigen::RowVectorXd &dObjAct, const Eigen::VectorXd &constrAct, Eigen::VectorXd &lagMultAct, Eigen::VectorXd &deltaDVAct, bool &isNan, bool &isLarge, bool attemptConstrOnly) |
Provides the infrastructure for calculation of the RSVS equations.
Definition in file SQPstep.hpp.
void ResizeLagrangianMultiplier | ( | const RSVScalc & | calcobj, |
Eigen::VectorXd & | lagMultAct, | ||
bool & | isLarge, | ||
bool & | isNan | ||
) |
Resizes the lagrangian multiplier LagMultAct based on whether any of its values are nan or too large.
This uses the RSVScalc object to guide the resizing operation if it is needed.
[in] | calcobj | The calculation object. |
[in,out] | lagMultAct | The std::vector of active lagrangian multipliers. |
[out] | isLarge | Returns if lagMultAct is too large. |
[out] | isNan | Returns if lagMultAct has Nan values. |
bool SQPsens_sparse | ( | const Eigen::MatrixXd & | sensMult, |
const MatrixXd_sparse & | sensInv, | ||
Eigen::MatrixXd & | sensRes | ||
) |
Calculation for the sparse SQP sensitivity.
sensMult | The sensitivity multiplier. | |
[in] | sensInv | The Matrix to invert. |
sensRes | The sensitivity result. |
T | The Eigen object type to use. Should take a RSVScalc::HLag as a constructor and support a solve method. |
Definition at line 257 of file SQPstep.hpp.
bool SQPstep | ( | const RSVScalc & | calcobj, |
const Eigen::MatrixXd & | dConstrAct, | ||
const Eigen::RowVectorXd & | dObjAct, | ||
const Eigen::VectorXd & | constrAct, | ||
Eigen::VectorXd & | lagMultAct, | ||
Eigen::VectorXd & | deltaDVAct, | ||
bool & | isNan, | ||
bool & | isLarge, | ||
bool | attemptConstrOnly = true |
||
) |
Template for calculation of an SQP step.
This template cannot be deduced and needs the developer to pass the required solver class when it is called.
Instantiation options: Eigen::HouseholderQR Eigen::ColPivHouseholderQR Eigen::LLT<Eigen::MatrixXd> (*) <- needs a full type to be defined (see below) Eigen::PartialPivLU
For stability info https://eigen.tuxfamily.org/dox/group__TutorialLinearAlgebra.html https://eigen.tuxfamily.org/dox/group__DenseDecompositionBenchmark.html
[in] | calcobj | The calculation object |
[in] | dConstrAct | Active constraint jacobian dh/dx |
[in] | dObjAct | Active objective jacobian dJ/dx |
[in] | constrAct | Active constraint vector |
lagMultAct | The active lagrangian multipliers | |
deltaDVAct | The active SQP step to take | |
[out] | isNan | Indicates if lagMult is nan |
[out] | isLarge | Indicates if lagMult is large |
[in] | attemptConstrOnly | Should the step algorithm attempt using only the constraint to step. |
T | The Eigen object template type to use. A full type will be defined using T<Eigen::MatrixXd>. |
This template cannot be deduced and needs the developer to pass the required solver class when it is called.
Instantiation options: Eigen::HouseholderQR<Eigen::MatrixXd> Eigen::ColPivHouseholderQR<Eigen::MatrixXd> Eigen::LLT<Eigen::MatrixXd> Eigen::PartialPivLU<Eigen::MatrixXd>
For stability info https://eigen.tuxfamily.org/dox/group__TutorialLinearAlgebra.html https://eigen.tuxfamily.org/dox/group__DenseDecompositionBenchmark.html
[in] | calcobj | The calculation object |
[in] | dConstrAct | Active constraint jacobian dh/dx |
[in] | dObjAct | Active objective jacobian dJ/dx |
[in] | constrAct | Active constraint vector |
lagMultAct | The active lagrangian multipliers | |
deltaDVAct | The active SQP step to take | |
[out] | isNan | Indicates if lagMult is nan |
[out] | isLarge | Indicates if lagMult is large |
[in] | attemptConstrOnly | Should the step algorithm attempt using only the constraint to step. |
T | The Eigen object type to use. Should take a RSVScalc::HLag as a constructor and support a solve method. |
Definition at line 142 of file SQPstep.hpp.