OptiPack Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
OptiPack_NonlinearCG_def.hpp
Go to the documentation of this file.
1 /*
2 // @HEADER
3 // ***********************************************************************
4 //
5 // OptiPack: Collection of simple Thyra-based Optimization ANAs
6 // Copyright (2009) Sandia Corporation
7 //
8 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
9 // license for use of this work by or on behalf of the U.S. Government.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov)
39 //
40 // ***********************************************************************
41 // @HEADER
42 */
43 
44 #ifndef OPTIPACK_NONLINEAR_CG_DEF_HPP
45 #define OPTIPACK_NONLINEAR_CG_DEF_HPP
46 
47 
50 #include "OptiPack_UnconstrainedOptMeritFunc1D.hpp"
51 #include "Thyra_ModelEvaluatorHelpers.hpp"
52 #include "Thyra_VectorStdOps.hpp"
53 #include "Teuchos_VerboseObjectParameterListHelpers.hpp"
54 #include "Teuchos_StandardParameterEntryValidators.hpp"
55 #include "Teuchos_Tuple.hpp"
56 
57 
58 namespace OptiPack {
59 
60 
61 template<typename Scalar>
62 RCP<Teuchos::ParameterEntryValidator>
63 NonlinearCG<Scalar>::solverType_validator_ = Teuchos::null;
64 
65 
66 // Constructor/Initializers/Accessors
67 
68 
69 template<typename Scalar>
71  : paramIndex_(-1),
72  responseIndex_(-1),
73  solverType_(NonlinearCGUtils::solverType_default_integral_val),
74  alpha_init_(NonlinearCGUtils::alpha_init_default),
75  alpha_reinit_(NonlinearCGUtils::alpha_reinit_default),
76  and_conv_tests_(NonlinearCGUtils::and_conv_tests_default),
77  minIters_(NonlinearCGUtils::minIters_default),
78  maxIters_(NonlinearCGUtils::maxIters_default),
79  g_reduct_tol_(NonlinearCGUtils::g_reduct_tol_default),
80  g_grad_tol_(NonlinearCGUtils::g_grad_tol_default),
81  g_mag_(NonlinearCGUtils::g_mag_default),
82  numIters_(0)
83 {}
84 
85 
86 template<typename Scalar>
88  const RCP<const Thyra::ModelEvaluator<Scalar> > &model,
89  const int paramIndex,
90  const int responseIndex,
91  const RCP<GlobiPack::LineSearchBase<Scalar> > &linesearch
92  )
93 {
94  // ToDo: Validate input objects!
95  model_ = model.assert_not_null();
96  paramIndex_ = paramIndex;
97  responseIndex_ = responseIndex;
98  linesearch_ = linesearch.assert_not_null();
99 }
100 
101 
102 template<typename Scalar>
104 {
105  return solverType_;
106 }
107 
108 
109 template<typename Scalar>
112 {
113  return alpha_init_;
114 }
115 
116 
117 template<typename Scalar>
119 {
120  return alpha_reinit_;
121 }
122 
123 
124 template<typename Scalar>
126 {
127  return and_conv_tests_;
128 }
129 
130 
131 template<typename Scalar>
133 {
134  return minIters_;
135 }
136 
137 
138 template<typename Scalar>
140 {
141  return maxIters_;
142 }
143 
144 
145 template<typename Scalar>
148 {
149  return g_reduct_tol_;
150 }
151 
152 
153 template<typename Scalar>
156 {
157  return g_grad_tol_;
158 }
159 
160 
161 template<typename Scalar>
164 {
165  return g_mag_;
166 }
167 
168 
169 // Overridden from ParameterListAcceptor (simple forwarding functions)
170 
171 
172 template<typename Scalar>
174 {
175  //typedef ScalarTraits<Scalar> ST; // unused
176  typedef ScalarTraits<ScalarMag> SMT;
177  namespace NCGU = NonlinearCGUtils;
178  using Teuchos::getParameter;
179  using Teuchos::getIntegralValue;
180  paramList->validateParametersAndSetDefaults(*this->getValidParameters());
181  solverType_ = getIntegralValue<NCGU::ESolverTypes>(*paramList, NCGU::solverType_name);
182  alpha_init_ = getParameter<double>(*paramList, NCGU::alpha_init_name);
183  alpha_reinit_ = getParameter<bool>(*paramList, NCGU::alpha_reinit_name);
184  and_conv_tests_ = getParameter<bool>(*paramList, NCGU::and_conv_tests_name);
185  minIters_ = getParameter<int>(*paramList, NCGU::minIters_name);
186  maxIters_ = getParameter<int>(*paramList, NCGU::maxIters_name);
187  g_reduct_tol_ = getParameter<double>(*paramList, NCGU::g_reduct_tol_name);
188  g_grad_tol_ = getParameter<double>(*paramList, NCGU::g_grad_tol_name);
189  g_mag_ = getParameter<double>(*paramList, NCGU::g_mag_name);
190  TEUCHOS_ASSERT_INEQUALITY( alpha_init_, >, SMT::zero() );
191  TEUCHOS_ASSERT_INEQUALITY( minIters_, >=, 0 );
192  TEUCHOS_ASSERT_INEQUALITY( minIters_, <, maxIters_ );
193  TEUCHOS_ASSERT_INEQUALITY( g_reduct_tol_, >=, SMT::zero() );
194  TEUCHOS_ASSERT_INEQUALITY( g_grad_tol_, >=, SMT::zero() );
195  TEUCHOS_ASSERT_INEQUALITY( g_mag_, >, SMT::zero() );
196  Teuchos::readVerboseObjectSublist(&*paramList, this);
197  setMyParamList(paramList);
198 }
199 
200 
201 template<typename Scalar>
204 {
205  using Teuchos::tuple;
206  namespace NCGU = NonlinearCGUtils;
207  static RCP<const ParameterList> validPL;
208  if (is_null(validPL)) {
211  solverType_validator_ =
212  Teuchos::stringToIntegralParameterEntryValidator<NCGU::ESolverTypes>(
213  tuple<std::string>(
214  "FR",
215  "PR+",
216  "FR-PR",
217  "HS"
218  ),
219  tuple<std::string>(
220  "Fletcher-Reeves Method",
221  "Polak-Ribiere Method",
222  "Fletcher-Reeves Polak-Ribiere Hybrid Method",
223  "Hestenes-Stiefel Method"
224  ),
225  tuple<NCGU::ESolverTypes>(
230  ),
232  );
234  "Set the type of nonlinear CG solver algorithm to use.",
235  solverType_validator_ );
239  pl->set( NCGU::minIters_name, NCGU::minIters_default );
240  pl->set( NCGU::maxIters_name, NCGU::maxIters_default );
244  Teuchos::setupVerboseObjectSublist(&*pl);
245  validPL = pl;
246  // ToDo: Add documentation for these parameters
247  }
248  return validPL;
249 }
250 
251 
252 // Solve
253 
254 
255 template<typename Scalar>
258  const Ptr<Thyra::VectorBase<Scalar> > &p_inout,
259  const Ptr<ScalarMag> &g_opt_out,
260  const Ptr<const ScalarMag> &g_reduct_tol_in,
261  const Ptr<const ScalarMag> &g_grad_tol_in,
262  const Ptr<const ScalarMag> &alpha_init_in,
263  const Ptr<int> &numIters_out
264  )
265 {
266 
267  typedef ScalarTraits<Scalar> ST;
268  typedef ScalarTraits<ScalarMag> SMT;
269 
270  using Teuchos::null;
271  using Teuchos::as;
272  using Teuchos::tuple;
273  using Teuchos::rcpFromPtr;
274  using Teuchos::optInArg;
275  using Teuchos::inOutArg;
276  using GlobiPack::computeValue;
278  using Thyra::VectorSpaceBase;
279  using Thyra::VectorBase;
280  using Thyra::MultiVectorBase;
281  using Thyra::scalarProd;
282  using Thyra::createMember;
283  using Thyra::createMembers;
284  using Thyra::get_ele;
285  using Thyra::norm;
286  using Thyra::V_StV;
287  using Thyra::Vt_S;
288  using Thyra::eval_g_DgDp;
289  typedef Thyra::Ordinal Ordinal;
290  typedef Thyra::ModelEvaluatorBase MEB;
291  namespace NCGU = NonlinearCGUtils;
292  using std::max;
293 
294  // Validate input
295 
296  g_opt_out.assert_not_null();
297 
298  // Set streams
299 
300  const RCP<Teuchos::FancyOStream> out = this->getOStream();
301  linesearch_->setOStream(out);
302 
303  // Determine what step constants will be computed
304 
305  const bool compute_beta_PR =
306  (
307  solverType_ == NCGU::NONLINEAR_CG_PR_PLUS
308  ||
309  solverType_ == NCGU::NONLINEAR_CG_FR_PR
310  );
311 
312  const bool compute_beta_HS = (solverType_ == NCGU::NONLINEAR_CG_HS);
313 
314  //
315  // A) Set up the storage for the algorithm
316  //
317 
319  pointEvaluator = defaultPolyLineSearchPointEvaluator<Scalar>();
320 
322  meritFunc = unconstrainedOptMeritFunc1D<Scalar>(
323  model_, paramIndex_, responseIndex_ );
324 
326  p_space = model_->get_p_space(paramIndex_),
327  g_space = model_->get_g_space(responseIndex_);
328 
329  // Stoarge for current iteration
331  p_k = rcpFromPtr(p_inout), // Current solution for p
332  p_kp1 = createMember(p_space), // Trial point for p (in line search)
333  g_vec = createMember(g_space), // Vector (size 1) form of objective g(p)
334  g_grad_k = createMember(p_space), // Gradient of g DgDp^T
335  d_k = createMember(p_space), // Search direction
336  g_grad_k_diff_km1 = null; // g_grad_k - g_grad_km1 (if needed)
337 
338  // Storage for previous iteration
340  g_grad_km1 = null, // Will allocate if we need it!
341  d_km1 = null; // Will allocate if we need it!
342  ScalarMag
343  alpha_km1 = SMT::zero(),
344  g_km1 = SMT::zero(),
345  g_grad_km1_inner_g_grad_km1 = SMT::zero();
346 
347  if (compute_beta_PR || compute_beta_HS) {
348  g_grad_km1 = createMember(p_space);
349  g_grad_k_diff_km1 = createMember(p_space);
350  }
351 
352  if (compute_beta_HS) {
353  d_km1 = createMember(p_space);
354  }
355 
356  //
357  // B) Do the nonlinear CG iterations
358  //
359 
360  *out << "\nStarting nonlinear CG iterations ...\n";
361 
362  if (and_conv_tests_) {
363  *out << "\nNOTE: Using AND of convergence tests!\n";
364  }
365  else {
366  *out << "\nNOTE: Using OR of convergence tests!\n";
367  }
368 
369  const Scalar alpha_init =
370  ( !is_null(alpha_init_in) ? *alpha_init_in : alpha_init_ );
371  const Scalar g_reduct_tol =
372  ( !is_null(g_reduct_tol_in) ? *g_reduct_tol_in : g_reduct_tol_ );
373  const Scalar g_grad_tol =
374  ( !is_null(g_grad_tol_in) ? *g_grad_tol_in : g_grad_tol_ );
375 
376  const Ordinal globalDim = p_space->dim();
377 
378  bool foundSolution = false;
379  bool fatalLinesearchFailure = false;
380  bool restart = true;
381  int numConsecutiveLineSearchFailures = 0;
382 
383  int numConsecutiveIters = 0;
384 
385  for (numIters_ = 0; numIters_ < maxIters_; ++numIters_, ++numConsecutiveIters) {
386 
387  Teuchos::OSTab tab(out);
388 
389  *out << "\nNonlinear CG Iteration k = " << numIters_ << "\n";
390 
391  Teuchos::OSTab tab2(out);
392 
393  //
394  // B.1) Evaluate the point (on first iteration)
395  //
396 
397  eval_g_DgDp(
398  *model_, paramIndex_, *p_k, responseIndex_,
399  numIters_ == 0 ? g_vec.ptr() : null, // Only on first iteration
400  MEB::Derivative<Scalar>(g_grad_k, MEB::DERIV_MV_GRADIENT_FORM) );
401 
402  const ScalarMag g_k = get_ele(*g_vec, 0);
403  // Above: If numIters_ > 0, then g_vec was updated in meritFunc->eval(...).
404 
405  //
406  // B.2) Check for convergence
407  //
408 
409  // B.2.a) ||g_k - g_km1|| |g_k + g_mag| <= g_reduct_tol
410 
411  bool g_reduct_converged = false;
412 
413  if (numIters_ > 0) {
414 
415  const ScalarMag g_reduct = g_k - g_km1;
416 
417  *out << "\ng_k - g_km1 = "<<g_reduct<<"\n";
418 
419  const ScalarMag g_reduct_err =
420  SMT::magnitude(g_reduct / SMT::magnitude(g_k + g_mag_));
421 
422  g_reduct_converged = (g_reduct_err <= g_reduct_tol);
423 
424  *out << "\nCheck convergence: |g_k - g_km1| / |g_k + g_mag| = "<<g_reduct_err
425  << (g_reduct_converged ? " <= " : " > ")
426  << "g_reduct_tol = "<<g_reduct_tol<<"\n";
427 
428  }
429 
430  // B.2.b) ||g_grad_k|| g_mag <= g_grad_tol
431 
432  const Scalar g_grad_k_inner_g_grad_k = scalarProd<Scalar>(*g_grad_k, *g_grad_k);
433  const ScalarMag norm_g_grad_k = ST::magnitude(ST::squareroot(g_grad_k_inner_g_grad_k));
434 
435  *out << "\n||g_grad_k|| = "<<norm_g_grad_k << "\n";
436 
437  const ScalarMag g_grad_err = norm_g_grad_k / g_mag_;
438 
439  const bool g_grad_converged = (g_grad_err <= g_grad_tol);
440 
441  *out << "\nCheck convergence: ||g_grad_k|| / g_mag = "<<g_grad_err
442  << (g_grad_converged ? " <= " : " > ")
443  << "g_grad_tol = "<<g_grad_tol<<"\n";
444 
445  // B.2.c) Convergence status
446 
447  bool isConverged = false;
448  if (and_conv_tests_) {
449  isConverged = g_reduct_converged && g_grad_converged;
450  }
451  else {
452  isConverged = g_reduct_converged || g_grad_converged;
453  }
454 
455  if (isConverged) {
456  if (numIters_ < minIters_) {
457  *out << "\nnumIters="<<numIters_<<" < minIters="<<minIters_
458  << ", continuing on!\n";
459  }
460  else {
461  *out << "\nFound solution, existing algorithm!\n";
462  foundSolution = true;
463  }
464  }
465  else {
466  *out << "\nNot converged!\n";
467  }
468 
469  if (foundSolution) {
470  break;
471  }
472 
473  //
474  // B.3) Compute the search direction d_k
475  //
476 
477  if (numConsecutiveIters == globalDim) {
478 
479  *out << "\nThe number of consecutive iterations exceeds the"
480  << " global dimension so restarting!\n";
481 
482  restart = true;
483 
484  }
485 
486  if (restart) {
487 
488  *out << "\nResetting search direction back to steppest descent!\n";
489 
490  // d_k = -g_grad_k
491  V_StV( d_k.ptr(), as<Scalar>(-1.0), *g_grad_k );
492 
493  restart = false;
494 
495  }
496  else {
497 
498  // g_grad_k - g_grad_km1
499  if (!is_null(g_grad_k_diff_km1)) {
500  V_VmV( g_grad_k_diff_km1.ptr(), *g_grad_k, *g_grad_km1 );
501  }
502 
503  // beta_FR = inner(g_grad_k, g_grad_k) / inner(g_grad_km1, g_grad_km1)
504  const Scalar beta_FR =
505  g_grad_k_inner_g_grad_k / g_grad_km1_inner_g_grad_km1;
506  *out << "\nbeta_FR = " << beta_FR << "\n";
507  // NOTE: Computing beta_FR is free so we might as well just do it!
508 
509  // beta_PR = inner(g_grad_k, g_grad_k - g_grad_km1) /
510  // inner(g_grad_km1, g_grad_km1)
511  Scalar beta_PR = ST::zero();
512  if (compute_beta_PR) {
513  beta_PR =
514  inner(*g_grad_k, *g_grad_k_diff_km1) / g_grad_km1_inner_g_grad_km1;
515  *out << "\nbeta_PR = " << beta_PR << "\n";
516  }
517 
518  // beta_HS = inner(g_grad_k, g_grad_k - g_grad_km1) /
519  // inner(g_grad_k - g_grad_km1, d_km1)
520  Scalar beta_HS = ST::zero();
521  if (compute_beta_HS) {
522  beta_HS =
523  inner(*g_grad_k, *g_grad_k_diff_km1) / inner(*g_grad_k_diff_km1, *d_km1);
524  *out << "\nbeta_HS = " << beta_HS << "\n";
525  }
526 
527  Scalar beta_k = ST::zero();
528  switch(solverType_) {
529  case NCGU::NONLINEAR_CG_FR: {
530  beta_k = beta_FR;
531  break;
532  }
534  beta_k = max(beta_PR, ST::zero());
535  break;
536  }
538  // NOTE: This does not seem to be working :-(
539  if (numConsecutiveIters < 2) {
540  beta_k = beta_PR;
541  }
542  else if (beta_PR < -beta_FR)
543  beta_k = -beta_FR;
544  else if (ST::magnitude(beta_PR) <= beta_FR)
545  beta_k = beta_PR;
546  else // beta_PR > beta_FR
547  beta_k = beta_FR;
548  break;
549  }
550  case NCGU::NONLINEAR_CG_HS: {
551  beta_k = beta_HS;
552  break;
553  }
554  default:
556  }
557  *out << "\nbeta_k = " << beta_k << "\n";
558 
559  // d_k = beta_k * d_last + -g_grad_k
560  if (!is_null(d_km1))
561  V_StV( d_k.ptr(), beta_k, *d_km1 );
562  else
563  Vt_S( d_k.ptr(), beta_k );
564  Vp_StV( d_k.ptr(), as<Scalar>(-1.0), *g_grad_k );
565 
566  }
567 
568  //
569  // B.4) Perform the line search
570  //
571 
572  // B.4.a) Compute the initial step length
573 
574  Scalar alpha_k = as<Scalar>(-1.0);
575 
576  if (numIters_ == 0) {
577  alpha_k = alpha_init;
578  }
579  else {
580  if (alpha_reinit_) {
581  alpha_k = alpha_init;
582  }
583  else {
584  alpha_k = alpha_km1;
585  // ToDo: Implement better logic from Nocedal and Wright for selecting
586  // this step length after first iteration!
587  }
588  }
589 
590  // B.4.b) Perform the linesearch (computing updated quantities in process)
591 
592  pointEvaluator->initialize(tuple<RCP<const VectorBase<Scalar> > >(p_k, d_k)());
593 
594  ScalarMag g_grad_k_inner_d_k = ST::zero();
595 
596  // Set up the merit function to only compute the value
597  meritFunc->setEvaluationQuantities(pointEvaluator, p_kp1, g_vec, null);
598 
599  PointEval1D<ScalarMag> point_k(ST::zero(), g_k);
600  if (linesearch_->requiresBaseDeriv()) {
601  g_grad_k_inner_d_k = scalarProd(*g_grad_k, *d_k);
602  point_k.Dphi = g_grad_k_inner_d_k;
603  }
604 
605  ScalarMag g_kp1 = computeValue(*meritFunc, alpha_k);
606  // NOTE: The above call updates p_kp1 and g_vec as well!
607 
608  PointEval1D<ScalarMag> point_kp1(alpha_k, g_kp1);
609 
610  const bool linesearchResult = linesearch_->doLineSearch(
611  *meritFunc, point_k, inOutArg(point_kp1), null );
612 
613  alpha_k = point_kp1.alpha;
614  g_kp1 = point_kp1.phi;
615 
616  if (linesearchResult) {
617  numConsecutiveLineSearchFailures = 0;
618  }
619  else {
620  if (numConsecutiveLineSearchFailures==0) {
621  *out << "\nLine search failure, resetting the search direction!\n";
622  restart = true;
623  }
624  if (numConsecutiveLineSearchFailures==1) {
625  *out << "\nLine search failure on last iteration also, terminating algorithm!\n";
626  fatalLinesearchFailure = true;
627  }
628  ++numConsecutiveLineSearchFailures;
629  }
630 
631  if (fatalLinesearchFailure) {
632  break;
633  }
634 
635  //
636  // B.5) Transition to the next iteration
637  //
638 
639  alpha_km1 = alpha_k;
640  g_km1 = g_k;
641  g_grad_km1_inner_g_grad_km1 = g_grad_k_inner_g_grad_k;
642  std::swap(p_k, p_kp1);
643  if (!is_null(g_grad_km1))
644  std::swap(g_grad_km1, g_grad_k);
645  if (!is_null(d_km1))
646  std::swap(d_k, d_km1);
647 
648 #ifdef TEUCHOS_DEBUG
649  // Make sure we compute these correctly before they are used!
650  V_S(g_grad_k.ptr(), ST::nan());
651  V_S(p_kp1.ptr(), ST::nan());
652 #endif
653 
654  }
655 
656  //
657  // C) Final clean up
658  //
659 
660  // Get the most current value of g(p)
661  *g_opt_out = get_ele(*g_vec, 0);
662 
663  // Make sure that the final value for p has been copied in!
664  V_V( p_inout, *p_k );
665 
666  if (!is_null(numIters_out)) {
667  *numIters_out = numIters_;
668  }
669 
670  if (numIters_ == maxIters_) {
671  *out << "\nMax nonlinear CG iterations exceeded!\n";
672  }
673 
674  if (foundSolution) {
676  }
677  else if(fatalLinesearchFailure) {
679  }
680 
681  // Else, the max number of iterations was exceeded
683 
684 }
685 
686 
687 } // namespace OptiPack
688 
689 
690 #endif // OPTIPACK_NONLINEAR_CG_DEF_HPP
bool is_null(const boost::shared_ptr< T > &p)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
#define TEUCHOS_ASSERT_INEQUALITY(val1, comp, val2)
RCP< const ParameterList > getValidParameters() const
const int maxIters_default
NonlinearCGUtils::ESolverTypes get_solverType() const
void setParameterList(RCP< ParameterList > const &paramList)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
NonlinearCG()
Construct with default parameters.
Ptr< T > ptr() const
void validateParametersAndSetDefaults(ParameterList const &validParamList, int const depth=1000)
const ESolverTypes solverType_default_integral_val
TypeTo as(const TypeFrom &t)
const Ptr< T > & assert_not_null() const
void initialize(const RCP< const Thyra::ModelEvaluator< Scalar > > &model, const int paramIndex, const int responseIndex, const RCP< GlobiPack::LineSearchBase< Scalar > > &linesearch)
Initialize.
const int minIters_default
Fletcher-Reeves Polak-Ribiere Hybrid Method.
Teuchos_Ordinal Ordinal
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
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.