Actual source code: asfls.c
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
14: Many other reformulations are possible and combinations of
15: feasible/infeasible and linesearch/trust region are possible.
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.
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.
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.
36: The theory for the feasible version follows from the feasible descent
37: algorithm framework.
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 407439, 1996.
45: . * - Ferris, Kanzow, Munson, "Feasible Descent Algorithms for Mixed
46: Complementarity Problems," Mathematical Programming, 86,
47: pages 475497, 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: */
55: static PetscErrorCode TaoSetUp_ASFLS(Tao tao)
56: {
57: TAO_SSLS *asls = (TAO_SSLS *)tao->data;
59: VecDuplicate(tao->solution,&tao->gradient);
60: VecDuplicate(tao->solution,&tao->stepdirection);
61: VecDuplicate(tao->solution,&asls->ff);
62: VecDuplicate(tao->solution,&asls->dpsi);
63: VecDuplicate(tao->solution,&asls->da);
64: VecDuplicate(tao->solution,&asls->db);
65: VecDuplicate(tao->solution,&asls->t1);
66: VecDuplicate(tao->solution,&asls->t2);
67: VecDuplicate(tao->solution, &asls->w);
68: asls->fixed = NULL;
69: asls->free = NULL;
70: asls->J_sub = NULL;
71: asls->Jpre_sub = NULL;
72: asls->r1 = NULL;
73: asls->r2 = NULL;
74: asls->r3 = NULL;
75: asls->dxfree = NULL;
76: return 0;
77: }
79: static PetscErrorCode Tao_ASLS_FunctionGradient(TaoLineSearch ls, Vec X, PetscReal *fcn, Vec G, void *ptr)
80: {
81: Tao tao = (Tao)ptr;
82: TAO_SSLS *asls = (TAO_SSLS *)tao->data;
84: TaoComputeConstraints(tao, X, tao->constraints);
85: VecFischer(X,tao->constraints,tao->XL,tao->XU,asls->ff);
86: VecNorm(asls->ff,NORM_2,&asls->merit);
87: *fcn = 0.5*asls->merit*asls->merit;
88: TaoComputeJacobian(tao,tao->solution,tao->jacobian,tao->jacobian_pre);
90: MatDFischer(tao->jacobian, tao->solution, tao->constraints,tao->XL, tao->XU, asls->t1, asls->t2,asls->da, asls->db);
91: VecPointwiseMult(asls->t1, asls->ff, asls->db);
92: MatMultTranspose(tao->jacobian,asls->t1,G);
93: VecPointwiseMult(asls->t1, asls->ff, asls->da);
94: VecAXPY(G,1.0,asls->t1);
95: return 0;
96: }
98: static PetscErrorCode TaoDestroy_ASFLS(Tao tao)
99: {
100: TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
102: VecDestroy(&ssls->ff);
103: VecDestroy(&ssls->dpsi);
104: VecDestroy(&ssls->da);
105: VecDestroy(&ssls->db);
106: VecDestroy(&ssls->w);
107: VecDestroy(&ssls->t1);
108: VecDestroy(&ssls->t2);
109: VecDestroy(&ssls->r1);
110: VecDestroy(&ssls->r2);
111: VecDestroy(&ssls->r3);
112: VecDestroy(&ssls->dxfree);
113: MatDestroy(&ssls->J_sub);
114: MatDestroy(&ssls->Jpre_sub);
115: ISDestroy(&ssls->fixed);
116: ISDestroy(&ssls->free);
117: PetscFree(tao->data);
118: tao->data = NULL;
119: return 0;
120: }
122: static PetscErrorCode TaoSolve_ASFLS(Tao tao)
123: {
124: TAO_SSLS *asls = (TAO_SSLS *)tao->data;
125: PetscReal psi,ndpsi, normd, innerd, t=0;
126: PetscInt nf;
127: TaoLineSearchConvergedReason ls_reason;
129: /* Assume that Setup has been called!
130: Set the structure for the Jacobian and create a linear solver. */
132: TaoComputeVariableBounds(tao);
133: TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_ASLS_FunctionGradient,tao);
134: TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);
135: TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);
137: VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);
139: /* Calculate the function value and fischer function value at the
140: current iterate */
141: TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,asls->dpsi);
142: VecNorm(asls->dpsi,NORM_2,&ndpsi);
144: tao->reason = TAO_CONTINUE_ITERATING;
145: while (1) {
146: /* Check the converged criteria */
147: PetscInfo(tao,"iter %D, merit: %g, ||dpsi||: %g\n",tao->niter,(double)asls->merit,(double)ndpsi);
148: TaoLogConvergenceHistory(tao,asls->merit,ndpsi,0.0,tao->ksp_its);
149: TaoMonitor(tao,tao->niter,asls->merit,ndpsi,0.0,t);
150: (*tao->ops->convergencetest)(tao,tao->cnvP);
151: if (TAO_CONTINUE_ITERATING != tao->reason) break;
153: /* Call general purpose update function */
154: if (tao->ops->update) {
155: (*tao->ops->update)(tao, tao->niter, tao->user_update);
156: }
157: tao->niter++;
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.
165: This rule gives superlinear asymptotic convergence
166: asls->atol = min(0.5, asls->merit*sqrt(asls->merit));
167: asls->rtol = 0.0;
169: This rule gives quadratic asymptotic convergence
170: asls->atol = min(0.5, asls->merit*asls->merit);
171: asls->rtol = 0.0;
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.
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: MatNorm(tao->jacobian,NORM_1,&asls->identifier);
182: asls->identifier = PetscMin(asls->merit, 1e-2) / (1 + asls->identifier);
184: VecSet(asls->t1,-asls->identifier);
185: VecSet(asls->t2, asls->identifier);
187: ISDestroy(&asls->fixed);
188: ISDestroy(&asls->free);
189: VecWhichBetweenOrEqual(asls->t1, asls->db, asls->t2, &asls->fixed);
190: ISComplementVec(asls->fixed,asls->t1, &asls->free);
192: ISGetSize(asls->fixed,&nf);
193: PetscInfo(tao,"Number of fixed variables: %D\n", nf);
195: /* We now have our partition. Now calculate the direction in the
196: fixed variable space. */
197: TaoVecGetSubVec(asls->ff, asls->fixed, tao->subset_type, 0.0, &asls->r1);
198: TaoVecGetSubVec(asls->da, asls->fixed, tao->subset_type, 1.0, &asls->r2);
199: VecPointwiseDivide(asls->r1,asls->r1,asls->r2);
200: VecSet(tao->stepdirection,0.0);
201: VecISAXPY(tao->stepdirection, asls->fixed, 1.0,asls->r1);
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. */
208: TaoVecGetSubVec(asls->da, asls->free, tao->subset_type, 0.0, &asls->r1);
209: TaoVecGetSubVec(asls->ff, asls->free, tao->subset_type, 0.0, &asls->r2);
210: TaoVecGetSubVec(asls->db, asls->free, tao->subset_type, 1.0, &asls->r3);
211: VecPointwiseDivide(asls->r1,asls->r1, asls->r3);
212: VecPointwiseDivide(asls->r2,asls->r2, asls->r3);
214: /* r1 is the diagonal perturbation
215: r2 is the right hand side
216: r3 is no longer needed
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. */
222: MatMult(tao->jacobian, tao->stepdirection, asls->t1);
223: TaoVecGetSubVec(asls->t1,asls->free,tao->subset_type,0.0,&asls->r3);
224: VecAXPY(asls->r2, -1.0, asls->r3);
226: /* Calculate the reduced problem matrix and the direction */
227: TaoMatGetSubMat(tao->jacobian, asls->free, asls->w, tao->subset_type,&asls->J_sub);
228: if (tao->jacobian != tao->jacobian_pre) {
229: TaoMatGetSubMat(tao->jacobian_pre, asls->free, asls->w, tao->subset_type, &asls->Jpre_sub);
230: } else {
231: MatDestroy(&asls->Jpre_sub);
232: asls->Jpre_sub = asls->J_sub;
233: PetscObjectReference((PetscObject)(asls->Jpre_sub));
234: }
235: MatDiagonalSet(asls->J_sub, asls->r1,ADD_VALUES);
236: TaoVecGetSubVec(tao->stepdirection, asls->free, tao->subset_type, 0.0, &asls->dxfree);
237: VecSet(asls->dxfree, 0.0);
239: /* Calculate the reduced direction. (Really negative of Newton
240: direction. Therefore, rest of the code uses -d.) */
241: KSPReset(tao->ksp);
242: KSPSetOperators(tao->ksp, asls->J_sub, asls->Jpre_sub);
243: KSPSolve(tao->ksp, asls->r2, asls->dxfree);
244: KSPGetIterationNumber(tao->ksp,&tao->ksp_its);
245: tao->ksp_tot_its+=tao->ksp_its;
247: /* Add the direction in the free variables back into the real direction. */
248: VecISAXPY(tao->stepdirection, asls->free, 1.0,asls->dxfree);
250: /* Check the projected real direction for descent and if not, use the negative
251: gradient direction. */
252: VecCopy(tao->stepdirection, asls->w);
253: VecScale(asls->w, -1.0);
254: VecBoundGradientProjection(asls->w, tao->solution, tao->XL, tao->XU, asls->w);
255: VecNorm(asls->w, NORM_2, &normd);
256: VecDot(asls->w, asls->dpsi, &innerd);
258: if (innerd >= -asls->delta*PetscPowReal(normd, asls->rho)) {
259: PetscInfo(tao,"Gradient direction: %5.4e.\n", (double)innerd);
260: PetscInfo(tao, "Iteration %D: newton direction not descent\n", tao->niter);
261: VecCopy(asls->dpsi, tao->stepdirection);
262: VecDot(asls->dpsi, tao->stepdirection, &innerd);
263: }
265: VecScale(tao->stepdirection, -1.0);
266: innerd = -innerd;
268: /* We now have a correct descent direction. Apply a linesearch to
269: find the new iterate. */
270: TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);
271: TaoLineSearchApply(tao->linesearch, tao->solution, &psi,asls->dpsi, tao->stepdirection, &t, &ls_reason);
272: VecNorm(asls->dpsi, NORM_2, &ndpsi);
273: }
274: return 0;
275: }
277: /* ---------------------------------------------------------- */
278: /*MC
279: TAOASFLS - Active-set feasible linesearch algorithm for solving
280: complementarity constraints
282: Options Database Keys:
283: + -tao_ssls_delta - descent test fraction
284: - -tao_ssls_rho - descent test power
286: Level: beginner
287: M*/
288: PETSC_EXTERN PetscErrorCode TaoCreate_ASFLS(Tao tao)
289: {
290: TAO_SSLS *asls;
291: const char *armijo_type = TAOLINESEARCHARMIJO;
293: PetscNewLog(tao,&asls);
294: tao->data = (void*)asls;
295: tao->ops->solve = TaoSolve_ASFLS;
296: tao->ops->setup = TaoSetUp_ASFLS;
297: tao->ops->view = TaoView_SSLS;
298: tao->ops->setfromoptions = TaoSetFromOptions_SSLS;
299: tao->ops->destroy = TaoDestroy_ASFLS;
300: tao->subset_type = TAO_SUBSET_SUBVEC;
301: asls->delta = 1e-10;
302: asls->rho = 2.1;
303: asls->fixed = NULL;
304: asls->free = NULL;
305: asls->J_sub = NULL;
306: asls->Jpre_sub = NULL;
307: asls->w = NULL;
308: asls->r1 = NULL;
309: asls->r2 = NULL;
310: asls->r3 = NULL;
311: asls->t1 = NULL;
312: asls->t2 = NULL;
313: asls->dxfree = NULL;
314: asls->identifier = 1e-5;
316: TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);
317: PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);
318: TaoLineSearchSetType(tao->linesearch, armijo_type);
319: TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);
320: TaoLineSearchSetFromOptions(tao->linesearch);
322: KSPCreate(((PetscObject)tao)->comm, &tao->ksp);
323: PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);
324: KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);
325: KSPSetFromOptions(tao->ksp);
327: /* Override default settings (unless already changed) */
328: if (!tao->max_it_changed) tao->max_it = 2000;
329: if (!tao->max_funcs_changed) tao->max_funcs = 4000;
330: if (!tao->gttol_changed) tao->gttol = 0;
331: if (!tao->grtol_changed) tao->grtol = 0;
332: #if defined(PETSC_USE_REAL_SINGLE)
333: if (!tao->gatol_changed) tao->gatol = 1.0e-6;
334: if (!tao->fmin_changed) tao->fmin = 1.0e-4;
335: #else
336: if (!tao->gatol_changed) tao->gatol = 1.0e-16;
337: if (!tao->fmin_changed) tao->fmin = 1.0e-8;
338: #endif
339: return 0;
340: }