|
NOX
Development
|
Generic object that provides constraints through model evaluator responses. More...
#include <LOCA_Tpetra_ConstraintModelEvaluator.hpp>
Public Member Functions | |
| ConstraintModelEvaluator (const Teuchos::RCP<::Thyra::ModelEvaluator< double >> &model, const LOCA::ParameterVector &pVec, const std::vector< std::string > &constraintResponseNames, const NOX::Abstract::Vector &cloneVec) | |
| ConstraintModelEvaluator (const LOCA::MultiContinuation::ConstraintModelEvaluator &cme, NOX::CopyType type=NOX::DeepCopy) | |
| void | copy (const ConstraintInterface &source) |
| Copy. | |
|
Teuchos::RCP < LOCA::MultiContinuation::ConstraintInterface > | clone (NOX::CopyType type=NOX::DeepCopy) const |
| Cloning function. | |
| int | numConstraints () const |
| Return number of constraints. | |
| void | setX (const NOX::Abstract::Vector &x) |
| Set the solution vector to x. | |
| void | setParam (int paramID, double val) |
| Set parameter value given a parameter indices corresponding to the LOCA::ParameterVector. | |
| void | setParams (const std::vector< int > ¶mIDs, const NOX::Abstract::MultiVector::DenseMatrix &vals) |
| Set parameter value given a parameter indices corresponding to the LOCA::ParameterVector. | |
| NOX::Abstract::Group::ReturnType | computeConstraints () |
| Compute constraint residuals. | |
| NOX::Abstract::Group::ReturnType | computeDX () |
| Compute derivative of constraints w.r.t. solution vector x. | |
| NOX::Abstract::Group::ReturnType | computeDP (const std::vector< int > ¶mIDs, NOX::Abstract::MultiVector::DenseMatrix &dgdp, bool isValidG) |
| Compute derivative of constraints w.r.t. supplied parameters. More... | |
| bool | isConstraints () const |
Return true if constraint residuals are valid. | |
| bool | isDX () const |
Return true if derivative of constraint w.r.t. x is valid. | |
|
const NOX::Abstract::MultiVector::DenseMatrix & | getConstraints () const |
| Return constraint residuals. | |
| bool | isDXZero () const |
Return true if solution component of constraint derivatives is zero. | |
| NOX::Abstract::MultiVector * | getDX () const |
| Return solution component of constraint derivatives. More... | |
| const LOCA::ParameterVector | getParams () const |
Public Member Functions inherited from LOCA::MultiContinuation::ConstraintInterfaceMVDX | |
| ConstraintInterfaceMVDX () | |
| Constructor. | |
| virtual | ~ConstraintInterfaceMVDX () |
| Destructor. | |
| virtual NOX::Abstract::Group::ReturnType | multiplyDX (double alpha, const NOX::Abstract::MultiVector &input_x, NOX::Abstract::MultiVector::DenseMatrix &result_p) const |
| Compute result_p = alpha * dg/dx * input_x. More... | |
| virtual NOX::Abstract::Group::ReturnType | addDX (Teuchos::ETransp transb, double alpha, const NOX::Abstract::MultiVector::DenseMatrix &b, double beta, NOX::Abstract::MultiVector &result_x) const |
| Compute result_x = alpha * dg/dx^T * op(b) + beta * result_x. More... | |
Public Member Functions inherited from LOCA::MultiContinuation::ConstraintInterface | |
| ConstraintInterface () | |
| Constructor. | |
| virtual | ~ConstraintInterface () |
| Destructor. | |
| virtual void | preProcessContinuationStep (LOCA::Abstract::Iterator::StepStatus) |
| Perform any preprocessing before a continuation step starts. More... | |
| virtual void | postProcessContinuationStep (LOCA::Abstract::Iterator::StepStatus) |
| Perform any postprocessing after a continuation step finishes. More... | |
Generic object that provides constraints through model evaluator responses.
| LOCA::MultiContinuation::ConstraintModelEvaluator::ConstraintModelEvaluator | ( | const Teuchos::RCP<::Thyra::ModelEvaluator< double >> & | model, |
| const LOCA::ParameterVector & | pVec, | ||
| const std::vector< std::string > & | constraintResponseNames, | ||
| const NOX::Abstract::Vector & | cloneVec | ||
| ) |
Constructor
| model | Model evaluator that provides constraints as responses. |
| pVec | The independent parameters for constraints. |
| constrantResponseNames | The names of the responses used as constraint equations. |
| cloneVec | NOX vector used to clone data structures with same space/map. |
References NOX::Abstract::Vector::clone(), NOX::Abstract::Vector::createMultiVector(), LOCA::ParameterVector::getNamesVector(), is_null(), LOCA::ParameterVector::length(), NOX::ShapeCopy, TEUCHOS_ASSERT, and TEUCHOS_TEST_FOR_EXCEPTION.
|
virtual |
Compute derivative of constraints w.r.t. supplied parameters.
The first column of dgdp should be filled with the constraint residuals
if isValidG is false. If isValidG is true, then the dgdp contains
on input.
Implements LOCA::MultiContinuation::ConstraintInterface.
References computeConstraints(), isConstraints(), LOCA::ParameterVector::length(), Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >::numCols(), Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >::numRows(), NOX::Abstract::Group::Ok, TEUCHOS_ASSERT, and Teuchos::VERB_EXTREME.
|
virtual |
Return solution component of constraint derivatives.
May return NULL if constraint derivative is zero
Implements LOCA::MultiContinuation::ConstraintInterfaceMVDX.
References Teuchos::RCP< T >::get().
1.8.5