45 #include "Teuchos_UnitTestHarness.hpp"
47 #include "Thyra_TestingTools.hpp"
49 #include "OptiPack_NonlinearCG.hpp"
50 #include "GlobiPack_ArmijoPolyInterpLineSearch.hpp"
51 #include "GlobiPack_BrentsLineSearch.hpp"
52 #include "Thyra_DiagonalQuadraticResponseOnlyModelEvaluator.hpp"
53 #include "Thyra_DefaultSpmdVectorSpace.hpp"
54 #include "Thyra_ModelEvaluatorHelpers.hpp"
55 #include "Thyra_VectorStdOps.hpp"
56 #include "Thyra_SpmdVectorSpaceBase.hpp"
57 #include "RTOpPack_RTOpTHelpers.hpp"
58 #include "Teuchos_DefaultComm.hpp"
73 using Teuchos::rcpFromRef;
74 using Teuchos::rcp_dynamic_cast;
76 using Teuchos::outArg;
78 using Teuchos::parameterList;
81 using Thyra::createMember;
86 using Thyra::VectorSpaceBase;
87 using Thyra::VectorBase;
89 using Thyra::DiagonalQuadraticResponseOnlyModelEvaluator;
90 using Thyra::diagonalQuadraticResponseOnlyModelEvaluator;
92 using GlobiPack::armijoQuadraticLineSearch;
94 using GlobiPack::brentsLineSearch;
96 using OptiPack::nonlinearCG;
97 namespace NCGU = OptiPack::NonlinearCGUtils;
100 std::string g_solverType =
"FR";
102 int g_globalDim = 16;
104 double g_solve_tol_scale = 10.0;
106 double g_error_tol_scale = 1000.0;
108 int g_nonlin_max_iters = 20;
110 double g_nonlin_term_factor = 1e-2;
112 double g_nonlin_solve_tol = 1e-4;
114 double g_nonlin_error_tol = 1e-3;
120 "solver-type", &g_solverType,
121 "Type type of nonlinear solver. Just pass in blah to see valid options." );
123 "global-dim", &g_globalDim,
124 "Number of global vector elements over all processes" );
126 "solve-tol-scale", &g_solve_tol_scale,
127 "Floating point tolerance for nonlinear CG solve for linear CG tests" );
129 "error-tol-scale", &g_error_tol_scale,
130 "Floating point tolerance for error checks for linear CG tests" );
132 "nonlin-max-iters", &g_nonlin_max_iters,
133 "Max nubmer of CG iterations for general nonlinear problem" );
135 "nonlin-term-factor", &g_nonlin_term_factor,
136 "Scale factor for cubic term in objective for general nonlinear problem" );
138 "nonlin-solve-tol", &g_nonlin_solve_tol,
139 "Floating point tolerance for general nonlinear CG solve" );
141 "nonlin-error-tol", &g_nonlin_error_tol,
142 "Floating point tolerance for error checks for general nonlinear CG solve" );
147 template<
class Scalar>
148 const RCP<DiagonalQuadraticResponseOnlyModelEvaluator<Scalar> >
155 const RCP<const Teuchos::Comm<Thyra::Ordinal> > comm =
158 const int numProcs = comm->getSize();
160 "Error, the number of processors can not be greater than the global"
161 " dimension of the vectors!." );
162 const int localDim = globalDim / numProcs;
163 const int localDimRemainder = globalDim % numProcs;
165 "Error, the number of processors must divide into the global number"
166 " of elements exactly for now!." );
168 const RCP<DiagonalQuadraticResponseOnlyModelEvaluator<Scalar> > model =
169 diagonalQuadraticResponseOnlyModelEvaluator<Scalar>(localDim);
170 const RCP<const VectorSpaceBase<Scalar> > p_space = model->get_p_space(0);
171 const RCP<VectorBase<Scalar> > ps = createMember(p_space);
172 const Scalar ps_val = 2.0;
173 V_S(ps.ptr(), ps_val);
174 model->setSolutionVector(ps);
175 model->setScalarOffset(g_offset);
182 template<
class Scalar>
183 const RCP<NonlinearCG<Scalar> >
184 createNonlinearCGSolver(
186 const RCP<FancyOStream> &out
192 const RCP<ArmijoPolyInterpLineSearch<Scalar> > linesearch =
193 armijoQuadraticLineSearch<Scalar>();
194 const RCP<ParameterList> lsPL = parameterList();
195 lsPL->set(
"Armijo Slope Fraction", 0.0);
196 lsPL->set(
"Min Backtrack Fraction", 0.0);
197 lsPL->set(
"Max Backtrack Fraction", 1e+50);
198 lsPL->set(
"Min Num Iterations", 1);
199 lsPL->set(
"Max Num Iterations", 2);
200 linesearch->setParameterList(lsPL);
202 const RCP<NonlinearCG<Scalar> > cgSolver =
203 nonlinearCG<Scalar>(model, 0, 0, linesearch);
205 const RCP<ParameterList> pl = parameterList();
206 pl->set(
"Solver Type", g_solverType);
207 cgSolver->setParameterList(pl);
209 cgSolver->setOStream(out);
220 template<
class Scalar>
232 const ArrayView<
const ConstSubVectorView<Scalar> > &sub_vecs,
233 const ArrayView<
const SubVectorView<Scalar> > &targ_sub_vecs,
234 const Ptr<ReductTarget> &reduct_obj
239 const SubVectorView<Scalar> &z = targ_sub_vecs[0];
240 const Ordinal z_global_offset = z.globalOffset();
241 const Ordinal z_sub_dim = z.subDim();
242 iter_t z_val = z.values().begin();
243 const ptrdiff_t z_val_s = z.stride();
245 for (
int i = 0; i < z_sub_dim; ++i, z_val += z_val_s ) {
246 const Ordinal global_i = z_global_offset + i;
247 if (range_.in_range(global_i)) {
248 *z_val = as<Scalar>(global_i + 1);
269 typedef ScalarTraits<Scalar> ST;
270 typedef typename ST::magnitudeType ScalarMag;
272 namespace NCGU = OptiPack::NonlinearCGUtils;
273 const RCP<NonlinearCG<Scalar> > cgSolver = nonlinearCG<Scalar>();
292 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( NonlinearCG, parseParamsDefaultParams, Scalar )
294 typedef ScalarTraits<Scalar> ST;
295 typedef typename ST::magnitudeType ScalarMag;
297 namespace NCGU = OptiPack::NonlinearCGUtils;
298 const RCP<NonlinearCG<Scalar> > cgSolver = nonlinearCG<Scalar>();
299 const RCP<ParameterList> pl = parameterList();
300 cgSolver->setParameterList(pl);
318 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( NonlinearCG, parseParams, Scalar )
320 typedef ScalarTraits<Scalar> ST;
321 typedef typename ST::magnitudeType ScalarMag;
323 namespace NCGU = OptiPack::NonlinearCGUtils;
324 const RCP<NonlinearCG<Scalar> > cgSolver = nonlinearCG<Scalar>();
326 const std::string solverTypeStrVal =
"PR+";
327 const double alpha_init = 0.9;
328 const bool alpha_reinit =
true;
329 const int minIters = 92;
330 const int maxIters = 99;
331 const double g_reduct_tol = 2.5;
332 const double g_grad_tol = 2.8;
333 const double g_mag = 3.1;
336 const RCP<ParameterList> pl = parameterList();
337 pl->set(
"Solver Type", solverTypeStrVal);
338 pl->set(
"Initial Linesearch Step Length", alpha_init);
339 pl->set(
"Reinitlaize Linesearch Step Length", alpha_reinit);
340 pl->set(
"Min Num Iterations", minIters);
341 pl->set(
"Max Num Iterations", maxIters);
342 pl->set(
"Objective Reduction Tol", g_reduct_tol);
343 pl->set(
"Objective Gradient Tol", g_grad_tol);
344 pl->set(
"Objective Magnitude", g_mag);
345 cgSolver->setParameterList(pl);
346 const ScalarMag tol = SMT::eps();
364 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( NonlinearCG, printValidParams, Scalar )
366 const RCP<NonlinearCG<Scalar> > cgSolver = nonlinearCG<Scalar>();
367 const RCP<const ParameterList> validPL = cgSolver->getValidParameters();
369 out <<
"\nvalidPL:\n";
370 validPL->print(out, PO().indent(2).showTypes(
true).showFlags(
true).showDoc(
true));
380 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( NonlinearCG, oneEigenVal, Scalar )
383 using Teuchos::optInArg;
384 typedef ScalarTraits<Scalar> ST;
385 typedef typename ST::magnitudeType ScalarMag;
388 const ScalarMag g_offset = as<ScalarMag>(5.0);
389 const RCP<DiagonalQuadraticResponseOnlyModelEvaluator<Scalar> > model =
390 createModel<Scalar>(g_globalDim, g_offset);
391 const RCP<const VectorSpaceBase<Scalar> > p_space = model->get_p_space(0);
392 const int dim = p_space->dim();
394 const RCP<NonlinearCG<Scalar> > cgSolver =
395 createNonlinearCGSolver<Scalar>(model, rcpFromRef(out));
397 const RCP<VectorBase<Scalar> > p = createMember(p_space);
398 V_S( p.ptr(), ST::zero() );
400 ScalarMag g_opt = -1.0;
401 const ScalarMag tol = as<Scalar>(g_solve_tol_scale * dim) * ST::eps();
402 const ScalarMag alpha_init = 10.0;
406 cgSolver->doSolve( p.ptr(), outArg(g_opt),
407 optInArg(tol), optInArg(tol), optInArg(alpha_init), outArg(numIters) );
411 const ScalarMag err_tol = as<Scalar>(g_error_tol_scale * dim) * ST::eps();
415 const bool result = Thyra::testRelNormDiffErr<Scalar>(
417 "ps", *model->getSolutionVector(),
419 "2*err_tol", as<ScalarMag>(2.0)*err_tol,
422 if (!result) success =
false;
433 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( NonlinearCG, partialEigenVal, Scalar )
437 using Teuchos::optInArg;
438 typedef ScalarTraits<Scalar> ST;
439 typedef typename ST::magnitudeType ScalarMag;
442 const ScalarMag g_offset = as<ScalarMag>(5.0);
443 const RCP<DiagonalQuadraticResponseOnlyModelEvaluator<Scalar> > model =
444 createModel<Scalar>(g_globalDim, g_offset);
446 const RCP<const VectorSpaceBase<Scalar> > p_space = model->get_p_space(0);
447 const int dim = p_space->dim();
449 const int numUniqueEigenVals = 3;
451 const RCP<VectorBase<Scalar> > diag = createMember(p_space);
452 V_S(diag.ptr(), ST::one());
453 applyOp<Scalar>( TOpAssignValToGlobalIndex<Scalar>(Range1D(0,numUniqueEigenVals-1)),
454 null, tuple(diag.ptr())(),
null );
455 out <<
"diag =\n" << *diag;
456 model->setDiagonalVector(diag);
459 const RCP<NonlinearCG<Scalar> > cgSolver =
460 createNonlinearCGSolver<Scalar>(model, rcpFromRef(out));
462 const RCP<ParameterList> pl = cgSolver->getNonconstParameterList();
463 const int minIters = numUniqueEigenVals;
464 const int maxIters = minIters + 1;
466 pl->set(
"Max Num Iterations", maxIters);
467 cgSolver->setParameterList(pl);
469 const RCP<VectorBase<Scalar> > p = createMember(p_space);
470 V_S( p.ptr(), ST::zero() );
472 ScalarMag g_opt = -1.0;
473 const ScalarMag tol = as<Scalar>(g_solve_tol_scale * dim) * ST::eps();
474 const ScalarMag alpha_init = 10.0;
477 cgSolver->doSolve( p.ptr(), outArg(g_opt),
478 optInArg(tol), optInArg(tol), optInArg(alpha_init), outArg(numIters) );
482 const ScalarMag err_tol = as<Scalar>(g_error_tol_scale * dim) * ST::eps();
486 const bool result = Thyra::testRelNormDiffErr<Scalar>(
488 "ps", *model->getSolutionVector(),
490 "2*err_tol", as<ScalarMag>(2.0)*err_tol,
493 if (!result) success =
false;
504 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( NonlinearCG, fullEigenVal, Scalar )
507 using Teuchos::optInArg;
508 typedef ScalarTraits<Scalar> ST;
509 typedef typename ST::magnitudeType ScalarMag;
512 const ScalarMag g_offset = as<ScalarMag>(5.0);
513 const RCP<DiagonalQuadraticResponseOnlyModelEvaluator<Scalar> > model =
514 createModel<Scalar>(g_globalDim, g_offset);
515 const RCP<const VectorSpaceBase<Scalar> > p_space = model->get_p_space(0);
516 const int dim = p_space->dim();
518 const RCP<VectorBase<Scalar> > diag = createMember(p_space);
519 applyOp<Scalar>( TOpAssignValToGlobalIndex<Scalar>(),
520 null, tuple(diag.ptr())(),
null );
521 out <<
"diag =\n" << *diag;
522 model->setDiagonalVector(diag);
525 const RCP<NonlinearCG<Scalar> > cgSolver =
526 createNonlinearCGSolver<Scalar>(model, rcpFromRef(out));
528 const RCP<ParameterList> pl = cgSolver->getNonconstParameterList();
529 const int minIters = dim;
530 const int maxIters = minIters+2;
532 pl->set(
"Max Num Iterations", maxIters);
533 cgSolver->setParameterList(pl);
535 const RCP<VectorBase<Scalar> > p = createMember(p_space);
536 V_S( p.ptr(), ST::zero() );
538 ScalarMag g_opt = -1.0;
539 const ScalarMag tol = as<Scalar>(g_solve_tol_scale * dim) * ST::eps();
540 const ScalarMag alpha_init = 10.0;
543 cgSolver->doSolve( p.ptr(), outArg(g_opt),
544 optInArg(tol), optInArg(tol), optInArg(alpha_init), outArg(numIters) );
548 const ScalarMag err_tol = as<Scalar>(g_error_tol_scale * dim) * ST::eps();
552 const bool result = Thyra::testRelNormDiffErr<Scalar>(
554 "ps", *model->getSolutionVector(),
556 "2*err_tol", as<ScalarMag>(2.0)*err_tol,
559 if (!result) success =
false;
572 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( NonlinearCG, fullEigenValScalarProd, Scalar )
576 using Teuchos::optInArg;
577 typedef ScalarTraits<Scalar> ST;
578 typedef typename ST::magnitudeType ScalarMag;
581 const ScalarMag g_offset = as<ScalarMag>(5.0);
582 const RCP<DiagonalQuadraticResponseOnlyModelEvaluator<Scalar> > model =
583 createModel<Scalar>(g_globalDim, g_offset);
585 const RCP<const VectorSpaceBase<Scalar> > p_space = model->get_p_space(0);
586 const int dim = p_space->dim();
589 const RCP<VectorBase<Scalar> > diag = createMember(p_space);
590 applyOp<Scalar>( TOpAssignValToGlobalIndex<Scalar>(),
591 null, tuple(diag.ptr())(),
null );
592 out <<
"diag =\n" << *diag;
593 model->setDiagonalVector(diag);
596 const int numUniqueEigenVals = 3;
598 const RCP<VectorBase<Scalar> > diag_bar = createMember(p_space);
599 V_S(diag_bar.ptr(), ST::one());
600 applyOp<Scalar>( TOpAssignValToGlobalIndex<Scalar>(Range1D(0, numUniqueEigenVals-1)),
601 null, tuple(diag_bar.ptr())(),
null );
604 out <<
"diag_bar =\n" << *diag_bar;
605 model->setDiagonalBarVector(diag_bar);
608 const RCP<NonlinearCG<Scalar> > cgSolver =
609 createNonlinearCGSolver<Scalar>(model, rcpFromRef(out));
611 const RCP<ParameterList> pl = cgSolver->getNonconstParameterList();
612 const int minIters = numUniqueEigenVals;
613 const int maxIters = minIters + 2;
614 pl->set(
"Max Num Iterations", maxIters);
615 cgSolver->setParameterList(pl);
617 const RCP<VectorBase<Scalar> > p = createMember(p_space);
618 V_S( p.ptr(), ST::zero() );
620 ScalarMag g_opt = -1.0;
621 const ScalarMag tol = as<Scalar>(g_solve_tol_scale * dim) * ST::eps();
622 const ScalarMag alpha_init = 10.0;
625 cgSolver->doSolve( p.ptr(), outArg(g_opt),
626 optInArg(tol), optInArg(tol), optInArg(alpha_init), outArg(numIters) );
630 const ScalarMag err_tol = as<Scalar>(g_error_tol_scale * dim) * ST::eps();
634 const bool result = Thyra::testRelNormDiffErr<Scalar>(
636 "ps", *model->getSolutionVector(),
638 "2*err_tol", as<ScalarMag>(2.0)*err_tol,
641 if (!result) success =
false;
652 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( NonlinearCG, generalNonlinearProblem, Scalar )
655 using Teuchos::optInArg;
656 typedef ScalarTraits<Scalar> ST;
657 typedef typename ST::magnitudeType ScalarMag;
660 const ScalarMag g_offset = as<ScalarMag>(5.0);
661 const RCP<DiagonalQuadraticResponseOnlyModelEvaluator<Scalar> > model =
662 createModel<Scalar>(g_globalDim, g_offset);
664 const RCP<const VectorSpaceBase<Scalar> > p_space = model->get_p_space(0);
667 const RCP<VectorBase<Scalar> > diag = createMember(p_space);
668 applyOp<Scalar>( TOpAssignValToGlobalIndex<Scalar>(),
669 null, tuple(diag.ptr())(),
null );
670 out <<
"diag =\n" << *diag;
671 model->setDiagonalVector(diag);
674 const ScalarMag nonlinearTermFactor = as<ScalarMag>(g_nonlin_term_factor);
675 model->setNonlinearTermFactor(nonlinearTermFactor);
677 RCP<BrentsLineSearch<Scalar> > linesearch = brentsLineSearch<Scalar>();
679 const RCP<NonlinearCG<Scalar> > cgSolver =
680 nonlinearCG<Scalar>(model, 0, 0, linesearch);
682 cgSolver->setOStream(rcpFromRef(out));
684 const RCP<ParameterList> pl = parameterList();
685 pl->set(
"Solver Type", g_solverType);
686 pl->set(
"Reinitlaize Linesearch Step Length",
false);
687 pl->set(
"Max Num Iterations", g_nonlin_max_iters);
688 cgSolver->setParameterList(pl);
690 const RCP<VectorBase<Scalar> > p = createMember(p_space);
691 V_S( p.ptr(), ST::zero() );
693 ScalarMag g_opt = -1.0;
694 const ScalarMag tol = as<ScalarMag>(g_nonlin_solve_tol);
695 const ScalarMag alpha_init = 5.0;
698 cgSolver->doSolve( p.ptr(), outArg(g_opt),
699 optInArg(tol), optInArg(tol), optInArg(alpha_init), outArg(numIters) );
703 const ScalarMag err_tol = as<ScalarMag>(g_nonlin_error_tol);
706 const bool result = Thyra::testRelNormDiffErr<Scalar>(
708 "ps", *model->getSolutionVector(),
710 "2*err_tol", as<ScalarMag>(2.0)*err_tol,
713 if (!result) success =
false;
718 generalNonlinearProblem )
726 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( NonlinearCG, generalNonlinearProblem_PL, Scalar )
729 using Teuchos::optInArg;
730 typedef ScalarTraits<Scalar> ST;
731 typedef typename ST::magnitudeType ScalarMag;
734 const ScalarMag g_offset = as<ScalarMag>(5.0);
735 const RCP<DiagonalQuadraticResponseOnlyModelEvaluator<Scalar> > model =
736 createModel<Scalar>(g_globalDim, g_offset);
738 const RCP<const VectorSpaceBase<Scalar> > p_space = model->get_p_space(0);
741 const RCP<VectorBase<Scalar> > diag = createMember(p_space);
742 applyOp<Scalar>( TOpAssignValToGlobalIndex<Scalar>(),
743 null, tuple(diag.ptr())(),
null );
744 out <<
"diag =\n" << *diag;
745 model->setDiagonalVector(diag);
748 const ScalarMag nonlinearTermFactor = as<ScalarMag>(g_nonlin_term_factor);
749 model->setNonlinearTermFactor(nonlinearTermFactor);
751 RCP<BrentsLineSearch<Scalar> > linesearch = brentsLineSearch<Scalar>();
753 const RCP<NonlinearCG<Scalar> > cgSolver =
754 nonlinearCG<Scalar>(model, 0, 0, linesearch);
756 cgSolver->setOStream(rcpFromRef(out));
758 const RCP<VectorBase<Scalar> > p = createMember(p_space);
759 V_S( p.ptr(), ST::zero() );
761 const double tol = as<double>(g_nonlin_solve_tol);
762 const double alpha_init = as<double>(5.0);
763 const RCP<ParameterList> pl = parameterList();
764 pl->set(
"Solver Type", g_solverType);
765 pl->set(
"Initial Linesearch Step Length", alpha_init);
766 pl->set(
"Reinitlaize Linesearch Step Length",
false);
767 pl->set(
"Max Num Iterations", g_nonlin_max_iters);
768 pl->set(
"Objective Reduction Tol", tol);
769 pl->set(
"Objective Gradient Tol", tol);
770 cgSolver->setParameterList(pl);
772 ScalarMag g_opt = -1.0;
775 cgSolver->doSolve( p.ptr(), outArg(g_opt),
776 null,
null, null, outArg(numIters) );
780 const ScalarMag err_tol = as<ScalarMag>(g_nonlin_error_tol);
783 const bool result = Thyra::testRelNormDiffErr<Scalar>(
785 "ps", *model->getSolutionVector(),
787 "2*err_tol", as<ScalarMag>(2.0)*err_tol,
790 if (!result) success =
false;
795 generalNonlinearProblem_PL )
const int maxIters_default
TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL(Teuchos_Conditions, NumberConditionSerialization, T)
#define TEST_INEQUALITY(v1, v2)
static CommandLineProcessor & getCLP()
const double g_mag_default
static Teuchos::RCP< const Comm< OrdinalType > > getComm()
Concrete class implementing several nonlinear CG algorithms.
const double g_grad_tol_default
const double g_reduct_tol_default
const double alpha_init_default
#define TEUCHOS_TEST_FOR_EXCEPT_MSG(throw_exception_test, msg)
void setOption(const char option_true[], const char option_false[], bool *option_val, const char documentation[]=NULL)
virtual void apply_op_impl(const ArrayView< const ConstSubVectorView< Scalar > > &sub_vecs, const ArrayView< const SubVectorView< Scalar > > &targ_sub_vecs, const Ptr< ReductTarget > &reduct_obj) const =0
virtual bool coord_invariant_impl() const
const ESolverTypes solverType_default_integral_val
TypeTo as(const TypeFrom &t)
#define TEST_FLOATING_EQUALITY(v1, v2, tol)
basic_FancyOStream< char > FancyOStream
#define TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT_REAL_SCALAR_TYPES(TEST_GROUP, TEST_NAME)
void validate_apply_op(const RTOpT< Scalar > &op, const int allowed_num_sub_vecs, const int allowed_num_targ_sub_vecs, const bool expect_reduct_obj, const ArrayView< const ConstSubVectorView< Scalar > > &sub_vecs, const ArrayView< const SubVectorView< Scalar > > &targ_sub_vecs, const Ptr< const ReductTarget > &reduct_obj)
#define TEST_EQUALITY(v1, v2)
const int minIters_default
const bool alpha_reinit_default
#define TEST_COMPARE(v1, comp, v2)