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 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: + -dmmg_snes_monitor
532: . -dmmg_jacobian_fd
533: . -dmmg_jacobian_ad
534: . -dmmg_jacobian_mf_fd_operator
535: . -dmmg_jacobian_mf_fd
536: . -dmmg_jacobian_mf_ad_operator
537: . -dmmg_jacobian_mf_ad
538: . -dmmg_iscoloring_type
539: - -dmmg_jacobian_period <p> - Indicates how often in the SNES solve the Jacobian is recomputed (on all levels)
540: as suggested by Florin Dobrian if p is -1 then Jacobian is computed only on first
541: SNES iteration (i.e. -1 is equivalent to infinity)
543: Level: advanced
545: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNESLocal()
547: @*/
548: PetscErrorCode DMMGSetSNES(DMMG *dmmg,PetscErrorCode (*function)(SNES,Vec,Vec,void*),PetscErrorCode (*jacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
549: {
550: PetscErrorCode ierr;
551: PetscInt i,nlevels = dmmg[0]->nlevels,period = 1;
552: PetscTruth snesmonitor,mffdoperator,mffd,fdjacobian;
553: PetscTruth useFAS = PETSC_FALSE, fasBlock, fasGMRES;
554: PetscTruth monitor, monitorAll;
555: PetscInt fasPresmooth = 1, fasPostsmooth = 1, fasCoarsesmooth = 1, fasMaxIter = 2;
556: PetscReal fasRtol = 1.0e-8, fasAbstol = 1.0e-50;
557: #if defined(PETSC_HAVE_ADIC)
558: PetscTruth mfadoperator,mfad,adjacobian;
559: #endif
560: PetscViewerASCIIMonitor ascii;
561: MPI_Comm comm;
562: PetscCookie cookie;
565: if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");
566: if (!jacobian) jacobian = DMMGComputeJacobianWithFD;
567: PetscObjectGetCookie((PetscObject) dmmg[0]->dm, &cookie);
569: PetscOptionsBegin(dmmg[0]->comm,dmmg[0]->prefix,"DMMG Options","SNES");
570: PetscOptionsName("-dmmg_monitor","Monitor DMMG iterations","DMMG",&monitor);
571: PetscOptionsName("-dmmg_monitor_all","Monitor DMMG iterations","DMMG",&monitorAll);
572: PetscOptionsTruth("-dmmg_fas","Use the Full Approximation Scheme","DMMGSetSNES",useFAS,&useFAS,PETSC_NULL);
573: PetscOptionsName("-dmmg_fas_block","Use point-block smoothing","DMMG",&fasBlock);
574: PetscOptionsName("-dmmg_fas_ngmres","Use Nonlinear GMRES","DMMG",&fasGMRES);
575: PetscOptionsInt("-dmmg_fas_presmooth","Number of downward smoother iterates","DMMG",fasPresmooth,&fasPresmooth,PETSC_NULL);
576: PetscOptionsInt("-dmmg_fas_postsmooth","Number of upward smoother iterates","DMMG",fasPostsmooth,&fasPostsmooth,PETSC_NULL);
577: PetscOptionsInt("-dmmg_fas_coarsesmooth","Number of coarse smoother iterates","DMMG",fasCoarsesmooth,&fasCoarsesmooth,PETSC_NULL);
578: PetscOptionsReal("-dmmg_fas_rtol","Relative tolerance for FAS","DMMG",fasRtol,&fasRtol,PETSC_NULL);
579: PetscOptionsReal("-dmmg_fas_atol","Absolute tolerance for FAS","DMMG",fasAbstol,&fasAbstol,PETSC_NULL);
580: PetscOptionsInt("-dmmg_fas_max_its","Maximum number of iterates per smoother","DMMG",fasMaxIter,&fasMaxIter,PETSC_NULL);
582: PetscOptionsName("-dmmg_snes_monitor","Monitor nonlinear convergence","SNESMonitorSet",&snesmonitor);
584: PetscOptionsName("-dmmg_jacobian_fd","Compute sparse Jacobian explicitly with finite differencing","DMMGSetSNES",&fdjacobian);
585: if (fdjacobian) jacobian = DMMGComputeJacobianWithFD;
586: #if defined(PETSC_HAVE_ADIC)
587: PetscOptionsName("-dmmg_jacobian_ad","Compute sparse Jacobian explicitly with ADIC (automatic differentiation)","DMMGSetSNES",&adjacobian);
588: if (adjacobian) jacobian = DMMGComputeJacobianWithAdic;
589: #endif
591: PetscOptionsTruthGroupBegin("-dmmg_jacobian_mf_fd_operator","Apply Jacobian via matrix free finite differencing","DMMGSetSNES",&mffdoperator);
592: PetscOptionsTruthGroupEnd("-dmmg_jacobian_mf_fd","Apply Jacobian via matrix free finite differencing even in computing preconditioner","DMMGSetSNES",&mffd);
593: if (mffd) mffdoperator = PETSC_TRUE;
594: #if defined(PETSC_HAVE_ADIC)
595: PetscOptionsTruthGroupBegin("-dmmg_jacobian_mf_ad_operator","Apply Jacobian via matrix free ADIC (automatic differentiation)","DMMGSetSNES",&mfadoperator);
596: PetscOptionsTruthGroupEnd("-dmmg_jacobian_mf_ad","Apply Jacobian via matrix free ADIC (automatic differentiation) even in computing preconditioner","DMMGSetSNES",&mfad);
597: if (mfad) mfadoperator = PETSC_TRUE;
598: #endif
599: PetscOptionsEnum("-dmmg_iscoloring_type","Type of ISColoring","None",ISColoringTypes,(PetscEnum)dmmg[0]->isctype,(PetscEnum*)&dmmg[0]->isctype,PETSC_NULL);
600:
601: PetscOptionsEnd();
603: /* create solvers for each level */
604: for (i=0; i<nlevels; i++) {
605: SNESCreate(dmmg[i]->comm,&dmmg[i]->snes);
606: SNESSetFunction(dmmg[i]->snes,dmmg[i]->b,function,dmmg[i]);
607: SNESSetOptionsPrefix(dmmg[i]->snes,dmmg[i]->prefix);
608: SNESGetKSP(dmmg[i]->snes,&dmmg[i]->ksp);
609: if (snesmonitor) {
610: PetscObjectGetComm((PetscObject)dmmg[i]->snes,&comm);
611: PetscViewerASCIIMonitorCreate(comm,"stdout",nlevels-i,&ascii);
612: SNESMonitorSet(dmmg[i]->snes,SNESMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);
613: }
615: if (mffdoperator) {
616: MatCreateSNESMF(dmmg[i]->snes,&dmmg[i]->J);
617: VecDuplicate(dmmg[i]->x,&dmmg[i]->work1);
618: VecDuplicate(dmmg[i]->x,&dmmg[i]->work2);
619: MatMFFDSetFunction(dmmg[i]->J,(PetscErrorCode (*)(void*, Vec,Vec))SNESComputeFunction,dmmg[i]->snes);
620: if (mffd) {
621: dmmg[i]->B = dmmg[i]->J;
622: jacobian = DMMGComputeJacobianWithMF;
623: PetscObjectReference((PetscObject) dmmg[i]->B);
624: }
625: #if defined(PETSC_HAVE_ADIC)
626: } else if (mfadoperator) {
627: MatRegisterDAAD();
628: MatCreateDAAD((DA)dmmg[i]->dm,&dmmg[i]->J);
629: MatDAADSetCtx(dmmg[i]->J,dmmg[i]->user);
630: if (mfad) {
631: dmmg[i]->B = dmmg[i]->J;
632: jacobian = DMMGComputeJacobianWithMF;
633: PetscObjectReference((PetscObject) dmmg[i]->B);
634: }
635: #endif
636: }
637:
638: if (!useFAS) {
639: if (!dmmg[i]->B) {
640: DMGetMatrix(dmmg[i]->dm,dmmg[i]->mtype,&dmmg[i]->B);
641: }
642: if (!dmmg[i]->J) {
643: dmmg[i]->J = dmmg[i]->B;
644: PetscObjectReference((PetscObject) dmmg[i]->B);
645: }
646: }
648: DMMGSetUpLevel(dmmg,dmmg[i]->ksp,i+1);
649:
650: /*
651: if the number of levels is > 1 then we want the coarse solve in the grid sequencing to use LU
652: when possible
653: */
654: if (nlevels > 1 && i == 0) {
655: PC pc;
656: KSP cksp;
657: PetscTruth flg1,flg2,flg3;
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: PCSetType(pc,PCLU);
667: }
668: }
670: dmmg[i]->computejacobian = jacobian;
671: dmmg[i]->computefunction = function;
672: if (useFAS) {
673: if (cookie == DA_COOKIE) {
674: #if defined(PETSC_HAVE_ADIC)
675: if (fasBlock) {
676: dmmg[i]->solve = DMMGSolveFASb;
677: } else if(fasGMRES) {
678: dmmg[i]->solve = DMMGSolveFAS_NGMRES;
679: } else {
680: dmmg[i]->solve = DMMGSolveFAS4;
681: }
682: #else
683: SETERRQ(PETSC_ERR_SUP, "Must use ADIC for structured FAS.");
684: #endif
685: } else {
686: #if defined(PETSC_HAVE_SIEVE)
687: dmmg[i]->solve = DMMGSolveFAS_Mesh;
688: #endif
689: }
690: } else {
691: dmmg[i]->solve = DMMGSolveSNES;
692: }
693: }
695: if (jacobian == DMMGComputeJacobianWithFD) {
696: ISColoring iscoloring;
698: for (i=0; i<nlevels; i++) {
699: DMGetColoring(dmmg[i]->dm,dmmg[0]->isctype,&iscoloring);
700: MatFDColoringCreate(dmmg[i]->B,iscoloring,&dmmg[i]->fdcoloring);
701: ISColoringDestroy(iscoloring);
702: if (function == DMMGFormFunction) function = DMMGFormFunctionFD;
703: MatFDColoringSetFunction(dmmg[i]->fdcoloring,(PetscErrorCode(*)(void))function,dmmg[i]);
704: MatFDColoringSetFromOptions(dmmg[i]->fdcoloring);
705: }
706: #if defined(PETSC_HAVE_ADIC)
707: } else if (jacobian == DMMGComputeJacobianWithAdic) {
708: for (i=0; i<nlevels; i++) {
709: ISColoring iscoloring;
710: DMGetColoring(dmmg[i]->dm,IS_COLORING_GHOSTED,&iscoloring);
711: MatSetColoring(dmmg[i]->B,iscoloring);
712: ISColoringDestroy(iscoloring);
713: }
714: #endif
715: }
716: if (!useFAS) {
717: for (i=0; i<nlevels; i++) {
718: SNESSetJacobian(dmmg[i]->snes,dmmg[i]->J,dmmg[i]->B,DMMGComputeJacobian_Multigrid,dmmg);
719: SNESSetFromOptions(dmmg[i]->snes);
720: }
722: PetscOptionsGetInt(PETSC_NULL,"-dmmg_jacobian_period",&period,PETSC_NULL);
723: for (i=0; i<nlevels; i++) {
724: dmmg[i]->updatejacobian = PETSC_TRUE;
725: dmmg[i]->updatejacobianperiod = period;
726: }
727: }
729: /* Create interpolation scaling */
730: for (i=1; i<nlevels; i++) {
731: DMGetInterpolationScale(dmmg[i-1]->dm,dmmg[i]->dm,dmmg[i]->R,&dmmg[i]->Rscale);
732: }
734: if (useFAS) {
735: PetscTruth flg;
737: for(i = 0; i < nlevels; i++) {
738: VecDuplicate(dmmg[i]->b,&dmmg[i]->w);
739: if (cookie == DA_COOKIE) {
740: #if defined(PETSC_HAVE_ADIC)
741: NLFCreate_DAAD(&dmmg[i]->nlf);
742: NLFDAADSetDA_DAAD(dmmg[i]->nlf,(DA)dmmg[i]->dm);
743: NLFDAADSetCtx_DAAD(dmmg[i]->nlf,dmmg[i]->user);
744: NLFDAADSetResidual_DAAD(dmmg[i]->nlf,dmmg[i]->r);
745: NLFDAADSetNewtonIterations_DAAD(dmmg[i]->nlf,fasMaxIter);
746: #endif
747: } else {
748: }
750: dmmg[i]->monitor = monitor;
751: dmmg[i]->monitorall = monitorAll;
752: dmmg[i]->presmooth = fasPresmooth;
753: dmmg[i]->postsmooth = fasPostsmooth;
754: dmmg[i]->coarsesmooth = fasCoarsesmooth;
755: dmmg[i]->rtol = fasRtol;
756: dmmg[i]->abstol = fasAbstol;
757: }
759: PetscOptionsHasName(0, "-dmmg_fas_view", &flg);
760: if (flg) {
761: for(i = 0; i < nlevels; i++) {
762: if (i == 0) {
763: PetscPrintf(dmmg[i]->comm,"FAS Solver Parameters\n");
764: PetscPrintf(dmmg[i]->comm," rtol %G atol %G\n",dmmg[i]->rtol,dmmg[i]->abstol);
765: PetscPrintf(dmmg[i]->comm," coarsesmooths %D\n",dmmg[i]->coarsesmooth);
766: PetscPrintf(dmmg[i]->comm," Newton iterations %D\n",fasMaxIter);
767: } else {
768: PetscPrintf(dmmg[i]->comm," level %D presmooths %D\n",i,dmmg[i]->presmooth);
769: PetscPrintf(dmmg[i]->comm," postsmooths %D\n",dmmg[i]->postsmooth);
770: PetscPrintf(dmmg[i]->comm," Newton iterations %D\n",fasMaxIter);
771: }
772: if (fasBlock) {
773: PetscPrintf(dmmg[i]->comm," using point-block smoothing\n");
774: } else if(fasGMRES) {
775: PetscPrintf(dmmg[i]->comm," using non-linear gmres\n");
776: }
777: }
778: }
779: }
780: return(0);
781: }
785: /*@
786: DMMGSetSNESLocalFD - Sets the local user function that is used to approximately compute the Jacobian
787: via finite differences.
789: Collective on DMMG
791: Input Parameter:
792: + dmmg - the context
793: - function - the function that defines the nonlinear system
795: Level: intermediate
797: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGSetSNESLocal()
799: @*/
800: PetscErrorCode DMMGSetSNESLocalFD(DMMG *dmmg,DALocalFunction1 function)
801: {
802: PetscInt i,nlevels = dmmg[0]->nlevels;
805: for (i=0; i<nlevels; i++) {
806: dmmg[i]->lfj = (PetscErrorCode (*)(void))function;
807: }
808: return(0);
809: }
814: /*@
815: DMMGGetSNESLocal - Returns the local functions for residual and Jacobian evaluation.
817: Collective on DMMG
819: Input Parameter:
820: . dmmg - the context
822: Output Parameters:
823: + function - the function that defines the nonlinear system
824: - jacobian - function defines the local part of the Jacobian
826: Level: intermediate
828: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGSetSNESLocal()
829: @*/
830: PetscErrorCode DMMGGetSNESLocal(DMMG *dmmg,DALocalFunction1 *function, DALocalFunction1 *jacobian)
831: {
832: PetscCookie cookie;
836: PetscObjectGetCookie((PetscObject) dmmg[0]->dm, &cookie);
837: if (cookie == DA_COOKIE) {
838: DAGetLocalFunction((DA) dmmg[0]->dm, function);
839: DAGetLocalJacobian((DA) dmmg[0]->dm, jacobian);
840: } else {
841: #ifdef PETSC_HAVE_SIEVE
842: MeshGetLocalFunction((Mesh) dmmg[0]->dm, (PetscErrorCode (**)(Mesh,SectionReal,SectionReal,void*)) function);
843: MeshGetLocalJacobian((Mesh) dmmg[0]->dm, (PetscErrorCode (**)(Mesh,SectionReal,Mat,void*)) jacobian);
844: #else
845: SETERRQ(PETSC_ERR_SUP, "Unstructured grids only supported when Sieve is enabled.\nReconfigure with --with-sieve.");
846: #endif
847: }
848: return(0);
849: }
851: /*M
852: DMMGSetSNESLocal - Sets the local user function that defines the nonlinear set of equations
853: that will use the grid hierarchy and (optionally) its derivative.
855: Collective on DMMG
857: Synopsis:
858: PetscErrorCode DMMGSetSNESLocal(DMMG *dmmg,DALocalFunction1 function, DALocalFunction1 jacobian,
859: DALocalFunction1 ad_function, DALocalFunction1 admf_function);
861: Input Parameter:
862: + dmmg - the context
863: . function - the function that defines the nonlinear system
864: . jacobian - function defines the local part of the Jacobian
865: . ad_function - the name of the function with an ad_ prefix. This is ignored if ADIC is
866: not installed
867: - admf_function - the name of the function with an ad_ prefix. This is ignored if ADIC is
868: not installed
870: Options Database Keys:
871: + -dmmg_snes_monitor
872: . -dmmg_jacobian_fd
873: . -dmmg_jacobian_ad
874: . -dmmg_jacobian_mf_fd_operator
875: . -dmmg_jacobian_mf_fd
876: . -dmmg_jacobian_mf_ad_operator
877: . -dmmg_jacobian_mf_ad
878: - -dmmg_jacobian_period <p> - Indicates how often in the SNES solve the Jacobian is recomputed (on all levels)
879: as suggested by Florin Dobrian if p is -1 then Jacobian is computed only on first
880: SNES iteration (i.e. -1 is equivalent to infinity)
883: Level: intermediate
885: Notes:
886: If ADIC or ADIFOR have been installed, this routine can use ADIC or ADIFOR to compute
887: the derivative; however, that function cannot call other functions except those in
888: standard C math libraries.
890: If ADIC/ADIFOR have not been installed and the Jacobian is not provided, this routine
891: uses finite differencing to approximate the Jacobian.
893: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES()
895: M*/
899: PetscErrorCode DMMGSetSNESLocal_Private(DMMG *dmmg,DALocalFunction1 function,DALocalFunction1 jacobian,DALocalFunction1 ad_function,DALocalFunction1 admf_function)
900: {
902: PetscInt i,nlevels = dmmg[0]->nlevels;
903: PetscCookie cookie;
904: PetscErrorCode (*computejacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*) = 0;
908: if (jacobian) computejacobian = DMMGComputeJacobian;
909: #if defined(PETSC_HAVE_ADIC)
910: else if (ad_function) computejacobian = DMMGComputeJacobianWithAdic;
911: #endif
912: CHKMEMQ;
913: PetscObjectGetCookie((PetscObject) dmmg[0]->dm,&cookie);
914: if (cookie == DA_COOKIE) {
915: PetscTruth flag;
916: PetscOptionsHasName(PETSC_NULL, "-dmmg_form_function_ghost", &flag);
917: if (flag) {
918: DMMGSetSNES(dmmg,DMMGFormFunctionGhost,computejacobian);
919: } else {
920: DMMGSetSNES(dmmg,DMMGFormFunction,computejacobian);
921: }
922: for (i=0; i<nlevels; i++) {
923: DASetLocalFunction((DA)dmmg[i]->dm,function);
924: dmmg[i]->lfj = (PetscErrorCode (*)(void))function;
925: DASetLocalJacobian((DA)dmmg[i]->dm,jacobian);
926: DASetLocalAdicFunction((DA)dmmg[i]->dm,ad_function);
927: DASetLocalAdicMFFunction((DA)dmmg[i]->dm,admf_function);
928: }
929: CHKMEMQ;
930: } else {
931: #ifdef PETSC_HAVE_SIEVE
932: DMMGSetSNES(dmmg, DMMGFormFunctionMesh, DMMGComputeJacobianMesh);
933: for (i=0; i<nlevels; i++) {
934: MeshSetLocalFunction((Mesh) dmmg[i]->dm, (PetscErrorCode (*)(Mesh,SectionReal,SectionReal,void*)) function);
935: dmmg[i]->lfj = (PetscErrorCode (*)(void)) function;
936: MeshSetLocalJacobian((Mesh) dmmg[i]->dm, (PetscErrorCode (*)(Mesh,SectionReal,Mat,void*)) jacobian);
937: // Setup a work section
938: SectionReal defaultSec, constantSec;
939: PetscTruth hasConstant;
941: MeshGetSectionReal((Mesh) dmmg[i]->dm, "default", &defaultSec);
942: MeshHasSectionReal((Mesh) dmmg[i]->dm, "constant", &hasConstant);
943: if (!hasConstant) {
944: SectionRealDuplicate(defaultSec, &constantSec);
945: PetscObjectSetName((PetscObject) constantSec, "constant");
946: MeshSetSectionReal((Mesh) dmmg[i]->dm, constantSec);
947: SectionRealDestroy(constantSec);
948: }
949: }
950: CHKMEMQ;
951: #else
952: SETERRQ(PETSC_ERR_SUP, "Unstructured grids only supported when Sieve is enabled.\nReconfigure with --with-sieve.");
953: #endif
954: }
955: CHKMEMQ;
956: return(0);
957: }
961: PetscErrorCode DMMGFunctioni(void* ctx,PetscInt i,Vec u,PetscScalar* r)
962: {
963: DMMG dmmg = (DMMG)ctx;
964: Vec U = dmmg->lwork1;
966: VecScatter gtol;
969: /* copy u into interior part of U */
970: DAGetScatter((DA)dmmg->dm,0,>ol,0);
971: VecScatterBegin(gtol,u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL);
972: VecScatterEnd(gtol,u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL);
973: DAFormFunctioni1((DA)dmmg->dm,i,U,r,dmmg->user);
974: return(0);
975: }
979: PetscErrorCode DMMGFunctionib(PetscInt i,Vec u,PetscScalar* r,void* ctx)
980: {
981: DMMG dmmg = (DMMG)ctx;
982: Vec U = dmmg->lwork1;
984: VecScatter gtol;
987: /* copy u into interior part of U */
988: DAGetScatter((DA)dmmg->dm,0,>ol,0);
989: VecScatterBegin(gtol,u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL);
990: VecScatterEnd(gtol,u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL);
991: DAFormFunctionib1((DA)dmmg->dm,i,U,r,dmmg->user);
992: return(0);
993: }
997: PetscErrorCode DMMGFunctioniBase(void* ctx,Vec u)
998: {
999: DMMG dmmg = (DMMG)ctx;
1000: Vec U = dmmg->lwork1;
1004: DAGlobalToLocalBegin((DA)dmmg->dm,u,INSERT_VALUES,U);
1005: DAGlobalToLocalEnd((DA)dmmg->dm,u,INSERT_VALUES,U);
1006: return(0);
1007: }
1011: 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*))
1012: {
1014: PetscInt i,nlevels = dmmg[0]->nlevels;
1017: for (i=0; i<nlevels; i++) {
1018: DASetLocalFunctioni((DA)dmmg[i]->dm,functioni);
1019: DASetLocalAdicFunctioni((DA)dmmg[i]->dm,adi);
1020: DASetLocalAdicMFFunctioni((DA)dmmg[i]->dm,adimf);
1021: MatMFFDSetFunctioni(dmmg[i]->J,DMMGFunctioni);
1022: MatMFFDSetFunctioniBase(dmmg[i]->J,DMMGFunctioniBase);
1023: if (!dmmg[i]->lwork1) {
1024: DACreateLocalVector((DA)dmmg[i]->dm,&dmmg[i]->lwork1);
1025: }
1026: }
1027: return(0);
1028: }
1032: 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*))
1033: {
1035: PetscInt i,nlevels = dmmg[0]->nlevels;
1038: for (i=0; i<nlevels; i++) {
1039: DASetLocalFunctionib((DA)dmmg[i]->dm,functioni);
1040: DASetLocalAdicFunctionib((DA)dmmg[i]->dm,adi);
1041: DASetLocalAdicMFFunctionib((DA)dmmg[i]->dm,adimf);
1042: if (!dmmg[i]->lwork1) {
1043: DACreateLocalVector((DA)dmmg[i]->dm,&dmmg[i]->lwork1);
1044: }
1045: }
1046: return(0);
1047: }
1049: static PetscErrorCode (*localfunc)(void) = 0;
1053: /*
1054: Uses the DM object to call the user provided function with the correct calling
1055: sequence.
1056: */
1057: PetscErrorCode DMMGInitialGuess_Local(DMMG dmmg,Vec x)
1058: {
1062: (*dmmg->dm->ops->forminitialguess)(dmmg->dm,localfunc,x,0);
1063: return(0);
1064: }
1068: /*@C
1069: DMMGSetInitialGuessLocal - sets code to compute the initial guess for each level
1071: Collective on DMMG
1073: Input Parameter:
1074: + dmmg - the context
1075: - localguess - the function that computes the initial guess
1077: Level: intermediate
1079: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGSetInitialGuess(), DMMGSetSNESLocal()
1081: @*/
1082: PetscErrorCode DMMGSetInitialGuessLocal(DMMG *dmmg,PetscErrorCode (*localguess)(void))
1083: {
1087: localfunc = localguess; /* stash into ugly static for now */
1089: DMMGSetInitialGuess(dmmg,DMMGInitialGuess_Local);
1090: return(0);
1091: }
1095: /*@C
1096: DMMGSetISColoringType - type of coloring used to compute Jacobian via finite differencing
1098: Collective on DMMG
1100: Input Parameter:
1101: + dmmg - the context
1102: - isctype - IS_COLORING_GHOSTED (default) or IS_COLORING_GLOBAL
1104: Options Database:
1105: . -dmmg_iscoloring_type <ghosted or global>
1107: Notes: ghosted coloring requires using DMMGSetSNESLocal()
1109: Level: intermediate
1111: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGSetInitialGuess(), DMMGSetSNESLocal()
1113: @*/
1114: PetscErrorCode DMMGSetISColoringType(DMMG *dmmg,ISColoringType isctype)
1115: {
1117: dmmg[0]->isctype = isctype;
1118: return(0);
1119: }