44 #ifndef OPTIPACK_NONLINEAR_CG_DECL_HPP 45 #define OPTIPACK_NONLINEAR_CG_DECL_HPP 49 #include "Thyra_ModelEvaluator.hpp" 51 #include "Teuchos_Describable.hpp" 52 #include "Teuchos_VerboseObject.hpp" 53 #include "Teuchos_ParameterListAcceptorDefaultBase.hpp" 54 #include "Teuchos_ParameterEntryValidator.hpp" 60 namespace NonlinearCGUtils {
87 template<
typename Scalar>
89 :
public Teuchos::Describable,
90 public Teuchos::VerboseObject<NonlinearCG<Scalar> >,
91 public Teuchos::ParameterListAcceptorDefaultBase
96 typedef typename ScalarTraits<Scalar>::magnitudeType
ScalarMag;
107 const int paramIndex,
108 const int responseIndex,
174 const Ptr<Thyra::VectorBase<Scalar> > &p,
175 const Ptr<ScalarMag> &g_opt,
176 const Ptr<const ScalarMag> &g_reduct_tol = Teuchos::null,
177 const Ptr<const ScalarMag> &g_grad_tol = Teuchos::null,
178 const Ptr<const ScalarMag> &alpha_init = Teuchos::null,
179 const Ptr<int> &numIters = Teuchos::null
189 RCP<const Thyra::ModelEvaluator<Scalar> >
model_;
206 static RCP<Teuchos::ParameterEntryValidator>
216 template<
typename Scalar>
217 const RCP<NonlinearCG<Scalar> >
228 template<
typename Scalar>
229 const RCP<NonlinearCG<Scalar> >
232 const int paramIndex,
233 const int responseIndex,
237 const RCP<NonlinearCG<Scalar> > solver =
239 solver->initialize(model, paramIndex, responseIndex, linesearch);
247 namespace NonlinearCGUtils {
294 #endif // OPTIPACK_NONLINEAR_CG_DECL_HPP bool get_and_conv_tests() const
const std::string and_conv_tests_name
const std::string alpha_init_name
const int maxIters_default
ScalarMag get_alpha_init() const
const RCP< NonlinearCG< Scalar > > nonlinearCG()
Nonmember constructor.
const bool and_conv_tests_default
RCP< const Thyra::ModelEvaluator< Scalar > > model_
const std::string maxIters_name
const double g_mag_default
const std::string solverType_name
Concrete class implementing several nonlinear CG algorithms.
bool get_alpha_reinit() const
const double g_grad_tol_default
RCP< GlobiPack::LineSearchBase< Scalar > > linesearch_
const std::string g_reduct_tol_name
const double g_reduct_tol_default
static RCP< Teuchos::ParameterEntryValidator > solverType_validator_
void setParameterList(RCP< ParameterList > const ¶mList)
const double alpha_init_default
RCP< const ParameterList > getValidParameters() const
NonlinearCG()
Construct with default parameters.
const std::string g_mag_name
const ESolverTypes solverType_default_integral_val
ScalarMag get_g_reduct_tol() const
void initialize(const RCP< const Thyra::ModelEvaluator< Scalar > > &model, const int paramIndex, const int responseIndex, const RCP< GlobiPack::LineSearchBase< Scalar > > &linesearch)
Initialize.
NonlinearCGUtils::ESolverTypes get_solverType() const
ScalarMag get_g_grad_tol() const
Fletcher-Reeves Polak-Ribiere Hybrid Method.
const RCP< NonlinearCG< Scalar > > nonlinearCG(const RCP< const Thyra::ModelEvaluator< Scalar > > &model, const int paramIndex, const int responseIndex, const RCP< GlobiPack::LineSearchBase< Scalar > > &linesearch)
Nonmember constructor.
ScalarMag get_g_mag() const
const std::string alpha_reinit_name
const int minIters_default
const std::string solverType_default
const std::string minIters_name
NonlinearCGUtils::ESolveReturn doSolve(const Ptr< Thyra::VectorBase< Scalar > > &p, const Ptr< ScalarMag > &g_opt, const Ptr< const ScalarMag > &g_reduct_tol=Teuchos::null, const Ptr< const ScalarMag > &g_grad_tol=Teuchos::null, const Ptr< const ScalarMag > &alpha_init=Teuchos::null, const Ptr< int > &numIters=Teuchos::null)
Perform a solve.
NonlinearCGUtils::ESolverTypes solverType_
const bool alpha_reinit_default
const std::string g_grad_tol_name
ScalarTraits< Scalar >::magnitudeType ScalarMag