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