Actual source code: damgsnes.c
1: #define PETSCSNES_DLL
2:
3: #include petscda.h
4: #include ../src/dm/da/daimpl.h
5: /* It appears that preprocessor directives are not respected by bfort */
6: #ifdef PETSC_HAVE_SIEVE
7: #include petscmesh.h
8: #endif
9: #include petscmg.h
10: #include petscdmmg.h
12: #if defined(PETSC_HAVE_ADIC)
19: #endif
20: #if defined(PETSC_HAVE_SIEVE)
22: #endif
25: EXTERN PetscErrorCode NLFCreate_DAAD(NLF*);
26: EXTERN PetscErrorCode NLFDAADSetDA_DAAD(NLF,DA);
27: EXTERN PetscErrorCode NLFDAADSetCtx_DAAD(NLF,void*);
28: EXTERN PetscErrorCode NLFDAADSetResidual_DAAD(NLF,Vec);
29: EXTERN PetscErrorCode NLFDAADSetNewtonIterations_DAAD(NLF,PetscInt);
32: /*
33: period of -1 indicates update only on zeroth iteration of SNES
34: */
35: #define ShouldUpdate(l,it) (((dmmg[l-1]->updatejacobianperiod == -1) && (it == 0)) || \
36: ((dmmg[l-1]->updatejacobianperiod > 0) && !(it % dmmg[l-1]->updatejacobianperiod)))
37: /*
38: Evaluates the Jacobian on all of the grids. It is used by DMMG to provide the
39: ComputeJacobian() function that SNESSetJacobian() requires.
40: */
43: PetscErrorCode DMMGComputeJacobian_Multigrid(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
44: {
45: DMMG *dmmg = (DMMG*)ptr;
47: PetscInt i,nlevels = dmmg[0]->nlevels,it;
48: KSP ksp,lksp;
49: PC pc;
50: PetscTruth ismg,galerkin = PETSC_FALSE;
51: Vec W;
52: MatStructure flg;
55: if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as user context which should contain DMMG");
56: SNESGetIterationNumber(snes,&it);
58: /* compute Jacobian on finest grid */
59: if (dmmg[nlevels-1]->updatejacobian && ShouldUpdate(nlevels,it)) {
60: (*DMMGGetFine(dmmg)->computejacobian)(snes,X,J,B,flag,DMMGGetFine(dmmg));
61: } else {
62: PetscInfo3(0,"Skipping Jacobian, SNES iteration %D frequence %D level %D\n",it,dmmg[nlevels-1]->updatejacobianperiod,nlevels-1);
63: *flag = SAME_PRECONDITIONER;
64: }
65: MatMFFDSetBase(DMMGGetFine(dmmg)->J,X,PETSC_NULL);
67: /* create coarser grid Jacobians for preconditioner if multigrid is the preconditioner */
68: SNESGetKSP(snes,&ksp);
69: KSPGetPC(ksp,&pc);
70: PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
71: if (ismg) {
72: PCMGGetGalerkin(pc,&galerkin);
73: }
75: if (!galerkin) {
76: for (i=nlevels-1; i>0; i--) {
77: if (!dmmg[i-1]->w) {
78: VecDuplicate(dmmg[i-1]->x,&dmmg[i-1]->w);
79: }
80: W = dmmg[i-1]->w;
81: /* restrict X to coarser grid */
82: MatRestrict(dmmg[i]->R,X,W);
83: X = W;
84: /* scale to "natural" scaling for that grid */
85: VecPointwiseMult(X,X,dmmg[i]->Rscale);
86: /* tell the base vector for matrix free multiplies */
87: MatMFFDSetBase(dmmg[i-1]->J,X,PETSC_NULL);
88: /* compute Jacobian on coarse grid */
89: if (dmmg[i-1]->updatejacobian && ShouldUpdate(i,it)) {
90: (*dmmg[i-1]->computejacobian)(snes,X,&dmmg[i-1]->J,&dmmg[i-1]->B,&flg,dmmg[i-1]);
91: flg = SAME_NONZERO_PATTERN;
92: } else {
93: PetscInfo3(0,"Skipping Jacobian, SNES iteration %D frequence %D level %D\n",it,dmmg[i-1]->updatejacobianperiod,i-1);
94: flg = SAME_PRECONDITIONER;
95: }
96: if (ismg) {
97: PCMGGetSmoother(pc,i-1,&lksp);
98: KSPSetOperators(lksp,dmmg[i-1]->J,dmmg[i-1]->B,flg);
99: }
100: }
101: }
102: return(0);
103: }
105: /* ---------------------------------------------------------------------------*/
110: /*
111: DMMGFormFunction - This is a universal global FormFunction used by the DMMG code
112: when the user provides a local function.
114: Input Parameters:
115: + snes - the SNES context
116: . X - input vector
117: - ptr - optional user-defined context, as set by SNESSetFunction()
119: Output Parameter:
120: . F - function vector
122: */
123: PetscErrorCode DMMGFormFunction(SNES snes,Vec X,Vec F,void *ptr)
124: {
125: DMMG dmmg = (DMMG)ptr;
127: Vec localX;
128: DA da = (DA)dmmg->dm;
131: DAGetLocalVector(da,&localX);
132: /*
133: Scatter ghost points to local vector, using the 2-step process
134: DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
135: */
136: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
137: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
138: DAFormFunction1(da,localX,F,dmmg->user);
139: DARestoreLocalVector(da,&localX);
140: return(0);
141: }
145: PetscErrorCode DMMGFormFunctionGhost(SNES snes,Vec X,Vec F,void *ptr)
146: {
147: DMMG dmmg = (DMMG)ptr;
149: Vec localX, localF;
150: DA da = (DA)dmmg->dm;
153: DAGetLocalVector(da,&localX);
154: DAGetLocalVector(da,&localF);
155: /*
156: Scatter ghost points to local vector, using the 2-step process
157: DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
158: */
159: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
160: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
161: VecSet(F, 0.0);
162: VecSet(localF, 0.0);
163: DAFormFunction1(da,localX,localF,dmmg->user);
164: DALocalToGlobalBegin(da,localF,F);
165: DALocalToGlobalEnd(da,localF,F);
166: DARestoreLocalVector(da,&localX);
167: DARestoreLocalVector(da,&localF);
168: return(0);
169: }
171: #ifdef PETSC_HAVE_SIEVE
174: /*
175: DMMGFormFunctionMesh - This is a universal global FormFunction used by the DMMG code
176: when the user provides a local function.
178: Input Parameters:
179: + snes - the SNES context
180: . X - input vector
181: - ptr - This is the DMMG object
183: Output Parameter:
184: . F - function vector
186: */
187: PetscErrorCode DMMGFormFunctionMesh(SNES snes, Vec X, Vec F, void *ptr)
188: {
189: DMMG dmmg = (DMMG) ptr;
190: Mesh mesh = (Mesh) dmmg->dm;
191: SectionReal sectionF, section;
195: MeshGetSectionReal(mesh, "default", §ion);
196: SectionRealDuplicate(section, §ionF);
197: SectionRealToVec(section, mesh, SCATTER_REVERSE, X);
198: MeshFormFunction(mesh, section, sectionF, dmmg->user);
199: SectionRealToVec(sectionF, mesh, SCATTER_FORWARD, F);
200: SectionRealDestroy(sectionF);
201: SectionRealDestroy(section);
202: return(0);
203: }
207: /*
208: DMMGComputeJacobianMesh - This is a universal global FormJacobian used by the DMMG code
209: when the user provides a local function.
211: Input Parameters:
212: + snes - the SNES context
213: . X - input vector
214: - ptr - This is the DMMG object
216: Output Parameter:
217: . F - function vector
219: */
220: PetscErrorCode DMMGComputeJacobianMesh(SNES snes, Vec X, Mat *J, Mat *B, MatStructure *flag, void *ptr)
221: {
222: DMMG dmmg = (DMMG) ptr;
223: Mesh mesh = (Mesh) dmmg->dm;
224: SectionReal sectionX;
228: MeshGetSectionReal(mesh, "default", §ionX);
229: SectionRealToVec(sectionX, mesh, SCATTER_REVERSE, X);
230: MeshFormJacobian(mesh, sectionX, *B, dmmg->user);
231: /* Assemble true Jacobian; if it is different */
232: if (*J != *B) {
233: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
234: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
235: }
236: MatSetOption(*B, MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
237: *flag = SAME_NONZERO_PATTERN;
238: SectionRealDestroy(sectionX);
239: return(0);
240: }
241: #endif
245: /*
246: DMMGFormFunctionFD - This is a universal global FormFunction used by the DMMG code
247: when the user provides a local function used to compute the Jacobian via FD.
249: Input Parameters:
250: + snes - the SNES context
251: . X - input vector
252: - ptr - optional user-defined context, as set by SNESSetFunction()
254: Output Parameter:
255: . F - function vector
257: */
258: PetscErrorCode DMMGFormFunctionFD(SNES snes,Vec X,Vec F,void *ptr)
259: {
260: DMMG dmmg = (DMMG)ptr;
262: Vec localX;
263: DA da = (DA)dmmg->dm;
264: PetscInt N,n;
265:
267: /* determine whether X=localX */
268: DAGetLocalVector(da,&localX);
269: VecGetSize(X,&N);
270: VecGetSize(localX,&n);
272: if (n != N){ /* X != localX */
273: /* Scatter ghost points to local vector, using the 2-step process
274: DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
275: */
276: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
277: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
278: } else {
279: DARestoreLocalVector(da,&localX);
280: localX = X;
281: }
282: DAFormFunction(da,dmmg->lfj,localX,F,dmmg->user);
283: if (n != N){
284: DARestoreLocalVector(da,&localX);
285: }
286: return(0);
287: }
291: /*@C
292: SNESDAFormFunction - This is a universal function evaluation routine that
293: may be used with SNESSetFunction() as long as the user context has a DA
294: as its first record and the user has called DASetLocalFunction().
296: Collective on SNES
298: Input Parameters:
299: + snes - the SNES context
300: . X - input vector
301: . F - function vector
302: - ptr - pointer to a structure that must have a DA as its first entry. For example this
303: could be a DMMG, this ptr must have been passed into SNESDAFormFunction() as the context
305: Level: intermediate
307: .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(),
308: SNESSetFunction(), SNESSetJacobian()
310: @*/
311: PetscErrorCode SNESDAFormFunction(SNES snes,Vec X,Vec F,void *ptr)
312: {
314: Vec localX;
315: DA da = *(DA*)ptr;
316: PetscInt N,n;
317:
319: if (!da) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Looks like you called SNESSetFromFuntion(snes,SNESDAFormFunction,) without the DA context");
321: /* determine whether X=localX */
322: DAGetLocalVector(da,&localX);
323: VecGetSize(X,&N);
324: VecGetSize(localX,&n);
325:
326:
327: if (n != N){ /* X != localX */
328: /* Scatter ghost points to local vector, using the 2-step process
329: DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
330: */
331: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
332: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
333: } else {
334: DARestoreLocalVector(da,&localX);
335: localX = X;
336: }
337: DAFormFunction1(da,localX,F,ptr);
338: if (n != N){
339: if (PetscExceptionValue(ierr)) {
340: PetscErrorCode pDARestoreLocalVector(da,&localX);CHKERRQ(pierr);
341: }
342:
343: DARestoreLocalVector(da,&localX);
344: }
345: return(0);
346: }
348: /* ------------------------------------------------------------------------------*/
349: #include private/matimpl.h
352: PetscErrorCode DMMGComputeJacobianWithFD(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
353: {
355: DMMG dmmg = (DMMG)ctx;
356: MatFDColoring color = (MatFDColoring)dmmg->fdcoloring;
357:
359: if (color->ctype == IS_COLORING_GHOSTED){
360: DA da=(DA)dmmg->dm;
361: Vec x1_loc;
362: DAGetLocalVector(da,&x1_loc);
363: DAGlobalToLocalBegin(da,x1,INSERT_VALUES,x1_loc);
364: DAGlobalToLocalEnd(da,x1,INSERT_VALUES,x1_loc);
365: SNESDefaultComputeJacobianColor(snes,x1_loc,J,B,flag,dmmg->fdcoloring);
366: DARestoreLocalVector(da,&x1_loc);
367: } else {
368: SNESDefaultComputeJacobianColor(snes,x1,J,B,flag,dmmg->fdcoloring);
369: }
370: return(0);
371: }
375: PetscErrorCode DMMGComputeJacobianWithMF(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
376: {
378:
380: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
381: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
382: return(0);
383: }
387: /*
388: DMMGComputeJacobian - Evaluates the Jacobian when the user has provided
389: a local function evaluation routine.
390: */
391: PetscErrorCode DMMGComputeJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
392: {
393: DMMG dmmg = (DMMG) ptr;
395: Vec localX;
396: DA da = (DA) dmmg->dm;
399: DAGetLocalVector(da,&localX);
400: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
401: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
402: DAComputeJacobian1(da,localX,*B,dmmg->user);
403: DARestoreLocalVector(da,&localX);
404: /* Assemble true Jacobian; if it is different */
405: if (*J != *B) {
406: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
407: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
408: }
409: MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
410: *flag = SAME_NONZERO_PATTERN;
411: return(0);
412: }
416: /*
417: SNESDAComputeJacobianWithAdifor - This is a universal Jacobian evaluation routine
418: that may be used with SNESSetJacobian() from Fortran as long as the user context has
419: a DA as its first record and DASetLocalAdiforFunction() has been called.
421: Collective on SNES
423: Input Parameters:
424: + snes - the SNES context
425: . X - input vector
426: . J - Jacobian
427: . B - Jacobian used in preconditioner (usally same as J)
428: . flag - indicates if the matrix changed its structure
429: - ptr - optional user-defined context, as set by SNESSetFunction()
431: Level: intermediate
433: .seealso: DASetLocalFunction(), DASetLocalAdicFunction(), SNESSetFunction(), SNESSetJacobian()
435: */
436: PetscErrorCode SNESDAComputeJacobianWithAdifor(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
437: {
438: DA da = *(DA*) ptr;
440: Vec localX;
443: DAGetLocalVector(da,&localX);
444: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
445: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
446: DAComputeJacobian1WithAdifor(da,localX,*B,ptr);
447: DARestoreLocalVector(da,&localX);
448: /* Assemble true Jacobian; if it is different */
449: if (*J != *B) {
450: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
451: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
452: }
453: MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
454: *flag = SAME_NONZERO_PATTERN;
455: return(0);
456: }
460: /*
461: SNESDAComputeJacobian - This is a universal Jacobian evaluation routine for a
462: locally provided Jacobian.
464: Collective on SNES
466: Input Parameters:
467: + snes - the SNES context
468: . X - input vector
469: . J - Jacobian
470: . B - Jacobian used in preconditioner (usally same as J)
471: . flag - indicates if the matrix changed its structure
472: - ptr - optional user-defined context, as set by SNESSetFunction()
474: Level: intermediate
476: .seealso: DASetLocalFunction(), DASetLocalJacobian(), SNESSetFunction(), SNESSetJacobian()
478: */
479: PetscErrorCode SNESDAComputeJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
480: {
481: DA da = *(DA*) ptr;
483: Vec localX;
486: DAGetLocalVector(da,&localX);
487: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
488: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
489: DAComputeJacobian1(da,localX,*B,ptr);
490: DARestoreLocalVector(da,&localX);
491: /* Assemble true Jacobian; if it is different */
492: if (*J != *B) {
493: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
494: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
495: }
496: MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
497: *flag = SAME_NONZERO_PATTERN;
498: return(0);
499: }
503: PetscErrorCode DMMGSolveSNES(DMMG *dmmg,PetscInt level)
504: {
506: PetscInt nlevels = dmmg[0]->nlevels;
509: dmmg[0]->nlevels = level+1;
510: SNESSolve(dmmg[level]->snes,PETSC_NULL,dmmg[level]->x);
511: dmmg[0]->nlevels = nlevels;
512: return(0);
513: }
515: /* ===========================================================================================================*/
519: /*@C
520: DMMGSetSNES - Sets the nonlinear function that defines the nonlinear set of equations
521: to be solved using the grid hierarchy.
523: Collective on DMMG
525: Input Parameter:
526: + dmmg - the context
527: . function - the function that defines the nonlinear system
528: - jacobian - optional function to compute Jacobian
530: Options Database Keys:
531: + -snes_monitor
532: . -dmmg_coloring_from_mat - use graph coloring on the actual matrix nonzero structure instead of getting the coloring from the DM
533: . -dmmg_jacobian_fd
534: . -dmmg_jacobian_ad
535: . -dmmg_jacobian_mf_fd_operator
536: . -dmmg_jacobian_mf_fd
537: . -dmmg_jacobian_mf_ad_operator
538: . -dmmg_jacobian_mf_ad
539: . -dmmg_iscoloring_type
540: -
541: The period at which the Jacobian is recomputed can be set differently for different levels
542: of the Jacobian (for example lag all Jacobians except on the finest level).
543: There is no user interface currently for setting a different period on the different levels, one must set the
544: fields dmmg[i]->updatejacobian and dmmg[i]->updatejacobianperiod directly in the DMMG data structure.
545:
547: Level: advanced
549: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNESLocal(), DMMGSetFromOptions()
551: @*/
552: PetscErrorCode DMMGSetSNES(DMMG *dmmg,PetscErrorCode (*function)(SNES,Vec,Vec,void*),PetscErrorCode (*jacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
553: {
554: PetscErrorCode ierr;
555: PetscInt i,nlevels = dmmg[0]->nlevels;
556: PetscTruth mffdoperator,mffd,fdjacobian;
557: PetscTruth useFAS = PETSC_FALSE, fasBlock, fasGMRES;
558: PetscTruth monitor, monitorAll;
559: PetscInt fasPresmooth = 1, fasPostsmooth = 1, fasCoarsesmooth = 1, fasMaxIter = 2;
560: PetscReal fasRtol = 1.0e-8, fasAbstol = 1.0e-50;
561: #if defined(PETSC_HAVE_ADIC)
562: PetscTruth mfadoperator,mfad,adjacobian;
563: #endif
564: PetscCookie cookie;
567: if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");
568: if (!jacobian) jacobian = DMMGComputeJacobianWithFD;
569: PetscObjectGetCookie((PetscObject) dmmg[0]->dm, &cookie);
571: PetscOptionsBegin(dmmg[0]->comm,dmmg[0]->prefix,"DMMG Options","SNES");
572: PetscOptionsName("-dmmg_monitor","Monitor DMMG iterations","DMMG",&monitor);
573: PetscOptionsName("-dmmg_monitor_all","Monitor DMMG iterations","DMMG",&monitorAll);
574: PetscOptionsTruth("-dmmg_fas","Use the Full Approximation Scheme","DMMGSetSNES",useFAS,&useFAS,PETSC_NULL);
575: PetscOptionsName("-dmmg_fas_block","Use point-block smoothing","DMMG",&fasBlock);
576: PetscOptionsName("-dmmg_fas_ngmres","Use Nonlinear GMRES","DMMG",&fasGMRES);
577: PetscOptionsInt("-dmmg_fas_presmooth","Number of downward smoother iterates","DMMG",fasPresmooth,&fasPresmooth,PETSC_NULL);
578: PetscOptionsInt("-dmmg_fas_postsmooth","Number of upward smoother iterates","DMMG",fasPostsmooth,&fasPostsmooth,PETSC_NULL);
579: PetscOptionsInt("-dmmg_fas_coarsesmooth","Number of coarse smoother iterates","DMMG",fasCoarsesmooth,&fasCoarsesmooth,PETSC_NULL);
580: PetscOptionsReal("-dmmg_fas_rtol","Relative tolerance for FAS","DMMG",fasRtol,&fasRtol,PETSC_NULL);
581: PetscOptionsReal("-dmmg_fas_atol","Absolute tolerance for FAS","DMMG",fasAbstol,&fasAbstol,PETSC_NULL);
582: PetscOptionsInt("-dmmg_fas_max_its","Maximum number of iterates per smoother","DMMG",fasMaxIter,&fasMaxIter,PETSC_NULL);
584: PetscOptionsTruth("-dmmg_coloring_from_mat","Compute the coloring directly from the matrix nonzero structure","DMMGSetSNES",dmmg[0]->getcoloringfrommat,&dmmg[0]->getcoloringfrommat,PETSC_NULL);
586: PetscOptionsName("-dmmg_jacobian_fd","Compute sparse Jacobian explicitly with finite differencing","DMMGSetSNES",&fdjacobian);
587: if (fdjacobian) jacobian = DMMGComputeJacobianWithFD;
588: #if defined(PETSC_HAVE_ADIC)
589: PetscOptionsName("-dmmg_jacobian_ad","Compute sparse Jacobian explicitly with ADIC (automatic differentiation)","DMMGSetSNES",&adjacobian);
590: if (adjacobian) jacobian = DMMGComputeJacobianWithAdic;
591: #endif
593: PetscOptionsTruthGroupBegin("-dmmg_jacobian_mf_fd_operator","Apply Jacobian via matrix free finite differencing","DMMGSetSNES",&mffdoperator);
594: PetscOptionsTruthGroupEnd("-dmmg_jacobian_mf_fd","Apply Jacobian via matrix free finite differencing even in computing preconditioner","DMMGSetSNES",&mffd);
595: if (mffd) mffdoperator = PETSC_TRUE;
596: #if defined(PETSC_HAVE_ADIC)
597: PetscOptionsTruthGroupBegin("-dmmg_jacobian_mf_ad_operator","Apply Jacobian via matrix free ADIC (automatic differentiation)","DMMGSetSNES",&mfadoperator);
598: PetscOptionsTruthGroupEnd("-dmmg_jacobian_mf_ad","Apply Jacobian via matrix free ADIC (automatic differentiation) even in computing preconditioner","DMMGSetSNES",&mfad);
599: if (mfad) mfadoperator = PETSC_TRUE;
600: #endif
601: PetscOptionsEnum("-dmmg_iscoloring_type","Type of ISColoring","None",ISColoringTypes,(PetscEnum)dmmg[0]->isctype,(PetscEnum*)&dmmg[0]->isctype,PETSC_NULL);
602:
603: PetscOptionsEnd();
605: /* create solvers for each level */
606: for (i=0; i<nlevels; i++) {
607: SNESCreate(dmmg[i]->comm,&dmmg[i]->snes);
608: PetscObjectIncrementTabLevel((PetscObject)dmmg[i]->snes,PETSC_NULL,nlevels - i - 1);
609: SNESSetFunction(dmmg[i]->snes,dmmg[i]->b,function,dmmg[i]);
610: SNESSetOptionsPrefix(dmmg[i]->snes,dmmg[i]->prefix);
611: SNESGetKSP(dmmg[i]->snes,&dmmg[i]->ksp);
613: if (mffdoperator) {
614: MatCreateSNESMF(dmmg[i]->snes,&dmmg[i]->J);
615: VecDuplicate(dmmg[i]->x,&dmmg[i]->work1);
616: VecDuplicate(dmmg[i]->x,&dmmg[i]->work2);
617: MatMFFDSetFunction(dmmg[i]->J,(PetscErrorCode (*)(void*, Vec,Vec))SNESComputeFunction,dmmg[i]->snes);
618: if (mffd) {
619: dmmg[i]->B = dmmg[i]->J;
620: jacobian = DMMGComputeJacobianWithMF;
621: PetscObjectReference((PetscObject) dmmg[i]->B);
622: }
623: #if defined(PETSC_HAVE_ADIC)
624: } else if (mfadoperator) {
625: MatRegisterDAAD();
626: MatCreateDAAD((DA)dmmg[i]->dm,&dmmg[i]->J);
627: MatDAADSetCtx(dmmg[i]->J,dmmg[i]->user);
628: if (mfad) {
629: dmmg[i]->B = dmmg[i]->J;
630: jacobian = DMMGComputeJacobianWithMF;
631: PetscObjectReference((PetscObject) dmmg[i]->B);
632: }
633: #endif
634: }
635:
636: if (!useFAS) {
637: if (!dmmg[i]->B) {
638: DMGetMatrix(dmmg[i]->dm,dmmg[i]->mtype,&dmmg[i]->B);
639: }
640: if (!dmmg[i]->J) {
641: dmmg[i]->J = dmmg[i]->B;
642: PetscObjectReference((PetscObject) dmmg[i]->B);
643: }
644: }
646: DMMGSetUpLevel(dmmg,dmmg[i]->ksp,i+1);
647:
648: /*
649: if the number of levels is > 1 then we want the coarse solve in the grid sequencing to use LU
650: when possible
651: */
652: if (nlevels > 1 && i == 0) {
653: PC pc;
654: PetscMPIInt size;
655: MPI_Comm comm;
656: PetscTruth flg1,flg2,flg3;
657: KSP cksp;
659: KSPGetPC(dmmg[i]->ksp,&pc);
660: PCMGGetCoarseSolve(pc,&cksp);
661: KSPGetPC(cksp,&pc);
662: PetscTypeCompare((PetscObject)pc,PCILU,&flg1);
663: PetscTypeCompare((PetscObject)pc,PCSOR,&flg2);
664: PetscTypeCompare((PetscObject)pc,PETSC_NULL,&flg3);
665: if (flg1 || flg2 || flg3) {
666: KSPSetType(dmmg[i]->ksp,KSPPREONLY);
667: PetscObjectGetComm((PetscObject) pc,&comm);
668: MPI_Comm_size(comm,&size);
669: if (size > 1) {
670: PCSetType(pc,PCREDUNDANT);
671: } else {
672: PCSetType(pc,PCLU);
673: }
674: }
675: }
677: dmmg[i]->computejacobian = jacobian;
678: dmmg[i]->computefunction = function;
679: if (useFAS) {
680: if (cookie == DM_COOKIE) {
681: #if defined(PETSC_HAVE_ADIC)
682: if (fasBlock) {
683: dmmg[i]->solve = DMMGSolveFASb;
684: } else if(fasGMRES) {
685: dmmg[i]->solve = DMMGSolveFAS_NGMRES;
686: } else {
687: dmmg[i]->solve = DMMGSolveFAS4;
688: }
689: #else
690: SETERRQ(PETSC_ERR_SUP, "Must use ADIC for structured FAS.");
691: #endif
692: } else {
693: #if defined(PETSC_HAVE_SIEVE)
694: dmmg[i]->solve = DMMGSolveFAS_Mesh;
695: #endif
696: }
697: } else {
698: dmmg[i]->solve = DMMGSolveSNES;
699: }
700: }
702: if (jacobian == DMMGComputeJacobianWithFD) {
703: ISColoring iscoloring;
705: for (i=0; i<nlevels; i++) {
706: if (dmmg[0]->getcoloringfrommat) {
707: MatGetColoring(dmmg[i]->B,(MatColoringType)MATCOLORING_SL,&iscoloring);
708: } else {
709: DMGetColoring(dmmg[i]->dm,dmmg[0]->isctype,&iscoloring);
710: }
711: MatFDColoringCreate(dmmg[i]->B,iscoloring,&dmmg[i]->fdcoloring);
712: ISColoringDestroy(iscoloring);
713: if (function == DMMGFormFunction) function = DMMGFormFunctionFD;
714: MatFDColoringSetFunction(dmmg[i]->fdcoloring,(PetscErrorCode(*)(void))function,dmmg[i]);
715: MatFDColoringSetFromOptions(dmmg[i]->fdcoloring);
716: }
717: #if defined(PETSC_HAVE_ADIC)
718: } else if (jacobian == DMMGComputeJacobianWithAdic) {
719: for (i=0; i<nlevels; i++) {
720: ISColoring iscoloring;
721: DMGetColoring(dmmg[i]->dm,IS_COLORING_GHOSTED,&iscoloring);
722: MatSetColoring(dmmg[i]->B,iscoloring);
723: ISColoringDestroy(iscoloring);
724: }
725: #endif
726: }
727: if (!useFAS) {
728: for (i=0; i<nlevels; i++) {
729: SNESSetJacobian(dmmg[i]->snes,dmmg[i]->J,dmmg[i]->B,DMMGComputeJacobian_Multigrid,dmmg);
730: }
731: }
733: /* Create interpolation scaling */
734: for (i=1; i<nlevels; i++) {
735: DMGetInterpolationScale(dmmg[i-1]->dm,dmmg[i]->dm,dmmg[i]->R,&dmmg[i]->Rscale);
736: }
738: if (useFAS) {
739: PetscTruth flg;
741: for(i = 0; i < nlevels; i++) {
742: VecDuplicate(dmmg[i]->b,&dmmg[i]->w);
743: if (cookie == DM_COOKIE) {
744: #if defined(PETSC_HAVE_ADIC)
745: NLFCreate_DAAD(&dmmg[i]->nlf);
746: NLFDAADSetDA_DAAD(dmmg[i]->nlf,(DA)dmmg[i]->dm);
747: NLFDAADSetCtx_DAAD(dmmg[i]->nlf,dmmg[i]->user);
748: NLFDAADSetResidual_DAAD(dmmg[i]->nlf,dmmg[i]->r);
749: NLFDAADSetNewtonIterations_DAAD(dmmg[i]->nlf,fasMaxIter);
750: #endif
751: } else {
752: }
754: dmmg[i]->monitor = monitor;
755: dmmg[i]->monitorall = monitorAll;
756: dmmg[i]->presmooth = fasPresmooth;
757: dmmg[i]->postsmooth = fasPostsmooth;
758: dmmg[i]->coarsesmooth = fasCoarsesmooth;
759: dmmg[i]->rtol = fasRtol;
760: dmmg[i]->abstol = fasAbstol;
761: }
763: PetscOptionsHasName(0, "-dmmg_fas_view", &flg);
764: if (flg) {
765: for(i = 0; i < nlevels; i++) {
766: if (i == 0) {
767: PetscPrintf(dmmg[i]->comm,"FAS Solver Parameters\n");
768: PetscPrintf(dmmg[i]->comm," rtol %G atol %G\n",dmmg[i]->rtol,dmmg[i]->abstol);
769: PetscPrintf(dmmg[i]->comm," coarsesmooths %D\n",dmmg[i]->coarsesmooth);
770: PetscPrintf(dmmg[i]->comm," Newton iterations %D\n",fasMaxIter);
771: } else {
772: PetscPrintf(dmmg[i]->comm," level %D presmooths %D\n",i,dmmg[i]->presmooth);
773: PetscPrintf(dmmg[i]->comm," postsmooths %D\n",dmmg[i]->postsmooth);
774: PetscPrintf(dmmg[i]->comm," Newton iterations %D\n",fasMaxIter);
775: }
776: if (fasBlock) {
777: PetscPrintf(dmmg[i]->comm," using point-block smoothing\n");
778: } else if(fasGMRES) {
779: PetscPrintf(dmmg[i]->comm," using non-linear gmres\n");
780: }
781: }
782: }
783: }
784: return(0);
785: }
789: /*@C
790: DMMGSetFromOptions - Sets various options associated with the DMMG object
792: Collective on DMMG
794: Input Parameter:
795: . dmmg - the context
798: Notes: this is currently only supported for use with DMMGSetSNES() NOT DMMGSetKSP()
800: Most options are checked in DMMGSetSNES(); this does call SNESSetFromOptions()
802: Level: advanced
804: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNESLocal(), DMMGSetSNES()
806: @*/
807: PetscErrorCode DMMGSetFromOptions(DMMG *dmmg)
808: {
809: PetscErrorCode ierr;
810: PetscInt i,nlevels = dmmg[0]->nlevels;
813: if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");
815: for (i=0; i<nlevels; i++) {
816: SNESSetFromOptions(dmmg[i]->snes);
817: }
818: return(0);
819: }
823: /*@C
824: DMMGSetSNESLocalFD - Sets the local user function that is used to approximately compute the Jacobian
825: via finite differences.
827: Collective on DMMG
829: Input Parameter:
830: + dmmg - the context
831: - function - the function that defines the nonlinear system
833: Level: intermediate
835: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGSetSNESLocal()
837: @*/
838: PetscErrorCode DMMGSetSNESLocalFD(DMMG *dmmg,DALocalFunction1 function)
839: {
840: PetscInt i,nlevels = dmmg[0]->nlevels;
843: for (i=0; i<nlevels; i++) {
844: dmmg[i]->lfj = (PetscErrorCode (*)(void))function;
845: }
846: return(0);
847: }
852: /*@C
853: DMMGGetSNESLocal - Returns the local functions for residual and Jacobian evaluation.
855: Collective on DMMG
857: Input Parameter:
858: . dmmg - the context
860: Output Parameters:
861: + function - the function that defines the nonlinear system
862: - jacobian - function defines the local part of the Jacobian
864: Level: intermediate
866: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGSetSNESLocal()
867: @*/
868: PetscErrorCode DMMGGetSNESLocal(DMMG *dmmg,DALocalFunction1 *function, DALocalFunction1 *jacobian)
869: {
870: PetscCookie cookie;
874: PetscObjectGetCookie((PetscObject) dmmg[0]->dm, &cookie);
875: if (cookie == DM_COOKIE) {
876: DAGetLocalFunction((DA) dmmg[0]->dm, function);
877: DAGetLocalJacobian((DA) dmmg[0]->dm, jacobian);
878: } else {
879: #ifdef PETSC_HAVE_SIEVE
880: MeshGetLocalFunction((Mesh) dmmg[0]->dm, (PetscErrorCode (**)(Mesh,SectionReal,SectionReal,void*)) function);
881: MeshGetLocalJacobian((Mesh) dmmg[0]->dm, (PetscErrorCode (**)(Mesh,SectionReal,Mat,void*)) jacobian);
882: #else
883: SETERRQ(PETSC_ERR_SUP, "Unstructured grids only supported when Sieve is enabled.\nReconfigure with --with-sieve.");
884: #endif
885: }
886: return(0);
887: }
889: /*MC
890: DMMGSetSNESLocal - Sets the local user function that defines the nonlinear set of equations
891: that will use the grid hierarchy and (optionally) its derivative.
893: Collective on DMMG
895: Synopsis:
896: PetscErrorCode DMMGSetSNESLocal(DMMG *dmmg,DALocalFunction1 function, DALocalFunction1 jacobian,
897: DALocalFunction1 ad_function, DALocalFunction1 admf_function);
899: Input Parameter:
900: + dmmg - the context
901: . function - the function that defines the nonlinear system
902: . jacobian - function defines the local part of the Jacobian
903: . ad_function - the name of the function with an ad_ prefix. This is ignored currently
904: - admf_function - the name of the function with an ad_ prefix. This is ignored currently
906: Options Database Keys:
907: + -dmmg_jacobian_fd
908: . -dmmg_jacobian_ad
909: . -dmmg_jacobian_mf_fd_operator
910: . -dmmg_jacobian_mf_fd
911: . -dmmg_jacobian_mf_ad_operator
912: - -dmmg_jacobian_mf_ad
915: Level: intermediate
917: Notes:
918: If the Jacobian is not provided, this routine
919: uses finite differencing to approximate the Jacobian.
921: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES()
923: M*/
927: PetscErrorCode DMMGSetSNESLocal_Private(DMMG *dmmg,DALocalFunction1 function,DALocalFunction1 jacobian,DALocalFunction1 ad_function,DALocalFunction1 admf_function)
928: {
930: PetscInt i,nlevels = dmmg[0]->nlevels;
931: PetscCookie cookie;
932: PetscErrorCode (*computejacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*) = 0;
936: if (jacobian) computejacobian = DMMGComputeJacobian;
937: #if defined(PETSC_HAVE_ADIC)
938: else if (ad_function) computejacobian = DMMGComputeJacobianWithAdic;
939: #endif
940: CHKMEMQ;
941: PetscObjectGetCookie((PetscObject) dmmg[0]->dm,&cookie);
942: if (cookie == DM_COOKIE) {
943: PetscTruth flag;
944: /* it makes no sense to use an option to decide on ghost, it depends on whether the
945: formfunctionlocal computes ghost values in F or not. */
946: PetscOptionsHasName(PETSC_NULL, "-dmmg_form_function_ghost", &flag);
947: if (flag) {
948: DMMGSetSNES(dmmg,DMMGFormFunctionGhost,computejacobian);
949: } else {
950: DMMGSetSNES(dmmg,DMMGFormFunction,computejacobian);
951: }
952: for (i=0; i<nlevels; i++) {
953: DASetLocalFunction((DA)dmmg[i]->dm,function);
954: dmmg[i]->lfj = (PetscErrorCode (*)(void))function;
955: DASetLocalJacobian((DA)dmmg[i]->dm,jacobian);
956: DASetLocalAdicFunction((DA)dmmg[i]->dm,ad_function);
957: DASetLocalAdicMFFunction((DA)dmmg[i]->dm,admf_function);
958: }
959: CHKMEMQ;
960: } else {
961: #ifdef PETSC_HAVE_SIEVE
962: DMMGSetSNES(dmmg, DMMGFormFunctionMesh, DMMGComputeJacobianMesh);
963: for (i=0; i<nlevels; i++) {
964: MeshSetLocalFunction((Mesh) dmmg[i]->dm, (PetscErrorCode (*)(Mesh,SectionReal,SectionReal,void*)) function);
965: dmmg[i]->lfj = (PetscErrorCode (*)(void)) function;
966: MeshSetLocalJacobian((Mesh) dmmg[i]->dm, (PetscErrorCode (*)(Mesh,SectionReal,Mat,void*)) jacobian);
967: // Setup a work section
968: SectionReal defaultSec, constantSec;
969: PetscTruth hasConstant;
971: MeshGetSectionReal((Mesh) dmmg[i]->dm, "default", &defaultSec);
972: MeshHasSectionReal((Mesh) dmmg[i]->dm, "constant", &hasConstant);
973: if (!hasConstant) {
974: SectionRealDuplicate(defaultSec, &constantSec);
975: PetscObjectSetName((PetscObject) constantSec, "constant");
976: MeshSetSectionReal((Mesh) dmmg[i]->dm, constantSec);
977: SectionRealDestroy(constantSec);
978: }
979: }
980: CHKMEMQ;
981: #else
982: SETERRQ(PETSC_ERR_SUP, "Unstructured grids only supported when Sieve is enabled.\nReconfigure with --with-sieve.");
983: #endif
984: }
985: CHKMEMQ;
986: return(0);
987: }
991: PetscErrorCode DMMGFunctioni(void* ctx,PetscInt i,Vec u,PetscScalar* r)
992: {
993: DMMG dmmg = (DMMG)ctx;
994: Vec U = dmmg->lwork1;
996: VecScatter gtol;
999: /* copy u into interior part of U */
1000: DAGetScatter((DA)dmmg->dm,0,>ol,0);
1001: VecScatterBegin(gtol,u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL);
1002: VecScatterEnd(gtol,u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL);
1003: DAFormFunctioni1((DA)dmmg->dm,i,U,r,dmmg->user);
1004: return(0);
1005: }
1009: PetscErrorCode DMMGFunctionib(PetscInt i,Vec u,PetscScalar* r,void* ctx)
1010: {
1011: DMMG dmmg = (DMMG)ctx;
1012: Vec U = dmmg->lwork1;
1014: VecScatter gtol;
1017: /* copy u into interior part of U */
1018: DAGetScatter((DA)dmmg->dm,0,>ol,0);
1019: VecScatterBegin(gtol,u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL);
1020: VecScatterEnd(gtol,u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL);
1021: DAFormFunctionib1((DA)dmmg->dm,i,U,r,dmmg->user);
1022: return(0);
1023: }
1027: PetscErrorCode DMMGFunctioniBase(void* ctx,Vec u)
1028: {
1029: DMMG dmmg = (DMMG)ctx;
1030: Vec U = dmmg->lwork1;
1034: DAGlobalToLocalBegin((DA)dmmg->dm,u,INSERT_VALUES,U);
1035: DAGlobalToLocalEnd((DA)dmmg->dm,u,INSERT_VALUES,U);
1036: return(0);
1037: }
1041: PetscErrorCode DMMGSetSNESLocali_Private(DMMG *dmmg,PetscErrorCode (*functioni)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*),PetscErrorCode (*adi)(DALocalInfo*,MatStencil*,void*,void*,void*),PetscErrorCode (*adimf)(DALocalInfo*,MatStencil*,void*,void*,void*))
1042: {
1044: PetscInt i,nlevels = dmmg[0]->nlevels;
1047: for (i=0; i<nlevels; i++) {
1048: DASetLocalFunctioni((DA)dmmg[i]->dm,functioni);
1049: DASetLocalAdicFunctioni((DA)dmmg[i]->dm,adi);
1050: DASetLocalAdicMFFunctioni((DA)dmmg[i]->dm,adimf);
1051: MatMFFDSetFunctioni(dmmg[i]->J,DMMGFunctioni);
1052: MatMFFDSetFunctioniBase(dmmg[i]->J,DMMGFunctioniBase);
1053: if (!dmmg[i]->lwork1) {
1054: DACreateLocalVector((DA)dmmg[i]->dm,&dmmg[i]->lwork1);
1055: }
1056: }
1057: return(0);
1058: }
1062: PetscErrorCode DMMGSetSNESLocalib_Private(DMMG *dmmg,PetscErrorCode (*functioni)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*),PetscErrorCode (*adi)(DALocalInfo*,MatStencil*,void*,void*,void*),PetscErrorCode (*adimf)(DALocalInfo*,MatStencil*,void*,void*,void*))
1063: {
1065: PetscInt i,nlevels = dmmg[0]->nlevels;
1068: for (i=0; i<nlevels; i++) {
1069: DASetLocalFunctionib((DA)dmmg[i]->dm,functioni);
1070: DASetLocalAdicFunctionib((DA)dmmg[i]->dm,adi);
1071: DASetLocalAdicMFFunctionib((DA)dmmg[i]->dm,adimf);
1072: if (!dmmg[i]->lwork1) {
1073: DACreateLocalVector((DA)dmmg[i]->dm,&dmmg[i]->lwork1);
1074: }
1075: }
1076: return(0);
1077: }
1079: static PetscErrorCode (*localfunc)(void) = 0;
1083: /*
1084: Uses the DM object to call the user provided function with the correct calling
1085: sequence.
1086: */
1087: PetscErrorCode DMMGInitialGuess_Local(DMMG dmmg,Vec x)
1088: {
1092: (*dmmg->dm->ops->forminitialguess)(dmmg->dm,localfunc,x,0);
1093: return(0);
1094: }
1098: /*@C
1099: DMMGSetInitialGuessLocal - sets code to compute the initial guess for each level
1101: Collective on DMMG
1103: Input Parameter:
1104: + dmmg - the context
1105: - localguess - the function that computes the initial guess
1107: Level: intermediate
1109: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGSetInitialGuess(), DMMGSetSNESLocal()
1111: @*/
1112: PetscErrorCode DMMGSetInitialGuessLocal(DMMG *dmmg,PetscErrorCode (*localguess)(void))
1113: {
1117: localfunc = localguess; /* stash into ugly static for now */
1119: DMMGSetInitialGuess(dmmg,DMMGInitialGuess_Local);
1120: return(0);
1121: }
1125: /*@C
1126: DMMGSetISColoringType - type of coloring used to compute Jacobian via finite differencing
1128: Collective on DMMG
1130: Input Parameter:
1131: + dmmg - the context
1132: - isctype - IS_COLORING_GHOSTED (default) or IS_COLORING_GLOBAL
1134: Options Database:
1135: . -dmmg_iscoloring_type <ghosted or global>
1137: Notes: ghosted coloring requires using DMMGSetSNESLocal()
1139: Level: intermediate
1141: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGSetInitialGuess(), DMMGSetSNESLocal()
1143: @*/
1144: PetscErrorCode DMMGSetISColoringType(DMMG *dmmg,ISColoringType isctype)
1145: {
1147: dmmg[0]->isctype = isctype;
1148: return(0);
1149: }
1154: /*@C
1155: DMMGSetUp - Called after DMMGSetSNES() and (optionally) DMMGSetFromOptions()
1157: Collective on DMMG
1159: Input Parameter:
1160: . dmmg - the context
1162: Notes: Currently this must be called by the user (if they want to). It checks to see if fieldsplit preconditioner
1163: is being used and manages it.
1165: Level: advanced
1167: .seealso DMMGCreate(), DMMGDestroy(), DMMG, DMMGSetSNES(), DMMGSetKSP(), DMMGSolve(), DMMGSetMatType()
1169: @*/
1170: PetscErrorCode DMMGSetUp(DMMG *dmmg)
1171: {
1173: PetscInt i,nDM;
1174: PetscTruth fieldsplit,dmcomposite;
1175: KSP ksp;
1176: PC pc;
1177: IS *fields;
1180: SNESGetKSP(DMMGGetSNES(dmmg),&ksp);
1181: KSPGetPC(ksp,&pc);
1182: PetscTypeCompare((PetscObject)pc,PCFIELDSPLIT,&fieldsplit);
1183: PetscTypeCompare((PetscObject)DMMGGetDM(dmmg),"DMComposite",&dmcomposite);
1184: if (fieldsplit && dmcomposite) {
1185: PetscInfo(ksp,"Setting up physics based fieldsplit preconditioner\n");
1186: DMCompositeGetNumberDM((DMComposite)DMMGGetDM(dmmg),&nDM);
1187: DMCompositeGetGlobalISs((DMComposite)DMMGGetDM(dmmg),&fields);
1188: for (i=0; i<nDM; i++) {
1189: PCFieldSplitSetIS(pc,fields[i]);
1190: ISDestroy(fields[i]);
1191: }
1192: PetscFree(fields);
1193: }
1195: return(0);
1196: }