Actual source code: asils.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  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: */


 56: static PetscErrorCode TaoSetUp_ASILS(Tao tao)
 57: {
 58:   TAO_SSLS       *asls = (TAO_SSLS *)tao->data;

 62:   VecDuplicate(tao->solution,&tao->gradient);
 63:   VecDuplicate(tao->solution,&tao->stepdirection);
 64:   VecDuplicate(tao->solution,&asls->ff);
 65:   VecDuplicate(tao->solution,&asls->dpsi);
 66:   VecDuplicate(tao->solution,&asls->da);
 67:   VecDuplicate(tao->solution,&asls->db);
 68:   VecDuplicate(tao->solution,&asls->t1);
 69:   VecDuplicate(tao->solution,&asls->t2);
 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:   return(0);
 80: }

 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;

 89:   TaoComputeConstraints(tao, X, tao->constraints);
 90:   VecFischer(X,tao->constraints,tao->XL,tao->XU,asls->ff);
 91:   VecNorm(asls->ff,NORM_2,&asls->merit);
 92:   *fcn = 0.5*asls->merit*asls->merit;

 94:   TaoComputeJacobian(tao,tao->solution,tao->jacobian,tao->jacobian_pre);
 95:   MatDFischer(tao->jacobian, tao->solution, tao->constraints,tao->XL, tao->XU, asls->t1, asls->t2,asls->da, asls->db);
 96:   VecPointwiseMult(asls->t1, asls->ff, asls->db);
 97:   MatMultTranspose(tao->jacobian,asls->t1,G);
 98:   VecPointwiseMult(asls->t1, asls->ff, asls->da);
 99:   VecAXPY(G,1.0,asls->t1);
100:   return(0);
101: }

103: static PetscErrorCode TaoDestroy_ASILS(Tao tao)
104: {
105:   TAO_SSLS       *ssls = (TAO_SSLS *)tao->data;

109:   VecDestroy(&ssls->ff);
110:   VecDestroy(&ssls->dpsi);
111:   VecDestroy(&ssls->da);
112:   VecDestroy(&ssls->db);
113:   VecDestroy(&ssls->w);
114:   VecDestroy(&ssls->t1);
115:   VecDestroy(&ssls->t2);
116:   VecDestroy(&ssls->r1);
117:   VecDestroy(&ssls->r2);
118:   VecDestroy(&ssls->r3);
119:   VecDestroy(&ssls->dxfree);
120:   MatDestroy(&ssls->J_sub);
121:   MatDestroy(&ssls->Jpre_sub);
122:   ISDestroy(&ssls->fixed);
123:   ISDestroy(&ssls->free);
124:   PetscFree(tao->data);
125:   return(0);
126: }

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;

137:   /* Assume that Setup has been called!
138:      Set the structure for the Jacobian and create a linear solver. */

140:   TaoComputeVariableBounds(tao);
141:   TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_ASLS_FunctionGradient,tao);
142:   TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);

144:   /* Calculate the function value and fischer function value at the
145:      current iterate */
146:   TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,asls->dpsi);
147:   VecNorm(asls->dpsi,NORM_2,&ndpsi);

149:   tao->reason = TAO_CONTINUE_ITERATING;
150:   while (1) {
151:     /* Check the termination criteria */
152:     PetscInfo3(tao,"iter %D, merit: %g, ||dpsi||: %g\n",tao->niter, (double)asls->merit,  (double)ndpsi);
153:     TaoLogConvergenceHistory(tao,asls->merit,ndpsi,0.0,tao->ksp_its);
154:     TaoMonitor(tao,tao->niter,asls->merit,ndpsi,0.0,t);
155:     (*tao->ops->convergencetest)(tao,tao->cnvP);
156:     if (TAO_CONTINUE_ITERATING != tao->reason) break;
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:     PetscInfo1(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:     if (!asls->w && (tao->subset_type == TAO_SUBSET_MASK || tao->subset_type == TAO_SUBSET_MATRIXFREE)) {
228:       VecDuplicate(tao->solution, &asls->w);
229:     }
230:     TaoMatGetSubMat(tao->jacobian, asls->free, asls->w, tao->subset_type,&asls->J_sub);
231:     if (tao->jacobian != tao->jacobian_pre) {
232:       TaoMatGetSubMat(tao->jacobian_pre, asls->free, asls->w, tao->subset_type, &asls->Jpre_sub);
233:     } else {
234:       MatDestroy(&asls->Jpre_sub);
235:       asls->Jpre_sub = asls->J_sub;
236:       PetscObjectReference((PetscObject)(asls->Jpre_sub));
237:     }
238:     MatDiagonalSet(asls->J_sub, asls->r1,ADD_VALUES);
239:     TaoVecGetSubVec(tao->stepdirection, asls->free, tao->subset_type, 0.0, &asls->dxfree);
240:     VecSet(asls->dxfree, 0.0);

242:     /* Calculate the reduced direction.  (Really negative of Newton
243:        direction.  Therefore, rest of the code uses -d.) */
244:     KSPReset(tao->ksp);
245:     KSPSetOperators(tao->ksp, asls->J_sub, asls->Jpre_sub);
246:     KSPSolve(tao->ksp, asls->r2, asls->dxfree);
247:     KSPGetIterationNumber(tao->ksp,&tao->ksp_its);
248:     tao->ksp_tot_its+=tao->ksp_its;

250:     /* Add the direction in the free variables back into the real direction. */
251:     VecISAXPY(tao->stepdirection, asls->free, 1.0,asls->dxfree);

253:     /* Check the real direction for descent and if not, use the negative
254:        gradient direction. */
255:     VecNorm(tao->stepdirection, NORM_2, &normd);
256:     VecDot(tao->stepdirection, asls->dpsi, &innerd);

258:     if (innerd <= asls->delta*PetscPowReal(normd, asls->rho)) {
259:       PetscInfo1(tao,"Gradient direction: %5.4e.\n", (double)innerd);
260:       PetscInfo1(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:    TAOASILS - Active-set infeasible 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_ASILS(Tao tao)
289: {
290:   TAO_SSLS       *asls;
292:   const char     *armijo_type = TAOLINESEARCHARMIJO;

295:   PetscNewLog(tao,&asls);
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;

317:   asls->identifier = 1e-5;

319:   TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);
320:   PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);
321:   TaoLineSearchSetType(tao->linesearch, armijo_type);
322:   TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);
323:   TaoLineSearchSetFromOptions(tao->linesearch);

325:   KSPCreate(((PetscObject)tao)->comm, &tao->ksp);
326:   PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);
327:   KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);
328:   KSPSetFromOptions(tao->ksp);

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:   return(0);
343: }