xref: /petsc/src/tao/complementarity/impls/ssls/ssfls.c (revision 9566063d113dddea24716c546802770db7481bc0)
1aaa7dc30SBarry Smith #include <../src/tao/complementarity/impls/ssls/ssls.h>
2a7e14dcfSSatish Balay 
3441846f8SBarry Smith PetscErrorCode TaoSetUp_SSFLS(Tao tao)
4a7e14dcfSSatish Balay {
5a7e14dcfSSatish Balay   TAO_SSLS       *ssls = (TAO_SSLS *)tao->data;
6a7e14dcfSSatish Balay 
7a7e14dcfSSatish Balay   PetscFunctionBegin;
8*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution,&tao->gradient));
9*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution,&tao->stepdirection));
10*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution,&ssls->w));
11*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution,&ssls->ff));
12*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution,&ssls->dpsi));
13*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution,&ssls->da));
14*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution,&ssls->db));
15*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution,&ssls->t1));
16*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution,&ssls->t2));
17a7e14dcfSSatish Balay   if (!tao->XL) {
18*9566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&tao->XL));
19a7e14dcfSSatish Balay   }
20a7e14dcfSSatish Balay   if (!tao->XU) {
21*9566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&tao->XU));
22a7e14dcfSSatish Balay   }
23*9566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU));
24a7e14dcfSSatish Balay   PetscFunctionReturn(0);
25a7e14dcfSSatish Balay }
26a7e14dcfSSatish Balay 
27441846f8SBarry Smith static PetscErrorCode TaoSolve_SSFLS(Tao tao)
28a7e14dcfSSatish Balay {
29a7e14dcfSSatish Balay   TAO_SSLS                     *ssls = (TAO_SSLS *)tao->data;
30a7e14dcfSSatish Balay   PetscReal                    psi, ndpsi, normd, innerd, t=0;
31a7e14dcfSSatish Balay   PetscReal                    delta, rho;
32e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason;
33a7e14dcfSSatish Balay 
34a7e14dcfSSatish Balay   PetscFunctionBegin;
35a7e14dcfSSatish Balay   /* Assume that Setup has been called!
36a7e14dcfSSatish Balay      Set the structure for the Jacobian and create a linear solver. */
37a7e14dcfSSatish Balay   delta = ssls->delta;
38a7e14dcfSSatish Balay   rho = ssls->rho;
39a7e14dcfSSatish Balay 
40*9566063dSJacob Faibussowitsch   PetscCall(TaoComputeVariableBounds(tao));
41a7e14dcfSSatish Balay   /* Project solution inside bounds */
42*9566063dSJacob Faibussowitsch   PetscCall(VecMedian(tao->XL,tao->solution,tao->XU,tao->solution));
43*9566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao));
44*9566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao));
45a7e14dcfSSatish Balay 
46a7e14dcfSSatish Balay   /* Calculate the function value and fischer function value at the
47a7e14dcfSSatish Balay      current iterate */
48*9566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi));
49*9566063dSJacob Faibussowitsch   PetscCall(VecNorm(ssls->dpsi,NORM_2,&ndpsi));
50a7e14dcfSSatish Balay 
51763847b4SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
52d90ca52eSBarry Smith   while (PETSC_TRUE) {
53*9566063dSJacob Faibussowitsch     PetscCall(PetscInfo(tao, "iter: %D, merit: %g, ndpsi: %g\n",tao->niter, (double)ssls->merit, (double)ndpsi));
54a7e14dcfSSatish Balay     /* Check the termination criteria */
55*9566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao,ssls->merit,ndpsi,0.0,tao->ksp_its));
56*9566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao,tao->niter,ssls->merit,ndpsi,0.0,t));
57*9566063dSJacob Faibussowitsch     PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP));
58763847b4SAlp Dener     if (tao->reason!=TAO_CONTINUE_ITERATING) break;
59e1e80dc8SAlp Dener 
60e1e80dc8SAlp Dener     /* Call general purpose update function */
61e1e80dc8SAlp Dener     if (tao->ops->update) {
62*9566063dSJacob Faibussowitsch       PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update));
63e1e80dc8SAlp Dener     }
64e6d4cb7fSJason Sarich     tao->niter++;
65a7e14dcfSSatish Balay 
66a7e14dcfSSatish Balay     /* Calculate direction.  (Really negative of newton direction.  Therefore,
67a7e14dcfSSatish Balay        rest of the code uses -d.) */
68*9566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre));
69*9566063dSJacob Faibussowitsch     PetscCall(KSPSolve(tao->ksp,ssls->ff,tao->stepdirection));
70*9566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(tao->ksp,&tao->ksp_its));
71b0026674SJason Sarich     tao->ksp_tot_its+=tao->ksp_its;
72a7e14dcfSSatish Balay 
73*9566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->stepdirection,ssls->w));
74*9566063dSJacob Faibussowitsch     PetscCall(VecScale(ssls->w,-1.0));
75*9566063dSJacob Faibussowitsch     PetscCall(VecBoundGradientProjection(ssls->w,tao->solution,tao->XL,tao->XU,ssls->w));
76a7e14dcfSSatish Balay 
77*9566063dSJacob Faibussowitsch     PetscCall(VecNorm(ssls->w,NORM_2,&normd));
78*9566063dSJacob Faibussowitsch     PetscCall(VecDot(ssls->w,ssls->dpsi,&innerd));
79a7e14dcfSSatish Balay 
80a7e14dcfSSatish Balay     /* Make sure that we have a descent direction */
81d90ca52eSBarry Smith     if (innerd >= -delta*PetscPowReal(normd, rho)) {
82*9566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "newton direction not descent\n"));
83*9566063dSJacob Faibussowitsch       PetscCall(VecCopy(ssls->dpsi,tao->stepdirection));
84*9566063dSJacob Faibussowitsch       PetscCall(VecDot(ssls->w,ssls->dpsi,&innerd));
85a7e14dcfSSatish Balay     }
86a7e14dcfSSatish Balay 
87*9566063dSJacob Faibussowitsch     PetscCall(VecScale(tao->stepdirection, -1.0));
88a7e14dcfSSatish Balay     innerd = -innerd;
89a7e14dcfSSatish Balay 
90*9566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,1.0));
91*9566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason));
92*9566063dSJacob Faibussowitsch     PetscCall(VecNorm(ssls->dpsi,NORM_2,&ndpsi));
93a7e14dcfSSatish Balay   }
94a7e14dcfSSatish Balay   PetscFunctionReturn(0);
95a7e14dcfSSatish Balay }
96a7e14dcfSSatish Balay 
97441846f8SBarry Smith PetscErrorCode TaoDestroy_SSFLS(Tao tao)
98a7e14dcfSSatish Balay {
99a7e14dcfSSatish Balay   TAO_SSLS       *ssls = (TAO_SSLS *)tao->data;
100a7e14dcfSSatish Balay 
101a7e14dcfSSatish Balay   PetscFunctionBegin;
102*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->ff));
103*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->w));
104*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->dpsi));
105*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->da));
106*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->db));
107*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->t1));
108*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ssls->t2));
109*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
110a7e14dcfSSatish Balay   PetscFunctionReturn(0);
111a7e14dcfSSatish Balay }
112a7e14dcfSSatish Balay 
113a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
1141522df2eSJason Sarich /*MC
1151522df2eSJason Sarich    TAOSSFLS - Semi-smooth feasible linesearch algorithm for solving
1161522df2eSJason Sarich        complementarity constraints
1171522df2eSJason Sarich 
1181522df2eSJason Sarich    Options Database Keys:
1191522df2eSJason Sarich + -tao_ssls_delta - descent test fraction
1201522df2eSJason Sarich - -tao_ssls_rho - descent test power
1211522df2eSJason Sarich 
1221eb8069cSJason Sarich    Level: beginner
1231522df2eSJason Sarich M*/
1241eb8069cSJason Sarich 
125728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_SSFLS(Tao tao)
126a7e14dcfSSatish Balay {
127a7e14dcfSSatish Balay   TAO_SSLS       *ssls;
1288caf6e8cSBarry Smith   const char     *armijo_type = TAOLINESEARCHARMIJO;
129a7e14dcfSSatish Balay 
130a7e14dcfSSatish Balay   PetscFunctionBegin;
131*9566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(tao,&ssls));
132a7e14dcfSSatish Balay   tao->data = (void*)ssls;
133a7e14dcfSSatish Balay   tao->ops->solve=TaoSolve_SSFLS;
134a7e14dcfSSatish Balay   tao->ops->setup=TaoSetUp_SSFLS;
135a7e14dcfSSatish Balay   tao->ops->view=TaoView_SSLS;
136a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_SSLS;
137a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_SSFLS;
138a7e14dcfSSatish Balay 
139a7e14dcfSSatish Balay   ssls->delta = 1e-10;
140a7e14dcfSSatish Balay   ssls->rho = 2.1;
141a7e14dcfSSatish Balay 
142*9566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch));
143*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
144*9566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch,armijo_type));
145*9566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix));
146*9566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
14747a47007SBarry Smith   /* Linesearch objective and objectivegradient routines are  set in solve routine */
148*9566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm,&tao->ksp));
149*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
150*9566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix));
151a7e14dcfSSatish Balay 
1526552cf8aSJason Sarich   /* Override default settings (unless already changed) */
1536552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 2000;
1546552cf8aSJason Sarich   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
1556552cf8aSJason Sarich   if (!tao->gttol_changed) tao->gttol = 0;
1566552cf8aSJason Sarich   if (!tao->grtol_changed) tao->grtol = 0;
1576f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
1586552cf8aSJason Sarich   if (!tao->gatol_changed) tao->gatol = 1.0e-6;
1596552cf8aSJason Sarich   if (!tao->fmin_changed)  tao->fmin = 1.0e-4;
1606f4723b1SBarry Smith #else
1616552cf8aSJason Sarich   if (!tao->gatol_changed) tao->gatol = 1.0e-16;
1626552cf8aSJason Sarich   if (!tao->fmin_changed)  tao->fmin = 1.0e-8;
1636f4723b1SBarry Smith #endif
164a7e14dcfSSatish Balay   PetscFunctionReturn(0);
165a7e14dcfSSatish Balay }
166