Actual source code: viss.c

petsc-dev 2014-02-02
Report Typos and Errors
  2: #include <../src/snes/impls/vi/ss/vissimpl.h> /*I "petscsnes.h" I*/
  3: #include <../include/petsc-private/kspimpl.h>
  4: #include <../include/petsc-private/matimpl.h>
  5: #include <../include/petsc-private/dmimpl.h>


  8: /*
  9:   SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem.

 11:   Input Parameter:
 12: . phi - the semismooth function

 14:   Output Parameter:
 15: . merit - the merit function
 16: . phinorm - ||phi||

 18:   Notes:
 19:   The merit function for the mixed complementarity problem is defined as
 20:      merit = 0.5*phi^T*phi
 21: */
 24: static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal *merit,PetscReal *phinorm)
 25: {

 29:   VecNormBegin(phi,NORM_2,phinorm);
 30:   VecNormEnd(phi,NORM_2,phinorm);

 32:   *merit = 0.5*(*phinorm)*(*phinorm);
 33:   return(0);
 34: }

 36: PETSC_STATIC_INLINE PetscScalar Phi(PetscScalar a,PetscScalar b)
 37: {
 38:   return a + b - PetscSqrtScalar(a*a + b*b);
 39: }

 41: PETSC_STATIC_INLINE PetscScalar DPhi(PetscScalar a,PetscScalar b)
 42: {
 43:   if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return 1.0 - a/ PetscSqrtScalar(a*a + b*b);
 44:   else return .5;
 45: }

 47: /*
 48:    SNESVIComputeFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form.

 50:    Input Parameters:
 51: .  snes - the SNES context
 52: .  x - current iterate
 53: .  functx - user defined function context

 55:    Output Parameters:
 56: .  phi - Semismooth function

 58: */
 61: static PetscErrorCode SNESVIComputeFunction(SNES snes,Vec X,Vec phi,void *functx)
 62: {
 63:   PetscErrorCode    ierr;
 64:   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS*)snes->data;
 65:   Vec               Xl  = snes->xl,Xu = snes->xu,F = snes->vec_func;
 66:   PetscScalar       *phi_arr,*x_arr,*f_arr,*l,*u;
 67:   PetscInt          i,nlocal;

 70:   (*vi->computeuserfunction)(snes,X,F,functx);
 71:   VecGetLocalSize(X,&nlocal);
 72:   VecGetArray(X,&x_arr);
 73:   VecGetArray(F,&f_arr);
 74:   VecGetArray(Xl,&l);
 75:   VecGetArray(Xu,&u);
 76:   VecGetArray(phi,&phi_arr);

 78:   for (i=0; i < nlocal; i++) {
 79:     if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */
 80:       phi_arr[i] = f_arr[i];
 81:     } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) {                      /* upper bound on variable only */
 82:       phi_arr[i] = -Phi(u[i] - x_arr[i],-f_arr[i]);
 83:     } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) {                       /* lower bound on variable only */
 84:       phi_arr[i] = Phi(x_arr[i] - l[i],f_arr[i]);
 85:     } else if (l[i] == u[i]) {
 86:       phi_arr[i] = l[i] - x_arr[i];
 87:     } else {                                                /* both bounds on variable */
 88:       phi_arr[i] = Phi(x_arr[i] - l[i],-Phi(u[i] - x_arr[i],-f_arr[i]));
 89:     }
 90:   }

 92:   VecRestoreArray(X,&x_arr);
 93:   VecRestoreArray(F,&f_arr);
 94:   VecRestoreArray(Xl,&l);
 95:   VecRestoreArray(Xu,&u);
 96:   VecRestoreArray(phi,&phi_arr);
 97:   return(0);
 98: }

100: /*
101:    SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the
102:                                           the semismooth jacobian.
103: */
106: PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes,Vec X,Vec F,Mat jac,Vec Da,Vec Db)
107: {
109:   PetscScalar    *l,*u,*x,*f,*da,*db,da1,da2,db1,db2;
110:   PetscInt       i,nlocal;

113:   VecGetArray(X,&x);
114:   VecGetArray(F,&f);
115:   VecGetArray(snes->xl,&l);
116:   VecGetArray(snes->xu,&u);
117:   VecGetArray(Da,&da);
118:   VecGetArray(Db,&db);
119:   VecGetLocalSize(X,&nlocal);

121:   for (i=0; i< nlocal; i++) {
122:     if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */
123:       da[i] = 0;
124:       db[i] = 1;
125:     } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) {                     /* upper bound on variable only */
126:       da[i] = DPhi(u[i] - x[i], -f[i]);
127:       db[i] = DPhi(-f[i],u[i] - x[i]);
128:     } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) {                      /* lower bound on variable only */
129:       da[i] = DPhi(x[i] - l[i], f[i]);
130:       db[i] = DPhi(f[i],x[i] - l[i]);
131:     } else if (l[i] == u[i]) {                              /* fixed variable */
132:       da[i] = 1;
133:       db[i] = 0;
134:     } else {                                                /* upper and lower bounds on variable */
135:       da1   = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i]));
136:       db1   = DPhi(-Phi(u[i] - x[i], -f[i]),x[i] - l[i]);
137:       da2   = DPhi(u[i] - x[i], -f[i]);
138:       db2   = DPhi(-f[i],u[i] - x[i]);
139:       da[i] = da1 + db1*da2;
140:       db[i] = db1*db2;
141:     }
142:   }

144:   VecRestoreArray(X,&x);
145:   VecRestoreArray(F,&f);
146:   VecRestoreArray(snes->xl,&l);
147:   VecRestoreArray(snes->xu,&u);
148:   VecRestoreArray(Da,&da);
149:   VecRestoreArray(Db,&db);
150:   return(0);
151: }

153: /*
154:    SNESVIComputeJacobian - Computes the jacobian of the semismooth function.The Jacobian for the semismooth function is an element of the B-subdifferential of the Fischer-Burmeister function for complementarity problems.

156:    Input Parameters:
157: .  Da       - Diagonal shift vector for the semismooth jacobian.
158: .  Db       - Row scaling vector for the semismooth jacobian.

160:    Output Parameters:
161: .  jac      - semismooth jacobian
162: .  jac_pre  - optional preconditioning matrix

164:    Notes:
165:    The semismooth jacobian matrix is given by
166:    jac = Da + Db*jacfun
167:    where Db is the row scaling matrix stored as a vector,
168:          Da is the diagonal perturbation matrix stored as a vector
169:    and   jacfun is the jacobian of the original nonlinear function.
170: */
173: PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre,Vec Da, Vec Db)
174: {

177:   /* Do row scaling  and add diagonal perturbation */
178:   MatDiagonalScale(jac,Db,NULL);
179:   MatDiagonalSet(jac,Da,ADD_VALUES);
180:   if (jac != jac_pre) { /* If jac and jac_pre are different */
181:     MatDiagonalScale(jac_pre,Db,NULL);
182:     MatDiagonalSet(jac_pre,Da,ADD_VALUES);
183:   }
184:   return(0);
185: }

187: /*
188:    SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi.

190:    Input Parameters:
191:    phi - semismooth function.
192:    H   - semismooth jacobian

194:    Output Parameters:
195:    dpsi - merit function gradient

197:    Notes:
198:   The merit function gradient is computed as follows
199:         dpsi = H^T*phi
200: */
203: PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi)
204: {

208:   MatMultTranspose(H,phi,dpsi);
209:   return(0);
210: }



214: /*
215:    SNESSolve_VINEWTONSSLS - Solves the complementarity problem with a semismooth Newton
216:    method using a line search.

218:    Input Parameters:
219: .  snes - the SNES context

221:    Output Parameter:
222: .  outits - number of iterations until termination

224:    Application Interface Routine: SNESSolve()

226:    Notes:
227:    This implements essentially a semismooth Newton method with a
228:    line search. The default line search does not do any line seach
229:    but rather takes a full newton step.

231:    Developer Note: the code in this file should be slightly modified so that this routine need not exist and the SNESSolve_NEWTONLS() routine is called directly with the appropriate wrapped function and Jacobian evaluations

233: */
236: PetscErrorCode SNESSolve_VINEWTONSSLS(SNES snes)
237: {
238:   SNES_VINEWTONSSLS  *vi = (SNES_VINEWTONSSLS*)snes->data;
239:   PetscErrorCode     ierr;
240:   PetscInt           maxits,i,lits;
241:   PetscBool          lssucceed;
242:   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
243:   PetscReal          gnorm,xnorm=0,ynorm;
244:   Vec                Y,X,F;
245:   KSPConvergedReason kspreason;
246:   DM                 dm;
247:   DMSNES             sdm;

250:   SNESGetDM(snes,&dm);
251:   DMGetDMSNES(dm,&sdm);

253:   vi->computeuserfunction   = sdm->ops->computefunction;
254:   sdm->ops->computefunction = SNESVIComputeFunction;

256:   snes->numFailures            = 0;
257:   snes->numLinearSolveFailures = 0;
258:   snes->reason                 = SNES_CONVERGED_ITERATING;

260:   maxits = snes->max_its;               /* maximum number of iterations */
261:   X      = snes->vec_sol;               /* solution vector */
262:   F      = snes->vec_func;              /* residual vector */
263:   Y      = snes->work[0];               /* work vectors */

265:   PetscObjectSAWsTakeAccess((PetscObject)snes);
266:   snes->iter = 0;
267:   snes->norm = 0.0;
268:   PetscObjectSAWsGrantAccess((PetscObject)snes);

270:   SNESVIProjectOntoBounds(snes,X);
271:   SNESComputeFunction(snes,X,vi->phi);
272:   if (snes->domainerror) {
273:     snes->reason              = SNES_DIVERGED_FUNCTION_DOMAIN;
274:     sdm->ops->computefunction = vi->computeuserfunction;
275:     return(0);
276:   }
277:   /* Compute Merit function */
278:   SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);

280:   VecNormBegin(X,NORM_2,&xnorm);        /* xnorm <- ||x||  */
281:   VecNormEnd(X,NORM_2,&xnorm);
282:   if (PetscIsInfOrNanReal(vi->merit)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");

284:   PetscObjectSAWsTakeAccess((PetscObject)snes);
285:   snes->norm = vi->phinorm;
286:   PetscObjectSAWsGrantAccess((PetscObject)snes);
287:   SNESLogConvergenceHistory(snes,vi->phinorm,0);
288:   SNESMonitor(snes,0,vi->phinorm);

290:   /* test convergence */
291:   (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);
292:   if (snes->reason) {
293:     sdm->ops->computefunction = vi->computeuserfunction;
294:     return(0);
295:   }

297:   for (i=0; i<maxits; i++) {

299:     /* Call general purpose update function */
300:     if (snes->ops->update) {
301:       (*snes->ops->update)(snes, snes->iter);
302:     }

304:     /* Solve J Y = Phi, where J is the semismooth jacobian */

306:     /* Get the jacobian -- note that the function must be the original function for snes_fd and snes_fd_color to work for this*/
307:     sdm->ops->computefunction = vi->computeuserfunction;
308:     SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
309:     sdm->ops->computefunction = SNESVIComputeFunction;

311:     /* Get the diagonal shift and row scaling vectors */
312:     SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);
313:     /* Compute the semismooth jacobian */
314:     SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);
315:     /* Compute the merit function gradient */
316:     SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);
317:     KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);
318:     KSPSolve(snes->ksp,vi->phi,Y);
319:     KSPGetConvergedReason(snes->ksp,&kspreason);

321:     if (kspreason < 0) {
322:       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
323:         PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);
324:         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
325:         break;
326:       }
327:     }
328:     KSPGetIterationNumber(snes->ksp,&lits);
329:     snes->linear_its += lits;
330:     PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);
331:     /*
332:     if (snes->ops->precheck) {
333:       PetscBool changed_y = PETSC_FALSE;
334:       (*snes->ops->precheck)(snes,X,Y,snes->precheck,&changed_y);
335:     }

337:     if (PetscLogPrintInfo) {
338:       SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);
339:     }
340:     */
341:     /* Compute a (scaled) negative update in the line search routine:
342:          Y <- X - lambda*Y
343:        and evaluate G = function(Y) (depends on the line search).
344:     */
345:     VecCopy(Y,snes->vec_sol_update);
346:     ynorm = 1; gnorm = vi->phinorm;
347:     SNESLineSearchApply(snes->linesearch, X, vi->phi, &gnorm, Y);
348:     SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);
349:     SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm);
350:     PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)vi->phinorm,(double)gnorm,(double)ynorm,(int)lssucceed);
351:     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
352:     if (snes->domainerror) {
353:       snes->reason              = SNES_DIVERGED_FUNCTION_DOMAIN;
354:       sdm->ops->computefunction = vi->computeuserfunction;
355:       return(0);
356:     }
357:     if (!lssucceed) {
358:       if (++snes->numFailures >= snes->maxFailures) {
359:         PetscBool ismin;
360:         snes->reason = SNES_DIVERGED_LINE_SEARCH;
361:         SNESVICheckLocalMin_Private(snes,snes->jacobian,vi->phi,X,gnorm,&ismin);
362:         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
363:         break;
364:       }
365:     }
366:     /* Update function and solution vectors */
367:     vi->phinorm = gnorm;
368:     vi->merit   = 0.5*vi->phinorm*vi->phinorm;
369:     /* Monitor convergence */
370:     PetscObjectSAWsTakeAccess((PetscObject)snes);
371:     snes->iter = i+1;
372:     snes->norm = vi->phinorm;
373:     PetscObjectSAWsGrantAccess((PetscObject)snes);
374:     SNESLogConvergenceHistory(snes,snes->norm,lits);
375:     SNESMonitor(snes,snes->iter,snes->norm);
376:     /* Test for convergence, xnorm = || X || */
377:     if (snes->ops->converged != SNESConvergedSkip) { VecNorm(X,NORM_2,&xnorm); }
378:     (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);
379:     if (snes->reason) break;
380:   }
381:   if (i == maxits) {
382:     PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);
383:     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
384:   }
385:   sdm->ops->computefunction = vi->computeuserfunction;
386:   return(0);
387: }

389: /* -------------------------------------------------------------------------- */
390: /*
391:    SNESSetUp_VINEWTONSSLS - Sets up the internal data structures for the later use
392:    of the SNES nonlinear solver.

394:    Input Parameter:
395: .  snes - the SNES context
396: .  x - the solution vector

398:    Application Interface Routine: SNESSetUp()

400:    Notes:
401:    For basic use of the SNES solvers, the user need not explicitly call
402:    SNESSetUp(), since these actions will automatically occur during
403:    the call to SNESSolve().
404:  */
407: PetscErrorCode SNESSetUp_VINEWTONSSLS(SNES snes)
408: {
409:   PetscErrorCode    ierr;
410:   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS*) snes->data;

413:   SNESSetUp_VI(snes);
414:   VecDuplicate(snes->vec_sol, &vi->dpsi);
415:   VecDuplicate(snes->vec_sol, &vi->phi);
416:   VecDuplicate(snes->vec_sol, &vi->Da);
417:   VecDuplicate(snes->vec_sol, &vi->Db);
418:   VecDuplicate(snes->vec_sol, &vi->z);
419:   VecDuplicate(snes->vec_sol, &vi->t);
420:   return(0);
421: }
422: /* -------------------------------------------------------------------------- */
425: PetscErrorCode SNESReset_VINEWTONSSLS(SNES snes)
426: {
427:   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS*) snes->data;
428:   PetscErrorCode    ierr;

431:   SNESReset_VI(snes);
432:   VecDestroy(&vi->dpsi);
433:   VecDestroy(&vi->phi);
434:   VecDestroy(&vi->Da);
435:   VecDestroy(&vi->Db);
436:   VecDestroy(&vi->z);
437:   VecDestroy(&vi->t);
438:   return(0);
439: }

441: /* -------------------------------------------------------------------------- */
442: /*
443:    SNESSetFromOptions_VINEWTONSSLS - Sets various parameters for the SNESVI method.

445:    Input Parameter:
446: .  snes - the SNES context

448:    Application Interface Routine: SNESSetFromOptions()
449: */
452: static PetscErrorCode SNESSetFromOptions_VINEWTONSSLS(SNES snes)
453: {
455:   SNESLineSearch linesearch;

458:   SNESSetFromOptions_VI(snes);
459:   PetscOptionsHead("SNES semismooth method options");
460:   PetscOptionsTail();
461:   /* set up the default line search */
462:   if (!snes->linesearch) {
463:     SNESGetLineSearch(snes, &linesearch);
464:     SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);
465:     SNESLineSearchBTSetAlpha(linesearch, 0.0);
466:   }
467:   return(0);
468: }


471: /* -------------------------------------------------------------------------- */
472: /*MC
473:       SNESVINEWTONSSLS - Semi-smooth solver for variational inequalities based on Newton's method

475:    Options Database:
476: +   -snes_vi_type <ss,rs,rsaug> a semi-smooth solver, a reduced space active set method, and a reduced space active set method that does not eliminate the active constraints from the Jacobian instead augments the Jacobian with additional variables that enforce the constraints
477: -   -snes_vi_monitor - prints the number of active constraints at each iteration.

479:    Level: beginner

481:    References:
482:    - T. S. Munson, F. Facchinei, M. C. Ferris, A. Fischer, and C. Kanzow. The semismooth
483:      algorithm for large scale complementarity problems. INFORMS Journal on Computing, 13 (2001).

485: .seealso:  SNESVISetVariableBounds(), SNESVISetComputeVariableBounds(), SNESCreate(), SNES, SNESSetType(), SNESVINEWTONRSLS, SNESVINEWTONSSLS, SNESNEWTONTR, SNESLineSearchSet(),
486:            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
487:            SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()

489: M*/
492: PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONSSLS(SNES snes)
493: {
494:   PetscErrorCode    ierr;
495:   SNES_VINEWTONSSLS *vi;

498:   snes->ops->reset          = SNESReset_VINEWTONSSLS;
499:   snes->ops->setup          = SNESSetUp_VINEWTONSSLS;
500:   snes->ops->solve          = SNESSolve_VINEWTONSSLS;
501:   snes->ops->destroy        = SNESDestroy_VI;
502:   snes->ops->setfromoptions = SNESSetFromOptions_VINEWTONSSLS;
503:   snes->ops->view           = NULL;

505:   snes->usesksp = PETSC_TRUE;
506:   snes->usespc  = PETSC_FALSE;

508:   PetscNewLog(snes,&vi);
509:   snes->data = (void*)vi;

511:   PetscObjectComposeFunction((PetscObject)snes,"SNESVISetVariableBounds_C",SNESVISetVariableBounds_VI);
512:   PetscObjectComposeFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",SNESVISetComputeVariableBounds_VI);
513:   return(0);
514: }