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