xref: /petsc/src/tao/complementarity/impls/asls/asils.c (revision 5b6bfdb9644f185dbf5e5a09b808ec241507e1e7)
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 
56 static PetscErrorCode TaoSetUp_ASILS(Tao tao)
57 {
58   TAO_SSLS       *asls = (TAO_SSLS *)tao->data;
59   PetscErrorCode ierr;
60 
61   PetscFunctionBegin;
62   ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
63   ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
64   ierr = VecDuplicate(tao->solution,&asls->ff);CHKERRQ(ierr);
65   ierr = VecDuplicate(tao->solution,&asls->dpsi);CHKERRQ(ierr);
66   ierr = VecDuplicate(tao->solution,&asls->da);CHKERRQ(ierr);
67   ierr = VecDuplicate(tao->solution,&asls->db);CHKERRQ(ierr);
68   ierr = VecDuplicate(tao->solution,&asls->t1);CHKERRQ(ierr);
69   ierr = VecDuplicate(tao->solution,&asls->t2);CHKERRQ(ierr);
70   asls->fixed = NULL;
71   asls->free = NULL;
72   asls->J_sub = NULL;
73   asls->Jpre_sub = NULL;
74   asls->w = NULL;
75   asls->r1 = NULL;
76   asls->r2 = NULL;
77   asls->r3 = NULL;
78   asls->dxfree = NULL;
79   PetscFunctionReturn(0);
80 }
81 
82 static PetscErrorCode Tao_ASLS_FunctionGradient(TaoLineSearch ls, Vec X, PetscReal *fcn,  Vec G, void *ptr)
83 {
84   Tao            tao = (Tao)ptr;
85   TAO_SSLS       *asls = (TAO_SSLS *)tao->data;
86   PetscErrorCode ierr;
87 
88   PetscFunctionBegin;
89   ierr = TaoComputeConstraints(tao, X, tao->constraints);CHKERRQ(ierr);
90   ierr = VecFischer(X,tao->constraints,tao->XL,tao->XU,asls->ff);CHKERRQ(ierr);
91   ierr = VecNorm(asls->ff,NORM_2,&asls->merit);CHKERRQ(ierr);
92   *fcn = 0.5*asls->merit*asls->merit;
93 
94   ierr = TaoComputeJacobian(tao,tao->solution,tao->jacobian,tao->jacobian_pre);CHKERRQ(ierr);
95   ierr = MatDFischer(tao->jacobian, tao->solution, tao->constraints,tao->XL, tao->XU, asls->t1, asls->t2,asls->da, asls->db);CHKERRQ(ierr);
96   ierr = VecPointwiseMult(asls->t1, asls->ff, asls->db);CHKERRQ(ierr);
97   ierr = MatMultTranspose(tao->jacobian,asls->t1,G);CHKERRQ(ierr);
98   ierr = VecPointwiseMult(asls->t1, asls->ff, asls->da);CHKERRQ(ierr);
99   ierr = VecAXPY(G,1.0,asls->t1);CHKERRQ(ierr);
100   PetscFunctionReturn(0);
101 }
102 
103 static PetscErrorCode TaoDestroy_ASILS(Tao tao)
104 {
105   TAO_SSLS       *ssls = (TAO_SSLS *)tao->data;
106   PetscErrorCode ierr;
107 
108   PetscFunctionBegin;
109   ierr = VecDestroy(&ssls->ff);CHKERRQ(ierr);
110   ierr = VecDestroy(&ssls->dpsi);CHKERRQ(ierr);
111   ierr = VecDestroy(&ssls->da);CHKERRQ(ierr);
112   ierr = VecDestroy(&ssls->db);CHKERRQ(ierr);
113   ierr = VecDestroy(&ssls->w);CHKERRQ(ierr);
114   ierr = VecDestroy(&ssls->t1);CHKERRQ(ierr);
115   ierr = VecDestroy(&ssls->t2);CHKERRQ(ierr);
116   ierr = VecDestroy(&ssls->r1);CHKERRQ(ierr);
117   ierr = VecDestroy(&ssls->r2);CHKERRQ(ierr);
118   ierr = VecDestroy(&ssls->r3);CHKERRQ(ierr);
119   ierr = VecDestroy(&ssls->dxfree);CHKERRQ(ierr);
120   ierr = MatDestroy(&ssls->J_sub);CHKERRQ(ierr);
121   ierr = MatDestroy(&ssls->Jpre_sub);CHKERRQ(ierr);
122   ierr = ISDestroy(&ssls->fixed);CHKERRQ(ierr);
123   ierr = ISDestroy(&ssls->free);CHKERRQ(ierr);
124   ierr = PetscFree(tao->data);CHKERRQ(ierr);
125   PetscFunctionReturn(0);
126 }
127 
128 static PetscErrorCode TaoSolve_ASILS(Tao tao)
129 {
130   TAO_SSLS                     *asls = (TAO_SSLS *)tao->data;
131   PetscReal                    psi,ndpsi, normd, innerd, t=0;
132   PetscInt                     nf;
133   PetscErrorCode               ierr;
134   TaoLineSearchConvergedReason ls_reason;
135 
136   PetscFunctionBegin;
137   /* Assume that Setup has been called!
138      Set the structure for the Jacobian and create a linear solver. */
139 
140   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
141   ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_ASLS_FunctionGradient,tao);CHKERRQ(ierr);
142   ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);CHKERRQ(ierr);
143 
144   /* Calculate the function value and fischer function value at the
145      current iterate */
146   ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,asls->dpsi);CHKERRQ(ierr);
147   ierr = VecNorm(asls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr);
148 
149   tao->reason = TAO_CONTINUE_ITERATING;
150   while (1) {
151     /* Check the termination criteria */
152     ierr = PetscInfo3(tao,"iter %D, merit: %g, ||dpsi||: %g\n",tao->niter, (double)asls->merit,  (double)ndpsi);CHKERRQ(ierr);
153     ierr = TaoLogConvergenceHistory(tao,asls->merit,ndpsi,0.0,tao->ksp_its);CHKERRQ(ierr);
154     ierr = TaoMonitor(tao,tao->niter,asls->merit,ndpsi,0.0,t);CHKERRQ(ierr);
155     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
156     if (TAO_CONTINUE_ITERATING != tao->reason) break;
157     tao->niter++;
158 
159     /* We are going to solve a linear system of equations.  We need to
160        set the tolerances for the solve so that we maintain an asymptotic
161        rate of convergence that is superlinear.
162        Note: these tolerances are for the reduced system.  We really need
163        to make sure that the full system satisfies the full-space conditions.
164 
165        This rule gives superlinear asymptotic convergence
166        asls->atol = min(0.5, asls->merit*sqrt(asls->merit));
167        asls->rtol = 0.0;
168 
169        This rule gives quadratic asymptotic convergence
170        asls->atol = min(0.5, asls->merit*asls->merit);
171        asls->rtol = 0.0;
172 
173        Calculate a free and fixed set of variables.  The fixed set of
174        variables are those for the d_b is approximately equal to zero.
175        The definition of approximately changes as we approach the solution
176        to the problem.
177 
178        No one rule is guaranteed to work in all cases.  The following
179        definition is based on the norm of the Jacobian matrix.  If the
180        norm is large, the tolerance becomes smaller. */
181     ierr = MatNorm(tao->jacobian,NORM_1,&asls->identifier);CHKERRQ(ierr);
182     asls->identifier = PetscMin(asls->merit, 1e-2) / (1 + asls->identifier);
183 
184     ierr = VecSet(asls->t1,-asls->identifier);CHKERRQ(ierr);
185     ierr = VecSet(asls->t2, asls->identifier);CHKERRQ(ierr);
186 
187     ierr = ISDestroy(&asls->fixed);CHKERRQ(ierr);
188     ierr = ISDestroy(&asls->free);CHKERRQ(ierr);
189     ierr = VecWhichBetweenOrEqual(asls->t1, asls->db, asls->t2, &asls->fixed);CHKERRQ(ierr);
190     ierr = ISComplementVec(asls->fixed,asls->t1, &asls->free);CHKERRQ(ierr);
191 
192     ierr = ISGetSize(asls->fixed,&nf);CHKERRQ(ierr);
193     ierr = PetscInfo1(tao,"Number of fixed variables: %D\n", nf);CHKERRQ(ierr);
194 
195     /* We now have our partition.  Now calculate the direction in the
196        fixed variable space. */
197     ierr = TaoVecGetSubVec(asls->ff, asls->fixed, tao->subset_type, 0.0, &asls->r1);CHKERRQ(ierr);
198     ierr = TaoVecGetSubVec(asls->da, asls->fixed, tao->subset_type, 1.0, &asls->r2);CHKERRQ(ierr);
199     ierr = VecPointwiseDivide(asls->r1,asls->r1,asls->r2);CHKERRQ(ierr);
200     ierr = VecSet(tao->stepdirection,0.0);CHKERRQ(ierr);
201     ierr = VecISAXPY(tao->stepdirection, asls->fixed,1.0,asls->r1);CHKERRQ(ierr);
202 
203     /* Our direction in the Fixed Variable Set is fixed.  Calculate the
204        information needed for the step in the Free Variable Set.  To
205        do this, we need to know the diagonal perturbation and the
206        right hand side. */
207 
208     ierr = TaoVecGetSubVec(asls->da, asls->free, tao->subset_type, 0.0, &asls->r1);CHKERRQ(ierr);
209     ierr = TaoVecGetSubVec(asls->ff, asls->free, tao->subset_type, 0.0, &asls->r2);CHKERRQ(ierr);
210     ierr = TaoVecGetSubVec(asls->db, asls->free, tao->subset_type, 1.0, &asls->r3);CHKERRQ(ierr);
211     ierr = VecPointwiseDivide(asls->r1,asls->r1, asls->r3);CHKERRQ(ierr);
212     ierr = VecPointwiseDivide(asls->r2,asls->r2, asls->r3);CHKERRQ(ierr);
213 
214     /* r1 is the diagonal perturbation
215        r2 is the right hand side
216        r3 is no longer needed
217 
218        Now need to modify r2 for our direction choice in the fixed
219        variable set:  calculate t1 = J*d, take the reduced vector
220        of t1 and modify r2. */
221 
222     ierr = MatMult(tao->jacobian, tao->stepdirection, asls->t1);CHKERRQ(ierr);
223     ierr = TaoVecGetSubVec(asls->t1,asls->free,tao->subset_type,0.0,&asls->r3);CHKERRQ(ierr);
224     ierr = VecAXPY(asls->r2, -1.0, asls->r3);CHKERRQ(ierr);
225 
226     /* Calculate the reduced problem matrix and the direction */
227     if (!asls->w && (tao->subset_type == TAO_SUBSET_MASK || tao->subset_type == TAO_SUBSET_MATRIXFREE)) {
228       ierr = VecDuplicate(tao->solution, &asls->w);CHKERRQ(ierr);
229     }
230     ierr = TaoMatGetSubMat(tao->jacobian, asls->free, asls->w, tao->subset_type,&asls->J_sub);CHKERRQ(ierr);
231     if (tao->jacobian != tao->jacobian_pre) {
232       ierr = TaoMatGetSubMat(tao->jacobian_pre, asls->free, asls->w, tao->subset_type, &asls->Jpre_sub);CHKERRQ(ierr);
233     } else {
234       ierr = MatDestroy(&asls->Jpre_sub);CHKERRQ(ierr);
235       asls->Jpre_sub = asls->J_sub;
236       ierr = PetscObjectReference((PetscObject)(asls->Jpre_sub));CHKERRQ(ierr);
237     }
238     ierr = MatDiagonalSet(asls->J_sub, asls->r1,ADD_VALUES);CHKERRQ(ierr);
239     ierr = TaoVecGetSubVec(tao->stepdirection, asls->free, tao->subset_type, 0.0, &asls->dxfree);CHKERRQ(ierr);
240     ierr = VecSet(asls->dxfree, 0.0);CHKERRQ(ierr);
241 
242     /* Calculate the reduced direction.  (Really negative of Newton
243        direction.  Therefore, rest of the code uses -d.) */
244     ierr = KSPReset(tao->ksp);CHKERRQ(ierr);
245     ierr = KSPSetOperators(tao->ksp, asls->J_sub, asls->Jpre_sub);CHKERRQ(ierr);
246     ierr = KSPSolve(tao->ksp, asls->r2, asls->dxfree);CHKERRQ(ierr);
247     ierr = KSPGetIterationNumber(tao->ksp,&tao->ksp_its);CHKERRQ(ierr);
248     tao->ksp_tot_its+=tao->ksp_its;
249 
250     /* Add the direction in the free variables back into the real direction. */
251     ierr = VecISAXPY(tao->stepdirection, asls->free, 1.0,asls->dxfree);CHKERRQ(ierr);
252 
253     /* Check the real direction for descent and if not, use the negative
254        gradient direction. */
255     ierr = VecNorm(tao->stepdirection, NORM_2, &normd);CHKERRQ(ierr);
256     ierr = VecDot(tao->stepdirection, asls->dpsi, &innerd);CHKERRQ(ierr);
257 
258     if (innerd <= asls->delta*PetscPowReal(normd, asls->rho)) {
259       ierr = PetscInfo1(tao,"Gradient direction: %5.4e.\n", (double)innerd);CHKERRQ(ierr);
260       ierr = PetscInfo1(tao, "Iteration %D: newton direction not descent\n", tao->niter);CHKERRQ(ierr);
261       ierr = VecCopy(asls->dpsi, tao->stepdirection);CHKERRQ(ierr);
262       ierr = VecDot(asls->dpsi, tao->stepdirection, &innerd);CHKERRQ(ierr);
263     }
264 
265     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
266     innerd = -innerd;
267 
268     /* We now have a correct descent direction.  Apply a linesearch to
269        find the new iterate. */
270     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr);
271     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &psi,asls->dpsi, tao->stepdirection, &t, &ls_reason);CHKERRQ(ierr);
272     ierr = VecNorm(asls->dpsi, NORM_2, &ndpsi);CHKERRQ(ierr);
273   }
274   PetscFunctionReturn(0);
275 }
276 
277 /* ---------------------------------------------------------- */
278 /*MC
279    TAOASILS - Active-set infeasible linesearch algorithm for solving
280        complementarity constraints
281 
282    Options Database Keys:
283 + -tao_ssls_delta - descent test fraction
284 - -tao_ssls_rho - descent test power
285 
286   Level: beginner
287 M*/
288 PETSC_EXTERN PetscErrorCode TaoCreate_ASILS(Tao tao)
289 {
290   TAO_SSLS       *asls;
291   PetscErrorCode ierr;
292   const char     *armijo_type = TAOLINESEARCHARMIJO;
293 
294   PetscFunctionBegin;
295   ierr = PetscNewLog(tao,&asls);CHKERRQ(ierr);
296   tao->data = (void*)asls;
297   tao->ops->solve = TaoSolve_ASILS;
298   tao->ops->setup = TaoSetUp_ASILS;
299   tao->ops->view = TaoView_SSLS;
300   tao->ops->setfromoptions = TaoSetFromOptions_SSLS;
301   tao->ops->destroy = TaoDestroy_ASILS;
302   tao->subset_type = TAO_SUBSET_SUBVEC;
303   asls->delta = 1e-10;
304   asls->rho = 2.1;
305   asls->fixed = NULL;
306   asls->free = NULL;
307   asls->J_sub = NULL;
308   asls->Jpre_sub = NULL;
309   asls->w = NULL;
310   asls->r1 = NULL;
311   asls->r2 = NULL;
312   asls->r3 = NULL;
313   asls->t1 = NULL;
314   asls->t2 = NULL;
315   asls->dxfree = NULL;
316 
317   asls->identifier = 1e-5;
318 
319   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
320   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
321   ierr = TaoLineSearchSetType(tao->linesearch, armijo_type);CHKERRQ(ierr);
322   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
323   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
324 
325   ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr);
326   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
327   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
328   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
329 
330   /* Override default settings (unless already changed) */
331   if (!tao->max_it_changed) tao->max_it = 2000;
332   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
333   if (!tao->gttol_changed) tao->gttol = 0;
334   if (!tao->grtol_changed) tao->grtol = 0;
335 #if defined(PETSC_USE_REAL_SINGLE)
336   if (!tao->gatol_changed) tao->gatol = 1.0e-6;
337   if (!tao->fmin_changed)  tao->fmin = 1.0e-4;
338 #else
339   if (!tao->gatol_changed) tao->gatol = 1.0e-16;
340   if (!tao->fmin_changed) tao->fmin = 1.0e-8;
341 #endif
342   PetscFunctionReturn(0);
343 }
344