Actual source code: virs.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: #include <../src/snes/impls/vi/rs/virsimpl.h> /*I "petscsnes.h" I*/
  3: #include <petsc/private/kspimpl.h>
  4: #include <petsc/private/matimpl.h>
  5: #include <petsc/private/dmimpl.h>
  6: #include <petsc/private/vecimpl.h>

 10: /*
 11:    SNESVIGetInactiveSet - Gets the global indices for the inactive set variables (these correspond to the degrees of freedom the linear
 12:      system is solved on)

 14:    Input parameter
 15: .  snes - the SNES context

 17:    Output parameter
 18: .  inact - inactive set index set

 20:  */
 21: PetscErrorCode SNESVIGetInactiveSet(SNES snes,IS *inact)
 22: {
 23:   SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS*)snes->data;

 26:   *inact = vi->IS_inact;
 27:   return(0);
 28: }

 30: /*
 31:     Provides a wrapper to a DM to allow it to be used to generated the interpolation/restriction from the DM for the smaller matrices and vectors
 32:   defined by the reduced space method.

 34:     Simple calls the regular DM interpolation and restricts it to operation on the variables not associated with active constraints.

 36: <*/
 37: typedef struct {
 38:   PetscInt n;                                              /* size of vectors in the reduced DM space */
 39:   IS       inactive;

 41:   PetscErrorCode (*createinterpolation)(DM,DM,Mat*,Vec*);  /* DM's original routines */
 42:   PetscErrorCode (*coarsen)(DM, MPI_Comm, DM*);
 43:   PetscErrorCode (*createglobalvector)(DM,Vec*);

 45:   DM dm;                                                  /* when destroying this object we need to reset the above function into the base DM */
 46: } DM_SNESVI;

 50: /*
 51:      DMCreateGlobalVector_SNESVI - Creates global vector of the size of the reduced space

 53: */
 54: PetscErrorCode  DMCreateGlobalVector_SNESVI(DM dm,Vec *vec)
 55: {
 57:   PetscContainer isnes;
 58:   DM_SNESVI      *dmsnesvi;

 61:   PetscObjectQuery((PetscObject)dm,"VI",(PetscObject*)&isnes);
 62:   if (!isnes) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Composed SNES is missing");
 63:   PetscContainerGetPointer(isnes,(void**)&dmsnesvi);
 64:   VecCreateMPI(PetscObjectComm((PetscObject)dm),dmsnesvi->n,PETSC_DETERMINE,vec);
 65:   return(0);
 66: }

 70: /*
 71:      DMCreateInterpolation_SNESVI - Modifieds the interpolation obtained from the DM by removing all rows and columns associated with active constraints.

 73: */
 74: PetscErrorCode  DMCreateInterpolation_SNESVI(DM dm1,DM dm2,Mat *mat,Vec *vec)
 75: {
 77:   PetscContainer isnes;
 78:   DM_SNESVI      *dmsnesvi1,*dmsnesvi2;
 79:   Mat            interp;

 82:   PetscObjectQuery((PetscObject)dm1,"VI",(PetscObject*)&isnes);
 83:   if (!isnes) SETERRQ(PetscObjectComm((PetscObject)dm1),PETSC_ERR_PLIB,"Composed VI data structure is missing");
 84:   PetscContainerGetPointer(isnes,(void**)&dmsnesvi1);
 85:   PetscObjectQuery((PetscObject)dm2,"VI",(PetscObject*)&isnes);
 86:   if (!isnes) SETERRQ(PetscObjectComm((PetscObject)dm2),PETSC_ERR_PLIB,"Composed VI data structure is missing");
 87:   PetscContainerGetPointer(isnes,(void**)&dmsnesvi2);

 89:   (*dmsnesvi1->createinterpolation)(dm1,dm2,&interp,NULL);
 90:   MatGetSubMatrix(interp,dmsnesvi2->inactive,dmsnesvi1->inactive,MAT_INITIAL_MATRIX,mat);
 91:   MatDestroy(&interp);
 92:   *vec = 0;
 93:   return(0);
 94: }

 96: extern PetscErrorCode  DMSetVI(DM,IS);

100: /*
101:      DMCoarsen_SNESVI - Computes the regular coarsened DM then computes additional information about its inactive set

103: */
104: PetscErrorCode  DMCoarsen_SNESVI(DM dm1,MPI_Comm comm,DM *dm2)
105: {
107:   PetscContainer isnes;
108:   DM_SNESVI      *dmsnesvi1;
109:   Vec            finemarked,coarsemarked;
110:   IS             inactive;
111:   Mat            inject;
112:   const PetscInt *index;
113:   PetscInt       n,k,cnt = 0,rstart,*coarseindex;
114:   PetscScalar    *marked;

117:   PetscObjectQuery((PetscObject)dm1,"VI",(PetscObject*)&isnes);
118:   if (!isnes) SETERRQ(PetscObjectComm((PetscObject)dm1),PETSC_ERR_PLIB,"Composed VI data structure is missing");
119:   PetscContainerGetPointer(isnes,(void**)&dmsnesvi1);

121:   /* get the original coarsen */
122:   (*dmsnesvi1->coarsen)(dm1,comm,dm2);

124:   /* not sure why this extra reference is needed, but without the dm2 disappears too early */
125:   /* Updating the KSPCreateVecs() to avoid using DMGetGlobalVector() when matrix is available removes the need for this reference? */
126:   /*  PetscObjectReference((PetscObject)*dm2);*/

128:   /* need to set back global vectors in order to use the original injection */
129:   DMClearGlobalVectors(dm1);

131:   dm1->ops->createglobalvector = dmsnesvi1->createglobalvector;

133:   DMCreateGlobalVector(dm1,&finemarked);
134:   DMCreateGlobalVector(*dm2,&coarsemarked);

136:   /*
137:      fill finemarked with locations of inactive points
138:   */
139:   ISGetIndices(dmsnesvi1->inactive,&index);
140:   ISGetLocalSize(dmsnesvi1->inactive,&n);
141:   VecSet(finemarked,0.0);
142:   for (k=0; k<n; k++) {
143:     VecSetValue(finemarked,index[k],1.0,INSERT_VALUES);
144:   }
145:   VecAssemblyBegin(finemarked);
146:   VecAssemblyEnd(finemarked);

148:   DMCreateInjection(*dm2,dm1,&inject);
149:   MatRestrict(inject,finemarked,coarsemarked);
150:   MatDestroy(&inject);

152:   /*
153:      create index set list of coarse inactive points from coarsemarked
154:   */
155:   VecGetLocalSize(coarsemarked,&n);
156:   VecGetOwnershipRange(coarsemarked,&rstart,NULL);
157:   VecGetArray(coarsemarked,&marked);
158:   for (k=0; k<n; k++) {
159:     if (marked[k] != 0.0) cnt++;
160:   }
161:   PetscMalloc1(cnt,&coarseindex);
162:   cnt  = 0;
163:   for (k=0; k<n; k++) {
164:     if (marked[k] != 0.0) coarseindex[cnt++] = k + rstart;
165:   }
166:   VecRestoreArray(coarsemarked,&marked);
167:   ISCreateGeneral(PetscObjectComm((PetscObject)coarsemarked),cnt,coarseindex,PETSC_OWN_POINTER,&inactive);

169:   DMClearGlobalVectors(dm1);

171:   dm1->ops->createglobalvector = DMCreateGlobalVector_SNESVI;

173:   DMSetVI(*dm2,inactive);

175:   VecDestroy(&finemarked);
176:   VecDestroy(&coarsemarked);
177:   ISDestroy(&inactive);
178:   return(0);
179: }

183: PetscErrorCode DMDestroy_SNESVI(DM_SNESVI *dmsnesvi)
184: {

188:   /* reset the base methods in the DM object that were changed when the DM_SNESVI was reset */
189:   dmsnesvi->dm->ops->createinterpolation = dmsnesvi->createinterpolation;
190:   dmsnesvi->dm->ops->coarsen             = dmsnesvi->coarsen;
191:   dmsnesvi->dm->ops->createglobalvector  = dmsnesvi->createglobalvector;
192:   /* need to clear out this vectors because some of them may not have a reference to the DM
193:     but they are counted as having references to the DM in DMDestroy() */
194:   DMClearGlobalVectors(dmsnesvi->dm);

196:   ISDestroy(&dmsnesvi->inactive);
197:   PetscFree(dmsnesvi);
198:   return(0);
199: }

203: /*
204:      DMSetVI - Marks a DM as associated with a VI problem. This causes the interpolation/restriction operators to
205:                be restricted to only those variables NOT associated with active constraints.

207: */
208: PetscErrorCode  DMSetVI(DM dm,IS inactive)
209: {
211:   PetscContainer isnes;
212:   DM_SNESVI      *dmsnesvi;

215:   if (!dm) return(0);

217:   PetscObjectReference((PetscObject)inactive);

219:   PetscObjectQuery((PetscObject)dm,"VI",(PetscObject*)&isnes);
220:   if (!isnes) {
221:     PetscContainerCreate(PetscObjectComm((PetscObject)dm),&isnes);
222:     PetscContainerSetUserDestroy(isnes,(PetscErrorCode (*)(void*))DMDestroy_SNESVI);
223:     PetscNew(&dmsnesvi);
224:     PetscContainerSetPointer(isnes,(void*)dmsnesvi);
225:     PetscObjectCompose((PetscObject)dm,"VI",(PetscObject)isnes);
226:     PetscContainerDestroy(&isnes);

228:     dmsnesvi->createinterpolation = dm->ops->createinterpolation;
229:     dm->ops->createinterpolation  = DMCreateInterpolation_SNESVI;
230:     dmsnesvi->coarsen             = dm->ops->coarsen;
231:     dm->ops->coarsen              = DMCoarsen_SNESVI;
232:     dmsnesvi->createglobalvector  = dm->ops->createglobalvector;
233:     dm->ops->createglobalvector   = DMCreateGlobalVector_SNESVI;
234:   } else {
235:     PetscContainerGetPointer(isnes,(void**)&dmsnesvi);
236:     ISDestroy(&dmsnesvi->inactive);
237:   }
238:   DMClearGlobalVectors(dm);
239:   ISGetLocalSize(inactive,&dmsnesvi->n);

241:   dmsnesvi->inactive = inactive;
242:   dmsnesvi->dm       = dm;
243:   return(0);
244: }

248: /*
249:      DMDestroyVI - Frees the DM_SNESVI object contained in the DM
250:          - also resets the function pointers in the DM for createinterpolation() etc to use the original DM
251: */
252: PetscErrorCode  DMDestroyVI(DM dm)
253: {

257:   if (!dm) return(0);
258:   PetscObjectCompose((PetscObject)dm,"VI",(PetscObject)NULL);
259:   return(0);
260: }

262: /* --------------------------------------------------------------------------------------------------------*/




269: PetscErrorCode SNESCreateIndexSets_VINEWTONRSLS(SNES snes,Vec X,Vec F,IS *ISact,IS *ISinact)
270: {

274:   SNESVIGetActiveSetIS(snes,X,F,ISact);
275:   ISComplement(*ISact,X->map->rstart,X->map->rend,ISinact);
276:   return(0);
277: }

279: /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */
282: PetscErrorCode SNESCreateSubVectors_VINEWTONRSLS(SNES snes,PetscInt n,Vec *newv)
283: {
285:   Vec            v;

288:   VecCreate(PetscObjectComm((PetscObject)snes),&v);
289:   VecSetSizes(v,n,PETSC_DECIDE);
290:   VecSetType(v,VECSTANDARD);
291:   *newv = v;
292:   return(0);
293: }

295: /* Resets the snes PC and KSP when the active set sizes change */
298: PetscErrorCode SNESVIResetPCandKSP(SNES snes,Mat Amat,Mat Pmat)
299: {
301:   KSP            snesksp;

304:   SNESGetKSP(snes,&snesksp);
305:   KSPReset(snesksp);

307:   /*
308:   KSP                    kspnew;
309:   PC                     pcnew;
310:   const MatSolverPackage stype;


313:   KSPCreate(PetscObjectComm((PetscObject)snes),&kspnew);
314:   kspnew->pc_side = snesksp->pc_side;
315:   kspnew->rtol    = snesksp->rtol;
316:   kspnew->abstol    = snesksp->abstol;
317:   kspnew->max_it  = snesksp->max_it;
318:   KSPSetType(kspnew,((PetscObject)snesksp)->type_name);
319:   KSPGetPC(kspnew,&pcnew);
320:   PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name);
321:   PCSetOperators(kspnew->pc,Amat,Pmat);
322:   PCFactorGetMatSolverPackage(snesksp->pc,&stype);
323:   PCFactorSetMatSolverPackage(kspnew->pc,stype);
324:   KSPDestroy(&snesksp);
325:   snes->ksp = kspnew;
326:   PetscLogObjectParent((PetscObject)snes,(PetscObject)kspnew);
327:    KSPSetFromOptions(kspnew);*/
328:   return(0);
329: }

331: /* Variational Inequality solver using reduce space method. No semismooth algorithm is
332:    implemented in this algorithm. It basically identifies the active constraints and does
333:    a linear solve on the other variables (those not associated with the active constraints). */
336: PetscErrorCode SNESSolve_VINEWTONRSLS(SNES snes)
337: {
338:   SNES_VINEWTONRSLS    *vi = (SNES_VINEWTONRSLS*)snes->data;
339:   PetscErrorCode       ierr;
340:   PetscInt             maxits,i,lits;
341:   SNESLineSearchReason lssucceed;
342:   PetscReal            fnorm,gnorm,xnorm=0,ynorm;
343:   Vec                  Y,X,F;
344:   KSPConvergedReason   kspreason;
345:   KSP                  ksp;
346:   PC                   pc;

349:   /* Multigrid must use Galerkin for coarse grids with active set/reduced space methods; cannot rediscretize on coarser grids*/
350:   SNESGetKSP(snes,&ksp);
351:   KSPGetPC(ksp,&pc);
352:   PCMGSetGalerkin(pc,PETSC_TRUE);

354:   snes->numFailures            = 0;
355:   snes->numLinearSolveFailures = 0;
356:   snes->reason                 = SNES_CONVERGED_ITERATING;

358:   maxits = snes->max_its;               /* maximum number of iterations */
359:   X      = snes->vec_sol;               /* solution vector */
360:   F      = snes->vec_func;              /* residual vector */
361:   Y      = snes->work[0];               /* work vectors */

363:   SNESLineSearchSetVIFunctions(snes->linesearch, SNESVIProjectOntoBounds, SNESVIComputeInactiveSetFnorm);
364:   SNESLineSearchSetVecs(snes->linesearch, X, NULL, NULL, NULL, NULL);
365:   SNESLineSearchSetUp(snes->linesearch);

367:   PetscObjectSAWsTakeAccess((PetscObject)snes);
368:   snes->iter = 0;
369:   snes->norm = 0.0;
370:   PetscObjectSAWsGrantAccess((PetscObject)snes);

372:   SNESVIProjectOntoBounds(snes,X);
373:   SNESComputeFunction(snes,X,F);
374:   SNESVIComputeInactiveSetFnorm(snes,F,X,&fnorm);
375:   VecNormBegin(X,NORM_2,&xnorm);        /* xnorm <- ||x||  */
376:   VecNormEnd(X,NORM_2,&xnorm);
377:   SNESCheckFunctionNorm(snes,fnorm);
378:   PetscObjectSAWsTakeAccess((PetscObject)snes);
379:   snes->norm = fnorm;
380:   PetscObjectSAWsGrantAccess((PetscObject)snes);
381:   SNESLogConvergenceHistory(snes,fnorm,0);
382:   SNESMonitor(snes,0,fnorm);

384:   /* test convergence */
385:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
386:   if (snes->reason) return(0);


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

391:     IS         IS_act; /* _act -> active set _inact -> inactive set */
392:     IS         IS_redact; /* redundant active set */
393:     VecScatter scat_act,scat_inact;
394:     PetscInt   nis_act,nis_inact;
395:     Vec        Y_act,Y_inact,F_inact;
396:     Mat        jac_inact_inact,prejac_inact_inact;
397:     PetscBool  isequal;

399:     /* Call general purpose update function */
400:     if (snes->ops->update) {
401:       (*snes->ops->update)(snes, snes->iter);
402:     }
403:     SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);


406:     /* Create active and inactive index sets */

408:     /*original
409:     SNESVICreateIndexSets_RS(snes,X,F,&IS_act,&vi->IS_inact);
410:      */
411:     SNESVIGetActiveSetIS(snes,X,F,&IS_act);

413:     if (vi->checkredundancy) {
414:       (*vi->checkredundancy)(snes,IS_act,&IS_redact,vi->ctxP);
415:       if (IS_redact) {
416:         ISSort(IS_redact);
417:         ISComplement(IS_redact,X->map->rstart,X->map->rend,&vi->IS_inact);
418:         ISDestroy(&IS_redact);
419:       } else {
420:         ISComplement(IS_act,X->map->rstart,X->map->rend,&vi->IS_inact);
421:       }
422:     } else {
423:       ISComplement(IS_act,X->map->rstart,X->map->rend,&vi->IS_inact);
424:     }


427:     /* Create inactive set submatrix */
428:     MatGetSubMatrix(snes->jacobian,vi->IS_inact,vi->IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);

430:     if (0) {                    /* Dead code (temporary developer hack) */
431:       IS keptrows;
432:       MatFindNonzeroRows(jac_inact_inact,&keptrows);
433:       if (keptrows) {
434:         PetscInt       cnt,*nrows,k;
435:         const PetscInt *krows,*inact;
436:         PetscInt       rstart=jac_inact_inact->rmap->rstart;

438:         MatDestroy(&jac_inact_inact);
439:         ISDestroy(&IS_act);

441:         ISGetLocalSize(keptrows,&cnt);
442:         ISGetIndices(keptrows,&krows);
443:         ISGetIndices(vi->IS_inact,&inact);
444:         PetscMalloc1(cnt,&nrows);
445:         for (k=0; k<cnt; k++) nrows[k] = inact[krows[k]-rstart];
446:         ISRestoreIndices(keptrows,&krows);
447:         ISRestoreIndices(vi->IS_inact,&inact);
448:         ISDestroy(&keptrows);
449:         ISDestroy(&vi->IS_inact);

451:         ISCreateGeneral(PetscObjectComm((PetscObject)snes),cnt,nrows,PETSC_OWN_POINTER,&vi->IS_inact);
452:         ISComplement(vi->IS_inact,F->map->rstart,F->map->rend,&IS_act);
453:         MatGetSubMatrix(snes->jacobian,vi->IS_inact,vi->IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);
454:       }
455:     }
456:     DMSetVI(snes->dm,vi->IS_inact);
457:     /* remove later */

459:     /*
460:     VecView(vi->xu,PETSC_VIEWER_BINARY_(((PetscObject)(vi->xu))->comm));
461:     VecView(vi->xl,PETSC_VIEWER_BINARY_(((PetscObject)(vi->xl))->comm));
462:     VecView(X,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)X)));
463:     VecView(F,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)F)));
464:     ISView(vi->IS_inact,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)vi->IS_inact)));
465:      */

467:     /* Get sizes of active and inactive sets */
468:     ISGetLocalSize(IS_act,&nis_act);
469:     ISGetLocalSize(vi->IS_inact,&nis_inact);

471:     /* Create active and inactive set vectors */
472:     SNESCreateSubVectors_VINEWTONRSLS(snes,nis_inact,&F_inact);
473:     SNESCreateSubVectors_VINEWTONRSLS(snes,nis_act,&Y_act);
474:     SNESCreateSubVectors_VINEWTONRSLS(snes,nis_inact,&Y_inact);

476:     /* Create scatter contexts */
477:     VecScatterCreate(Y,IS_act,Y_act,NULL,&scat_act);
478:     VecScatterCreate(Y,vi->IS_inact,Y_inact,NULL,&scat_inact);

480:     /* Do a vec scatter to active and inactive set vectors */
481:     VecScatterBegin(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);
482:     VecScatterEnd(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);

484:     VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);
485:     VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);

487:     VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);
488:     VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);

490:     /* Active set direction = 0 */
491:     VecSet(Y_act,0);
492:     if (snes->jacobian != snes->jacobian_pre) {
493:       MatGetSubMatrix(snes->jacobian_pre,vi->IS_inact,vi->IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);
494:     } else prejac_inact_inact = jac_inact_inact;

496:     ISEqual(vi->IS_inact_prev,vi->IS_inact,&isequal);
497:     if (!isequal) {
498:       SNESVIResetPCandKSP(snes,jac_inact_inact,prejac_inact_inact);
499:       PCFieldSplitRestrictIS(pc,vi->IS_inact);
500:     }

502:     /*      ISView(vi->IS_inact,0); */
503:     /*      ISView(IS_act,0);*/
504:     /*      MatView(snes->jacobian_pre,0); */



508:     KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact);
509:     KSPSetUp(snes->ksp);
510:     {
511:       PC        pc;
512:       PetscBool flg;
513:       KSPGetPC(snes->ksp,&pc);
514:       PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&flg);
515:       if (flg) {
516:         KSP *subksps;
517:         PCFieldSplitGetSubKSP(pc,NULL,&subksps);
518:         KSPGetPC(subksps[0],&pc);
519:         PetscFree(subksps);
520:         PetscObjectTypeCompare((PetscObject)pc,PCBJACOBI,&flg);
521:         if (flg) {
522:           PetscInt       n,N = 101*101,j,cnts[3] = {0,0,0};
523:           const PetscInt *ii;

525:           ISGetSize(vi->IS_inact,&n);
526:           ISGetIndices(vi->IS_inact,&ii);
527:           for (j=0; j<n; j++) {
528:             if (ii[j] < N) cnts[0]++;
529:             else if (ii[j] < 2*N) cnts[1]++;
530:             else if (ii[j] < 3*N) cnts[2]++;
531:           }
532:           ISRestoreIndices(vi->IS_inact,&ii);

534:           PCBJacobiSetTotalBlocks(pc,3,cnts);
535:         }
536:       }
537:     }

539:     KSPSolve(snes->ksp,F_inact,Y_inact);
540:     VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);
541:     VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);
542:     VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);
543:     VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);

545:     VecDestroy(&F_inact);
546:     VecDestroy(&Y_act);
547:     VecDestroy(&Y_inact);
548:     VecScatterDestroy(&scat_act);
549:     VecScatterDestroy(&scat_inact);
550:     ISDestroy(&IS_act);
551:     if (!isequal) {
552:       ISDestroy(&vi->IS_inact_prev);
553:       ISDuplicate(vi->IS_inact,&vi->IS_inact_prev);
554:     }
555:     ISDestroy(&vi->IS_inact);
556:     MatDestroy(&jac_inact_inact);
557:     if (snes->jacobian != snes->jacobian_pre) {
558:       MatDestroy(&prejac_inact_inact);
559:     }

561:     KSPGetConvergedReason(snes->ksp,&kspreason);
562:     if (kspreason < 0) {
563:       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
564:         PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);
565:         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
566:         break;
567:       }
568:     }

570:     KSPGetIterationNumber(snes->ksp,&lits);
571:     snes->linear_its += lits;
572:     PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);
573:     /*
574:     if (snes->ops->precheck) {
575:       PetscBool changed_y = PETSC_FALSE;
576:       (*snes->ops->precheck)(snes,X,Y,snes->precheck,&changed_y);
577:     }

579:     if (PetscLogPrintInfo) {
580:       SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);
581:     }
582:     */
583:     /* Compute a (scaled) negative update in the line search routine:
584:          Y <- X - lambda*Y
585:        and evaluate G = function(Y) (depends on the line search).
586:     */
587:     VecCopy(Y,snes->vec_sol_update);
588:     ynorm = 1; gnorm = fnorm;
589:     SNESLineSearchApply(snes->linesearch, X, F, &gnorm, Y);
590:     SNESLineSearchGetReason(snes->linesearch, &lssucceed);
591:     SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm);
592:     PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);
593:     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
594:     if (snes->domainerror) {
595:       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
596:       DMDestroyVI(snes->dm);
597:       return(0);
598:     }
599:     if (lssucceed) {
600:       if (++snes->numFailures >= snes->maxFailures) {
601:         PetscBool ismin;
602:         snes->reason = SNES_DIVERGED_LINE_SEARCH;
603:         SNESVICheckLocalMin_Private(snes,snes->jacobian,F,X,gnorm,&ismin);
604:         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
605:         break;
606:       }
607:    }
608:    DMDestroyVI(snes->dm);
609:     /* Update function and solution vectors */
610:     fnorm = gnorm;
611:     /* Monitor convergence */
612:     PetscObjectSAWsTakeAccess((PetscObject)snes);
613:     snes->iter = i+1;
614:     snes->norm = fnorm;
615:     PetscObjectSAWsGrantAccess((PetscObject)snes);
616:     SNESLogConvergenceHistory(snes,snes->norm,lits);
617:     SNESMonitor(snes,snes->iter,snes->norm);
618:     /* Test for convergence, xnorm = || X || */
619:     if (snes->ops->converged != SNESConvergedSkip) { VecNorm(X,NORM_2,&xnorm); }
620:     (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
621:     if (snes->reason) break;
622:   }
623:   /* make sure that the VI information attached to the DM is removed if the for loop above was broken early due to some exceptional conditional */
624:   DMDestroyVI(snes->dm);
625:   if (i == maxits) {
626:     PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);
627:     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
628:   }
629:   return(0);
630: }

634: PetscErrorCode SNESVISetRedundancyCheck(SNES snes,PetscErrorCode (*func)(SNES,IS,IS*,void*),void *ctx)
635: {
636:   SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS*)snes->data;

640:   vi->checkredundancy = func;
641:   vi->ctxP            = ctx;
642:   return(0);
643: }

645: #if defined(PETSC_HAVE_MATLAB_ENGINE)
646: #include <engine.h>
647: #include <mex.h>
648: typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext;

652: PetscErrorCode SNESVIRedundancyCheck_Matlab(SNES snes,IS is_act,IS *is_redact,void *ctx)
653: {
654:   PetscErrorCode    ierr;
655:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
656:   int               nlhs  = 1, nrhs = 5;
657:   mxArray           *plhs[1], *prhs[5];
658:   long long int     l1      = 0, l2 = 0, ls = 0;
659:   PetscInt          *indices=NULL;


667:   /* Create IS for reduced active set of size 0, its size and indices will
668:    bet set by the Matlab function */
669:   ISCreateGeneral(PetscObjectComm((PetscObject)snes),0,indices,PETSC_OWN_POINTER,is_redact);
670:   /* call Matlab function in ctx */
671:   PetscMemcpy(&ls,&snes,sizeof(snes));
672:   PetscMemcpy(&l1,&is_act,sizeof(is_act));
673:   PetscMemcpy(&l2,is_redact,sizeof(is_act));
674:   prhs[0] = mxCreateDoubleScalar((double)ls);
675:   prhs[1] = mxCreateDoubleScalar((double)l1);
676:   prhs[2] = mxCreateDoubleScalar((double)l2);
677:   prhs[3] = mxCreateString(sctx->funcname);
678:   prhs[4] = sctx->ctx;
679:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESVIRedundancyCheckInternal");
680:   mxGetScalar(plhs[0]);
681:   mxDestroyArray(prhs[0]);
682:   mxDestroyArray(prhs[1]);
683:   mxDestroyArray(prhs[2]);
684:   mxDestroyArray(prhs[3]);
685:   mxDestroyArray(plhs[0]);
686:   return(0);
687: }

691: PetscErrorCode SNESVISetRedundancyCheckMatlab(SNES snes,const char *func,mxArray *ctx)
692: {
693:   PetscErrorCode    ierr;
694:   SNESMatlabContext *sctx;

697:   /* currently sctx is memory bleed */
698:   PetscNew(&sctx);
699:   PetscStrallocpy(func,&sctx->funcname);
700:   sctx->ctx = mxDuplicateArray(ctx);
701:   SNESVISetRedundancyCheck(snes,SNESVIRedundancyCheck_Matlab,sctx);
702:   return(0);
703: }

705: #endif

707: /* -------------------------------------------------------------------------- */
708: /*
709:    SNESSetUp_VINEWTONRSLS - Sets up the internal data structures for the later use
710:    of the SNESVI nonlinear solver.

712:    Input Parameter:
713: .  snes - the SNES context

715:    Application Interface Routine: SNESSetUp()

717:    Notes:
718:    For basic use of the SNES solvers, the user need not explicitly call
719:    SNESSetUp(), since these actions will automatically occur during
720:    the call to SNESSolve().
721:  */
724: PetscErrorCode SNESSetUp_VINEWTONRSLS(SNES snes)
725: {
726:   PetscErrorCode    ierr;
727:   SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS*) snes->data;
728:   PetscInt          *indices;
729:   PetscInt          i,n,rstart,rend;
730:   SNESLineSearch    linesearch;

733:   SNESSetUp_VI(snes);

735:   /* Set up previous active index set for the first snes solve
736:    vi->IS_inact_prev = 0,1,2,....N */

738:   VecGetOwnershipRange(snes->vec_sol,&rstart,&rend);
739:   VecGetLocalSize(snes->vec_sol,&n);
740:   PetscMalloc1(n,&indices);
741:   for (i=0; i < n; i++) indices[i] = rstart + i;
742:   ISCreateGeneral(PetscObjectComm((PetscObject)snes),n,indices,PETSC_OWN_POINTER,&vi->IS_inact_prev);

744:   /* set the line search functions */
745:   if (!snes->linesearch) {
746:     SNESGetLineSearch(snes, &linesearch);
747:     SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);
748:   }
749:   return(0);
750: }
751: /* -------------------------------------------------------------------------- */
754: PetscErrorCode SNESReset_VINEWTONRSLS(SNES snes)
755: {
756:   SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS*) snes->data;
757:   PetscErrorCode    ierr;

760:   SNESReset_VI(snes);
761:   ISDestroy(&vi->IS_inact_prev);
762:   return(0);
763: }

765: /* -------------------------------------------------------------------------- */
766: /*MC
767:       SNESVINEWTONRSLS - Reduced space active set solvers for variational inequalities based on Newton's method

769:    Options Database:
770: +   -snes_type <vinewtonssls,vinewtonrsls> a semi-smooth solver, a reduced space active set method
771: -   -snes_vi_monitor - prints the number of active constraints at each iteration.

773:    Level: beginner

775:    References:
776: .  1. - T. S. Munson, and S. Benson. Flexible Complementarity Solvers for Large Scale
777:      Applications, Optimization Methods and Software, 21 (2006).

779: .seealso:  SNESVISetVariableBounds(), SNESVISetComputeVariableBounds(), SNESCreate(), SNES, SNESSetType(), SNESVINEWTONSSLS, SNESNEWTONTR, SNESLineSearchSet(),SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()

781: M*/
784: PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONRSLS(SNES snes)
785: {
786:   PetscErrorCode    ierr;
787:   SNES_VINEWTONRSLS *vi;

790:   snes->ops->reset          = SNESReset_VINEWTONRSLS;
791:   snes->ops->setup          = SNESSetUp_VINEWTONRSLS;
792:   snes->ops->solve          = SNESSolve_VINEWTONRSLS;
793:   snes->ops->destroy        = SNESDestroy_VI;
794:   snes->ops->setfromoptions = SNESSetFromOptions_VI;
795:   snes->ops->view           = NULL;
796:   snes->ops->converged      = SNESConvergedDefault_VI;

798:   snes->usesksp = PETSC_TRUE;
799:   snes->usespc  = PETSC_FALSE;

801:   PetscNewLog(snes,&vi);
802:   snes->data          = (void*)vi;
803:   vi->checkredundancy = NULL;

805:   PetscObjectComposeFunction((PetscObject)snes,"SNESVISetVariableBounds_C",SNESVISetVariableBounds_VI);
806:   PetscObjectComposeFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",SNESVISetComputeVariableBounds_VI);
807:   return(0);
808: }