[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
P4 preconditioner from Biros & Ghattas 2005. More...
#include <kkt_birosghattas_preconditioners.hpp>
Public Member Functions | |
KKT_P4_Preconditioner (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, const ROL::Ptr< ROL::Secant< Real > > secant, const bool use_approximate_preconditioner=false) | |
< Parallel std::cout that only outputs on mpi_rank==0 More... | |
void | vmult (dealiiSolverVectorWrappingROL< Real > &dst, const dealiiSolverVectorWrappingROL< Real > &src) const override |
Application of KKT preconditionner on vector src outputted into dst. | |
![]() | |
BirosGhattasPreconditioner (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, const ROL::Ptr< ROL::Secant< Real > > secant, const bool use_approximate_preconditioner=false) | |
Constructor. | |
virtual void | Tvmult (dealiiSolverVectorWrappingROL< Real > &dst, const dealiiSolverVectorWrappingROL< Real > &src) const |
Application of transposed KKT preconditioner on vector src outputted into dst. More... | |
Additional Inherited Members | |
![]() | |
const ROL::Ptr< ROL::Objective_SimOpt< Real > > | objective_ |
Objective function. | |
const ROL::Ptr< PHiLiP::FlowConstraints< PHILIP_DIM > > | equal_constraints_ |
Equality constraints. | |
const ROL::Ptr< const ROL::Vector_SimOpt< Real > > | design_variables_ |
Design variables. | |
const ROL::Ptr< const ROL::Vector< Real > > | lagrange_mult_ |
Lagrange multipliers. | |
const ROL::Ptr< const ROL::Vector< Real > > | simulation_variables_ |
Simulation design variables. | |
const ROL::Ptr< const ROL::Vector< Real > > | control_variables_ |
Control design variables. | |
const ROL::Ptr< ROL::Secant< Real > > | secant_ |
Secant method used to precondition the reduced Hessian. | |
const bool | use_approximate_preconditioner_ |
const unsigned int | mpi_rank |
MPI rank used to reset the deallog depth. | |
dealii::ConditionalOStream | pcout |
Parallel std::cout that only outputs on mpi_rank==0. | |
P4 preconditioner from Biros & Ghattas 2005.
Second order terms of the Lagrangian Hessian are ignored and exact inverses of the Jacobian (transpose) are used.
Definition at line 221 of file kkt_birosghattas_preconditioners.hpp.
|
inline |
< Parallel std::cout that only outputs on mpi_rank==0
Constructor.
Definition at line 253 of file kkt_birosghattas_preconditioners.hpp.