44 #ifndef OPTIPACK_NONLINEAR_CG_DECL_HPP
45 #define OPTIPACK_NONLINEAR_CG_DECL_HPP
48 #include "OptiPack_Types.hpp"
49 #include "Thyra_ModelEvaluator.hpp"
50 #include "GlobiPack_LineSearchBase.hpp"
51 #include "Teuchos_Describable.hpp"
52 #include "Teuchos_VerboseObject.hpp"
53 #include "Teuchos_ParameterListAcceptorDefaultBase.hpp"
54 #include "Teuchos_ParameterEntryValidator.hpp"
60 namespace NonlinearCGUtils {
66 SOLVE_LINSEARCH_FAILURE,
67 SOLVE_MAX_ITERS_EXCEEDED
87 template<
typename Scalar>
107 const int paramIndex,
108 const int responseIndex,
172 NonlinearCGUtils::ESolveReturn
179 const Ptr<int> &numIters = Teuchos::null
194 NonlinearCGUtils::ESolverTypes solverType_;
197 bool and_conv_tests_;
204 mutable int numIters_;
207 solverType_validator_;
216 template<
typename Scalar>
228 template<
typename Scalar>
232 const int paramIndex,
233 const int responseIndex,
239 solver->initialize(model, paramIndex, responseIndex, linesearch);
247 namespace NonlinearCGUtils {
249 const std::string solverType_name =
"Solver Type";
250 const std::string solverType_default =
"FR";
251 const ESolverTypes solverType_default_integral_val = NONLINEAR_CG_FR;
253 const std::string alpha_init_name =
"Initial Linesearch Step Length";
254 const double alpha_init_default = 1.0;
256 const std::string alpha_reinit_name =
"Reinitlaize Linesearch Step Length";
257 const bool alpha_reinit_default =
false;
259 const std::string and_conv_tests_name =
"AND Convergence Tests";
260 const bool and_conv_tests_default =
false;
262 const std::string minIters_name =
"Min Num Iterations";
263 const int minIters_default = 0;
265 const std::string maxIters_name =
"Max Num Iterations";
266 const int maxIters_default = 20;
268 const std::string g_reduct_tol_name =
"Objective Reduction Tol";
269 const double g_reduct_tol_default = 1e-5;
271 const std::string g_grad_tol_name =
"Objective Gradient Tol";
272 const double g_grad_tol_default = 1e-5;
274 const std::string g_mag_name =
"Objective Magnitude";
275 const double g_mag_default = 1.0;
294 #endif // OPTIPACK_NONLINEAR_CG_DECL_HPP
bool get_and_conv_tests() const
const RCP< NonlinearCG< Scalar > > nonlinearCG()
Nonmember constructor.
bool get_alpha_reinit() const
RCP< const ParameterList > getValidParameters() const
Concrete class implementing several nonlinear CG algorithms.
NonlinearCGUtils::ESolverTypes get_solverType() const
ScalarMag get_alpha_init() const
void setParameterList(RCP< ParameterList > const ¶mList)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
NonlinearCG()
Construct with default parameters.
ScalarMag get_g_reduct_tol() const
ScalarMag get_g_grad_tol() const
ScalarMag get_g_mag() const
void initialize(const RCP< const Thyra::ModelEvaluator< Scalar > > &model, const int paramIndex, const int responseIndex, const RCP< GlobiPack::LineSearchBase< Scalar > > &linesearch)
Initialize.
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.
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.
ScalarTraits< Scalar >::magnitudeType ScalarMag