xref: /petsc/src/tao/complementarity/impls/ssls/ssfls.c (revision 03047865b8d8757cf1cf9cda45785c1537b01dc1)
1aaa7dc30SBarry Smith #include <../src/tao/complementarity/impls/ssls/ssls.h>
2a7e14dcfSSatish Balay 
TaoSetUp_SSFLS(Tao tao)366976f2fSJacob Faibussowitsch static PetscErrorCode TaoSetUp_SSFLS(Tao tao)
4d71ae5a4SJacob Faibussowitsch {
5a7e14dcfSSatish Balay   TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
6a7e14dcfSSatish Balay 
7a7e14dcfSSatish Balay   PetscFunctionBegin;
89566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tao->gradient));
99566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
109566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &ssls->w));
119566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &ssls->ff));
129566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &ssls->dpsi));
139566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &ssls->da));
149566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &ssls->db));
159566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &ssls->t1));
169566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &ssls->t2));
179566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU));
183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19a7e14dcfSSatish Balay }
20a7e14dcfSSatish Balay 
TaoSolve_SSFLS(Tao tao)21d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_SSFLS(Tao tao)
22d71ae5a4SJacob Faibussowitsch {
23a7e14dcfSSatish Balay   TAO_SSLS                    *ssls = (TAO_SSLS *)tao->data;
24a7e14dcfSSatish Balay   PetscReal                    psi, ndpsi, normd, innerd, t = 0;
25a7e14dcfSSatish Balay   PetscReal                    delta, rho;
26e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason;
27a7e14dcfSSatish Balay 
28a7e14dcfSSatish Balay   PetscFunctionBegin;
29a7e14dcfSSatish Balay   /* Assume that Setup has been called!
30a7e14dcfSSatish Balay      Set the structure for the Jacobian and create a linear solver. */
31a7e14dcfSSatish Balay   delta = ssls->delta;
32a7e14dcfSSatish Balay   rho   = ssls->rho;
33a7e14dcfSSatish Balay 
349566063dSJacob Faibussowitsch   PetscCall(TaoComputeVariableBounds(tao));
35a7e14dcfSSatish Balay   /* Project solution inside bounds */
369566063dSJacob Faibussowitsch   PetscCall(VecMedian(tao->XL, tao->solution, tao->XU, tao->solution));
379566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch, Tao_SSLS_FunctionGradient, tao));
389566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetObjectiveRoutine(tao->linesearch, Tao_SSLS_Function, tao));
39a7e14dcfSSatish Balay 
40a7e14dcfSSatish Balay   /* Calculate the function value and fischer function value at the
41a7e14dcfSSatish Balay      current iterate */
429566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchComputeObjectiveAndGradient(tao->linesearch, tao->solution, &psi, ssls->dpsi));
439566063dSJacob Faibussowitsch   PetscCall(VecNorm(ssls->dpsi, NORM_2, &ndpsi));
44a7e14dcfSSatish Balay 
45763847b4SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
46d90ca52eSBarry Smith   while (PETSC_TRUE) {
4763a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(tao, "iter: %" PetscInt_FMT ", merit: %g, ndpsi: %g\n", tao->niter, (double)ssls->merit, (double)ndpsi));
48a7e14dcfSSatish Balay     /* Check the termination criteria */
499566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, ssls->merit, ndpsi, 0.0, tao->ksp_its));
509566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, ssls->merit, ndpsi, 0.0, t));
51dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
52763847b4SAlp Dener     if (tao->reason != TAO_CONTINUE_ITERATING) break;
53e1e80dc8SAlp Dener 
54e1e80dc8SAlp Dener     /* Call general purpose update function */
55dbbe0bcdSBarry Smith     PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
56e6d4cb7fSJason Sarich     tao->niter++;
57a7e14dcfSSatish Balay 
58a7e14dcfSSatish Balay     /* Calculate direction.  (Really negative of newton direction.  Therefore,
59a7e14dcfSSatish Balay        rest of the code uses -d.) */
609566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(tao->ksp, tao->jacobian, tao->jacobian_pre));
619566063dSJacob Faibussowitsch     PetscCall(KSPSolve(tao->ksp, ssls->ff, tao->stepdirection));
629566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(tao->ksp, &tao->ksp_its));
63b0026674SJason Sarich     tao->ksp_tot_its += tao->ksp_its;
64a7e14dcfSSatish Balay 
659566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->stepdirection, ssls->w));
669566063dSJacob Faibussowitsch     PetscCall(VecScale(ssls->w, -1.0));
679566063dSJacob Faibussowitsch     PetscCall(VecBoundGradientProjection(ssls->w, tao->solution, tao->XL, tao->XU, ssls->w));
68a7e14dcfSSatish Balay 
699566063dSJacob Faibussowitsch     PetscCall(VecNorm(ssls->w, NORM_2, &normd));
709566063dSJacob Faibussowitsch     PetscCall(VecDot(ssls->w, ssls->dpsi, &innerd));
71a7e14dcfSSatish Balay 
72a7e14dcfSSatish Balay     /* Make sure that we have a descent direction */
73d90ca52eSBarry Smith     if (innerd >= -delta * PetscPowReal(normd, rho)) {
749566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "newton direction not descent\n"));
759566063dSJacob Faibussowitsch       PetscCall(VecCopy(ssls->dpsi, tao->stepdirection));
769566063dSJacob Faibussowitsch       PetscCall(VecDot(ssls->w, ssls->dpsi, &innerd));
77a7e14dcfSSatish Balay     }
78a7e14dcfSSatish Balay 
799566063dSJacob Faibussowitsch     PetscCall(VecScale(tao->stepdirection, -1.0));
80a7e14dcfSSatish Balay     innerd = -innerd;
81a7e14dcfSSatish Balay 
829566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
839566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &psi, ssls->dpsi, tao->stepdirection, &t, &ls_reason));
849566063dSJacob Faibussowitsch     PetscCall(VecNorm(ssls->dpsi, NORM_2, &ndpsi));
85a7e14dcfSSatish Balay   }
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87a7e14dcfSSatish Balay }
88a7e14dcfSSatish Balay 
TaoDestroy_SSFLS(Tao tao)8966976f2fSJacob Faibussowitsch static PetscErrorCode TaoDestroy_SSFLS(Tao tao)
90d71ae5a4SJacob Faibussowitsch {
91a7e14dcfSSatish Balay   TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
92a7e14dcfSSatish Balay 
93a7e14dcfSSatish Balay   PetscFunctionBegin;
949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->ff));
959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->w));
969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->dpsi));
979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->da));
989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->db));
999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->t1));
1009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->t2));
101a958fbfcSStefano Zampini   PetscCall(KSPDestroy(&tao->ksp));
1029566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
1033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104a7e14dcfSSatish Balay }
105a7e14dcfSSatish Balay 
1061522df2eSJason Sarich /*MC
1071522df2eSJason Sarich    TAOSSFLS - Semi-smooth feasible linesearch algorithm for solving
1081522df2eSJason Sarich        complementarity constraints
1091522df2eSJason Sarich 
1101522df2eSJason Sarich    Options Database Keys:
1111522df2eSJason Sarich + -tao_ssls_delta - descent test fraction
1121522df2eSJason Sarich - -tao_ssls_rho - descent test power
1131522df2eSJason Sarich 
1141eb8069cSJason Sarich    Level: beginner
1151522df2eSJason Sarich M*/
1161eb8069cSJason Sarich 
TaoCreate_SSFLS(Tao tao)117d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_SSFLS(Tao tao)
118d71ae5a4SJacob Faibussowitsch {
119a7e14dcfSSatish Balay   TAO_SSLS   *ssls;
1208caf6e8cSBarry Smith   const char *armijo_type = TAOLINESEARCHARMIJO;
121a7e14dcfSSatish Balay 
122a7e14dcfSSatish Balay   PetscFunctionBegin;
1234dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ssls));
124a7e14dcfSSatish Balay   tao->data                = (void *)ssls;
125a7e14dcfSSatish Balay   tao->ops->solve          = TaoSolve_SSFLS;
126a7e14dcfSSatish Balay   tao->ops->setup          = TaoSetUp_SSFLS;
127a7e14dcfSSatish Balay   tao->ops->view           = TaoView_SSLS;
128a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_SSLS;
129a7e14dcfSSatish Balay   tao->ops->destroy        = TaoDestroy_SSFLS;
130a7e14dcfSSatish Balay 
131a7e14dcfSSatish Balay   ssls->delta = 1e-10;
132a7e14dcfSSatish Balay   ssls->rho   = 2.1;
133a7e14dcfSSatish Balay 
1349566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
1359566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
1369566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, armijo_type));
1379566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
1389566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
13947a47007SBarry Smith   /* Linesearch objective and objectivegradient routines are  set in solve routine */
1409566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
1419566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
1429566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix));
143a7e14dcfSSatish Balay 
1446552cf8aSJason Sarich   /* Override default settings (unless already changed) */
145*606f75f6SBarry Smith   PetscCall(TaoParametersInitialize(tao));
146*606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, max_it, 2000);
147*606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, max_funcs, 4000);
148*606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, gttol, 0);
149*606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, grtol, 0);
150*606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, gatol, PetscDefined(USE_REAL_SINGLE) ? 1.0e-6 : 1.0e-16);
151*606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, fmin, PetscDefined(USE_REAL_SINGLE) ? 1.0e-4 : 1.0e-8);
1523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
153a7e14dcfSSatish Balay }
154