Actual source code: damgsnes.c
1: /*$Id: damgsnes.c,v 1.51 2001/08/06 21:18:06 bsmith Exp $*/
2:
3: #include petscda.h
4: #include petscmg.h
7: /*
8: period of -1 indicates update only on zeroth iteration of SNES
9: */
10: #define ShouldUpdate(l,it) (((dmmg[l-1]->updatejacobianperiod == -1) && (it == 0)) || \
11: ((dmmg[l-1]->updatejacobianperiod > 0) && !(it % dmmg[l-1]->updatejacobianperiod)))
12: /*
13: Evaluates the Jacobian on all of the grids. It is used by DMMG to provide the
14: ComputeJacobian() function that SNESSetJacobian() requires.
15: */
18: int DMMGComputeJacobian_Multigrid(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
19: {
20: DMMG *dmmg = (DMMG*)ptr;
21: int ierr,i,nlevels = dmmg[0]->nlevels,it;
22: SLES sles,lsles;
23: PC pc;
24: PetscTruth ismg;
25: Vec W;
26: MatStructure flg;
29: if (!dmmg) SETERRQ(1,"Passing null as user context which should contain DMMG");
30: SNESGetIterationNumber(snes,&it);
32: /* compute Jacobian on finest grid */
33: if (dmmg[nlevels-1]->updatejacobian && ShouldUpdate(nlevels,it)) {
34: (*DMMGGetFine(dmmg)->computejacobian)(snes,X,J,B,flag,DMMGGetFine(dmmg));
35: } else {
36: PetscLogInfo(0,"DMMGComputeJacobian_Multigrid:Skipping Jacobian, SNES iteration %d frequence %d level %d\n",it,dmmg[nlevels-1]->updatejacobianperiod,nlevels-1);
37: *flag = SAME_PRECONDITIONER;
38: }
39: MatSNESMFSetBase(DMMGGetFine(dmmg)->J,X);
41: /* create coarser grid Jacobians for preconditioner if multigrid is the preconditioner */
42: SNESGetSLES(snes,&sles);
43: SLESGetPC(sles,&pc);
44: PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
45: if (ismg) {
47: MGGetSmoother(pc,nlevels-1,&lsles);
48: SLESSetOperators(lsles,DMMGGetFine(dmmg)->J,DMMGGetFine(dmmg)->B,*flag);
50: for (i=nlevels-1; i>0; i--) {
52: if (!dmmg[i-1]->w) {
53: VecDuplicate(dmmg[i-1]->x,&dmmg[i-1]->w);
54: }
56: W = dmmg[i-1]->w;
57: /* restrict X to coarser grid */
58: MatRestrict(dmmg[i]->R,X,W);
59: X = W;
61: /* scale to "natural" scaling for that grid */
62: VecPointwiseMult(dmmg[i]->Rscale,X,X);
64: /* tell the base vector for matrix free multiplies */
65: MatSNESMFSetBase(dmmg[i-1]->J,X);
67: /* compute Jacobian on coarse grid */
68: if (dmmg[i-1]->updatejacobian && ShouldUpdate(i,it)) {
69: (*dmmg[i-1]->computejacobian)(snes,X,&dmmg[i-1]->J,&dmmg[i-1]->B,&flg,dmmg[i-1]);
70: } else {
71: PetscLogInfo(0,"DMMGComputeJacobian_Multigrid:Skipping Jacobian, SNES iteration %d frequence %d level %d\n",it,dmmg[i-1]->updatejacobianperiod,i-1);
72: flg = SAME_PRECONDITIONER;
73: }
75: MGGetSmoother(pc,i-1,&lsles);
76: SLESSetOperators(lsles,dmmg[i-1]->J,dmmg[i-1]->B,flg);
77: }
78: }
79: return(0);
80: }
82: /* ---------------------------------------------------------------------------*/
86: /*
87: DMMGFormFunction - This is a universal global FormFunction used by the DMMG code
88: when the user provides a local function.
90: Input Parameters:
91: + snes - the SNES context
92: . X - input vector
93: - ptr - optional user-defined context, as set by SNESSetFunction()
95: Output Parameter:
96: . F - function vector
98: */
99: int DMMGFormFunction(SNES snes,Vec X,Vec F,void *ptr)
100: {
101: DMMG dmmg = (DMMG)ptr;
102: int ierr;
103: Vec localX;
104: DA da = (DA)dmmg->dm;
107: DAGetLocalVector(da,&localX);
108: /*
109: Scatter ghost points to local vector, using the 2-step process
110: DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
111: */
112: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
113: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
114: DAFormFunction1(da,localX,F,dmmg->user);
115: DARestoreLocalVector(da,&localX);
116: return(0);
117: }
121: /*@C
122: SNESDAFormFunction - This is a universal function evaluation routine that
123: may be used with SNESSetFunction() as long as the user context has a DA
124: as its first record and the user has called DASetLocalFunction().
126: Collective on SNES
128: Input Parameters:
129: + snes - the SNES context
130: . X - input vector
131: . F - function vector
132: - ptr - pointer to a structure that must have a DA as its first entry. For example this
133: could be a DMMG
135: Level: intermediate
137: .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(),
138: SNESSetFunction(), SNESSetJacobian()
140: @*/
141: int SNESDAFormFunction(SNES snes,Vec X,Vec F,void *ptr)
142: {
143: int ierr;
144: Vec localX;
145: DA da = *(DA*)ptr;
148: DAGetLocalVector(da,&localX);
149: /*
150: Scatter ghost points to local vector, using the 2-step process
151: DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
152: */
153: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
154: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
155: DAFormFunction1(da,localX,F,ptr);
156: DARestoreLocalVector(da,&localX);
157: return(0);
158: }
160: /* ---------------------------------------------------------------------------------------------------------------------------*/
164: int DMMGComputeJacobianWithFD(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
165: {
166: int ierr;
167: DMMG dmmg = (DMMG)ctx;
168:
170: SNESDefaultComputeJacobianColor(snes,x1,J,B,flag,dmmg->fdcoloring);
171: return(0);
172: }
176: int DMMGComputeJacobianWithMF(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
177: {
178: int ierr;
179:
181: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
182: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
183: return(0);
184: }
188: /*
189: DMMGComputeJacobianWithAdic - Evaluates the Jacobian via Adic when the user has provided
190: a local function evaluation routine.
191: */
192: int DMMGComputeJacobianWithAdic(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
193: {
194: DMMG dmmg = (DMMG) ptr;
195: int ierr;
196: Vec localX;
197: DA da = (DA) dmmg->dm;
200: DAGetLocalVector(da,&localX);
201: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
202: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
203: DAComputeJacobian1WithAdic(da,localX,*B,dmmg->user);
204: DARestoreLocalVector(da,&localX);
205: /* Assemble true Jacobian; if it is different */
206: if (*J != *B) {
207: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
208: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
209: }
210: MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR);
211: *flag = SAME_NONZERO_PATTERN;
212: return(0);
213: }
217: /*@
218: SNESDAComputeJacobianWithAdic - This is a universal Jacobian evaluation routine
219: that may be used with SNESSetJacobian() as long as the user context has a DA as
220: its first record and DASetLocalAdicFunction() has been called.
222: Collective on SNES
224: Input Parameters:
225: + snes - the SNES context
226: . X - input vector
227: . J - Jacobian
228: . B - Jacobian used in preconditioner (usally same as J)
229: . flag - indicates if the matrix changed its structure
230: - ptr - optional user-defined context, as set by SNESSetFunction()
232: Level: intermediate
234: .seealso: DASetLocalFunction(), DASetLocalAdicFunction(), SNESSetFunction(), SNESSetJacobian()
236: @*/
237: int SNESDAComputeJacobianWithAdic(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
238: {
239: DA da = *(DA*) ptr;
240: int ierr;
241: Vec localX;
244: DAGetLocalVector(da,&localX);
245: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
246: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
247: DAComputeJacobian1WithAdic(da,localX,*B,ptr);
248: DARestoreLocalVector(da,&localX);
249: /* Assemble true Jacobian; if it is different */
250: if (*J != *B) {
251: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
252: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
253: }
254: MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR);
255: *flag = SAME_NONZERO_PATTERN;
256: return(0);
257: }
261: /*
262: SNESDAComputeJacobianWithAdifor - This is a universal Jacobian evaluation routine
263: that may be used with SNESSetJacobian() from Fortran as long as the user context has
264: a DA as its first record and DASetLocalAdiforFunction() has been called.
266: Collective on SNES
268: Input Parameters:
269: + snes - the SNES context
270: . X - input vector
271: . J - Jacobian
272: . B - Jacobian used in preconditioner (usally same as J)
273: . flag - indicates if the matrix changed its structure
274: - ptr - optional user-defined context, as set by SNESSetFunction()
276: Level: intermediate
278: .seealso: DASetLocalFunction(), DASetLocalAdicFunction(), SNESSetFunction(), SNESSetJacobian()
280: */
281: int SNESDAComputeJacobianWithAdifor(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
282: {
283: DA da = *(DA*) ptr;
284: int ierr;
285: Vec localX;
288: DAGetLocalVector(da,&localX);
289: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
290: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
291: DAComputeJacobian1WithAdifor(da,localX,*B,ptr);
292: DARestoreLocalVector(da,&localX);
293: /* Assemble true Jacobian; if it is different */
294: if (*J != *B) {
295: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
296: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
297: }
298: MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR);
299: *flag = SAME_NONZERO_PATTERN;
300: return(0);
301: }
305: /*
306: SNESDAComputeJacobian - This is a universal Jacobian evaluation routine for a
307: locally provided Jacobian.
309: Collective on SNES
311: Input Parameters:
312: + snes - the SNES context
313: . X - input vector
314: . J - Jacobian
315: . B - Jacobian used in preconditioner (usally same as J)
316: . flag - indicates if the matrix changed its structure
317: - ptr - optional user-defined context, as set by SNESSetFunction()
319: Level: intermediate
321: .seealso: DASetLocalFunction(), DASetLocalJacobian(), SNESSetFunction(), SNESSetJacobian()
323: */
324: int SNESDAComputeJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
325: {
326: DA da = *(DA*) ptr;
327: int ierr;
328: Vec localX;
331: DAGetLocalVector(da,&localX);
332: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
333: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
334: DAComputeJacobian1(da,localX,*B,ptr);
335: DARestoreLocalVector(da,&localX);
336: /* Assemble true Jacobian; if it is different */
337: if (*J != *B) {
338: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
339: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
340: }
341: MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR);
342: *flag = SAME_NONZERO_PATTERN;
343: return(0);
344: }
348: int DMMGSolveSNES(DMMG *dmmg,int level)
349: {
350: int ierr,nlevels = dmmg[0]->nlevels,its;
353: dmmg[0]->nlevels = level+1;
354: SNESSolve(dmmg[level]->snes,dmmg[level]->x,&its);
355: dmmg[0]->nlevels = nlevels;
356: return(0);
357: }
359: EXTERN_C_BEGIN
360: extern int NLFCreate_DAAD(NLF*);
361: extern int NLFRelax_DAAD(NLF,MatSORType,int,Vec);
362: extern int NLFDAADSetDA_DAAD(NLF,DA);
363: extern int NLFDAADSetCtx_DAAD(NLF,void*);
364: extern int NLFDAADSetResidual_DAAD(NLF,Vec);
365: extern int NLFDAADSetNewtonIterations_DAAD(NLF,int);
366: EXTERN_C_END
368: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
369: #include src/sles/pc/impls/mg/mgimpl.h
370: /*
371: This is pre-beta FAS code. It's design should not be taken seriously!
372: */
375: int DMMGSolveFAS(DMMG *dmmg,int level)
376: {
377: int ierr,i,j,k;
378: PetscReal norm;
379: PetscScalar zero = 0.0,mone = -1.0,one = 1.0;
380: MG *mg;
381: PC pc;
382: SLES sles;
385: VecSet(&zero,dmmg[level]->r);
386: for (j=1; j<=level; j++) {
387: if (!dmmg[j]->inject) {
388: DMGetInjection(dmmg[j-1]->dm,dmmg[j]->dm,&dmmg[j]->inject);
389: }
390: }
392: SNESGetSLES(dmmg[level]->snes,&sles);
393: SLESGetPC(sles,&pc);
394: mg = ((MG*)pc->data);
396: for (i=0; i<100; i++) {
398: for (j=level; j>0; j--) {
400: /* Relax residual_fine - F(x_fine) = 0 */
401: for (k=0; k<dmmg[j]->presmooth; k++) {
402: NLFRelax_DAAD(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);
403: }
405: /* R*(residual_fine - F(x_fine)) */
406: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
407: VecAYPX(&mone,dmmg[j]->r,dmmg[j]->w);
409: if (j == level || dmmg[j]->monitorall) {
410: /* norm( residual_fine - f(x_fine) ) */
411: VecNorm(dmmg[j]->w,NORM_2,&norm);
412: if (j == level) {
413: if (norm < dmmg[level]->atol) goto theend;
414: if (i == 0) {
415: dmmg[level]->rrtol = norm*dmmg[level]->rtol;
416: } else {
417: if (norm < dmmg[level]->rrtol) goto theend;
418: }
419: }
420: }
422: if (dmmg[j]->monitorall) {
423: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
424: PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
425: }
426: MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);
427:
428: /* F(R*x_fine) */
429: VecScatterBegin(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);
430: VecScatterEnd(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);
431: DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);
433: /* residual_coarse = F(R*x_fine) + R*(residual_fine - F(x_fine)) */
434: VecAYPX(&one,dmmg[j-1]->w,dmmg[j-1]->r);
436: /* save R*x_fine into b (needed when interpolating compute x back up */
437: VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);
438: }
440: for (j=0; j<dmmg[0]->presmooth; j++) {
441: NLFRelax_DAAD(dmmg[0]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[0]->x);
442: }
443: if (dmmg[0]->monitorall){
444: DMMGFormFunction(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);
445: VecAXPY(&mone,dmmg[0]->r,dmmg[0]->w);
446: VecNorm(dmmg[0]->w,NORM_2,&norm);
447: for (k=0; k<level+1; k++) {PetscPrintf(dmmg[0]->comm," ");}
448: PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %g\n",norm);
449: }
451: for (j=1; j<=level; j++) {
452: /* x_fine = x_fine + R'*(x_coarse - R*x_fine) */
453: VecAXPY(&mone,dmmg[j-1]->b,dmmg[j-1]->x);
454: MatInterpolateAdd(mg[j]->restrct,dmmg[j-1]->x,dmmg[j]->x,dmmg[j]->x);
456: if (dmmg[j]->monitorall) {
457: /* norm( F(x_fine) - residual_fine ) */
458: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
459: VecAXPY(&mone,dmmg[j]->r,dmmg[j]->w);
460: VecNorm(dmmg[j]->w,NORM_2,&norm);
461: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
462: PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
463: }
465: /* Relax residual_fine - F(x_fine) = 0 */
466: for (k=0; k<dmmg[j]->postsmooth; k++) {
467: NLFRelax_DAAD(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);
468: }
470: if (dmmg[j]->monitorall) {
471: /* norm( F(x_fine) - residual_fine ) */
472: DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
473: VecAXPY(&mone,dmmg[j]->r,dmmg[j]->w);
474: VecNorm(dmmg[j]->w,NORM_2,&norm);
475: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
476: PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
477: }
478: }
480: if (dmmg[level]->monitor){
481: DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
482: VecNorm(dmmg[level]->w,NORM_2,&norm);
483: PetscPrintf(dmmg[level]->comm,"%d FAS function norm %g\n",i,norm);
484: }
485: }
486: theend:
487: return(0);
488: }
489: #endif
491: /* ===========================================================================================================*/
495: /*@C
496: DMMGSetSNES - Sets the nonlinear function that defines the nonlinear set of equations
497: to be solved using the grid hierarchy.
499: Collective on DMMG
501: Input Parameter:
502: + dmmg - the context
503: . function - the function that defines the nonlinear system
504: - jacobian - optional function to compute Jacobian
506: Options Database Keys:
507: + -dmmg_snes_monitor
508: . -dmmg_jacobian_fd
509: . -dmmg_jacobian_ad
510: . -dmmg_jacobian_mf_fd_operator
511: . -dmmg_jacobian_mf_fd
512: . -dmmg_jacobian_mf_ad_operator
513: . -dmmg_jacobian_mf_ad
514: - -dmmg_jacobian_period <p> - Indicates how often in the SNES solve the Jacobian is recomputed (on all levels)
515: as suggested by Florin Dobrian if p is -1 then Jacobian is computed only on first
516: SNES iteration (i.e. -1 is equivalent to infinity)
518: Level: advanced
520: .seealso DMMGCreate(), DMMGDestroy, DMMGSetSLES(), DMMGSetSNESLocal()
522: @*/
523: int DMMGSetSNES(DMMG *dmmg,int (*function)(SNES,Vec,Vec,void*),int (*jacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
524: {
525: int ierr,size,i,nlevels = dmmg[0]->nlevels,period = 1;
526: PetscTruth snesmonitor,mffdoperator,mffd,fdjacobian;
527: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
528: PetscTruth mfadoperator,mfad,adjacobian;
529: #endif
530: SLES sles;
531: PetscViewer ascii;
532: MPI_Comm comm;
535: if (!dmmg) SETERRQ(1,"Passing null as DMMG");
536: if (!jacobian) jacobian = DMMGComputeJacobianWithFD;
538: PetscOptionsBegin(dmmg[0]->comm,PETSC_NULL,"DMMG Options","SNES");
539: PetscOptionsName("-dmmg_snes_monitor","Monitor nonlinear convergence","SNESSetMonitor",&snesmonitor);
542: PetscOptionsName("-dmmg_jacobian_fd","Compute sparse Jacobian explicitly with finite differencing","DMMGSetSNES",&fdjacobian);
543: if (fdjacobian) jacobian = DMMGComputeJacobianWithFD;
544: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
545: PetscOptionsName("-dmmg_jacobian_ad","Compute sparse Jacobian explicitly with ADIC (automatic differentiation)","DMMGSetSNES",&adjacobian);
546: if (adjacobian) jacobian = DMMGComputeJacobianWithAdic;
547: #endif
549: PetscOptionsLogicalGroupBegin("-dmmg_jacobian_mf_fd_operator","Apply Jacobian via matrix free finite differencing","DMMGSetSNES",&mffdoperator);
550: PetscOptionsLogicalGroupEnd("-dmmg_jacobian_mf_fd","Apply Jacobian via matrix free finite differencing even in computing preconditioner","DMMGSetSNES",&mffd);
551: if (mffd) mffdoperator = PETSC_TRUE;
552: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
553: PetscOptionsLogicalGroupBegin("-dmmg_jacobian_mf_ad_operator","Apply Jacobian via matrix free ADIC (automatic differentiation)","DMMGSetSNES",&mfadoperator);
554: PetscOptionsLogicalGroupEnd("-dmmg_jacobian_mf_ad","Apply Jacobian via matrix free ADIC (automatic differentiation) even in computing preconditioner","DMMGSetSNES",&mfad);
555: if (mfad) mfadoperator = PETSC_TRUE;
556: #endif
557: PetscOptionsEnd();
559: /* create solvers for each level */
560: for (i=0; i<nlevels; i++) {
561: SNESCreate(dmmg[i]->comm,&dmmg[i]->snes);
562: if (snesmonitor) {
563: PetscObjectGetComm((PetscObject)dmmg[i]->snes,&comm);
564: PetscViewerASCIIOpen(comm,"stdout",&ascii);
565: PetscViewerASCIISetTab(ascii,nlevels-i);
566: SNESSetMonitor(dmmg[i]->snes,SNESDefaultMonitor,ascii,(int(*)(void*))PetscViewerDestroy);
567: }
569: if (mffdoperator) {
570: MatCreateSNESMF(dmmg[i]->snes,dmmg[i]->x,&dmmg[i]->J);
571: VecDuplicate(dmmg[i]->x,&dmmg[i]->work1);
572: VecDuplicate(dmmg[i]->x,&dmmg[i]->work2);
573: MatSNESMFSetFunction(dmmg[i]->J,dmmg[i]->work1,function,dmmg[i]);
574: if (mffd) {
575: dmmg[i]->B = dmmg[i]->J;
576: jacobian = DMMGComputeJacobianWithMF;
577: }
578: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
579: } else if (mfadoperator) {
580: MatRegisterDAAD();
581: MatCreateDAAD((DA)dmmg[i]->dm,&dmmg[i]->J);
582: MatDAADSetCtx(dmmg[i]->J,dmmg[i]->user);
583: if (mfad) {
584: dmmg[i]->B = dmmg[i]->J;
585: jacobian = DMMGComputeJacobianWithMF;
586: }
587: #endif
588: }
589:
590: if (!dmmg[i]->B) {
591: MPI_Comm_size(dmmg[i]->comm,&size);
592: DMGetMatrix(dmmg[i]->dm,MATAIJ,&dmmg[i]->B);
593: }
594: if (!dmmg[i]->J) {
595: dmmg[i]->J = dmmg[i]->B;
596: }
598: SNESGetSLES(dmmg[i]->snes,&sles);
599: DMMGSetUpLevel(dmmg,sles,i+1);
600:
601: /*
602: if the number of levels is > 1 then we want the coarse solve in the grid sequencing to use LU
603: when possible
604: */
605: if (nlevels > 1 && i == 0) {
606: PC pc;
607: SLES csles;
608: PetscTruth flg1,flg2,flg3;
610: SLESGetPC(sles,&pc);
611: MGGetCoarseSolve(pc,&csles);
612: SLESGetPC(csles,&pc);
613: PetscTypeCompare((PetscObject)pc,PCILU,&flg1);
614: PetscTypeCompare((PetscObject)pc,PCSOR,&flg2);
615: PetscTypeCompare((PetscObject)pc,PETSC_NULL,&flg3);
616: if (flg1 || flg2 || flg3) {
617: PCSetType(pc,PCLU);
618: }
619: }
621: SNESSetFromOptions(dmmg[i]->snes);
622: dmmg[i]->solve = DMMGSolveSNES;
623: dmmg[i]->computejacobian = jacobian;
624: dmmg[i]->computefunction = function;
625: }
628: if (jacobian == DMMGComputeJacobianWithFD) {
629: ISColoring iscoloring;
630: for (i=0; i<nlevels; i++) {
631: DMGetColoring(dmmg[i]->dm,IS_COLORING_LOCAL,&iscoloring);
632: MatFDColoringCreate(dmmg[i]->B,iscoloring,&dmmg[i]->fdcoloring);
633: ISColoringDestroy(iscoloring);
634: MatFDColoringSetFunction(dmmg[i]->fdcoloring,(int(*)(void))function,dmmg[i]);
635: MatFDColoringSetFromOptions(dmmg[i]->fdcoloring);
636: }
637: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
638: } else if (jacobian == DMMGComputeJacobianWithAdic) {
639: for (i=0; i<nlevels; i++) {
640: ISColoring iscoloring;
641: DMGetColoring(dmmg[i]->dm,IS_COLORING_GHOSTED,&iscoloring);
642: MatSetColoring(dmmg[i]->B,iscoloring);
643: ISColoringDestroy(iscoloring);
644: }
645: #endif
646: }
648: for (i=0; i<nlevels; i++) {
649: SNESSetJacobian(dmmg[i]->snes,dmmg[i]->J,dmmg[i]->B,DMMGComputeJacobian_Multigrid,dmmg);
650: SNESSetFunction(dmmg[i]->snes,dmmg[i]->b,function,dmmg[i]);
651: }
653: /* Create interpolation scaling */
654: for (i=1; i<nlevels; i++) {
655: DMGetInterpolationScale(dmmg[i-1]->dm,dmmg[i]->dm,dmmg[i]->R,&dmmg[i]->Rscale);
656: }
658: PetscOptionsGetInt(PETSC_NULL,"-dmmg_jacobian_period",&period,PETSC_NULL);
659: for (i=0; i<nlevels; i++) {
660: dmmg[i]->updatejacobian = PETSC_TRUE;
661: dmmg[i]->updatejacobianperiod = period;
662: }
664: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
665: {
666: PetscTruth flg;
667: PetscOptionsHasName(PETSC_NULL,"-dmmg_fas",&flg);
668: if (flg) {
669: int newton_its;
670: PetscOptionsHasName(0,"-fas_view",&flg);
671: for (i=0; i<nlevels; i++) {
672: NLFCreate_DAAD(&dmmg[i]->nlf);
673: NLFDAADSetDA_DAAD(dmmg[i]->nlf,(DA)dmmg[i]->dm);
674: NLFDAADSetCtx_DAAD(dmmg[i]->nlf,dmmg[i]->user);
675: NLFDAADSetResidual_DAAD(dmmg[i]->nlf,dmmg[i]->r);
676: VecDuplicate(dmmg[i]->b,&dmmg[i]->w);
678: dmmg[i]->monitor = PETSC_FALSE;
679: PetscOptionsHasName(0,"-dmmg_fas_monitor",&dmmg[i]->monitor);
680: dmmg[i]->monitorall = PETSC_FALSE;
681: PetscOptionsHasName(0,"-dmmg_fas_monitor_all",&dmmg[i]->monitorall);
682: dmmg[i]->presmooth = 2;
683: PetscOptionsGetInt(0,"-dmmg_fas_presmooth",&dmmg[i]->presmooth,0);
684: dmmg[i]->postsmooth = 2;
685: PetscOptionsGetInt(0,"-dmmg_fas_postsmooth",&dmmg[i]->postsmooth,0);
686: dmmg[i]->coarsesmooth = 2;
687: PetscOptionsGetInt(0,"-dmmg_fas_coarsesmooth",&dmmg[i]->coarsesmooth,0);
689: dmmg[i]->rtol = 1.e-8;
690: PetscOptionsGetReal(0,"-dmmg_fas_rtol",&dmmg[i]->rtol,0);
691: dmmg[i]->atol = 1.e-50;
692: PetscOptionsGetReal(0,"-dmmg_fas_atol",&dmmg[i]->atol,0);
694: newton_its = 2;
695: PetscOptionsGetInt(0,"-dmmg_fas_newton_its",&newton_its,0);
696: NLFDAADSetNewtonIterations_DAAD(dmmg[i]->nlf,newton_its);
698: if (flg) {
699: if (i == 0) {
700: PetscPrintf(dmmg[i]->comm,"FAS Solver Parameters\n");
701: PetscPrintf(dmmg[i]->comm," rtol %g atol %g\n",dmmg[i]->rtol,dmmg[i]->atol);
702: PetscPrintf(dmmg[i]->comm," coarsesmooths %d\n",dmmg[i]->coarsesmooth);
703: PetscPrintf(dmmg[i]->comm," Newton iterations %d\n",newton_its);
704: } else {
705: PetscPrintf(dmmg[i]->comm," level %d presmooths %d\n",i,dmmg[i]->presmooth);
706: PetscPrintf(dmmg[i]->comm," postsmooths %d\n",dmmg[i]->postsmooth);
707: PetscPrintf(dmmg[i]->comm," Newton iterations %d\n",newton_its);
708: }
709: }
710: dmmg[i]->solve = DMMGSolveFAS;
711: }
712: }
713: }
714: #endif
715:
716: return(0);
717: }
721: /*@C
722: DMMGSetInitialGuess - Sets the function that computes an initial guess, if not given
723: uses 0.
725: Collective on DMMG and SNES
727: Input Parameter:
728: + dmmg - the context
729: - guess - the function
731: Level: advanced
733: .seealso DMMGCreate(), DMMGDestroy, DMMGSetSLES()
735: @*/
736: int DMMGSetInitialGuess(DMMG *dmmg,int (*guess)(SNES,Vec,void*))
737: {
738: int i,nlevels = dmmg[0]->nlevels;
741: for (i=0; i<nlevels; i++) {
742: dmmg[i]->initialguess = guess;
743: }
744: return(0);
745: }
747: /*M
748: DMMGSetSNESLocal - Sets the local user function that defines the nonlinear set of equations
749: that will use the grid hierarchy and (optionally) its derivative.
751: Collective on DMMG
753: Synopsis:
754: int DMMGSetSNESLocal(DMMG *dmmg,DALocalFunction1 function, DALocalFunction1 jacobian,
755: DALocalFunction1 ad_function, DALocalFunction1 admf_function);
757: Input Parameter:
758: + dmmg - the context
759: . function - the function that defines the nonlinear system
760: . jacobian - function defines the local part of the Jacobian (not currently supported)
761: . ad_function - the name of the function with an ad_ prefix. This is ignored if ADIC is
762: not installed
763: - admf_function - the name of the function with an ad_ prefix. This is ignored if ADIC is
764: not installed
766: Options Database Keys:
767: + -dmmg_snes_monitor
768: . -dmmg_jacobian_fd
769: . -dmmg_jacobian_ad
770: . -dmmg_jacobian_mf_fd_operator
771: . -dmmg_jacobian_mf_fd
772: . -dmmg_jacobian_mf_ad_operator
773: . -dmmg_jacobian_mf_ad
774: - -dmmg_jacobian_period <p> - Indicates how often in the SNES solve the Jacobian is recomputed (on all levels)
775: as suggested by Florin Dobrian if p is -1 then Jacobian is computed only on first
776: SNES iteration (i.e. -1 is equivalent to infinity)
779: Level: intermediate
781: Notes:
782: If ADIC or ADIFOR have been installed, this routine can use ADIC or ADIFOR to compute
783: the derivative; however, that function cannot call other functions except those in
784: standard C math libraries.
786: If ADIC/ADIFOR have not been installed and the Jacobian is not provided, this routine
787: uses finite differencing to approximate the Jacobian.
789: .seealso DMMGCreate(), DMMGDestroy, DMMGSetSLES(), DMMGSetSNES()
791: M*/
795: int DMMGSetSNESLocal_Private(DMMG *dmmg,DALocalFunction1 function,DALocalFunction1 jacobian,DALocalFunction1 ad_function,DALocalFunction1 admf_function)
796: {
797: int ierr,i,nlevels = dmmg[0]->nlevels;
798: int (*computejacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*) = 0;
802: if (jacobian) computejacobian = SNESDAComputeJacobian;
803: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
804: else if (ad_function) computejacobian = DMMGComputeJacobianWithAdic;
805: #endif
807: DMMGSetSNES(dmmg,DMMGFormFunction,computejacobian);
808: for (i=0; i<nlevels; i++) {
809: DASetLocalFunction((DA)dmmg[i]->dm,function);
810: DASetLocalJacobian((DA)dmmg[i]->dm,jacobian);
811: DASetLocalAdicFunction((DA)dmmg[i]->dm,ad_function);
812: DASetLocalAdicMFFunction((DA)dmmg[i]->dm,admf_function);
813: }
814: return(0);
815: }
819: static int DMMGFunctioni(int i,Vec u,PetscScalar* r,void* ctx)
820: {
821: DMMG dmmg = (DMMG)ctx;
822: Vec U = dmmg->lwork1;
823: int ierr;
824: VecScatter gtol;
827: /* copy u into interior part of U */
828: DAGetScatter((DA)dmmg->dm,0,>ol,0);
829: VecScatterBegin(u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL,gtol);
830: VecScatterEnd(u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL,gtol);
832: DAFormFunctioni1((DA)dmmg->dm,i,U,r,dmmg->user);
833: return(0);
834: }
838: static int DMMGFunctioniBase(Vec u,void* ctx)
839: {
840: DMMG dmmg = (DMMG)ctx;
841: Vec U = dmmg->lwork1;
842: int ierr;
845: DAGlobalToLocalBegin((DA)dmmg->dm,u,INSERT_VALUES,U);
846: DAGlobalToLocalEnd((DA)dmmg->dm,u,INSERT_VALUES,U);
847: return(0);
848: }
852: int DMMGSetSNESLocali_Private(DMMG *dmmg,int (*functioni)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*),int (*adi)(DALocalInfo*,MatStencil*,void*,void*,void*),int (*adimf)(DALocalInfo*,MatStencil*,void*,void*,void*))
853: {
854: int ierr,i,nlevels = dmmg[0]->nlevels;
857: for (i=0; i<nlevels; i++) {
858: DASetLocalFunctioni((DA)dmmg[i]->dm,functioni);
859: DASetLocalAdicFunctioni((DA)dmmg[i]->dm,adi);
860: DASetLocalAdicMFFunctioni((DA)dmmg[i]->dm,adimf);
861: MatSNESMFSetFunctioni(dmmg[i]->J,DMMGFunctioni);
862: MatSNESMFSetFunctioniBase(dmmg[i]->J,DMMGFunctioniBase);
863: DACreateLocalVector((DA)dmmg[i]->dm,&dmmg[i]->lwork1);
864: }
865: return(0);
866: }
869: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
870: EXTERN_C_BEGIN
871: #include "adic/ad_utils.h"
872: EXTERN_C_END
876: int PetscADView(int N,int nc,double *ptr,PetscViewer viewer)
877: {
878: int i,j,nlen = PetscADGetDerivTypeSize();
879: char *cptr = (char*)ptr;
880: double *values;
883: for (i=0; i<N; i++) {
884: printf("Element %d value %g derivatives: ",i,*(double*)cptr);
885: values = PetscADGetGradArray(cptr);
886: for (j=0; j<nc; j++) {
887: printf("%g ",*values++);
888: }
889: printf("\n");
890: cptr += nlen;
891: }
893: return(0);
894: }
896: #endif