rsvs3D  0.0.0
Codes for the c++ implementation of the 3D RSVS
SQPstep.hpp File Reference

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)
 

Detailed Description

Provides the infrastructure for calculation of the RSVS equations.

Definition in file SQPstep.hpp.

Function Documentation

◆ ResizeLagrangianMultiplier()

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.

Parameters
[in]calcobjThe calculation object.
[in,out]lagMultActThe std::vector of active lagrangian multipliers.
[out]isLargeReturns if lagMultAct is too large.
[out]isNanReturns if lagMultAct has Nan values.

◆ SQPsens_sparse()

template<class T >
bool SQPsens_sparse ( const Eigen::MatrixXd &  sensMult,
const MatrixXd_sparse &  sensInv,
Eigen::MatrixXd &  sensRes 
)

Calculation for the sparse SQP sensitivity.

Parameters
sensMultThe sensitivity multiplier.
[in]sensInvThe Matrix to invert.
sensResThe sensitivity result.
Template Parameters
TThe Eigen object type to use. Should take a RSVScalc::HLag as a constructor and support a solve method.
Returns
true.

Definition at line 257 of file SQPstep.hpp.

◆ SQPstep()

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.

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

Parameters
[in]calcobjThe calculation object
[in]dConstrActActive constraint jacobian dh/dx
[in]dObjActActive objective jacobian dJ/dx
[in]constrActActive constraint vector
lagMultActThe active lagrangian multipliers
deltaDVActThe active SQP step to take
[out]isNanIndicates if lagMult is nan
[out]isLargeIndicates if lagMult is large
[in]attemptConstrOnlyShould the step algorithm attempt using only the constraint to step.
Template Parameters
TThe Eigen object template type to use. A full type will be defined using T<Eigen::MatrixXd>.
Returns
(isLarge || isNan), if true some form of failure was detected.

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

Parameters
[in]calcobjThe calculation object
[in]dConstrActActive constraint jacobian dh/dx
[in]dObjActActive objective jacobian dJ/dx
[in]constrActActive constraint vector
lagMultActThe active lagrangian multipliers
deltaDVActThe active SQP step to take
[out]isNanIndicates if lagMult is nan
[out]isLargeIndicates if lagMult is large
[in]attemptConstrOnlyShould the step algorithm attempt using only the constraint to step.
Template Parameters
TThe Eigen object type to use. Should take a RSVScalc::HLag as a constructor and support a solve method.
Returns
(isLarge || isNan), if true some form of failure was detected.

Definition at line 142 of file SQPstep.hpp.