xref: /petsc/src/tao/complementarity/impls/ssls/ssils.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1aaa7dc30SBarry Smith #include <../src/tao/complementarity/impls/ssls/ssls.h>
2a7e14dcfSSatish Balay 
39371c9d4SSatish Balay PetscErrorCode TaoSetUp_SSILS(Tao tao) {
4a7e14dcfSSatish Balay   TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
5a7e14dcfSSatish Balay 
6a7e14dcfSSatish Balay   PetscFunctionBegin;
79566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tao->gradient));
89566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
99566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &ssls->ff));
109566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &ssls->dpsi));
119566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &ssls->da));
129566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &ssls->db));
139566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &ssls->t1));
149566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &ssls->t2));
15a7e14dcfSSatish Balay   PetscFunctionReturn(0);
16a7e14dcfSSatish Balay }
17a7e14dcfSSatish Balay 
189371c9d4SSatish Balay PetscErrorCode TaoDestroy_SSILS(Tao tao) {
19a7e14dcfSSatish Balay   TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
20a7e14dcfSSatish Balay 
21a7e14dcfSSatish Balay   PetscFunctionBegin;
229566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->ff));
239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->dpsi));
249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->da));
259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->db));
269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->t1));
279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->t2));
28a958fbfcSStefano Zampini   PetscCall(KSPDestroy(&tao->ksp));
299566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
30a7e14dcfSSatish Balay   PetscFunctionReturn(0);
31a7e14dcfSSatish Balay }
32a7e14dcfSSatish Balay 
339371c9d4SSatish Balay static PetscErrorCode TaoSolve_SSILS(Tao tao) {
34a7e14dcfSSatish Balay   TAO_SSLS                    *ssls = (TAO_SSLS *)tao->data;
35a7e14dcfSSatish Balay   PetscReal                    psi, ndpsi, normd, innerd, t = 0;
36a7e14dcfSSatish Balay   PetscReal                    delta, rho;
37e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason;
38a7e14dcfSSatish Balay 
39a7e14dcfSSatish Balay   PetscFunctionBegin;
40a7e14dcfSSatish Balay   /* Assume that Setup has been called!
41a7e14dcfSSatish Balay      Set the structure for the Jacobian and create a linear solver. */
42a7e14dcfSSatish Balay   delta = ssls->delta;
43a7e14dcfSSatish Balay   rho   = ssls->rho;
44a7e14dcfSSatish Balay 
459566063dSJacob Faibussowitsch   PetscCall(TaoComputeVariableBounds(tao));
469566063dSJacob Faibussowitsch   PetscCall(VecMedian(tao->XL, tao->solution, tao->XU, tao->solution));
479566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch, Tao_SSLS_FunctionGradient, tao));
489566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetObjectiveRoutine(tao->linesearch, Tao_SSLS_Function, tao));
49a7e14dcfSSatish Balay 
50a7e14dcfSSatish Balay   /* Calculate the function value and fischer function value at the
51a7e14dcfSSatish Balay      current iterate */
529566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchComputeObjectiveAndGradient(tao->linesearch, tao->solution, &psi, ssls->dpsi));
539566063dSJacob Faibussowitsch   PetscCall(VecNorm(ssls->dpsi, NORM_2, &ndpsi));
54a7e14dcfSSatish Balay 
55763847b4SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
56d90ca52eSBarry Smith   while (PETSC_TRUE) {
5763a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(tao, "iter: %" PetscInt_FMT ", merit: %g, ndpsi: %g\n", tao->niter, (double)ssls->merit, (double)ndpsi));
58a7e14dcfSSatish Balay     /* Check the termination criteria */
599566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, ssls->merit, ndpsi, 0.0, tao->ksp_its));
609566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, ssls->merit, ndpsi, 0.0, t));
61dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
62763847b4SAlp Dener     if (tao->reason != TAO_CONTINUE_ITERATING) break;
63e1e80dc8SAlp Dener 
64e1e80dc8SAlp Dener     /* Call general purpose update function */
65dbbe0bcdSBarry Smith     PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
66e6d4cb7fSJason Sarich     tao->niter++;
67a7e14dcfSSatish Balay 
68a7e14dcfSSatish Balay     /* Calculate direction.  (Really negative of newton direction.  Therefore,
69a7e14dcfSSatish Balay        rest of the code uses -d.) */
709566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(tao->ksp, tao->jacobian, tao->jacobian_pre));
719566063dSJacob Faibussowitsch     PetscCall(KSPSolve(tao->ksp, ssls->ff, tao->stepdirection));
729566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(tao->ksp, &tao->ksp_its));
73b0026674SJason Sarich     tao->ksp_tot_its += tao->ksp_its;
749566063dSJacob Faibussowitsch     PetscCall(VecNorm(tao->stepdirection, NORM_2, &normd));
759566063dSJacob Faibussowitsch     PetscCall(VecDot(tao->stepdirection, ssls->dpsi, &innerd));
76a7e14dcfSSatish Balay 
77a7e14dcfSSatish Balay     /* Make sure that we have a descent direction */
781118d4bcSLisandro Dalcin     if (innerd <= delta * PetscPowReal(normd, rho)) {
799566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "newton direction not descent\n"));
809566063dSJacob Faibussowitsch       PetscCall(VecCopy(ssls->dpsi, tao->stepdirection));
819566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->stepdirection, ssls->dpsi, &innerd));
82a7e14dcfSSatish Balay     }
83a7e14dcfSSatish Balay 
849566063dSJacob Faibussowitsch     PetscCall(VecScale(tao->stepdirection, -1.0));
85a7e14dcfSSatish Balay     innerd = -innerd;
86a7e14dcfSSatish Balay 
879566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
889566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &psi, ssls->dpsi, tao->stepdirection, &t, &ls_reason));
899566063dSJacob Faibussowitsch     PetscCall(VecNorm(ssls->dpsi, NORM_2, &ndpsi));
90a7e14dcfSSatish Balay   }
91a7e14dcfSSatish Balay   PetscFunctionReturn(0);
92a7e14dcfSSatish Balay }
93a7e14dcfSSatish Balay 
94a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
951522df2eSJason Sarich /*MC
961522df2eSJason Sarich    TAOSSILS - semi-smooth infeasible linesearch algorithm for solving
971522df2eSJason Sarich        complementarity constraints
981522df2eSJason Sarich 
991522df2eSJason Sarich    Options Database Keys:
1001522df2eSJason Sarich + -tao_ssls_delta - descent test fraction
1011522df2eSJason Sarich - -tao_ssls_rho - descent test power
1021522df2eSJason Sarich 
1031eb8069cSJason Sarich    Level: beginner
1041522df2eSJason Sarich M*/
1059371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TaoCreate_SSILS(Tao tao) {
106a7e14dcfSSatish Balay   TAO_SSLS   *ssls;
1078caf6e8cSBarry Smith   const char *armijo_type = TAOLINESEARCHARMIJO;
108a7e14dcfSSatish Balay 
109a7e14dcfSSatish Balay   PetscFunctionBegin;
110*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ssls));
111a7e14dcfSSatish Balay   tao->data                = (void *)ssls;
112a7e14dcfSSatish Balay   tao->ops->solve          = TaoSolve_SSILS;
113a7e14dcfSSatish Balay   tao->ops->setup          = TaoSetUp_SSILS;
114a7e14dcfSSatish Balay   tao->ops->view           = TaoView_SSLS;
115a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_SSLS;
116a7e14dcfSSatish Balay   tao->ops->destroy        = TaoDestroy_SSILS;
117a7e14dcfSSatish Balay 
118a7e14dcfSSatish Balay   ssls->delta = 1e-10;
119a7e14dcfSSatish Balay   ssls->rho   = 2.1;
120a7e14dcfSSatish Balay 
1219566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
1229566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
1239566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, armijo_type));
1249566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
1259566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
126a7e14dcfSSatish Balay   /* Note: linesearch objective and objectivegradient routines are set in solve routine */
1279566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
1289566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
1299566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix));
130a7e14dcfSSatish Balay 
1316552cf8aSJason Sarich   /* Override default settings (unless already changed) */
1326552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 2000;
1336552cf8aSJason Sarich   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
1346552cf8aSJason Sarich   if (!tao->gttol_changed) tao->gttol = 0;
1356552cf8aSJason Sarich   if (!tao->grtol_changed) tao->grtol = 0;
1366f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
1376552cf8aSJason Sarich   if (!tao->gatol_changed) tao->gatol = 1.0e-6;
1386552cf8aSJason Sarich   if (!tao->fmin_changed) tao->fmin = 1.0e-4;
1396f4723b1SBarry Smith #else
1406552cf8aSJason Sarich   if (!tao->gatol_changed) tao->gatol = 1.0e-16;
1416552cf8aSJason Sarich   if (!tao->fmin_changed) tao->fmin = 1.0e-8;
1426f4723b1SBarry Smith #endif
143a7e14dcfSSatish Balay   PetscFunctionReturn(0);
144a7e14dcfSSatish Balay }
145