Actual source code: damg.c
1: #define PETSCSNES_DLL
2:
3: #include petscda.h
4: #include petscksp.h
5: #include petscmg.h
6: #include petscdmmg.h
7: #include private/pcimpl.h
9: /*
10: Code for almost fully managing multigrid/multi-level linear solvers for DA grids
11: */
15: /*@C
16: DMMGCreate - Creates a DA based multigrid solver object. This allows one to
17: easily implement MG methods on regular grids.
19: Collective on MPI_Comm
21: Input Parameter:
22: + comm - the processors that will share the grids and solution process
23: . nlevels - number of multigrid levels (if this is negative it CANNOT be reset with -dmmg_nlevels
24: - user - an optional user context
26: Output Parameters:
27: . - the context
29: Options Database:
30: + -dmmg_nlevels <levels> - number of levels to use
31: . -pc_mg_galerkin - use Galerkin approach to compute coarser matrices
32: - -dmmg_mat_type <type> - matrix type that DMMG should create, defaults to MATAIJ
34: Notes:
35: To provide a different user context for each level call DMMGSetUser() after calling
36: this routine
38: Level: advanced
40: .seealso DMMGDestroy(), DMMGSetUser(), DMMGGetUser(), DMMGSetMatType(), DMMGSetNullSpace(), DMMGSetInitialGuess(),
41: DMMGSetISColoringType()
43: @*/
44: PetscErrorCode DMMGCreate(MPI_Comm comm,PetscInt nlevels,void *user,DMMG **dmmg)
45: {
47: PetscInt i;
48: DMMG *p;
49: PetscTruth ftype;
50: char mtype[256];
53: if (nlevels < 0) {
54: nlevels = -nlevels;
55: } else {
56: PetscOptionsGetInt(0,"-dmmg_nlevels",&nlevels,PETSC_IGNORE);
57: }
58: if (nlevels < 1) SETERRQ(PETSC_ERR_USER,"Cannot set levels less than 1");
60: PetscMalloc(nlevels*sizeof(DMMG),&p);
61: for (i=0; i<nlevels; i++) {
62: PetscNew(struct _n_DMMG,&p[i]);
63: p[i]->nlevels = nlevels - i;
64: p[i]->comm = comm;
65: p[i]->user = user;
66: p[i]->updatejacobianperiod = 1;
67: p[i]->updatejacobian = PETSC_TRUE;
68: p[i]->isctype = IS_COLORING_GLOBAL;
69: PetscStrallocpy(MATAIJ,&p[i]->mtype);
70: }
71: *dmmg = p;
73: PetscOptionsGetString(PETSC_NULL,"-dmmg_mat_type",mtype,256,&ftype);
74: if (ftype) {
75: DMMGSetMatType(*dmmg,mtype);
76: }
77: return(0);
78: }
82: /*@C
83: DMMGSetMatType - Sets the type of matrices that DMMG will create for its solvers.
85: Collective on MPI_Comm
87: Input Parameters:
88: + dmmg - the DMMG object created with DMMGCreate()
89: - mtype - the matrix type, defaults to MATAIJ
91: Level: intermediate
93: .seealso DMMGDestroy(), DMMGSetUser(), DMMGGetUser(), DMMGCreate(), DMMGSetNullSpace()
95: @*/
96: PetscErrorCode DMMGSetMatType(DMMG *dmmg,const MatType mtype)
97: {
98: PetscInt i;
102: for (i=0; i<dmmg[0]->nlevels; i++) {
103: PetscFree(dmmg[i]->mtype);
104: PetscStrallocpy(mtype,&dmmg[i]->mtype);
105: }
106: return(0);
107: }
111: /*@C
112: DMMGSetOptionsPrefix - Sets the prefix used for the solvers inside a DMMG
114: Collective on MPI_Comm
116: Input Parameters:
117: + dmmg - the DMMG object created with DMMGCreate()
118: - prefix - the prefix string
120: Level: intermediate
122: .seealso DMMGDestroy(), DMMGSetUser(), DMMGGetUser(), DMMGCreate(), DMMGSetNullSpace()
124: @*/
125: PetscErrorCode DMMGSetOptionsPrefix(DMMG *dmmg,const char prefix[])
126: {
127: PetscInt i;
129:
131: for (i=0; i<dmmg[0]->nlevels; i++) {
132: PetscStrallocpy(prefix,&dmmg[i]->prefix);
133: }
134: return(0);
135: }
139: /*@C
140: DMMGDestroy - Destroys a DA based multigrid solver object.
142: Collective on DMMG
144: Input Parameter:
145: . - the context
147: Level: advanced
149: .seealso DMMGCreate()
151: @*/
152: PetscErrorCode DMMGDestroy(DMMG *dmmg)
153: {
155: PetscInt i,nlevels = dmmg[0]->nlevels;
158: if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");
160: for (i=1; i<nlevels; i++) {
161: if (dmmg[i]->R) {MatDestroy(dmmg[i]->R);}
162: }
163: for (i=0; i<nlevels; i++) {
164: PetscStrfree(dmmg[i]->prefix);
165: PetscStrfree(dmmg[i]->mtype);
166: if (dmmg[i]->dm) {DMDestroy(dmmg[i]->dm);}
167: if (dmmg[i]->x) {VecDestroy(dmmg[i]->x);}
168: if (dmmg[i]->b) {VecDestroy(dmmg[i]->b);}
169: if (dmmg[i]->r) {VecDestroy(dmmg[i]->r);}
170: if (dmmg[i]->work1) {VecDestroy(dmmg[i]->work1);}
171: if (dmmg[i]->w) {VecDestroy(dmmg[i]->w);}
172: if (dmmg[i]->work2) {VecDestroy(dmmg[i]->work2);}
173: if (dmmg[i]->lwork1) {VecDestroy(dmmg[i]->lwork1);}
174: if (dmmg[i]->B) {MatDestroy(dmmg[i]->B);}
175: if (dmmg[i]->J) {MatDestroy(dmmg[i]->J);}
176: if (dmmg[i]->Rscale) {VecDestroy(dmmg[i]->Rscale);}
177: if (dmmg[i]->fdcoloring){MatFDColoringDestroy(dmmg[i]->fdcoloring);}
178: if (dmmg[i]->ksp && !dmmg[i]->snes) {KSPDestroy(dmmg[i]->ksp);}
179: if (dmmg[i]->snes) {PetscObjectDestroy((PetscObject)dmmg[i]->snes);}
180: if (dmmg[i]->inject) {VecScatterDestroy(dmmg[i]->inject);}
181: PetscFree(dmmg[i]);
182: }
183: PetscFree(dmmg);
184: return(0);
185: }
189: /*@C
190: DMMGSetDM - Sets the coarse grid information for the grids
192: Collective on DMMG
194: Input Parameter:
195: + dmmg - the context
196: - dm - the DA or DMComposite object
198: Options Database Keys:
199: . -dmmg_refine: Use the input problem as the coarse level and refine. Otherwise, use it as the fine level and coarsen.
201: Level: advanced
203: .seealso DMMGCreate(), DMMGDestroy(), DMMGSetMatType()
205: @*/
206: PetscErrorCode DMMGSetDM(DMMG *dmmg, DM dm)
207: {
208: PetscInt nlevels = dmmg[0]->nlevels;
209: PetscTruth doRefine = PETSC_TRUE;
210: PetscInt i;
211: DM *hierarchy;
215: if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");
217: /* Create DM data structure for all the levels */
218: PetscOptionsGetTruth(PETSC_NULL, "-dmmg_refine", &doRefine, PETSC_IGNORE);
219: PetscObjectReference((PetscObject) dm);
220: PetscMalloc(nlevels*sizeof(DM),&hierarchy);
221: if (doRefine) {
222: DMRefineHierarchy(dm,nlevels-1,hierarchy);
223: dmmg[0]->dm = dm;
224: for(i=1; i<nlevels; ++i) {
225: dmmg[i]->dm = hierarchy[i-1];
226: }
227: } else {
228: dmmg[nlevels-1]->dm = dm;
229: DMCoarsenHierarchy(dm,nlevels-1,hierarchy);
230: for(i=0; i<nlevels-1; ++i) {
231: dmmg[nlevels-2-i]->dm = hierarchy[i];
232: }
233: }
234: PetscFree(hierarchy);
235: /* Cleanup old structures (should use some private Destroy() instead) */
236: for(i = 0; i < nlevels; ++i) {
237: if (dmmg[i]->B) {MatDestroy(dmmg[i]->B); dmmg[i]->B = PETSC_NULL;}
238: if (dmmg[i]->J) {MatDestroy(dmmg[i]->J); dmmg[i]->J = PETSC_NULL;}
239: }
241: /* Create work vectors and matrix for each level */
242: for (i=0; i<nlevels; i++) {
243: DMCreateGlobalVector(dmmg[i]->dm,&dmmg[i]->x);
244: VecDuplicate(dmmg[i]->x,&dmmg[i]->b);
245: VecDuplicate(dmmg[i]->x,&dmmg[i]->r);
246: }
248: /* Create interpolation/restriction between levels */
249: for (i=1; i<nlevels; i++) {
250: DMGetInterpolation(dmmg[i-1]->dm,dmmg[i]->dm,&dmmg[i]->R,PETSC_NULL);
251: }
252: return(0);
253: }
258: /*@C
259: DMMGSolve - Actually solves the (non)linear system defined with the DMMG
261: Collective on DMMG
263: Input Parameter:
264: . dmmg - the context
266: Level: advanced
268: Options Database:
269: + -dmmg_grid_sequence - use grid sequencing to get the initial solution for each level from the previous
270: - -dmmg_monitor_solution - display the solution at each iteration
272: Notes: For linear (KSP) problems may be called more than once, uses the same
273: matrices but recomputes the right hand side for each new solve. Call DMMGSetKSP()
274: to generate new matrices.
275:
276: .seealso DMMGCreate(), DMMGDestroy(), DMMG, DMMGSetSNES(), DMMGSetKSP(), DMMGSetUp(), DMMGSetMatType()
278: @*/
279: PetscErrorCode DMMGSolve(DMMG *dmmg)
280: {
282: PetscInt i,nlevels = dmmg[0]->nlevels;
283: PetscTruth gridseq = PETSC_FALSE,vecmonitor = PETSC_FALSE,flg;
286: PetscOptionsGetTruth(0,"-dmmg_grid_sequence",&gridseq,PETSC_NULL);
287: PetscOptionsGetTruth(0,"-dmmg_monitor_solution",&vecmonitor,PETSC_NULL);
288: if (gridseq) {
289: if (dmmg[0]->initialguess) {
290: (*dmmg[0]->initialguess)(dmmg[0],dmmg[0]->x);
291: if (dmmg[0]->ksp && !dmmg[0]->snes) {
292: KSPSetInitialGuessNonzero(dmmg[0]->ksp,PETSC_TRUE);
293: }
294: }
295: for (i=0; i<nlevels-1; i++) {
296: (*dmmg[i]->solve)(dmmg,i);
297: if (vecmonitor) {
298: VecView(dmmg[i]->x,PETSC_VIEWER_DRAW_(dmmg[i]->comm));
299: }
300: MatInterpolate(dmmg[i+1]->R,dmmg[i]->x,dmmg[i+1]->x);
301: if (dmmg[i+1]->ksp && !dmmg[i+1]->snes) {
302: KSPSetInitialGuessNonzero(dmmg[i+1]->ksp,PETSC_TRUE);
303: }
304: }
305: } else {
306: if (dmmg[nlevels-1]->initialguess) {
307: (*dmmg[nlevels-1]->initialguess)(dmmg[nlevels-1],dmmg[nlevels-1]->x);
308: }
309: }
311: /*VecView(dmmg[nlevels-1]->x,PETSC_VIEWER_DRAW_WORLD);*/
313: (*DMMGGetFine(dmmg)->solve)(dmmg,nlevels-1);
314: if (vecmonitor) {
315: VecView(dmmg[nlevels-1]->x,PETSC_VIEWER_DRAW_(dmmg[nlevels-1]->comm));
316: }
318: flg = PETSC_FALSE;
319: PetscOptionsGetTruth(PETSC_NULL,"-dmmg_view",&flg,PETSC_NULL);
320: if (flg && !PetscPreLoadingOn) {
321: PetscViewer viewer;
322: PetscViewerASCIIGetStdout(dmmg[0]->comm,&viewer);
323: DMMGView(dmmg,viewer);
324: }
325: flg = PETSC_FALSE;
326: PetscOptionsGetTruth(PETSC_NULL,"-dmmg_view_binary",&flg,PETSC_NULL);
327: if (flg && !PetscPreLoadingOn) {
328: DMMGView(dmmg,PETSC_VIEWER_BINARY_(dmmg[0]->comm));
329: }
330: return(0);
331: }
335: PetscErrorCode DMMGSolveKSP(DMMG *dmmg,PetscInt level)
336: {
340: if (dmmg[level]->rhs) {
341: CHKMEMQ;
342: (*dmmg[level]->rhs)(dmmg[level],dmmg[level]->b);
343: CHKMEMQ;
344: }
345: KSPSolve(dmmg[level]->ksp,dmmg[level]->b,dmmg[level]->x);
346: return(0);
347: }
349: /*
350: For each level (of grid sequencing) this sets the interpolation/restriction and
351: work vectors needed by the multigrid preconditioner within the KSP
352: (for nonlinear problems the KSP inside the SNES) of that level.
354: Also sets the KSP monitoring on all the levels if requested by user.
356: */
359: PetscErrorCode DMMGSetUpLevel(DMMG *dmmg,KSP ksp,PetscInt nlevels)
360: {
361: PetscErrorCode ierr;
362: PetscInt i;
363: PC pc;
364: PetscTruth ismg,ismf,isshell,ismffd;
365: KSP lksp; /* solver internal to the multigrid preconditioner */
366: MPI_Comm *comms;
369: if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");
371: /* use fgmres on outer iteration by default */
372: KSPSetType(ksp,KSPFGMRES);
373: KSPGetPC(ksp,&pc);
374: PCSetType(pc,PCMG);
375: PetscMalloc(nlevels*sizeof(MPI_Comm),&comms);
376: for (i=0; i<nlevels; i++) {
377: comms[i] = dmmg[i]->comm;
378: }
379: PCMGSetLevels(pc,nlevels,comms);
380: PetscFree(comms);
381: PCMGSetType(pc,PC_MG_FULL);
383: PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
384: if (ismg) {
385: /* set solvers for each level */
386: for (i=0; i<nlevels; i++) {
387: if (i < nlevels-1) { /* don't set for finest level, they are set in PCApply_MG()*/
388: PCMGSetX(pc,i,dmmg[i]->x);
389: PCMGSetRhs(pc,i,dmmg[i]->b);
390: }
391: if (i > 0) {
392: PCMGSetR(pc,i,dmmg[i]->r);
393: }
394: /* If using a matrix free multiply and did not provide an explicit matrix to build
395: the preconditioner then must use no preconditioner
396: */
397: PetscTypeCompare((PetscObject)dmmg[i]->B,MATSHELL,&isshell);
398: PetscTypeCompare((PetscObject)dmmg[i]->B,MATDAAD,&ismf);
399: PetscTypeCompare((PetscObject)dmmg[i]->B,MATMFFD,&ismffd);
400: if (isshell || ismf || ismffd) {
401: PC lpc;
402: PCMGGetSmoother(pc,i,&lksp);
403: KSPGetPC(lksp,&lpc);
404: PCSetType(lpc,PCNONE);
405: }
406: }
408: /* Set interpolation/restriction between levels */
409: for (i=1; i<nlevels; i++) {
410: PCMGSetInterpolation(pc,i,dmmg[i]->R);
411: PCMGSetRestriction(pc,i,dmmg[i]->R);
412: }
413: }
414: return(0);
415: }
419: /*@C
420: DMMGSetKSP - Sets the linear solver object that will use the grid hierarchy
422: Collective on DMMG
424: Input Parameter:
425: + dmmg - the context
426: . func - function to compute linear system matrix on each grid level
427: - rhs - function to compute right hand side on each level (need only work on the finest grid
428: if you do not use grid sequencing)
430: Level: advanced
432: Notes: For linear problems my be called more than once, reevaluates the matrices if it is called more
433: than once. Call DMMGSolve() directly several times to solve with the same matrix but different
434: right hand sides.
435:
436: .seealso DMMGCreate(), DMMGDestroy, DMMGSetDM(), DMMGSolve(), DMMGSetMatType()
438: @*/
439: PetscErrorCode DMMGSetKSP(DMMG *dmmg,PetscErrorCode (*rhs)(DMMG,Vec),PetscErrorCode (*func)(DMMG,Mat,Mat))
440: {
442: PetscInt i,nlevels = dmmg[0]->nlevels,level;
443: PetscTruth ismg,galerkin=PETSC_FALSE;
444: PC pc;
445: KSP lksp;
446:
448: if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");
450: if (!dmmg[0]->ksp) {
451: /* create solvers for each level if they don't already exist*/
452: for (i=0; i<nlevels; i++) {
454: KSPCreate(dmmg[i]->comm,&dmmg[i]->ksp);
455: PetscObjectIncrementTabLevel((PetscObject)dmmg[i]->ksp,PETSC_NULL,nlevels-i);
456: KSPSetOptionsPrefix(dmmg[i]->ksp,dmmg[i]->prefix);
457: DMMGSetUpLevel(dmmg,dmmg[i]->ksp,i+1);
458: KSPSetFromOptions(dmmg[i]->ksp);
460: /* if the multigrid is being run with Galerkin then these matrices do not need to be created except on the finest level
461: we do not take advantage of this because it may be that Galerkin has not yet been selected for the KSP object
462: These are also used if grid sequencing is selected for the linear problem. We should probably turn off grid sequencing
463: for the linear problem */
464: if (!dmmg[i]->B) {
465: DMGetMatrix(dmmg[i]->dm,dmmg[nlevels-1]->mtype,&dmmg[i]->B);
466: }
467: if (!dmmg[i]->J) {
468: dmmg[i]->J = dmmg[i]->B;
469: PetscObjectReference((PetscObject) dmmg[i]->J);
470: }
472: dmmg[i]->solve = DMMGSolveKSP;
473: dmmg[i]->rhs = rhs;
474: }
475: }
477: /* evalute matrix on each level */
478: KSPGetPC(dmmg[nlevels-1]->ksp,&pc);
479: PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
480: if (ismg) {
481: PCMGGetGalerkin(pc,&galerkin);
482: }
483: if (func) {
484: if (galerkin) {
485: (*func)(dmmg[nlevels-1],dmmg[nlevels-1]->J,dmmg[nlevels-1]->B);
486: } else {
487: for (i=0; i<nlevels; i++) {
488: (*func)(dmmg[i],dmmg[i]->J,dmmg[i]->B);
489: }
490: }
491: }
493: for (i=0; i<nlevels-1; i++) {
494: KSPSetOptionsPrefix(dmmg[i]->ksp,"dmmg_");
495: }
497: for (level=0; level<nlevels; level++) {
498: KSPSetOperators(dmmg[level]->ksp,dmmg[level]->J,dmmg[level]->B,SAME_NONZERO_PATTERN);
499: KSPGetPC(dmmg[level]->ksp,&pc);
500: if (ismg) {
501: for (i=0; i<=level; i++) {
502: PCMGGetSmoother(pc,i,&lksp);
503: KSPSetOperators(lksp,dmmg[i]->J,dmmg[i]->B,SAME_NONZERO_PATTERN);
504: }
505: }
506: }
508: return(0);
509: }
513: /*@C
514: DMMGView - prints information on a DA based multi-level preconditioner
516: Collective on DMMG and PetscViewer
518: Input Parameter:
519: + dmmg - the context
520: - viewer - the viewer
522: Level: advanced
524: .seealso DMMGCreate(), DMMGDestroy(), DMMGSetMatType()
526: @*/
527: PetscErrorCode DMMGView(DMMG *dmmg,PetscViewer viewer)
528: {
530: PetscInt i,nlevels = dmmg[0]->nlevels;
531: PetscMPIInt flag;
532: MPI_Comm comm;
533: PetscTruth iascii,isbinary;
538: PetscObjectGetComm((PetscObject)viewer,&comm);
539: MPI_Comm_compare(comm,dmmg[0]->comm,&flag);
540: if (flag != MPI_CONGRUENT && flag != MPI_IDENT) {
541: SETERRQ(PETSC_ERR_ARG_NOTSAMECOMM,"Different communicators in the DMMG and the PetscViewer");
542: }
544: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
545: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
546: if (isbinary) {
547: for (i=0; i<nlevels; i++) {
548: MatView(dmmg[i]->J,viewer);
549: }
550: for (i=1; i<nlevels; i++) {
551: MatView(dmmg[i]->R,viewer);
552: }
553: } else {
554: if (iascii) {
555: PetscViewerASCIIPrintf(viewer,"DMMG Object with %D levels\n",nlevels);
556: if (dmmg[0]->isctype == IS_COLORING_GLOBAL) {
557: PetscViewerASCIIPrintf(viewer,"Using global (nonghosted) Jacobian coloring computation\n");
558: } else {
559: PetscViewerASCIIPrintf(viewer,"Using ghosted Jacobian coloring computation\n");
560: }
561: }
562: for (i=0; i<nlevels; i++) {
563: PetscViewerASCIIPushTab(viewer);
564: DMView(dmmg[i]->dm,viewer);
565: PetscViewerASCIIPopTab(viewer);
566: }
567: if (iascii) {
568: PetscViewerASCIIPrintf(viewer,"Using matrix type %s\n",dmmg[nlevels-1]->mtype);
569: }
570: if (DMMGGetKSP(dmmg)) {
571: KSPView(DMMGGetKSP(dmmg),viewer);
572: } else if (DMMGGetSNES(dmmg)) {
573: SNESView(DMMGGetSNES(dmmg),viewer);
574: } else if (iascii) {
575: PetscViewerASCIIPrintf(viewer,"DMMG does not have a SNES or KSP set\n");
576: }
577: }
578: return(0);
579: }
583: /*@C
584: DMMGSetNullSpace - Indicates the null space in the linear operator (this is needed by the linear solver)
586: Collective on DMMG
588: Input Parameter:
589: + dmmg - the context
590: . has_cnst - is the constant vector in the null space
591: . n - number of null vectors (excluding the possible constant vector)
592: - func - a function that fills an array of vectors with the null vectors (must be orthonormal), may be PETSC_NULL
594: Level: advanced
596: .seealso DMMGCreate(), DMMGDestroy, DMMGSetDM(), DMMGSolve(), MatNullSpaceCreate(), KSPSetNullSpace(), DMMGSetMatType()
598: @*/
599: PetscErrorCode DMMGSetNullSpace(DMMG *dmmg,PetscTruth has_cnst,PetscInt n,PetscErrorCode (*func)(DMMG,Vec[]))
600: {
602: PetscInt i,j,nlevels = dmmg[0]->nlevels;
603: Vec *nulls = 0;
604: MatNullSpace nullsp;
605: KSP iksp;
606: PC pc,ipc;
607: PetscTruth ismg,isred;
610: if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");
611: if (!dmmg[0]->ksp) SETERRQ(PETSC_ERR_ORDER,"Must call AFTER DMMGSetKSP() or DMMGSetSNES()");
612: if ((n && !func) || (!n && func)) SETERRQ(PETSC_ERR_ARG_INCOMP,"Both n and func() must be set together");
613: if (n < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative number of vectors in null space n = %D",n)
615: for (i=0; i<nlevels; i++) {
616: if (n) {
617: VecDuplicateVecs(dmmg[i]->b,n,&nulls);
618: (*func)(dmmg[i],nulls);
619: }
620: MatNullSpaceCreate(dmmg[i]->comm,has_cnst,n,nulls,&nullsp);
621: KSPSetNullSpace(dmmg[i]->ksp,nullsp);
622: for (j=i; j<nlevels; j++) {
623: KSPGetPC(dmmg[j]->ksp,&pc);
624: PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
625: if (ismg) {
626: PCMGGetSmoother(pc,i,&iksp);
627: KSPSetNullSpace(iksp, nullsp);
628: }
629: }
630: MatNullSpaceDestroy(nullsp);
631: if (n) {
632: VecDestroyVecs(nulls,n);
633: }
634: }
635: /* make all the coarse grid solvers have LU shift since they are singular */
636: for (i=0; i<nlevels; i++) {
637: KSPGetPC(dmmg[i]->ksp,&pc);
638: PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
639: if (ismg) {
640: PCMGGetSmoother(pc,0,&iksp);
641: KSPGetPC(iksp,&ipc);
642: PetscTypeCompare((PetscObject)ipc,PCREDUNDANT,&isred);
643: if (isred) {
644: PCRedundantGetPC(ipc,&ipc);
645: }
646: PCFactorSetShiftType(ipc,MAT_SHIFT_POSITIVE_DEFINITE);
647: }
648: }
649: return(0);
650: }
654: /*@C
655: DMMGInitialGuessCurrent - Use with DMMGSetInitialGuess() to use the current value in the
656: solution vector (obtainable with DMMGGetx()) as the initial guess. Otherwise for linear
657: problems zero is used for the initial guess (unless grid sequencing is used). For nonlinear
658: problems this is not needed; it always uses the previous solution as the initial guess.
660: Collective on DMMG
662: Input Parameter:
663: + dmmg - the context
664: - vec - dummy argument
666: Level: intermediate
668: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGSetInitialGuess()
670: @*/
671: PetscErrorCode DMMGInitialGuessCurrent(DMMG dmmg,Vec vec)
672: {
674: return(0);
675: }
679: /*@C
680: DMMGSetInitialGuess - Sets the function that computes an initial guess.
682: Collective on DMMG
684: Input Parameter:
685: + dmmg - the context
686: - guess - the function
688: Notes: For nonlinear problems, if this is not set, then the current value in the
689: solution vector (obtained with DMMGGetX()) is used. Thus is if you doing 'time
690: stepping' it will use your current solution as the guess for the next timestep.
691: If grid sequencing is used (via -dmmg_grid_sequence) then the "guess" function
692: is used only on the coarsest grid.
693: For linear problems, if this is not set, then 0 is used as an initial guess.
694: If you would like the linear solver to also (like the nonlinear solver) use
695: the current solution vector as the initial guess then use DMMGInitialGuessCurrent()
696: as the function you pass in
698: Level: intermediate
701: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGInitialGuessCurrent(), DMMGSetGalekin(), DMMGSetMatType(), DMMGSetNullSpace()
703: @*/
704: PetscErrorCode DMMGSetInitialGuess(DMMG *dmmg,PetscErrorCode (*guess)(DMMG,Vec))
705: {
706: PetscInt i,nlevels = dmmg[0]->nlevels;
710: for (i=0; i<nlevels; i++) {
711: if (dmmg[i]->ksp && !dmmg[i]->snes) {
712: KSPSetInitialGuessNonzero(dmmg[i]->ksp,PETSC_TRUE);
713: }
714: dmmg[i]->initialguess = guess;
715: }
716: return(0);
717: }