xref: /petsc/src/tao/complementarity/impls/ssls/ssils.c (revision 60a53c4cce06db480fbcf932f4e3989079773c3a)
1 #include <../src/tao/complementarity/impls/ssls/ssls.h>
2 
3 PetscErrorCode TaoSetUp_SSILS(Tao tao)
4 {
5   TAO_SSLS       *ssls = (TAO_SSLS *)tao->data;
6 
7   PetscFunctionBegin;
8   PetscCall(VecDuplicate(tao->solution,&tao->gradient));
9   PetscCall(VecDuplicate(tao->solution,&tao->stepdirection));
10   PetscCall(VecDuplicate(tao->solution,&ssls->ff));
11   PetscCall(VecDuplicate(tao->solution,&ssls->dpsi));
12   PetscCall(VecDuplicate(tao->solution,&ssls->da));
13   PetscCall(VecDuplicate(tao->solution,&ssls->db));
14   PetscCall(VecDuplicate(tao->solution,&ssls->t1));
15   PetscCall(VecDuplicate(tao->solution,&ssls->t2));
16   PetscFunctionReturn(0);
17 }
18 
19 PetscErrorCode TaoDestroy_SSILS(Tao tao)
20 {
21   TAO_SSLS       *ssls = (TAO_SSLS *)tao->data;
22 
23   PetscFunctionBegin;
24   PetscCall(VecDestroy(&ssls->ff));
25   PetscCall(VecDestroy(&ssls->dpsi));
26   PetscCall(VecDestroy(&ssls->da));
27   PetscCall(VecDestroy(&ssls->db));
28   PetscCall(VecDestroy(&ssls->t1));
29   PetscCall(VecDestroy(&ssls->t2));
30   PetscCall(KSPDestroy(&tao->ksp));
31   PetscCall(PetscFree(tao->data));
32   PetscFunctionReturn(0);
33 }
34 
35 static PetscErrorCode TaoSolve_SSILS(Tao tao)
36 {
37   TAO_SSLS                     *ssls = (TAO_SSLS *)tao->data;
38   PetscReal                    psi, ndpsi, normd, innerd, t=0;
39   PetscReal                    delta, rho;
40   TaoLineSearchConvergedReason ls_reason;
41 
42   PetscFunctionBegin;
43   /* Assume that Setup has been called!
44      Set the structure for the Jacobian and create a linear solver. */
45   delta = ssls->delta;
46   rho = ssls->rho;
47 
48   PetscCall(TaoComputeVariableBounds(tao));
49   PetscCall(VecMedian(tao->XL,tao->solution,tao->XU,tao->solution));
50   PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao));
51   PetscCall(TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao));
52 
53   /* Calculate the function value and fischer function value at the
54      current iterate */
55   PetscCall(TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi));
56   PetscCall(VecNorm(ssls->dpsi,NORM_2,&ndpsi));
57 
58   tao->reason = TAO_CONTINUE_ITERATING;
59   while (PETSC_TRUE) {
60     PetscCall(PetscInfo(tao, "iter: %" PetscInt_FMT ", merit: %g, ndpsi: %g\n",tao->niter, (double)ssls->merit, (double)ndpsi));
61     /* Check the termination criteria */
62     PetscCall(TaoLogConvergenceHistory(tao,ssls->merit,ndpsi,0.0,tao->ksp_its));
63     PetscCall(TaoMonitor(tao,tao->niter,ssls->merit,ndpsi,0.0,t));
64     PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP));
65     if (tao->reason!=TAO_CONTINUE_ITERATING) break;
66 
67     /* Call general purpose update function */
68     if (tao->ops->update) {
69       PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update));
70     }
71     tao->niter++;
72 
73     /* Calculate direction.  (Really negative of newton direction.  Therefore,
74        rest of the code uses -d.) */
75     PetscCall(KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre));
76     PetscCall(KSPSolve(tao->ksp,ssls->ff,tao->stepdirection));
77     PetscCall(KSPGetIterationNumber(tao->ksp,&tao->ksp_its));
78     tao->ksp_tot_its+=tao->ksp_its;
79     PetscCall(VecNorm(tao->stepdirection,NORM_2,&normd));
80     PetscCall(VecDot(tao->stepdirection,ssls->dpsi,&innerd));
81 
82     /* Make sure that we have a descent direction */
83     if (innerd <= delta*PetscPowReal(normd, rho)) {
84       PetscCall(PetscInfo(tao, "newton direction not descent\n"));
85       PetscCall(VecCopy(ssls->dpsi,tao->stepdirection));
86       PetscCall(VecDot(tao->stepdirection,ssls->dpsi,&innerd));
87     }
88 
89     PetscCall(VecScale(tao->stepdirection, -1.0));
90     innerd = -innerd;
91 
92     PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,1.0));
93     PetscCall(TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason));
94     PetscCall(VecNorm(ssls->dpsi,NORM_2,&ndpsi));
95   }
96   PetscFunctionReturn(0);
97 }
98 
99 /* ---------------------------------------------------------- */
100 /*MC
101    TAOSSILS - semi-smooth infeasible linesearch algorithm for solving
102        complementarity constraints
103 
104    Options Database Keys:
105 + -tao_ssls_delta - descent test fraction
106 - -tao_ssls_rho - descent test power
107 
108    Level: beginner
109 M*/
110 PETSC_EXTERN PetscErrorCode TaoCreate_SSILS(Tao tao)
111 {
112   TAO_SSLS       *ssls;
113   const char     *armijo_type = TAOLINESEARCHARMIJO;
114 
115   PetscFunctionBegin;
116   PetscCall(PetscNewLog(tao,&ssls));
117   tao->data = (void*)ssls;
118   tao->ops->solve=TaoSolve_SSILS;
119   tao->ops->setup=TaoSetUp_SSILS;
120   tao->ops->view=TaoView_SSLS;
121   tao->ops->setfromoptions=TaoSetFromOptions_SSLS;
122   tao->ops->destroy = TaoDestroy_SSILS;
123 
124   ssls->delta = 1e-10;
125   ssls->rho = 2.1;
126 
127   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch));
128   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
129   PetscCall(TaoLineSearchSetType(tao->linesearch,armijo_type));
130   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix));
131   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
132   /* Note: linesearch objective and objectivegradient routines are set in solve routine */
133   PetscCall(KSPCreate(((PetscObject)tao)->comm,&tao->ksp));
134   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
135   PetscCall(KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix));
136 
137   /* Override default settings (unless already changed) */
138   if (!tao->max_it_changed) tao->max_it = 2000;
139   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
140   if (!tao->gttol_changed) tao->gttol = 0;
141   if (!tao->grtol_changed) tao->grtol = 0;
142 #if defined(PETSC_USE_REAL_SINGLE)
143   if (!tao->gatol_changed) tao->gatol = 1.0e-6;
144   if (!tao->fmin_changed)  tao->fmin = 1.0e-4;
145 #else
146   if (!tao->gatol_changed) tao->gatol = 1.0e-16;
147   if (!tao->fmin_changed)  tao->fmin = 1.0e-8;
148 #endif
149   PetscFunctionReturn(0);
150 }
151