Actual source code: virs.c
petsc-dev 2014-02-02
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: . ISact - active 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_prev;
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: VecScatter 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: PetscObjectReference((PetscObject)*dm2);
127: /* need to set back global vectors in order to use the original injection */
128: DMClearGlobalVectors(dm1);
130: dm1->ops->createglobalvector = dmsnesvi1->createglobalvector;
132: DMCreateGlobalVector(dm1,&finemarked);
133: DMCreateGlobalVector(*dm2,&coarsemarked);
135: /*
136: fill finemarked with locations of inactive points
137: */
138: ISGetIndices(dmsnesvi1->inactive,&index);
139: ISGetLocalSize(dmsnesvi1->inactive,&n);
140: VecSet(finemarked,0.0);
141: for (k=0; k<n; k++) {
142: VecSetValue(finemarked,index[k],1.0,INSERT_VALUES);
143: }
144: VecAssemblyBegin(finemarked);
145: VecAssemblyEnd(finemarked);
147: DMCreateInjection(*dm2,dm1,&inject);
148: VecScatterBegin(inject,finemarked,coarsemarked,INSERT_VALUES,SCATTER_FORWARD);
149: VecScatterEnd(inject,finemarked,coarsemarked,INSERT_VALUES,SCATTER_FORWARD);
150: VecScatterDestroy(&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,DIFFERENT_NONZERO_PATTERN);
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: PetscBool lssucceed;
342: MatStructure flg = DIFFERENT_NONZERO_PATTERN;
343: PetscReal fnorm,gnorm,xnorm=0,ynorm;
344: Vec Y,X,F;
345: KSPConvergedReason kspreason;
348: snes->numFailures = 0;
349: snes->numLinearSolveFailures = 0;
350: snes->reason = SNES_CONVERGED_ITERATING;
352: maxits = snes->max_its; /* maximum number of iterations */
353: X = snes->vec_sol; /* solution vector */
354: F = snes->vec_func; /* residual vector */
355: Y = snes->work[0]; /* work vectors */
357: SNESLineSearchSetVIFunctions(snes->linesearch, SNESVIProjectOntoBounds, SNESVIComputeInactiveSetFnorm);
358: SNESLineSearchSetVecs(snes->linesearch, X, NULL, NULL, NULL, NULL);
359: SNESLineSearchSetUp(snes->linesearch);
361: PetscObjectSAWsTakeAccess((PetscObject)snes);
362: snes->iter = 0;
363: snes->norm = 0.0;
364: PetscObjectSAWsGrantAccess((PetscObject)snes);
366: SNESVIProjectOntoBounds(snes,X);
367: SNESComputeFunction(snes,X,F);
368: if (snes->domainerror) {
369: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
370: return(0);
371: }
372: SNESVIComputeInactiveSetFnorm(snes,F,X,&fnorm);
373: VecNormBegin(X,NORM_2,&xnorm); /* xnorm <- ||x|| */
374: VecNormEnd(X,NORM_2,&xnorm);
375: if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
377: PetscObjectSAWsTakeAccess((PetscObject)snes);
378: snes->norm = fnorm;
379: PetscObjectSAWsGrantAccess((PetscObject)snes);
380: SNESLogConvergenceHistory(snes,fnorm,0);
381: SNESMonitor(snes,0,fnorm);
383: /* test convergence */
384: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
385: if (snes->reason) return(0);
388: for (i=0; i<maxits; i++) {
390: IS IS_act,IS_inact; /* _act -> active set _inact -> inactive set */
391: IS IS_redact; /* redundant active set */
392: VecScatter scat_act,scat_inact;
393: PetscInt nis_act,nis_inact;
394: Vec Y_act,Y_inact,F_inact;
395: Mat jac_inact_inact,prejac_inact_inact;
396: PetscBool isequal;
398: /* Call general purpose update function */
399: if (snes->ops->update) {
400: (*snes->ops->update)(snes, snes->iter);
401: }
402: SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
405: /* Create active and inactive index sets */
407: /*original
408: SNESVICreateIndexSets_RS(snes,X,F,&IS_act,&IS_inact);
409: */
410: SNESVIGetActiveSetIS(snes,X,F,&IS_act);
412: if (vi->checkredundancy) {
413: (*vi->checkredundancy)(snes,IS_act,&IS_redact,vi->ctxP);
414: if (IS_redact) {
415: ISSort(IS_redact);
416: ISComplement(IS_redact,X->map->rstart,X->map->rend,&IS_inact);
417: ISDestroy(&IS_redact);
418: } else {
419: ISComplement(IS_act,X->map->rstart,X->map->rend,&IS_inact);
420: }
421: } else {
422: ISComplement(IS_act,X->map->rstart,X->map->rend,&IS_inact);
423: }
426: /* Create inactive set submatrix */
427: MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);
429: if (0) { /* Dead code (temporary developer hack) */
430: IS keptrows;
431: MatFindNonzeroRows(jac_inact_inact,&keptrows);
432: if (keptrows) {
433: PetscInt cnt,*nrows,k;
434: const PetscInt *krows,*inact;
435: PetscInt rstart=jac_inact_inact->rmap->rstart;
437: MatDestroy(&jac_inact_inact);
438: ISDestroy(&IS_act);
440: ISGetLocalSize(keptrows,&cnt);
441: ISGetIndices(keptrows,&krows);
442: ISGetIndices(IS_inact,&inact);
443: PetscMalloc1(cnt,&nrows);
444: for (k=0; k<cnt; k++) nrows[k] = inact[krows[k]-rstart];
445: ISRestoreIndices(keptrows,&krows);
446: ISRestoreIndices(IS_inact,&inact);
447: ISDestroy(&keptrows);
448: ISDestroy(&IS_inact);
450: ISCreateGeneral(PetscObjectComm((PetscObject)snes),cnt,nrows,PETSC_OWN_POINTER,&IS_inact);
451: ISComplement(IS_inact,F->map->rstart,F->map->rend,&IS_act);
452: MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);
453: }
454: }
455: DMSetVI(snes->dm,IS_inact);
456: /* remove later */
458: /*
459: VecView(vi->xu,PETSC_VIEWER_BINARY_(((PetscObject)(vi->xu))->comm));
460: VecView(vi->xl,PETSC_VIEWER_BINARY_(((PetscObject)(vi->xl))->comm));
461: VecView(X,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)X)));
462: VecView(F,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)F)));
463: ISView(IS_inact,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)IS_inact)));
464: */
466: /* Get sizes of active and inactive sets */
467: ISGetLocalSize(IS_act,&nis_act);
468: ISGetLocalSize(IS_inact,&nis_inact);
470: /* Create active and inactive set vectors */
471: SNESCreateSubVectors_VINEWTONRSLS(snes,nis_inact,&F_inact);
472: SNESCreateSubVectors_VINEWTONRSLS(snes,nis_act,&Y_act);
473: SNESCreateSubVectors_VINEWTONRSLS(snes,nis_inact,&Y_inact);
475: /* Create scatter contexts */
476: VecScatterCreate(Y,IS_act,Y_act,NULL,&scat_act);
477: VecScatterCreate(Y,IS_inact,Y_inact,NULL,&scat_inact);
479: /* Do a vec scatter to active and inactive set vectors */
480: VecScatterBegin(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);
481: VecScatterEnd(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);
483: VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);
484: VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);
486: VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);
487: VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);
489: /* Active set direction = 0 */
490: VecSet(Y_act,0);
491: if (snes->jacobian != snes->jacobian_pre) {
492: MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);
493: } else prejac_inact_inact = jac_inact_inact;
495: ISEqual(vi->IS_inact_prev,IS_inact,&isequal);
496: if (!isequal) {
497: SNESVIResetPCandKSP(snes,jac_inact_inact,prejac_inact_inact);
498: flg = DIFFERENT_NONZERO_PATTERN;
499: }
501: /* ISView(IS_inact,0); */
502: /* ISView(IS_act,0);*/
503: /* MatView(snes->jacobian_pre,0); */
507: KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);
508: KSPSetUp(snes->ksp);
509: {
510: PC pc;
511: PetscBool flg;
512: KSPGetPC(snes->ksp,&pc);
513: PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&flg);
514: if (flg) {
515: KSP *subksps;
516: PCFieldSplitGetSubKSP(pc,NULL,&subksps);
517: KSPGetPC(subksps[0],&pc);
518: PetscFree(subksps);
519: PetscObjectTypeCompare((PetscObject)pc,PCBJACOBI,&flg);
520: if (flg) {
521: PetscInt n,N = 101*101,j,cnts[3] = {0,0,0};
522: const PetscInt *ii;
524: ISGetSize(IS_inact,&n);
525: ISGetIndices(IS_inact,&ii);
526: for (j=0; j<n; j++) {
527: if (ii[j] < N) cnts[0]++;
528: else if (ii[j] < 2*N) cnts[1]++;
529: else if (ii[j] < 3*N) cnts[2]++;
530: }
531: ISRestoreIndices(IS_inact,&ii);
533: PCBJacobiSetTotalBlocks(pc,3,cnts);
534: }
535: }
536: }
538: KSPSolve(snes->ksp,F_inact,Y_inact);
539: KSPGetConvergedReason(snes->ksp,&kspreason);
540: if (kspreason < 0) {
541: if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
542: PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);
543: snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
544: break;
545: }
546: }
548: VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);
549: VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);
550: VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);
551: VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);
553: VecDestroy(&F_inact);
554: VecDestroy(&Y_act);
555: VecDestroy(&Y_inact);
556: VecScatterDestroy(&scat_act);
557: VecScatterDestroy(&scat_inact);
558: ISDestroy(&IS_act);
559: if (!isequal) {
560: ISDestroy(&vi->IS_inact_prev);
561: ISDuplicate(IS_inact,&vi->IS_inact_prev);
562: }
563: ISDestroy(&IS_inact);
564: MatDestroy(&jac_inact_inact);
565: if (snes->jacobian != snes->jacobian_pre) {
566: MatDestroy(&prejac_inact_inact);
567: }
568: KSPGetIterationNumber(snes->ksp,&lits);
569: snes->linear_its += lits;
570: PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);
571: /*
572: if (snes->ops->precheck) {
573: PetscBool changed_y = PETSC_FALSE;
574: (*snes->ops->precheck)(snes,X,Y,snes->precheck,&changed_y);
575: }
577: if (PetscLogPrintInfo) {
578: SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);
579: }
580: */
581: /* Compute a (scaled) negative update in the line search routine:
582: Y <- X - lambda*Y
583: and evaluate G = function(Y) (depends on the line search).
584: */
585: VecCopy(Y,snes->vec_sol_update);
586: ynorm = 1; gnorm = fnorm;
587: SNESLineSearchApply(snes->linesearch, X, F, &gnorm, Y);
588: SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);
589: SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm);
590: PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);
591: if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
592: if (snes->domainerror) {
593: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
594: DMDestroyVI(snes->dm);
595: return(0);
596: }
597: if (!lssucceed) {
598: if (++snes->numFailures >= snes->maxFailures) {
599: PetscBool ismin;
600: snes->reason = SNES_DIVERGED_LINE_SEARCH;
601: SNESVICheckLocalMin_Private(snes,snes->jacobian,F,X,gnorm,&ismin);
602: if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
603: break;
604: }
605: }
606: /* Update function and solution vectors */
607: fnorm = gnorm;
608: /* Monitor convergence */
609: PetscObjectSAWsTakeAccess((PetscObject)snes);
610: snes->iter = i+1;
611: snes->norm = fnorm;
612: PetscObjectSAWsGrantAccess((PetscObject)snes);
613: SNESLogConvergenceHistory(snes,snes->norm,lits);
614: SNESMonitor(snes,snes->iter,snes->norm);
615: /* Test for convergence, xnorm = || X || */
616: if (snes->ops->converged != SNESConvergedSkip) { VecNorm(X,NORM_2,&xnorm); }
617: (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
618: if (snes->reason) break;
619: }
620: DMDestroyVI(snes->dm);
621: if (i == maxits) {
622: PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);
623: if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
624: }
625: return(0);
626: }
630: PetscErrorCode SNESVISetRedundancyCheck(SNES snes,PetscErrorCode (*func)(SNES,IS,IS*,void*),void *ctx)
631: {
632: SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS*)snes->data;
636: vi->checkredundancy = func;
637: vi->ctxP = ctx;
638: return(0);
639: }
641: #if defined(PETSC_HAVE_MATLAB_ENGINE)
642: #include <engine.h>
643: #include <mex.h>
644: typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext;
648: PetscErrorCode SNESVIRedundancyCheck_Matlab(SNES snes,IS is_act,IS *is_redact,void *ctx)
649: {
650: PetscErrorCode ierr;
651: SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
652: int nlhs = 1, nrhs = 5;
653: mxArray *plhs[1], *prhs[5];
654: long long int l1 = 0, l2 = 0, ls = 0;
655: PetscInt *indices=NULL;
663: /* Create IS for reduced active set of size 0, its size and indices will
664: bet set by the Matlab function */
665: ISCreateGeneral(PetscObjectComm((PetscObject)snes),0,indices,PETSC_OWN_POINTER,is_redact);
666: /* call Matlab function in ctx */
667: PetscMemcpy(&ls,&snes,sizeof(snes));
668: PetscMemcpy(&l1,&is_act,sizeof(is_act));
669: PetscMemcpy(&l2,is_redact,sizeof(is_act));
670: prhs[0] = mxCreateDoubleScalar((double)ls);
671: prhs[1] = mxCreateDoubleScalar((double)l1);
672: prhs[2] = mxCreateDoubleScalar((double)l2);
673: prhs[3] = mxCreateString(sctx->funcname);
674: prhs[4] = sctx->ctx;
675: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESVIRedundancyCheckInternal");
676: mxGetScalar(plhs[0]);
677: mxDestroyArray(prhs[0]);
678: mxDestroyArray(prhs[1]);
679: mxDestroyArray(prhs[2]);
680: mxDestroyArray(prhs[3]);
681: mxDestroyArray(plhs[0]);
682: return(0);
683: }
687: PetscErrorCode SNESVISetRedundancyCheckMatlab(SNES snes,const char *func,mxArray *ctx)
688: {
689: PetscErrorCode ierr;
690: SNESMatlabContext *sctx;
693: /* currently sctx is memory bleed */
694: PetscMalloc(sizeof(SNESMatlabContext),&sctx);
695: PetscStrallocpy(func,&sctx->funcname);
696: sctx->ctx = mxDuplicateArray(ctx);
697: SNESVISetRedundancyCheck(snes,SNESVIRedundancyCheck_Matlab,sctx);
698: return(0);
699: }
701: #endif
703: /* -------------------------------------------------------------------------- */
704: /*
705: SNESSetUp_VINEWTONRSLS - Sets up the internal data structures for the later use
706: of the SNESVI nonlinear solver.
708: Input Parameter:
709: . snes - the SNES context
710: . x - the solution vector
712: Application Interface Routine: SNESSetUp()
714: Notes:
715: For basic use of the SNES solvers, the user need not explicitly call
716: SNESSetUp(), since these actions will automatically occur during
717: the call to SNESSolve().
718: */
721: PetscErrorCode SNESSetUp_VINEWTONRSLS(SNES snes)
722: {
723: PetscErrorCode ierr;
724: SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS*) snes->data;
725: PetscInt *indices;
726: PetscInt i,n,rstart,rend;
727: SNESLineSearch linesearch;
730: SNESSetUp_VI(snes);
732: /* Set up previous active index set for the first snes solve
733: vi->IS_inact_prev = 0,1,2,....N */
735: VecGetOwnershipRange(snes->vec_sol,&rstart,&rend);
736: VecGetLocalSize(snes->vec_sol,&n);
737: PetscMalloc1(n,&indices);
738: for (i=0; i < n; i++) indices[i] = rstart + i;
739: ISCreateGeneral(PetscObjectComm((PetscObject)snes),n,indices,PETSC_OWN_POINTER,&vi->IS_inact_prev);
741: /* set the line search functions */
742: if (!snes->linesearch) {
743: SNESGetLineSearch(snes, &linesearch);
744: SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);
745: }
746: return(0);
747: }
748: /* -------------------------------------------------------------------------- */
751: PetscErrorCode SNESReset_VINEWTONRSLS(SNES snes)
752: {
753: SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS*) snes->data;
754: PetscErrorCode ierr;
757: SNESReset_VI(snes);
758: ISDestroy(&vi->IS_inact_prev);
759: return(0);
760: }
762: /* -------------------------------------------------------------------------- */
763: /*MC
764: SNESVINEWTONRSLS - Reduced space active set solvers for variational inequalities based on Newton's method
766: Options Database:
767: + -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
768: - -snes_vi_monitor - prints the number of active constraints at each iteration.
770: Level: beginner
772: References:
773: - T. S. Munson, and S. Benson. Flexible Complementarity Solvers for Large-Scale
774: Applications, Optimization Methods and Software, 21 (2006).
776: .seealso: SNESVISetVariableBounds(), SNESVISetComputeVariableBounds(), SNESCreate(), SNES, SNESSetType(), SNESVINEWTONRSLS, SNESVINEWTONSSLS, SNESNEWTONTR, SNESLineSearchSet(),
777: SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
778: SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
780: M*/
783: PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONRSLS(SNES snes)
784: {
785: PetscErrorCode ierr;
786: SNES_VINEWTONRSLS *vi;
789: snes->ops->reset = SNESReset_VINEWTONRSLS;
790: snes->ops->setup = SNESSetUp_VINEWTONRSLS;
791: snes->ops->solve = SNESSolve_VINEWTONRSLS;
792: snes->ops->destroy = SNESDestroy_VI;
793: snes->ops->setfromoptions = SNESSetFromOptions_VI;
794: snes->ops->view = NULL;
795: snes->ops->converged = SNESConvergedDefault_VI;
797: snes->usesksp = PETSC_TRUE;
798: snes->usespc = PETSC_FALSE;
800: PetscNewLog(snes,&vi);
801: snes->data = (void*)vi;
802: vi->checkredundancy = NULL;
804: PetscObjectComposeFunction((PetscObject)snes,"SNESVISetVariableBounds_C",SNESVISetVariableBounds_VI);
805: PetscObjectComposeFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",SNESVISetComputeVariableBounds_VI);
806: return(0);
807: }