Actual source code: virs.c
petsc-3.3-p5 2012-12-01
2: #include <../src/snes/impls/vi/rs/virsimpl.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>
9: /*
10: SNESVIGetInactiveSet - Gets the global indices for the inactive set variables (these correspond to the degrees of freedom the linear
11: system is solved on)
13: Input parameter
14: . snes - the SNES context
16: Output parameter
17: . ISact - active set index set
19: */
20: PetscErrorCode SNESVIGetInactiveSet(SNES snes,IS* inact)
21: {
22: SNES_VIRS *vi = (SNES_VIRS*)snes->data;
24: *inact = vi->IS_inact_prev;
25: return(0);
26: }
28: /*
29: 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
30: defined by the reduced space method.
32: Simple calls the regular DM interpolation and restricts it to operation on the variables not associated with active constraints.
34: <*/
35: typedef struct {
36: PetscInt n; /* size of vectors in the reduced DM space */
37: IS inactive;
38: PetscErrorCode (*createinterpolation)(DM,DM,Mat*,Vec*); /* DM's original routines */
39: PetscErrorCode (*coarsen)(DM, MPI_Comm, DM*);
40: PetscErrorCode (*createglobalvector)(DM,Vec*);
41: DM dm; /* when destroying this object we need to reset the above function into the base DM */
42: } DM_SNESVI;
46: /*
47: DMCreateGlobalVector_SNESVI - Creates global vector of the size of the reduced space
49: */
50: PetscErrorCode DMCreateGlobalVector_SNESVI(DM dm,Vec *vec)
51: {
52: PetscErrorCode ierr;
53: PetscContainer isnes;
54: DM_SNESVI *dmsnesvi;
57: PetscObjectQuery((PetscObject)dm,"VI",(PetscObject *)&isnes);
58: if (!isnes) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_PLIB,"Composed SNES is missing");
59: PetscContainerGetPointer(isnes,(void**)&dmsnesvi);
60: VecCreateMPI(((PetscObject)dm)->comm,dmsnesvi->n,PETSC_DETERMINE,vec);
61: return(0);
62: }
66: /*
67: DMCreateInterpolation_SNESVI - Modifieds the interpolation obtained from the DM by removing all rows and columns associated with active constraints.
69: */
70: PetscErrorCode DMCreateInterpolation_SNESVI(DM dm1,DM dm2,Mat *mat,Vec *vec)
71: {
72: PetscErrorCode ierr;
73: PetscContainer isnes;
74: DM_SNESVI *dmsnesvi1,*dmsnesvi2;
75: Mat interp;
78: PetscObjectQuery((PetscObject)dm1,"VI",(PetscObject *)&isnes);
79: if (!isnes) SETERRQ(((PetscObject)dm1)->comm,PETSC_ERR_PLIB,"Composed VI data structure is missing");
80: PetscContainerGetPointer(isnes,(void**)&dmsnesvi1);
81: PetscObjectQuery((PetscObject)dm2,"VI",(PetscObject *)&isnes);
82: if (!isnes) SETERRQ(((PetscObject)dm2)->comm,PETSC_ERR_PLIB,"Composed VI data structure is missing");
83: PetscContainerGetPointer(isnes,(void**)&dmsnesvi2);
84:
85: (*dmsnesvi1->createinterpolation)(dm1,dm2,&interp,PETSC_NULL);
86: MatGetSubMatrix(interp,dmsnesvi2->inactive,dmsnesvi1->inactive,MAT_INITIAL_MATRIX,mat);
87: MatDestroy(&interp);
88: *vec = 0;
89: return(0);
90: }
92: extern PetscErrorCode DMSetVI(DM,IS);
96: /*
97: DMCoarsen_SNESVI - Computes the regular coarsened DM then computes additional information about its inactive set
99: */
100: PetscErrorCode DMCoarsen_SNESVI(DM dm1,MPI_Comm comm,DM *dm2)
101: {
102: PetscErrorCode ierr;
103: PetscContainer isnes;
104: DM_SNESVI *dmsnesvi1;
105: Vec finemarked,coarsemarked;
106: IS inactive;
107: VecScatter inject;
108: const PetscInt *index;
109: PetscInt n,k,cnt = 0,rstart,*coarseindex;
110: PetscScalar *marked;
113: PetscObjectQuery((PetscObject)dm1,"VI",(PetscObject *)&isnes);
114: if (!isnes) SETERRQ(((PetscObject)dm1)->comm,PETSC_ERR_PLIB,"Composed VI data structure is missing");
115: PetscContainerGetPointer(isnes,(void**)&dmsnesvi1);
116:
117: /* get the original coarsen */
118: (*dmsnesvi1->coarsen)(dm1,comm,dm2);
120: /* not sure why this extra reference is needed, but without the dm2 disappears too early */
121: PetscObjectReference((PetscObject)*dm2);
123: /* need to set back global vectors in order to use the original injection */
124: DMClearGlobalVectors(dm1);
125: dm1->ops->createglobalvector = dmsnesvi1->createglobalvector;
126: DMCreateGlobalVector(dm1,&finemarked);
127: DMCreateGlobalVector(*dm2,&coarsemarked);
129: /*
130: fill finemarked with locations of inactive points
131: */
132: ISGetIndices(dmsnesvi1->inactive,&index);
133: ISGetLocalSize(dmsnesvi1->inactive,&n);
134: VecSet(finemarked,0.0);
135: for (k=0;k<n;k++){
136: VecSetValue(finemarked,index[k],1.0,INSERT_VALUES);
137: }
138: VecAssemblyBegin(finemarked);
139: VecAssemblyEnd(finemarked);
141: DMCreateInjection(*dm2,dm1,&inject);
142: VecScatterBegin(inject,finemarked,coarsemarked,INSERT_VALUES,SCATTER_FORWARD);
143: VecScatterEnd(inject,finemarked,coarsemarked,INSERT_VALUES,SCATTER_FORWARD);
144: VecScatterDestroy(&inject);
146: /*
147: create index set list of coarse inactive points from coarsemarked
148: */
149: VecGetLocalSize(coarsemarked,&n);
150: VecGetOwnershipRange(coarsemarked,&rstart,PETSC_NULL);
151: VecGetArray(coarsemarked,&marked);
152: for (k=0; k<n; k++) {
153: if (marked[k] != 0.0) cnt++;
154: }
155: PetscMalloc(cnt*sizeof(PetscInt),&coarseindex);
156: cnt = 0;
157: for (k=0; k<n; k++) {
158: if (marked[k] != 0.0) coarseindex[cnt++] = k + rstart;
159: }
160: VecRestoreArray(coarsemarked,&marked);
161: ISCreateGeneral(((PetscObject)coarsemarked)->comm,cnt,coarseindex,PETSC_OWN_POINTER,&inactive);
163: DMClearGlobalVectors(dm1);
164: dm1->ops->createglobalvector = DMCreateGlobalVector_SNESVI;
165: DMSetVI(*dm2,inactive);
167: VecDestroy(&finemarked);
168: VecDestroy(&coarsemarked);
169: ISDestroy(&inactive);
170: return(0);
171: }
175: PetscErrorCode DMDestroy_SNESVI(DM_SNESVI *dmsnesvi)
176: {
178:
180: /* reset the base methods in the DM object that were changed when the DM_SNESVI was reset */
181: dmsnesvi->dm->ops->createinterpolation = dmsnesvi->createinterpolation;
182: dmsnesvi->dm->ops->coarsen = dmsnesvi->coarsen;
183: dmsnesvi->dm->ops->createglobalvector = dmsnesvi->createglobalvector;
184: /* need to clear out this vectors because some of them may not have a reference to the DM
185: but they are counted as having references to the DM in DMDestroy() */
186: DMClearGlobalVectors(dmsnesvi->dm);
188: ISDestroy(&dmsnesvi->inactive);
189: PetscFree(dmsnesvi);
190: return(0);
191: }
195: /*
196: DMSetVI - Marks a DM as associated with a VI problem. This causes the interpolation/restriction operators to
197: be restricted to only those variables NOT associated with active constraints.
199: */
200: PetscErrorCode DMSetVI(DM dm,IS inactive)
201: {
202: PetscErrorCode ierr;
203: PetscContainer isnes;
204: DM_SNESVI *dmsnesvi;
207: if (!dm) return(0);
209: PetscObjectReference((PetscObject)inactive);
211: PetscObjectQuery((PetscObject)dm,"VI",(PetscObject *)&isnes);
212: if (!isnes) {
213: PetscContainerCreate(((PetscObject)dm)->comm,&isnes);
214: PetscContainerSetUserDestroy(isnes,(PetscErrorCode (*)(void*))DMDestroy_SNESVI);
215: PetscNew(DM_SNESVI,&dmsnesvi);
216: PetscContainerSetPointer(isnes,(void*)dmsnesvi);
217: PetscObjectCompose((PetscObject)dm,"VI",(PetscObject)isnes);
218: PetscContainerDestroy(&isnes);
219: dmsnesvi->createinterpolation = dm->ops->createinterpolation;
220: dm->ops->createinterpolation = DMCreateInterpolation_SNESVI;
221: dmsnesvi->coarsen = dm->ops->coarsen;
222: dm->ops->coarsen = DMCoarsen_SNESVI;
223: dmsnesvi->createglobalvector = dm->ops->createglobalvector;
224: dm->ops->createglobalvector = DMCreateGlobalVector_SNESVI;
225: } else {
226: PetscContainerGetPointer(isnes,(void**)&dmsnesvi);
227: ISDestroy(&dmsnesvi->inactive);
228: }
229: DMClearGlobalVectors(dm);
230: ISGetLocalSize(inactive,&dmsnesvi->n);
231: dmsnesvi->inactive = inactive;
232: dmsnesvi->dm = dm;
233: return(0);
234: }
238: /*
239: DMDestroyVI - Frees the DM_SNESVI object contained in the DM
240: - also resets the function pointers in the DM for createinterpolation() etc to use the original DM
241: */
242: PetscErrorCode DMDestroyVI(DM dm)
243: {
244: PetscErrorCode ierr;
247: if (!dm) return(0);
248: PetscObjectCompose((PetscObject)dm,"VI",(PetscObject)PETSC_NULL);
249: return(0);
250: }
252: /* --------------------------------------------------------------------------------------------------------*/
259: PetscErrorCode SNESCreateIndexSets_VIRS(SNES snes,Vec X,Vec F,IS* ISact,IS* ISinact)
260: {
261: PetscErrorCode ierr;
264: SNESVIGetActiveSetIS(snes,X,F,ISact);
265: ISComplement(*ISact,X->map->rstart,X->map->rend,ISinact);
266: return(0);
267: }
269: /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */
272: PetscErrorCode SNESCreateSubVectors_VIRS(SNES snes,PetscInt n,Vec* newv)
273: {
275: Vec v;
278: VecCreate(((PetscObject)snes)->comm,&v);
279: VecSetSizes(v,n,PETSC_DECIDE);
280: VecSetFromOptions(v);
281: *newv = v;
283: return(0);
284: }
286: /* Resets the snes PC and KSP when the active set sizes change */
289: PetscErrorCode SNESVIResetPCandKSP(SNES snes,Mat Amat,Mat Pmat)
290: {
291: PetscErrorCode ierr;
292: KSP snesksp;
295: SNESGetKSP(snes,&snesksp);
296: KSPReset(snesksp);
298: /*
299: KSP kspnew;
300: PC pcnew;
301: const MatSolverPackage stype;
304: KSPCreate(((PetscObject)snes)->comm,&kspnew);
305: kspnew->pc_side = snesksp->pc_side;
306: kspnew->rtol = snesksp->rtol;
307: kspnew->abstol = snesksp->abstol;
308: kspnew->max_it = snesksp->max_it;
309: KSPSetType(kspnew,((PetscObject)snesksp)->type_name);
310: KSPGetPC(kspnew,&pcnew);
311: PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name);
312: PCSetOperators(kspnew->pc,Amat,Pmat,DIFFERENT_NONZERO_PATTERN);
313: PCFactorGetMatSolverPackage(snesksp->pc,&stype);
314: PCFactorSetMatSolverPackage(kspnew->pc,stype);
315: KSPDestroy(&snesksp);
316: snes->ksp = kspnew;
317: PetscLogObjectParent(snes,kspnew);
318: KSPSetFromOptions(kspnew);*/
319: return(0);
320: }
322: /* Variational Inequality solver using reduce space method. No semismooth algorithm is
323: implemented in this algorithm. It basically identifies the active constraints and does
324: a linear solve on the other variables (those not associated with the active constraints). */
327: PetscErrorCode SNESSolve_VIRS(SNES snes)
328: {
329: SNES_VIRS *vi = (SNES_VIRS*)snes->data;
330: PetscErrorCode ierr;
331: PetscInt maxits,i,lits;
332: PetscBool lssucceed;
333: MatStructure flg = DIFFERENT_NONZERO_PATTERN;
334: PetscReal fnorm,gnorm,xnorm=0,ynorm;
335: Vec Y,X,F;
336: KSPConvergedReason kspreason;
340: snes->numFailures = 0;
341: snes->numLinearSolveFailures = 0;
342: snes->reason = SNES_CONVERGED_ITERATING;
344: maxits = snes->max_its; /* maximum number of iterations */
345: X = snes->vec_sol; /* solution vector */
346: F = snes->vec_func; /* residual vector */
347: Y = snes->work[0]; /* work vectors */
349: SNESLineSearchSetVIFunctions(snes->linesearch, SNESVIProjectOntoBounds, SNESVIComputeInactiveSetFnorm);
350: SNESLineSearchSetVecs(snes->linesearch, X, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL);
351: SNESLineSearchSetUp(snes->linesearch);
353: PetscObjectTakeAccess(snes);
354: snes->iter = 0;
355: snes->norm = 0.0;
356: PetscObjectGrantAccess(snes);
358: SNESVIProjectOntoBounds(snes,X);
359: SNESComputeFunction(snes,X,F);
360: if (snes->domainerror) {
361: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
362: return(0);
363: }
364: SNESVIComputeInactiveSetFnorm(snes,F,X,&fnorm);
365: VecNormBegin(X,NORM_2,&xnorm); /* xnorm <- ||x|| */
366: VecNormEnd(X,NORM_2,&xnorm);
367: if (PetscIsInfOrNanReal(fnorm)) SETERRQ(((PetscObject)X)->comm,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
369: PetscObjectTakeAccess(snes);
370: snes->norm = fnorm;
371: PetscObjectGrantAccess(snes);
372: SNESLogConvHistory(snes,fnorm,0);
373: SNESMonitor(snes,0,fnorm);
375: /* set parameter for default relative tolerance convergence test */
376: snes->ttol = fnorm*snes->rtol;
377: /* test convergence */
378: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
379: if (snes->reason) return(0);
382: for (i=0; i<maxits; i++) {
384: IS IS_act,IS_inact; /* _act -> active set _inact -> inactive set */
385: IS IS_redact; /* redundant active set */
386: VecScatter scat_act,scat_inact;
387: PetscInt nis_act,nis_inact;
388: Vec Y_act,Y_inact,F_inact;
389: Mat jac_inact_inact,prejac_inact_inact;
390: PetscBool isequal;
392: /* Call general purpose update function */
393: if (snes->ops->update) {
394: (*snes->ops->update)(snes, snes->iter);
395: }
396: SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
399: /* Create active and inactive index sets */
400:
401: /*original
402: SNESVICreateIndexSets_RS(snes,X,F,&IS_act,&IS_inact);
403: */
404: SNESVIGetActiveSetIS(snes,X,F,&IS_act);
405:
406: if (vi->checkredundancy) {
407: (*vi->checkredundancy)(snes,IS_act,&IS_redact,vi->ctxP);
408: if (IS_redact){
409: ISSort(IS_redact);
410: ISComplement(IS_redact,X->map->rstart,X->map->rend,&IS_inact);
411: ISDestroy(&IS_redact);
412: }
413: else {
414: ISComplement(IS_act,X->map->rstart,X->map->rend,&IS_inact);
415: }
416: } else {
417: ISComplement(IS_act,X->map->rstart,X->map->rend,&IS_inact);
418: }
419:
421: /* Create inactive set submatrix */
422: MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);
423:
424: if (0) { /* Dead code (temporary developer hack) */
425: IS keptrows;
426: MatFindNonzeroRows(jac_inact_inact,&keptrows);
427: if (keptrows) {
428: PetscInt cnt,*nrows,k;
429: const PetscInt *krows,*inact;
430: PetscInt rstart=jac_inact_inact->rmap->rstart;
432: MatDestroy(&jac_inact_inact);
433: ISDestroy(&IS_act);
435: ISGetLocalSize(keptrows,&cnt);
436: ISGetIndices(keptrows,&krows);
437: ISGetIndices(IS_inact,&inact);
438: PetscMalloc(cnt*sizeof(PetscInt),&nrows);
439: for (k=0; k<cnt; k++) {
440: nrows[k] = inact[krows[k]-rstart];
441: }
442: ISRestoreIndices(keptrows,&krows);
443: ISRestoreIndices(IS_inact,&inact);
444: ISDestroy(&keptrows);
445: ISDestroy(&IS_inact);
447: ISCreateGeneral(((PetscObject)snes)->comm,cnt,nrows,PETSC_OWN_POINTER,&IS_inact);
448: ISComplement(IS_inact,F->map->rstart,F->map->rend,&IS_act);
449: MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);
450: }
451: }
452: DMSetVI(snes->dm,IS_inact);
453: /* remove later */
455: /*
456: VecView(vi->xu,PETSC_VIEWER_BINARY_(((PetscObject)(vi->xu))->comm));
457: VecView(vi->xl,PETSC_VIEWER_BINARY_(((PetscObject)(vi->xl))->comm));
458: VecView(X,PETSC_VIEWER_BINARY_(((PetscObject)X)->comm));
459: VecView(F,PETSC_VIEWER_BINARY_(((PetscObject)F)->comm));
460: ISView(IS_inact,PETSC_VIEWER_BINARY_(((PetscObject)IS_inact)->comm));
461: */
463: /* Get sizes of active and inactive sets */
464: ISGetLocalSize(IS_act,&nis_act);
465: ISGetLocalSize(IS_inact,&nis_inact);
467: /* Create active and inactive set vectors */
468: SNESCreateSubVectors_VIRS(snes,nis_inact,&F_inact);
469: SNESCreateSubVectors_VIRS(snes,nis_act,&Y_act);
470: SNESCreateSubVectors_VIRS(snes,nis_inact,&Y_inact);
472: /* Create scatter contexts */
473: VecScatterCreate(Y,IS_act,Y_act,PETSC_NULL,&scat_act);
474: VecScatterCreate(Y,IS_inact,Y_inact,PETSC_NULL,&scat_inact);
476: /* Do a vec scatter to active and inactive set vectors */
477: VecScatterBegin(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);
478: VecScatterEnd(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);
480: VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);
481: VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);
483: VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);
484: VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);
485:
486: /* Active set direction = 0 */
487: VecSet(Y_act,0);
488: if (snes->jacobian != snes->jacobian_pre) {
489: MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);
490: } else prejac_inact_inact = jac_inact_inact;
492: ISEqual(vi->IS_inact_prev,IS_inact,&isequal);
493: if (!isequal) {
494: SNESVIResetPCandKSP(snes,jac_inact_inact,prejac_inact_inact);
495: flg = DIFFERENT_NONZERO_PATTERN;
496: }
497:
498: /* ISView(IS_inact,0); */
499: /* ISView(IS_act,0);*/
500: /* MatView(snes->jacobian_pre,0); */
502:
503:
504: KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);
505: KSPSetUp(snes->ksp);
506: {
507: PC pc;
508: PetscBool flg;
509: KSPGetPC(snes->ksp,&pc);
510: PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&flg);
511: if (flg) {
512: KSP *subksps;
513: PCFieldSplitGetSubKSP(pc,PETSC_NULL,&subksps);
514: KSPGetPC(subksps[0],&pc);
515: PetscFree(subksps);
516: PetscObjectTypeCompare((PetscObject)pc,PCBJACOBI,&flg);
517: if (flg) {
518: PetscInt n,N = 101*101,j,cnts[3] = {0,0,0};
519: const PetscInt *ii;
521: ISGetSize(IS_inact,&n);
522: ISGetIndices(IS_inact,&ii);
523: for (j=0; j<n; j++) {
524: if (ii[j] < N) cnts[0]++;
525: else if (ii[j] < 2*N) cnts[1]++;
526: else if (ii[j] < 3*N) cnts[2]++;
527: }
528: ISRestoreIndices(IS_inact,&ii);
530: PCBJacobiSetTotalBlocks(pc,3,cnts);
531: }
532: }
533: }
535: SNES_KSPSolve(snes,snes->ksp,F_inact,Y_inact);
536: KSPGetConvergedReason(snes->ksp,&kspreason);
537: if (kspreason < 0) {
538: if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
539: PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);
540: snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
541: break;
542: }
543: }
545: VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);
546: VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);
547: VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);
548: VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);
550: VecDestroy(&F_inact);
551: VecDestroy(&Y_act);
552: VecDestroy(&Y_inact);
553: VecScatterDestroy(&scat_act);
554: VecScatterDestroy(&scat_inact);
555: ISDestroy(&IS_act);
556: if (!isequal) {
557: ISDestroy(&vi->IS_inact_prev);
558: ISDuplicate(IS_inact,&vi->IS_inact_prev);
559: }
560: ISDestroy(&IS_inact);
561: MatDestroy(&jac_inact_inact);
562: if (snes->jacobian != snes->jacobian_pre) {
563: MatDestroy(&prejac_inact_inact);
564: }
565: KSPGetIterationNumber(snes->ksp,&lits);
566: snes->linear_its += lits;
567: PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);
568: /*
569: if (snes->ops->precheckstep) {
570: PetscBool changed_y = PETSC_FALSE;
571: (*snes->ops->precheckstep)(snes,X,Y,snes->precheck,&changed_y);
572: }
574: if (PetscLogPrintInfo){
575: SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);
576: }
577: */
578: /* Compute a (scaled) negative update in the line search routine:
579: Y <- X - lambda*Y
580: and evaluate G = function(Y) (depends on the line search).
581: */
582: VecCopy(Y,snes->vec_sol_update);
583: ynorm = 1; gnorm = fnorm;
584: SNESLineSearchApply(snes->linesearch, X, F, &gnorm, Y);
585: SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm);
586: PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);
587: if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
588: if (snes->domainerror) {
589: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
590: DMDestroyVI(snes->dm);
591: return(0);
592: }
593: SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);
595: if (!lssucceed) {
596: if (++snes->numFailures >= snes->maxFailures) {
597: PetscBool ismin;
598: snes->reason = SNES_DIVERGED_LINE_SEARCH;
599: SNESVICheckLocalMin_Private(snes,snes->jacobian,F,X,gnorm,&ismin);
600: if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
601: break;
602: }
603: }
604: /* Update function and solution vectors */
605: fnorm = gnorm;
606: /* Monitor convergence */
607: PetscObjectTakeAccess(snes);
608: snes->iter = i+1;
609: snes->norm = fnorm;
610: PetscObjectGrantAccess(snes);
611: SNESLogConvHistory(snes,snes->norm,lits);
612: SNESMonitor(snes,snes->iter,snes->norm);
613: /* Test for convergence, xnorm = || X || */
614: if (snes->ops->converged != SNESSkipConverged) { VecNorm(X,NORM_2,&xnorm); }
615: (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
616: if (snes->reason) break;
617: }
618: DMDestroyVI(snes->dm);
619: if (i == maxits) {
620: PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);
621: if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
622: }
623: return(0);
624: }
628: PetscErrorCode SNESVISetRedundancyCheck(SNES snes,PetscErrorCode (*func)(SNES,IS,IS*,void*),void *ctx)
629: {
630: SNES_VIRS *vi = (SNES_VIRS*)snes->data;
634: vi->checkredundancy = func;
635: vi->ctxP = ctx;
636: return(0);
637: }
639: #if defined(PETSC_HAVE_MATLAB_ENGINE)
640: #include <engine.h>
641: #include <mex.h>
642: typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext;
646: PetscErrorCode SNESVIRedundancyCheck_Matlab(SNES snes,IS is_act,IS* is_redact,void* ctx)
647: {
648: PetscErrorCode ierr;
649: SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
650: int nlhs = 1, nrhs = 5;
651: mxArray *plhs[1], *prhs[5];
652: long long int l1 = 0, l2 = 0, ls = 0;
653: PetscInt *indices=PETSC_NULL;
661: /* Create IS for reduced active set of size 0, its size and indices will
662: bet set by the Matlab function */
663: ISCreateGeneral(((PetscObject)snes)->comm,0,indices,PETSC_OWN_POINTER,is_redact);
664: /* call Matlab function in ctx */
665: PetscMemcpy(&ls,&snes,sizeof(snes));
666: PetscMemcpy(&l1,&is_act,sizeof(is_act));
667: PetscMemcpy(&l2,is_redact,sizeof(is_act));
668: prhs[0] = mxCreateDoubleScalar((double)ls);
669: prhs[1] = mxCreateDoubleScalar((double)l1);
670: prhs[2] = mxCreateDoubleScalar((double)l2);
671: prhs[3] = mxCreateString(sctx->funcname);
672: prhs[4] = sctx->ctx;
673: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESVIRedundancyCheckInternal");
674: mxGetScalar(plhs[0]);
675: mxDestroyArray(prhs[0]);
676: mxDestroyArray(prhs[1]);
677: mxDestroyArray(prhs[2]);
678: mxDestroyArray(prhs[3]);
679: mxDestroyArray(plhs[0]);
680: return(0);
681: }
685: PetscErrorCode SNESVISetRedundancyCheckMatlab(SNES snes,const char* func,mxArray* ctx)
686: {
687: PetscErrorCode ierr;
688: SNESMatlabContext *sctx;
691: /* currently sctx is memory bleed */
692: PetscMalloc(sizeof(SNESMatlabContext),&sctx);
693: PetscStrallocpy(func,&sctx->funcname);
694: sctx->ctx = mxDuplicateArray(ctx);
695: SNESVISetRedundancyCheck(snes,SNESVIRedundancyCheck_Matlab,sctx);
696: return(0);
697: }
698:
699: #endif
701: /* -------------------------------------------------------------------------- */
702: /*
703: SNESSetUp_VIRS - Sets up the internal data structures for the later use
704: of the SNESVI nonlinear solver.
706: Input Parameter:
707: . snes - the SNES context
708: . x - the solution vector
710: Application Interface Routine: SNESSetUp()
712: Notes:
713: For basic use of the SNES solvers, the user need not explicitly call
714: SNESSetUp(), since these actions will automatically occur during
715: the call to SNESSolve().
716: */
719: PetscErrorCode SNESSetUp_VIRS(SNES snes)
720: {
722: SNES_VIRS *vi = (SNES_VIRS*) snes->data;
723: PetscInt *indices;
724: PetscInt i,n,rstart,rend;
725: SNESLineSearch linesearch;
728: SNESSetUp_VI(snes);
730: /* Set up previous active index set for the first snes solve
731: vi->IS_inact_prev = 0,1,2,....N */
733: VecGetOwnershipRange(snes->vec_sol,&rstart,&rend);
734: VecGetLocalSize(snes->vec_sol,&n);
735: PetscMalloc(n*sizeof(PetscInt),&indices);
736: for(i=0;i < n; i++) indices[i] = rstart + i;
737: ISCreateGeneral(((PetscObject)snes)->comm,n,indices,PETSC_OWN_POINTER,&vi->IS_inact_prev);
739: /* set the line search functions */
740: if (!snes->linesearch) {
741: SNESGetSNESLineSearch(snes, &linesearch);
742: SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);
743: }
744: return(0);
745: }
746: /* -------------------------------------------------------------------------- */
749: PetscErrorCode SNESReset_VIRS(SNES snes)
750: {
751: SNES_VIRS *vi = (SNES_VIRS*) snes->data;
755: SNESReset_VI(snes);
756: ISDestroy(&vi->IS_inact_prev);
757: return(0);
758: }
760: /* -------------------------------------------------------------------------- */
761: /*MC
762: SNESVIRS - Reduced space active set solvers for variational inequalities based on Newton's method
764: Options Database:
765: + -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
766: - -snes_vi_monitor - prints the number of active constraints at each iteration.
768: Level: beginner
770: References:
771: - T. S. Munson, and S. Benson. Flexible Complementarity Solvers for Large-Scale
772: Applications, Optimization Methods and Software, 21 (2006).
774: .seealso: SNESVISetVariableBounds(), SNESVISetComputeVariableBounds(), SNESCreate(), SNES, SNESSetType(), SNESVIRS, SNESVISS, SNESTR, SNESLineSearchSet(),
775: SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
776: SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
778: M*/
779: EXTERN_C_BEGIN
782: PetscErrorCode SNESCreate_VIRS(SNES snes)
783: {
785: SNES_VIRS *vi;
788: snes->ops->reset = SNESReset_VIRS;
789: snes->ops->setup = SNESSetUp_VIRS;
790: snes->ops->solve = SNESSolve_VIRS;
791: snes->ops->destroy = SNESDestroy_VI;
792: snes->ops->setfromoptions = SNESSetFromOptions_VI;
793: snes->ops->view = PETSC_NULL;
794: snes->ops->converged = SNESDefaultConverged_VI;
796: snes->usesksp = PETSC_TRUE;
797: snes->usespc = PETSC_FALSE;
799: PetscNewLog(snes,SNES_VIRS,&vi);
800: snes->data = (void*)vi;
801: vi->checkredundancy = PETSC_NULL;
803: PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESVISetVariableBounds_C","SNESVISetVariableBounds_VI",SNESVISetVariableBounds_VI);
804: PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESVISetComputeVariableBounds_C","SNESVISetComputeVariableBounds_VI",SNESVISetComputeVariableBounds_VI);
805: return(0);
806: }
807: EXTERN_C_END