Actual source code: damgsnes.c
1: /*$Id: damgsnes.c,v 1.14 2001/04/10 19:37:08 bsmith Exp bsmith $*/
2:
3: #include petscda.h
4: #include petscmg.h
6: /*
7: These evaluate the Jacobian on all of the grids. It is used by DMMG to "replace"
8: the user provided Jacobian function. In fact, it calls the user provided one at each level.
9: */
10: /*
11: Version for matrix-free Jacobian
12: */
13: int DMMGComputeJacobian_MF(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
14: {
15: DMMG *dmmg = (DMMG*)ptr;
16: int ierr,i,nlevels = dmmg[0]->nlevels;
17: SLES sles,lsles;
18: PC pc;
19: PetscTruth ismg;
22: if (!dmmg) SETERRQ(1,"Passing null as user context which should contain DMMG");
24: /* The finest level matrix is "shared" by the corresponding SNES object so we need
25: only call MatAssemblyXXX() on it to indicate it is being used in a new solve */
26: MatAssemblyBegin(dmmg[nlevels-1]->J,MAT_FINAL_ASSEMBLY);
27: MatAssemblyEnd(dmmg[nlevels-1]->J,MAT_FINAL_ASSEMBLY);
29: /*
30: The other levels MUST be told the vector from which we are doing the differencing
31: */
32: SNESGetSLES(snes,&sles);
33: SLESGetPC(sles,&pc);
34: PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
35: if (ismg) {
37: MGGetSmoother(pc,nlevels-1,&lsles);
38: SLESSetOperators(lsles,DMMGGetFine(dmmg)->J,DMMGGetFine(dmmg)->J,*flag);
40: for (i=nlevels-1; i>0; i--) {
42: /* restrict X to coarse grid */
43: MatRestrict(dmmg[i]->R,X,dmmg[i-1]->work2);
44: X = dmmg[i-1]->work2;
46: /* scale to "natural" scaling for that grid */
47: VecPointwiseMult(dmmg[i]->Rscale,X,X);
49: MatSNESMFSetBase(dmmg[i-1]->J,X);
50: MatAssemblyBegin(dmmg[i-1]->J,MAT_FINAL_ASSEMBLY);
51: MatAssemblyEnd(dmmg[i-1]->J,MAT_FINAL_ASSEMBLY);
53: MGGetSmoother(pc,i-1,&lsles);
54: SLESSetOperators(lsles,dmmg[i-1]->J,dmmg[i-1]->B,*flag);
55: }
56: }
57: return(0);
58: }
60: /*
61: Version for user provided Jacobian
62: */
63: int DMMGComputeJacobian_User(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
64: {
65: DMMG *dmmg = (DMMG*)ptr;
66: int ierr,i,nlevels = dmmg[0]->nlevels;
67: SLES sles,lsles;
68: PC pc;
69: PetscTruth ismg;
72: if (!dmmg) SETERRQ(1,"Passing null as user context which should contain DMMG");
74: (*DMMGGetFine(dmmg)->computejacobian)(snes,X,J,B,flag,DMMGGetFine(dmmg));
76: /* create coarse grid jacobian for preconditioner if multigrid is the preconditioner */
77: SNESGetSLES(snes,&sles);
78: SLESGetPC(sles,&pc);
79: PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
80: if (ismg) {
82: MGGetSmoother(pc,nlevels-1,&lsles);
83: SLESSetOperators(lsles,DMMGGetFine(dmmg)->J,DMMGGetFine(dmmg)->J,*flag);
85: for (i=nlevels-1; i>0; i--) {
87: /* restrict X to coarse grid */
88: MatRestrict(dmmg[i]->R,X,dmmg[i-1]->x);
89: X = dmmg[i-1]->x;
91: /* scale to "natural" scaling for that grid */
92: VecPointwiseMult(dmmg[i]->Rscale,X,X);
94: /* form Jacobian on coarse grid */
95: (*dmmg[i-1]->computejacobian)(snes,X,&dmmg[i-1]->J,&dmmg[i-1]->B,flag,dmmg[i-1]);
97: MGGetSmoother(pc,i-1,&lsles);
98: SLESSetOperators(lsles,dmmg[i-1]->J,dmmg[i-1]->B,*flag);
99: }
100: }
101: return(0);
102: }
103: /*
104: Version for Jacobian computed via PETSc finite differencing. This is the same
105: as DMMGComputeJacobian_User() except passes in the fdcoloring as the private context
106: */
107: int DMMGComputeJacobian_FD(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
108: {
109: DMMG *dmmg = (DMMG*)ptr;
110: int ierr,i,nlevels = dmmg[0]->nlevels;
111: SLES sles,lsles;
112: PC pc;
113: PetscTruth ismg;
116: if (!dmmg) SETERRQ(1,"Passing null as user context which should contain DMMG");
118: (*DMMGGetFine(dmmg)->computejacobian)(snes,X,J,B,flag,DMMGGetFine(dmmg)->fdcoloring);
120: /* create coarse grid jacobian for preconditioner if multigrid is the preconditioner */
121: SNESGetSLES(snes,&sles);
122: SLESGetPC(sles,&pc);
123: PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
124: if (ismg) {
126: MGGetSmoother(pc,nlevels-1,&lsles);
127: SLESSetOperators(lsles,DMMGGetFine(dmmg)->J,DMMGGetFine(dmmg)->J,*flag);
129: for (i=nlevels-1; i>0; i--) {
131: /* restrict X to coarse grid */
132: MatRestrict(dmmg[i]->R,X,dmmg[i-1]->x);
133: X = dmmg[i-1]->x;
135: /* scale to "natural" scaling for that grid */
136: VecPointwiseMult(dmmg[i]->Rscale,X,X);
138: /* form Jacobian on coarse grid */
139: (*dmmg[i-1]->computejacobian)(snes,X,&dmmg[i-1]->J,&dmmg[i-1]->B,flag,dmmg[i-1]->fdcoloring);
141: MGGetSmoother(pc,i-1,&lsles);
142: SLESSetOperators(lsles,dmmg[i-1]->J,dmmg[i-1]->B,*flag);
143: }
144: }
145: return(0);
146: }
148: int DMMGSolveSNES(DMMG *dmmg,int level)
149: {
150: int ierr,nlevels = dmmg[0]->nlevels,its;
153: dmmg[0]->nlevels = level+1;
154: SNESSolve(dmmg[level]->snes,dmmg[level]->x,&its);
155: dmmg[0]->nlevels = nlevels;
156: return(0);
157: }
159: EXTERN int DMMGSetUpLevel(DMMG*,SLES,int);
161: /*@C
162: DMMGSetSNES - Sets the nonlinear function that defines the nonlinear set of equations
163: to be solved will use the grid hierarchy
165: Collective on DMMG
167: Input Parameter:
168: + dmmg - the context
169: . function - the function that defines the nonlinear system
170: - jacobian - optional function to compute Jacobian
172: Level: advanced
174: .seealso DMMGCreate(), DMMGDestroy, DMMGSetSLES()
176: @*/
177: int DMMGSetSNES(DMMG *dmmg,int (*function)(SNES,Vec,Vec,void*),int (*jacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
178: {
179: int ierr,i,nlevels = dmmg[0]->nlevels;
180: PetscTruth usefd,snesmonitor;
181: SLES sles;
182: PetscViewer ascii;
183: MPI_Comm comm;
186: if (!dmmg) SETERRQ(1,"Passing null as DMMG");
188: PetscOptionsHasName(PETSC_NULL,"-dmmg_snes_monitor",&snesmonitor);
189: /* create solvers for each level */
190: for (i=0; i<nlevels; i++) {
191: SNESCreate(dmmg[i]->comm,SNES_NONLINEAR_EQUATIONS,&dmmg[i]->snes);
192: if (snesmonitor) {
193: PetscObjectGetComm((PetscObject)dmmg[i]->snes,&comm);
194: PetscViewerASCIIOpen(comm,"stdout",&ascii);
195: PetscViewerASCIISetTab(ascii,nlevels-i);
196: SNESSetMonitor(dmmg[i]->snes,SNESDefaultMonitor,ascii,(int(*)(void*))PetscViewerDestroy);
197: }
198: if (dmmg[0]->matrixfree) {
199: MatCreateSNESMF(dmmg[i]->snes,dmmg[i]->x,&dmmg[i]->J);
200: if (!dmmg[i]->B) dmmg[i]->B = dmmg[i]->J;
201: if (i != nlevels-1) {
202: VecDuplicate(dmmg[i]->x,&dmmg[i]->work1);
203: VecDuplicate(dmmg[i]->x,&dmmg[i]->work2);
204: MatSNESMFSetFunction(dmmg[i]->J,dmmg[i]->work1,function,dmmg[i]);
205: }
206: }
208: SNESGetSLES(dmmg[i]->snes,&sles);
209: DMMGSetUpLevel(dmmg,sles,i+1);
210: SNESSetFromOptions(dmmg[i]->snes);
211: dmmg[i]->solve = DMMGSolveSNES;
212: dmmg[i]->computejacobian = jacobian;
213: dmmg[i]->computefunction = function;
214: }
216: PetscOptionsHasName(PETSC_NULL,"-dmmg_fd",&usefd);
217: if ((!jacobian && !dmmg[0]->matrixfree) || usefd) {
218: ISColoring iscoloring;
219: for (i=0; i<nlevels; i++) {
220: DMGetColoring(dmmg[i]->dm,IS_COLORING_GLOBAL,MATMPIAIJ,&iscoloring,PETSC_NULL);
221: MatFDColoringCreate(dmmg[i]->J,iscoloring,&dmmg[i]->fdcoloring);
222: ISColoringDestroy(iscoloring);
223: MatFDColoringSetFunction(dmmg[i]->fdcoloring,(int(*)(void))function,dmmg[i]);
224: MatFDColoringSetFromOptions(dmmg[i]->fdcoloring);
225: dmmg[i]->computejacobian = SNESDefaultComputeJacobianColor;
226: }
227: }
229: for (i=0; i<nlevels; i++) {
230: if (dmmg[i]->matrixfree) {
231: SNESSetJacobian(dmmg[i]->snes,dmmg[i]->J,dmmg[i]->B,DMMGComputeJacobian_MF,dmmg);
232: } else if (dmmg[i]->computejacobian == SNESDefaultComputeJacobianColor) {
233: SNESSetJacobian(dmmg[i]->snes,dmmg[i]->J,dmmg[i]->B,DMMGComputeJacobian_FD,dmmg);
234: } else {
235: SNESSetJacobian(dmmg[i]->snes,dmmg[i]->J,dmmg[i]->B,DMMGComputeJacobian_User,dmmg);
236: }
237: SNESSetFunction(dmmg[i]->snes,dmmg[i]->b,function,dmmg[i]);
238: }
240: /* Create interpolation scaling */
241: for (i=1; i<nlevels; i++) {
242: DMGetInterpolationScale(dmmg[i-1]->dm,dmmg[i]->dm,dmmg[i]->R,&dmmg[i]->Rscale);
243: }
245: for (i=0; i<nlevels-1; i++) {
246: SNESSetOptionsPrefix(dmmg[i]->snes,"dmmg_levels_");
247: }
249: return(0);
250: }
252: /*@C
253: DMMGSetInitialGuess - Sets the function that computes an initial guess, if not given
254: uses 0.
256: Collective on DMMG and SNES
258: Input Parameter:
259: + dmmg - the context
260: - guess - the function
262: Level: advanced
264: .seealso DMMGCreate(), DMMGDestroy, DMMGSetSLES()
266: @*/
267: int DMMGSetInitialGuess(DMMG *dmmg,int (*guess)(SNES,Vec,void*))
268: {
269: int i,nlevels = dmmg[0]->nlevels;
272: for (i=0; i<nlevels; i++) {
273: dmmg[i]->initialguess = guess;
274: }
275: return(0);
276: }
279: /*
280: DMMGFormFunction - This is a universal global FormFunction used by the DMMG code
281: when the user provides a local function.
283: Input Parameters:
284: . snes - the SNES context
285: . X - input vector
286: . ptr - optional user-defined context, as set by SNESSetFunction()
288: Output Parameter:
289: . F - function vector
291: */
292: int DMMGFormFunction(SNES snes,Vec X,Vec F,void *ptr)
293: {
294: DMMG dmmg = (DMMG)ptr;
295: int ierr;
296: Scalar **x,**f;
297: Vec localX;
298: DA da = (DA)dmmg->dm;
299: DALocalInfo info;
302: DAGetLocalVector((DA)dmmg->dm,&localX);
303: DAGetLocalInfo(da,&info);
305: /*
306: Scatter ghost points to local vector, using the 2-step process
307: DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
308: */
309: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
310: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
312: /*
313: Get pointers to vector data
314: */
315: DAVecGetArray((DA)dmmg->dm,localX,(void**)&x);
316: DAVecGetArray((DA)dmmg->dm,F,(void**)&f);
318: /*
319: Compute function over the locally owned part of the grid
320: */
321: (*dmmg->computefunctionlocal)(x,f,&info,dmmg->user);
323: /*
324: Restore vectors
325: */
326: DAVecRestoreArray((DA)dmmg->dm,localX,(void**)&x);
327: DAVecRestoreArray((DA)dmmg->dm,F,(void**)&f);
329: DARestoreLocalVector((DA)dmmg->dm,&localX);
331: return 0;
332: }
334: /*@C
335: DMMGSetSNESLocal - Sets the local user function that defines the nonlinear set of equations
336: that will use the grid hierarchy
338: Collective on DMMG
340: Input Parameter:
341: + dmmg - the context
342: . function - the function that defines the nonlinear system
343: - jacobian - optional function to compute Jacobian
346: Level: advanced
348: .seealso DMMGCreate(), DMMGDestroy, DMMGSetSLES()
350: @*/
351: int DMMGSetSNESLocal(DMMG *dmmg,int (*function)(Scalar **,Scalar**,DALocalInfo*,void*),int (*jacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
352: {
353: int ierr,i,nlevels = dmmg[0]->nlevels;
356: DMMGSetSNES(dmmg,DMMGFormFunction,jacobian);
357: for (i=0; i<nlevels; i++) {
358: dmmg[i]->computefunctionlocal = function;
359: }
360: return(0);
361: }
364: /* ---------------------------------------------------------------------------------------------------------------------------*/
366: typedef struct {
367: double value;
368: double grad[20];
369: } DERIV_TYPE;
371: #define DERIV_val(a) ((a).value)
372: #define DERIV_grad(a) ((a).grad)
374: /* ------------------------------------------------------------------- */
375: /*
376: DMMGFormJacobianWithAD - Evaluates the Jacobian via AD when the user has provide
377: a local form function
378: */
379: int DMMGFormJacobianWithAD(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
380: {
381: DMMG dmmg = (DMMG) ptr;
382: int ierr,i,j,k,l,row,col[25],color,*colors = dmmg->iscoloring->colors;
383: int xs,ys,xm,ym;
384: int gxs,gys,gxm,gym;
385: int mx,my;
386: Vec localX;
387: Scalar **x;
388: Mat jac = *B; /* Jacobian matrix */
389: DALocalInfo info;
391: DERIV_TYPE **ad_x, **ad_f, *ad_ptr, dummy;
392: Scalar av[25];
393: int dim,dof,stencil_size,stencil_width;
394: DAStencilType stencil_type;
395: DAPeriodicType periodicity;
397: DA da= (DA) dmmg->dm;
400: DAGetLocalInfo(da,&info);
403: DAGetLocalVector(da,&localX);
404: DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
405: DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
407: DAGetInfo(da,&dim,&mx,&my,0,0,0,0,&dof,&stencil_width,&periodicity,&stencil_type);
409: /* Verify that this DA type is supported */
410: if ((dim != 2) || (stencil_width != 1) || (stencil_type != DA_STENCIL_STAR)
411: || (periodicity != DA_NONPERIODIC)) {
412: SETERRQ(0,"This DA type is not yet supported. Sorry.n");
413: }
415: stencil_size = 5;
417: /*
418: Get local grid boundaries
419: */
420: DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
421: DAGetGhostCorners(da,&gxs,&gys,PETSC_NULL,&gxm,&gym,PETSC_NULL);
423: /* printf("myid = %d x[%d,%d], y[%d,%d], gx[%d,%d], gy[%d,%d]n",myid); */
427: /*
428: Get pointer to vector data
429: */
430: DAVecGetArray(da,localX,(void**)&x);
432: /* allocate space for derivative objects. Should be able to use something
433: like:
435: DAVecGetStructArray(da,localX,(void **)&ad_x,sizeof(DERIV_TYPE));
436: DAVecGetStructArray(da,F,(void **)&ad_f,sizeof(DERIV_TYPE));
438: instead. */
440: PetscMalloc((gym)*sizeof(DERIV_TYPE *),(void **)&ad_x);
441: ad_x -= gys;
442: for(j=gys;j<gys+gym;j++) {
443: PetscMalloc((dof*gxm)*sizeof(DERIV_TYPE),(void **)&ad_ptr);
444: ad_x[j] = ad_ptr - dof*gxs;
445: for(i=dof*gxs;i<dof*(gxs+gxm);i++) DERIV_val(ad_x[j][i]) = x[j][i];
446: }
448: PetscMalloc((ym)*sizeof(DERIV_TYPE *),(void**)&ad_f);
449: ad_f -= ys;
450: for(j=ys;j<ys+ym;j++) {
451: PetscMalloc((dof*xm)*sizeof(DERIV_TYPE),(void**)&ad_ptr);
452: ad_f[j] = ad_ptr - dof*xs;
453: }
455: ad_AD_Init(stencil_size*dof);
457: /* Fake ADIC into thinking there are stencil_size*dof independent variables - fix this */
458: for(i=0;i<stencil_size*dof;i++) ad_AD_SetIndep(dummy);
459: ad_AD_SetIndepDone();
461: /* Initialize seed matrix, using coloring */
462: for(j=gys;j<gys+gym;j++) {
463: for(i=gxs;i<gxs+gxm;i++) {
464: for(k=0;k<dof;k++) {
465: color = colors[(j-gys)*gxm+i-gxs];
466: ad_AD_ClearGrad(DERIV_grad(ad_x[j][i*dof+k]));
467: DERIV_grad(ad_x[j][i*dof+k])[color] = 1.0;
468: }
469: }
470: }
472: /*
473: Compute entries for the locally owned part of the Jacobian.
474: */
475: (*dmmg->ad_computefunctionlocal)((Scalar**)ad_x,(Scalar**)ad_f,&info,dmmg->user);
477: DAVecRestoreArray(da,localX,(void**)&x);
478: DARestoreLocalVector(da,&localX);
480: /* stick the values into the matrix */
481: MatADSetColoring_SeqAIJ(B,dmmg->iscoloring);
482: MatADSetValues_SeqAIJ(B,(Scalar**)ad_f,20);
484: /* Assemble true Jacobian; if it is different */
485: ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
486: ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
487: ierr = MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR);
488: *flag = SAME_NONZERO_PATTERN;
489: return(0);
490: }
493: /*@C
494: DMMGSetSNESLocalWithAD - Sets the local user function that defines the nonlinear set of equations
495: that will use the grid hierarchy and its AD derivate function
497: Collective on DMMG
499: Input Parameter:
500: + dmmg - the context
501: . function - the function that defines the nonlinear system
502: - jacobian - AD function to compute Jacobian
504: Level: advanced
506: .seealso DMMGCreate(), DMMGDestroy, DMMGSetSLES()
508: @*/
509: int DMMGSetSNESLocalWithAD(DMMG *dmmg,int (*function)(Scalar **,Scalar**,DALocalInfo*,void*),int (*ad_function)(Scalar **,Scalar**,DALocalInfo*,void*))
510: {
511: int ierr,i,nlevels = dmmg[0]->nlevels;
514: DMMGSetSNES(dmmg,DMMGFormFunction,DMMGFormJacobianWithAD);
515: for (i=0; i<nlevels; i++) {
516: dmmg[i]->computefunctionlocal = function;
517: dmmg[i]->ad_computefunctionlocal = ad_function;
518: }
519: return(0);
520: }