| [P]arallel [Hi]gh-order [Li]brary for [P]DEs
    Latest
    Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods | 
KKT_Operator to be used with dealii::SolverBase class. More...
#include <kkt_operator.hpp>
| Public Member Functions | |
| KKT_Operator (const ROL::Ptr< ROL::Objective< Real >> objective, const ROL::Ptr< ROL::Constraint< Real >> equal_constraints, const ROL::Ptr< const ROL::Vector< Real >> design_variables, const ROL::Ptr< const ROL::Vector< Real >> lagrange_mult) | |
| Constructor. | |
| unsigned int | size () | 
| Returns the size of the KKT system. | |
| void | vmult (dealiiSolverVectorWrappingROL< Real > &dst, const dealiiSolverVectorWrappingROL< Real > &src) const | 
| Application of KKT matrix on vector src outputted into dst. | |
| void | Tvmult (dealiiSolverVectorWrappingROL< Real > &dst, const dealiiSolverVectorWrappingROL< Real > &src) const | 
| Application of transposed KKT matrix on vector src outputted into dst.  More... | |
| void | print (const dealiiSolverVectorWrappingROL< Real > &vector_format) | 
| Print the KKT system if the program is run with 1 MPI process.  More... | |
| Protected Attributes | |
| const ROL::Ptr< ROL::Objective< Real > > | objective_ | 
| Objective function. | |
| const ROL::Ptr< ROL::Constraint< Real > > | equal_constraints_ | 
| Equality constraints. | |
| const ROL::Ptr< const ROL::Vector< Real > > | design_variables_ | 
| Design variables. | |
| const ROL::Ptr< const ROL::Vector< Real > > | lagrange_mult_ | 
| Lagrange multipliers. | |
| Private Attributes | |
| const ROL::Ptr< ROL::Vector< Real > > | temp_design_variables_size_vector_ | 
| Used to perform the Lagrangian Hessian in two steps. | |
KKT_Operator to be used with dealii::SolverBase class.
Definition at line 11 of file kkt_operator.hpp.
| 
 | inline | 
Print the KKT system if the program is run with 1 MPI process.
If more than 1 MPI process is used, we can't print out the matrix since the information is distributed
Definition at line 114 of file kkt_operator.hpp.
| 
 | inline | 
Application of transposed KKT matrix on vector src outputted into dst.
Same as vmult since KKT matrix is symmetric.
Definition at line 104 of file kkt_operator.hpp.