xref: /petsc/src/tao/complementarity/impls/asls/asils.c (revision 030f984af8d8bb4c203755d35bded3c05b3d83ce)
1 #include <../src/tao/complementarity/impls/ssls/ssls.h>
2 /*
3    Context for ASXLS
4      -- active-set      - reduced matrices formed
5                           - inherit properties of original system
6      -- semismooth (S)  - function not differentiable
7                         - merit function continuously differentiable
8                         - Fischer-Burmeister reformulation of complementarity
9                           - Billups composition for two finite bounds
10      -- infeasible (I)  - iterates not guaranteed to remain within bounds
11      -- feasible (F)    - iterates guaranteed to remain within bounds
12      -- linesearch (LS) - Armijo rule on direction
13 
14    Many other reformulations are possible and combinations of
15    feasible/infeasible and linesearch/trust region are possible.
16 
17    Basic theory
18      Fischer-Burmeister reformulation is semismooth with a continuously
19      differentiable merit function and strongly semismooth if the F has
20      lipschitz continuous derivatives.
21 
22      Every accumulation point generated by the algorithm is a stationary
23      point for the merit function.  Stationary points of the merit function
24      are solutions of the complementarity problem if
25        a.  the stationary point has a BD-regular subdifferential, or
26        b.  the Schur complement F'/F'_ff is a P_0-matrix where ff is the
27            index set corresponding to the free variables.
28 
29      If one of the accumulation points has a BD-regular subdifferential then
30        a.  the entire sequence converges to this accumulation point at
31            a local q-superlinear rate
32        b.  if in addition the reformulation is strongly semismooth near
33            this accumulation point, then the algorithm converges at a
34            local q-quadratic rate.
35 
36    The theory for the feasible version follows from the feasible descent
37    algorithm framework.
38 
39    References:
40      Billups, "Algorithms for Complementarity Problems and Generalized
41        Equations," Ph.D thesis, University of Wisconsin  Madison, 1995.
42      De Luca, Facchinei, Kanzow, "A Semismooth Equation Approach to the
43        Solution of Nonlinear Complementarity Problems," Mathematical
44        Programming, 75, 1996.
45      Ferris, Kanzow, Munson, "Feasible Descent Algorithms for Mixed
46        Complementarity Problems," Mathematical Programming, 86,
47        1999.
48      Fischer, "A Special Newton type Optimization Method," Optimization,
49        24, 1992
50      Munson, Facchinei, Ferris, Fischer, Kanzow, "The Semismooth Algorithm
51        for Large Scale Complementarity Problems," Technical Report,
52        University of Wisconsin  Madison, 1999.
53 */
54 
55 static PetscErrorCode TaoSetUp_ASILS(Tao tao)
56 {
57   TAO_SSLS       *asls = (TAO_SSLS *)tao->data;
58   PetscErrorCode ierr;
59 
60   PetscFunctionBegin;
61   ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
62   ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
63   ierr = VecDuplicate(tao->solution,&asls->ff);CHKERRQ(ierr);
64   ierr = VecDuplicate(tao->solution,&asls->dpsi);CHKERRQ(ierr);
65   ierr = VecDuplicate(tao->solution,&asls->da);CHKERRQ(ierr);
66   ierr = VecDuplicate(tao->solution,&asls->db);CHKERRQ(ierr);
67   ierr = VecDuplicate(tao->solution,&asls->t1);CHKERRQ(ierr);
68   ierr = VecDuplicate(tao->solution,&asls->t2);CHKERRQ(ierr);
69   asls->fixed = NULL;
70   asls->free = NULL;
71   asls->J_sub = NULL;
72   asls->Jpre_sub = NULL;
73   asls->w = NULL;
74   asls->r1 = NULL;
75   asls->r2 = NULL;
76   asls->r3 = NULL;
77   asls->dxfree = NULL;
78   PetscFunctionReturn(0);
79 }
80 
81 static PetscErrorCode Tao_ASLS_FunctionGradient(TaoLineSearch ls, Vec X, PetscReal *fcn,  Vec G, void *ptr)
82 {
83   Tao            tao = (Tao)ptr;
84   TAO_SSLS       *asls = (TAO_SSLS *)tao->data;
85   PetscErrorCode ierr;
86 
87   PetscFunctionBegin;
88   ierr = TaoComputeConstraints(tao, X, tao->constraints);CHKERRQ(ierr);
89   ierr = VecFischer(X,tao->constraints,tao->XL,tao->XU,asls->ff);CHKERRQ(ierr);
90   ierr = VecNorm(asls->ff,NORM_2,&asls->merit);CHKERRQ(ierr);
91   *fcn = 0.5*asls->merit*asls->merit;
92 
93   ierr = TaoComputeJacobian(tao,tao->solution,tao->jacobian,tao->jacobian_pre);CHKERRQ(ierr);
94   ierr = MatDFischer(tao->jacobian, tao->solution, tao->constraints,tao->XL, tao->XU, asls->t1, asls->t2,asls->da, asls->db);CHKERRQ(ierr);
95   ierr = VecPointwiseMult(asls->t1, asls->ff, asls->db);CHKERRQ(ierr);
96   ierr = MatMultTranspose(tao->jacobian,asls->t1,G);CHKERRQ(ierr);
97   ierr = VecPointwiseMult(asls->t1, asls->ff, asls->da);CHKERRQ(ierr);
98   ierr = VecAXPY(G,1.0,asls->t1);CHKERRQ(ierr);
99   PetscFunctionReturn(0);
100 }
101 
102 static PetscErrorCode TaoDestroy_ASILS(Tao tao)
103 {
104   TAO_SSLS       *ssls = (TAO_SSLS *)tao->data;
105   PetscErrorCode ierr;
106 
107   PetscFunctionBegin;
108   ierr = VecDestroy(&ssls->ff);CHKERRQ(ierr);
109   ierr = VecDestroy(&ssls->dpsi);CHKERRQ(ierr);
110   ierr = VecDestroy(&ssls->da);CHKERRQ(ierr);
111   ierr = VecDestroy(&ssls->db);CHKERRQ(ierr);
112   ierr = VecDestroy(&ssls->w);CHKERRQ(ierr);
113   ierr = VecDestroy(&ssls->t1);CHKERRQ(ierr);
114   ierr = VecDestroy(&ssls->t2);CHKERRQ(ierr);
115   ierr = VecDestroy(&ssls->r1);CHKERRQ(ierr);
116   ierr = VecDestroy(&ssls->r2);CHKERRQ(ierr);
117   ierr = VecDestroy(&ssls->r3);CHKERRQ(ierr);
118   ierr = VecDestroy(&ssls->dxfree);CHKERRQ(ierr);
119   ierr = MatDestroy(&ssls->J_sub);CHKERRQ(ierr);
120   ierr = MatDestroy(&ssls->Jpre_sub);CHKERRQ(ierr);
121   ierr = ISDestroy(&ssls->fixed);CHKERRQ(ierr);
122   ierr = ISDestroy(&ssls->free);CHKERRQ(ierr);
123   ierr = PetscFree(tao->data);CHKERRQ(ierr);
124   PetscFunctionReturn(0);
125 }
126 
127 static PetscErrorCode TaoSolve_ASILS(Tao tao)
128 {
129   TAO_SSLS                     *asls = (TAO_SSLS *)tao->data;
130   PetscReal                    psi,ndpsi, normd, innerd, t=0;
131   PetscInt                     nf;
132   PetscErrorCode               ierr;
133   TaoLineSearchConvergedReason ls_reason;
134 
135   PetscFunctionBegin;
136   /* Assume that Setup has been called!
137      Set the structure for the Jacobian and create a linear solver. */
138 
139   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
140   ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_ASLS_FunctionGradient,tao);CHKERRQ(ierr);
141   ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);CHKERRQ(ierr);
142 
143   /* Calculate the function value and fischer function value at the
144      current iterate */
145   ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,asls->dpsi);CHKERRQ(ierr);
146   ierr = VecNorm(asls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr);
147 
148   tao->reason = TAO_CONTINUE_ITERATING;
149   while (1) {
150     /* Check the termination criteria */
151     ierr = PetscInfo3(tao,"iter %D, merit: %g, ||dpsi||: %g\n",tao->niter, (double)asls->merit,  (double)ndpsi);CHKERRQ(ierr);
152     ierr = TaoLogConvergenceHistory(tao,asls->merit,ndpsi,0.0,tao->ksp_its);CHKERRQ(ierr);
153     ierr = TaoMonitor(tao,tao->niter,asls->merit,ndpsi,0.0,t);CHKERRQ(ierr);
154     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
155     if (TAO_CONTINUE_ITERATING != tao->reason) break;
156 
157     /* Call general purpose update function */
158     if (tao->ops->update) {
159       ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr);
160     }
161     tao->niter++;
162 
163     /* We are going to solve a linear system of equations.  We need to
164        set the tolerances for the solve so that we maintain an asymptotic
165        rate of convergence that is superlinear.
166        Note: these tolerances are for the reduced system.  We really need
167        to make sure that the full system satisfies the full-space conditions.
168 
169        This rule gives superlinear asymptotic convergence
170        asls->atol = min(0.5, asls->merit*sqrt(asls->merit));
171        asls->rtol = 0.0;
172 
173        This rule gives quadratic asymptotic convergence
174        asls->atol = min(0.5, asls->merit*asls->merit);
175        asls->rtol = 0.0;
176 
177        Calculate a free and fixed set of variables.  The fixed set of
178        variables are those for the d_b is approximately equal to zero.
179        The definition of approximately changes as we approach the solution
180        to the problem.
181 
182        No one rule is guaranteed to work in all cases.  The following
183        definition is based on the norm of the Jacobian matrix.  If the
184        norm is large, the tolerance becomes smaller. */
185     ierr = MatNorm(tao->jacobian,NORM_1,&asls->identifier);CHKERRQ(ierr);
186     asls->identifier = PetscMin(asls->merit, 1e-2) / (1 + asls->identifier);
187 
188     ierr = VecSet(asls->t1,-asls->identifier);CHKERRQ(ierr);
189     ierr = VecSet(asls->t2, asls->identifier);CHKERRQ(ierr);
190 
191     ierr = ISDestroy(&asls->fixed);CHKERRQ(ierr);
192     ierr = ISDestroy(&asls->free);CHKERRQ(ierr);
193     ierr = VecWhichBetweenOrEqual(asls->t1, asls->db, asls->t2, &asls->fixed);CHKERRQ(ierr);
194     ierr = ISComplementVec(asls->fixed,asls->t1, &asls->free);CHKERRQ(ierr);
195 
196     ierr = ISGetSize(asls->fixed,&nf);CHKERRQ(ierr);
197     ierr = PetscInfo1(tao,"Number of fixed variables: %D\n", nf);CHKERRQ(ierr);
198 
199     /* We now have our partition.  Now calculate the direction in the
200        fixed variable space. */
201     ierr = TaoVecGetSubVec(asls->ff, asls->fixed, tao->subset_type, 0.0, &asls->r1);CHKERRQ(ierr);
202     ierr = TaoVecGetSubVec(asls->da, asls->fixed, tao->subset_type, 1.0, &asls->r2);CHKERRQ(ierr);
203     ierr = VecPointwiseDivide(asls->r1,asls->r1,asls->r2);CHKERRQ(ierr);
204     ierr = VecSet(tao->stepdirection,0.0);CHKERRQ(ierr);
205     ierr = VecISAXPY(tao->stepdirection, asls->fixed,1.0,asls->r1);CHKERRQ(ierr);
206 
207     /* Our direction in the Fixed Variable Set is fixed.  Calculate the
208        information needed for the step in the Free Variable Set.  To
209        do this, we need to know the diagonal perturbation and the
210        right hand side. */
211 
212     ierr = TaoVecGetSubVec(asls->da, asls->free, tao->subset_type, 0.0, &asls->r1);CHKERRQ(ierr);
213     ierr = TaoVecGetSubVec(asls->ff, asls->free, tao->subset_type, 0.0, &asls->r2);CHKERRQ(ierr);
214     ierr = TaoVecGetSubVec(asls->db, asls->free, tao->subset_type, 1.0, &asls->r3);CHKERRQ(ierr);
215     ierr = VecPointwiseDivide(asls->r1,asls->r1, asls->r3);CHKERRQ(ierr);
216     ierr = VecPointwiseDivide(asls->r2,asls->r2, asls->r3);CHKERRQ(ierr);
217 
218     /* r1 is the diagonal perturbation
219        r2 is the right hand side
220        r3 is no longer needed
221 
222        Now need to modify r2 for our direction choice in the fixed
223        variable set:  calculate t1 = J*d, take the reduced vector
224        of t1 and modify r2. */
225 
226     ierr = MatMult(tao->jacobian, tao->stepdirection, asls->t1);CHKERRQ(ierr);
227     ierr = TaoVecGetSubVec(asls->t1,asls->free,tao->subset_type,0.0,&asls->r3);CHKERRQ(ierr);
228     ierr = VecAXPY(asls->r2, -1.0, asls->r3);CHKERRQ(ierr);
229 
230     /* Calculate the reduced problem matrix and the direction */
231     if (!asls->w && (tao->subset_type == TAO_SUBSET_MASK || tao->subset_type == TAO_SUBSET_MATRIXFREE)) {
232       ierr = VecDuplicate(tao->solution, &asls->w);CHKERRQ(ierr);
233     }
234     ierr = TaoMatGetSubMat(tao->jacobian, asls->free, asls->w, tao->subset_type,&asls->J_sub);CHKERRQ(ierr);
235     if (tao->jacobian != tao->jacobian_pre) {
236       ierr = TaoMatGetSubMat(tao->jacobian_pre, asls->free, asls->w, tao->subset_type, &asls->Jpre_sub);CHKERRQ(ierr);
237     } else {
238       ierr = MatDestroy(&asls->Jpre_sub);CHKERRQ(ierr);
239       asls->Jpre_sub = asls->J_sub;
240       ierr = PetscObjectReference((PetscObject)(asls->Jpre_sub));CHKERRQ(ierr);
241     }
242     ierr = MatDiagonalSet(asls->J_sub, asls->r1,ADD_VALUES);CHKERRQ(ierr);
243     ierr = TaoVecGetSubVec(tao->stepdirection, asls->free, tao->subset_type, 0.0, &asls->dxfree);CHKERRQ(ierr);
244     ierr = VecSet(asls->dxfree, 0.0);CHKERRQ(ierr);
245 
246     /* Calculate the reduced direction.  (Really negative of Newton
247        direction.  Therefore, rest of the code uses -d.) */
248     ierr = KSPReset(tao->ksp);CHKERRQ(ierr);
249     ierr = KSPSetOperators(tao->ksp, asls->J_sub, asls->Jpre_sub);CHKERRQ(ierr);
250     ierr = KSPSolve(tao->ksp, asls->r2, asls->dxfree);CHKERRQ(ierr);
251     ierr = KSPGetIterationNumber(tao->ksp,&tao->ksp_its);CHKERRQ(ierr);
252     tao->ksp_tot_its+=tao->ksp_its;
253 
254     /* Add the direction in the free variables back into the real direction. */
255     ierr = VecISAXPY(tao->stepdirection, asls->free, 1.0,asls->dxfree);CHKERRQ(ierr);
256 
257     /* Check the real direction for descent and if not, use the negative
258        gradient direction. */
259     ierr = VecNorm(tao->stepdirection, NORM_2, &normd);CHKERRQ(ierr);
260     ierr = VecDot(tao->stepdirection, asls->dpsi, &innerd);CHKERRQ(ierr);
261 
262     if (innerd <= asls->delta*PetscPowReal(normd, asls->rho)) {
263       ierr = PetscInfo1(tao,"Gradient direction: %5.4e.\n", (double)innerd);CHKERRQ(ierr);
264       ierr = PetscInfo1(tao, "Iteration %D: newton direction not descent\n", tao->niter);CHKERRQ(ierr);
265       ierr = VecCopy(asls->dpsi, tao->stepdirection);CHKERRQ(ierr);
266       ierr = VecDot(asls->dpsi, tao->stepdirection, &innerd);CHKERRQ(ierr);
267     }
268 
269     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
270     innerd = -innerd;
271 
272     /* We now have a correct descent direction.  Apply a linesearch to
273        find the new iterate. */
274     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr);
275     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &psi,asls->dpsi, tao->stepdirection, &t, &ls_reason);CHKERRQ(ierr);
276     ierr = VecNorm(asls->dpsi, NORM_2, &ndpsi);CHKERRQ(ierr);
277   }
278   PetscFunctionReturn(0);
279 }
280 
281 /* ---------------------------------------------------------- */
282 /*MC
283    TAOASILS - Active-set infeasible linesearch algorithm for solving
284        complementarity constraints
285 
286    Options Database Keys:
287 + -tao_ssls_delta - descent test fraction
288 - -tao_ssls_rho - descent test power
289 
290   Level: beginner
291 M*/
292 PETSC_EXTERN PetscErrorCode TaoCreate_ASILS(Tao tao)
293 {
294   TAO_SSLS       *asls;
295   PetscErrorCode ierr;
296   const char     *armijo_type = TAOLINESEARCHARMIJO;
297 
298   PetscFunctionBegin;
299   ierr = PetscNewLog(tao,&asls);CHKERRQ(ierr);
300   tao->data = (void*)asls;
301   tao->ops->solve = TaoSolve_ASILS;
302   tao->ops->setup = TaoSetUp_ASILS;
303   tao->ops->view = TaoView_SSLS;
304   tao->ops->setfromoptions = TaoSetFromOptions_SSLS;
305   tao->ops->destroy = TaoDestroy_ASILS;
306   tao->subset_type = TAO_SUBSET_SUBVEC;
307   asls->delta = 1e-10;
308   asls->rho = 2.1;
309   asls->fixed = NULL;
310   asls->free = NULL;
311   asls->J_sub = NULL;
312   asls->Jpre_sub = NULL;
313   asls->w = NULL;
314   asls->r1 = NULL;
315   asls->r2 = NULL;
316   asls->r3 = NULL;
317   asls->t1 = NULL;
318   asls->t2 = NULL;
319   asls->dxfree = NULL;
320 
321   asls->identifier = 1e-5;
322 
323   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
324   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
325   ierr = TaoLineSearchSetType(tao->linesearch, armijo_type);CHKERRQ(ierr);
326   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
327   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
328 
329   ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr);
330   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
331   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
332   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
333 
334   /* Override default settings (unless already changed) */
335   if (!tao->max_it_changed) tao->max_it = 2000;
336   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
337   if (!tao->gttol_changed) tao->gttol = 0;
338   if (!tao->grtol_changed) tao->grtol = 0;
339 #if defined(PETSC_USE_REAL_SINGLE)
340   if (!tao->gatol_changed) tao->gatol = 1.0e-6;
341   if (!tao->fmin_changed)  tao->fmin = 1.0e-4;
342 #else
343   if (!tao->gatol_changed) tao->gatol = 1.0e-16;
344   if (!tao->fmin_changed) tao->fmin = 1.0e-8;
345 #endif
346   PetscFunctionReturn(0);
347 }
348