Actual source code: dmdats.c
petsc-dev 2014-02-02
1: #include <petscdmda.h> /*I "petscdmda.h" I*/
2: #include <petsc-private/dmimpl.h>
3: #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/
4: #include <petscdraw.h>
6: /* This structure holds the user-provided DMDA callbacks */
7: typedef struct {
8: PetscErrorCode (*ifunctionlocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*);
9: PetscErrorCode (*rhsfunctionlocal)(DMDALocalInfo*,PetscReal,void*,void*,void*);
10: PetscErrorCode (*ijacobianlocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*);
11: PetscErrorCode (*rhsjacobianlocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*);
12: void *ifunctionlocalctx;
13: void *ijacobianlocalctx;
14: void *rhsfunctionlocalctx;
15: void *rhsjacobianlocalctx;
16: InsertMode ifunctionlocalimode;
17: InsertMode rhsfunctionlocalimode;
18: } DMTS_DA;
22: static PetscErrorCode DMTSDestroy_DMDA(DMTS sdm)
23: {
27: PetscFree(sdm->data);
28: return(0);
29: }
33: static PetscErrorCode DMTSDuplicate_DMDA(DMTS oldsdm,DMTS sdm)
34: {
38: PetscNewLog(sdm,(DMTS_DA**)&sdm->data);
39: if (oldsdm->data) {PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMTS_DA));}
40: return(0);
41: }
45: static PetscErrorCode DMDATSGetContext(DM dm,DMTS sdm,DMTS_DA **dmdats)
46: {
50: *dmdats = NULL;
51: if (!sdm->data) {
52: PetscNewLog(dm,(DMTS_DA**)&sdm->data);
53: sdm->ops->destroy = DMTSDestroy_DMDA;
54: sdm->ops->duplicate = DMTSDuplicate_DMDA;
55: }
56: *dmdats = (DMTS_DA*)sdm->data;
57: return(0);
58: }
62: /*
63: This function should eventually replace:
64: DMDAComputeFunction() and DMDAComputeFunction1()
65: */
66: static PetscErrorCode TSComputeIFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,Vec F,void *ctx)
67: {
69: DM dm;
70: DMTS_DA *dmdats = (DMTS_DA*)ctx;
71: DMDALocalInfo info;
72: Vec Xloc;
73: void *x,*f,*xdot;
79: if (!dmdats->ifunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context");
80: TSGetDM(ts,&dm);
81: DMGetLocalVector(dm,&Xloc);
82: DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);
83: DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);
84: DMDAGetLocalInfo(dm,&info);
85: DMDAVecGetArray(dm,Xloc,&x);
86: DMDAVecGetArray(dm,Xdot,&xdot);
87: switch (dmdats->ifunctionlocalimode) {
88: case INSERT_VALUES: {
89: DMDAVecGetArray(dm,F,&f);
90: CHKMEMQ;
91: (*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx);
92: CHKMEMQ;
93: DMDAVecRestoreArray(dm,F,&f);
94: } break;
95: case ADD_VALUES: {
96: Vec Floc;
97: DMGetLocalVector(dm,&Floc);
98: VecZeroEntries(Floc);
99: DMDAVecGetArray(dm,Floc,&f);
100: CHKMEMQ;
101: (*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx);
102: CHKMEMQ;
103: DMDAVecRestoreArray(dm,Floc,&f);
104: VecZeroEntries(F);
105: DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);
106: DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);
107: DMRestoreLocalVector(dm,&Floc);
108: } break;
109: default: SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->ifunctionlocalimode);
110: }
111: DMDAVecRestoreArray(dm,Xloc,&x);
112: DMRestoreLocalVector(dm,&Xloc);
113: DMDAVecRestoreArray(dm,Xdot,&xdot);
114: return(0);
115: }
119: static PetscErrorCode TSComputeIJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
120: {
122: DM dm;
123: DMTS_DA *dmdats = (DMTS_DA*)ctx;
124: DMDALocalInfo info;
125: Vec Xloc;
126: void *x,*xdot;
129: if (!dmdats->ifunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context");
130: TSGetDM(ts,&dm);
132: if (dmdats->ijacobianlocal) {
133: DMGetLocalVector(dm,&Xloc);
134: DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);
135: DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);
136: DMDAGetLocalInfo(dm,&info);
137: DMDAVecGetArray(dm,Xloc,&x);
138: DMDAVecGetArray(dm,Xdot,&xdot);
139: CHKMEMQ;
140: (*dmdats->ijacobianlocal)(&info,ptime,x,xdot,shift,*A,*B,mstr,dmdats->ijacobianlocalctx);
141: CHKMEMQ;
142: DMDAVecRestoreArray(dm,Xloc,&x);
143: DMDAVecRestoreArray(dm,Xdot,&xdot);
144: DMRestoreLocalVector(dm,&Xloc);
145: } else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"TSComputeIJacobian_DMDA() called without calling DMDATSSetIJacobian()");
146: /* This will be redundant if the user called both, but it's too common to forget. */
147: if (*A != *B) {
148: MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
149: MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
150: }
151: return(0);
152: }
156: /*
157: This function should eventually replace:
158: DMDAComputeFunction() and DMDAComputeFunction1()
159: */
160: static PetscErrorCode TSComputeRHSFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec F,void *ctx)
161: {
163: DM dm;
164: DMTS_DA *dmdats = (DMTS_DA*)ctx;
165: DMDALocalInfo info;
166: Vec Xloc;
167: void *x,*f;
173: if (!dmdats->rhsfunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context");
174: TSGetDM(ts,&dm);
175: DMGetLocalVector(dm,&Xloc);
176: DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);
177: DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);
178: DMDAGetLocalInfo(dm,&info);
179: DMDAVecGetArray(dm,Xloc,&x);
180: switch (dmdats->rhsfunctionlocalimode) {
181: case INSERT_VALUES: {
182: DMDAVecGetArray(dm,F,&f);
183: CHKMEMQ;
184: (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);
185: CHKMEMQ;
186: DMDAVecRestoreArray(dm,F,&f);
187: } break;
188: case ADD_VALUES: {
189: Vec Floc;
190: DMGetLocalVector(dm,&Floc);
191: VecZeroEntries(Floc);
192: DMDAVecGetArray(dm,Floc,&f);
193: CHKMEMQ;
194: (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);
195: CHKMEMQ;
196: DMDAVecRestoreArray(dm,Floc,&f);
197: VecZeroEntries(F);
198: DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);
199: DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);
200: DMRestoreLocalVector(dm,&Floc);
201: } break;
202: default: SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->rhsfunctionlocalimode);
203: }
204: DMDAVecRestoreArray(dm,Xloc,&x);
205: DMRestoreLocalVector(dm,&Xloc);
206: return(0);
207: }
211: static PetscErrorCode TSComputeRHSJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
212: {
214: DM dm;
215: DMTS_DA *dmdats = (DMTS_DA*)ctx;
216: DMDALocalInfo info;
217: Vec Xloc;
218: void *x;
221: if (!dmdats->rhsfunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context");
222: TSGetDM(ts,&dm);
224: if (dmdats->rhsjacobianlocal) {
225: DMGetLocalVector(dm,&Xloc);
226: DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);
227: DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);
228: DMDAGetLocalInfo(dm,&info);
229: DMDAVecGetArray(dm,Xloc,&x);
230: CHKMEMQ;
231: (*dmdats->rhsjacobianlocal)(&info,ptime,x,*A,*B,mstr,dmdats->rhsjacobianlocalctx);
232: CHKMEMQ;
233: DMDAVecRestoreArray(dm,Xloc,&x);
234: DMRestoreLocalVector(dm,&Xloc);
235: } else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"TSComputeRHSJacobian_DMDA() called without calling DMDATSSetRHSJacobian()");
236: /* This will be redundant if the user called both, but it's too common to forget. */
237: if (*A != *B) {
238: MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
239: MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
240: }
241: return(0);
242: }
247: /*@C
248: DMDATSSetRHSFunctionLocal - set a local residual evaluation function
250: Logically Collective
252: Input Arguments:
253: + dm - DM to associate callback with
254: . imode - insert mode for the residual
255: . func - local residual evaluation
256: - ctx - optional context for local residual evaluation
258: Calling sequence for func:
260: $ func(DMDALocalInfo info,PetscReal t,void *x,void *f,void *ctx)
262: + info - DMDALocalInfo defining the subdomain to evaluate the residual on
263: . t - time at which to evaluate residual
264: . x - array of local state information
265: . f - output array of local residual information
266: - ctx - optional user context
268: Level: beginner
270: .seealso: DMTSSetRHSFunction(), DMDATSSetRHSJacobianLocal(), DMDASNESSetFunctionLocal()
271: @*/
272: PetscErrorCode DMDATSSetRHSFunctionLocal(DM dm,InsertMode imode,DMDATSRHSFunctionLocal func,void *ctx)
273: {
275: DMTS sdm;
276: DMTS_DA *dmdats;
280: DMGetDMTSWrite(dm,&sdm);
281: DMDATSGetContext(dm,sdm,&dmdats);
282: dmdats->rhsfunctionlocalimode = imode;
283: dmdats->rhsfunctionlocal = func;
284: dmdats->rhsfunctionlocalctx = ctx;
285: DMTSSetRHSFunction(dm,TSComputeRHSFunction_DMDA,dmdats);
286: return(0);
287: }
291: /*@C
292: DMDATSSetRHSJacobianLocal - set a local residual evaluation function
294: Logically Collective
296: Input Arguments:
297: + dm - DM to associate callback with
298: . func - local RHS Jacobian evaluation routine
299: - ctx - optional context for local jacobian evaluation
301: Calling sequence for func:
303: $ func(DMDALocalInfo* info,PetscReal t,void* x,Mat J,Mat B,MatStructure *flg,void *ctx);
305: + info - DMDALocalInfo defining the subdomain to evaluate the residual on
306: . t - time at which to evaluate residual
307: . x - array of local state information
308: . J - Jacobian matrix
309: . B - preconditioner matrix; often same as J
310: . flg - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators())
311: - ctx - optional context passed above
313: Level: beginner
315: .seealso: DMTSSetRHSJacobian(), DMDATSSetRHSFunctionLocal(), DMDASNESSetJacobianLocal()
316: @*/
317: PetscErrorCode DMDATSSetRHSJacobianLocal(DM dm,DMDATSRHSJacobianLocal func,void *ctx)
318: {
320: DMTS sdm;
321: DMTS_DA *dmdats;
325: DMGetDMTSWrite(dm,&sdm);
326: DMDATSGetContext(dm,sdm,&dmdats);
327: dmdats->rhsjacobianlocal = func;
328: dmdats->rhsjacobianlocalctx = ctx;
329: DMTSSetRHSJacobian(dm,TSComputeRHSJacobian_DMDA,dmdats);
330: return(0);
331: }
336: /*@C
337: DMDATSSetIFunctionLocal - set a local residual evaluation function
339: Logically Collective
341: Input Arguments:
342: + dm - DM to associate callback with
343: . func - local residual evaluation
344: - ctx - optional context for local residual evaluation
346: Calling sequence for func:
347: + info - DMDALocalInfo defining the subdomain to evaluate the residual on
348: . t - time at which to evaluate residual
349: . x - array of local state information
350: . xdot - array of local time derivative information
351: . f - output array of local function evaluation information
352: - ctx - optional context passed above
354: Level: beginner
356: .seealso: DMTSSetIFunction(), DMDATSSetIJacobianLocal(), DMDASNESSetFunctionLocal()
357: @*/
358: PetscErrorCode DMDATSSetIFunctionLocal(DM dm,InsertMode imode,DMDATSIFunctionLocal func,void *ctx)
359: {
361: DMTS sdm;
362: DMTS_DA *dmdats;
366: DMGetDMTSWrite(dm,&sdm);
367: DMDATSGetContext(dm,sdm,&dmdats);
368: dmdats->ifunctionlocalimode = imode;
369: dmdats->ifunctionlocal = func;
370: dmdats->ifunctionlocalctx = ctx;
371: DMTSSetIFunction(dm,TSComputeIFunction_DMDA,dmdats);
372: return(0);
373: }
377: /*@C
378: DMDATSSetIJacobianLocal - set a local residual evaluation function
380: Logically Collective
382: Input Arguments:
383: + dm - DM to associate callback with
384: . func - local residual evaluation
385: - ctx - optional context for local residual evaluation
387: Calling sequence for func:
389: $ func(DMDALocalInfo* info,PetscReal t,void* x,void *xdot,Mat J,Mat B,MatStructure *flg,void *ctx);
391: + info - DMDALocalInfo defining the subdomain to evaluate the residual on
392: . t - time at which to evaluate the jacobian
393: . x - array of local state information
394: . xdot - time derivative at this state
395: . J - Jacobian matrix
396: . B - preconditioner matrix; often same as J
397: . flg - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators())
398: - ctx - optional context passed above
400: Level: beginner
402: .seealso: DMTSSetJacobian(), DMDATSSetIFunctionLocal(), DMDASNESSetJacobianLocal()
403: @*/
404: PetscErrorCode DMDATSSetIJacobianLocal(DM dm,DMDATSIJacobianLocal func,void *ctx)
405: {
407: DMTS sdm;
408: DMTS_DA *dmdats;
412: DMGetDMTSWrite(dm,&sdm);
413: DMDATSGetContext(dm,sdm,&dmdats);
414: dmdats->ijacobianlocal = func;
415: dmdats->ijacobianlocalctx = ctx;
416: DMTSSetIJacobian(dm,TSComputeIJacobian_DMDA,dmdats);
417: return(0);
418: }
422: PetscErrorCode TSMonitorDMDARayDestroy(void **mctx)
423: {
424: TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *) *mctx;
425: PetscErrorCode ierr;
428: if (rayctx->lgctx) {TSMonitorLGCtxDestroy(&rayctx->lgctx);}
429: VecDestroy(&rayctx->ray);
430: VecScatterDestroy(&rayctx->scatter);
431: PetscViewerDestroy(&rayctx->viewer);
432: PetscFree(rayctx);
433: return(0);
434: }
438: PetscErrorCode TSMonitorDMDARay(TS ts,PetscInt steps,PetscReal time,Vec u,void *mctx)
439: {
440: TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx*)mctx;
441: Vec solution;
442: PetscErrorCode ierr;
445: TSGetSolution(ts,&solution);
446: VecScatterBegin(rayctx->scatter,solution,rayctx->ray,INSERT_VALUES,SCATTER_FORWARD);
447: VecScatterEnd(rayctx->scatter,solution,rayctx->ray,INSERT_VALUES,SCATTER_FORWARD);
448: if (rayctx->viewer) {
449: VecView(rayctx->ray,rayctx->viewer);
450: }
451: return(0);
452: }
456: PetscErrorCode TSMonitorLGDMDARay(TS ts, PetscInt step, PetscReal ptime, Vec u, void *ctx)
457: {
458: TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *) ctx;
459: TSMonitorLGCtx lgctx = (TSMonitorLGCtx) rayctx->lgctx;
460: Vec v = rayctx->ray;
461: const PetscScalar *a;
462: PetscInt dim;
463: PetscErrorCode ierr;
466: VecScatterBegin(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD);
467: VecScatterEnd(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD);
468: if (!step) {
469: PetscDrawAxis axis;
471: PetscDrawLGGetAxis(lgctx->lg, &axis);
472: PetscDrawAxisSetLabels(axis, "Solution Ray as function of time", "Time", "Solution");
473: VecGetLocalSize(rayctx->ray, &dim);
474: PetscDrawLGSetDimension(lgctx->lg, dim);
475: PetscDrawLGReset(lgctx->lg);
476: }
477: VecGetArrayRead(v, &a);
478: #if defined(PETSC_USE_COMPLEX)
479: {
480: PetscReal *areal;
481: PetscInt i,n;
482: VecGetLocalSize(v, &n);
483: PetscMalloc1(n, &areal);
484: for (i = 0; i < n; ++i) areal[i] = PetscRealPart(a[i]);
485: PetscDrawLGAddCommonPoint(lgctx->lg, ptime, areal);
486: PetscFree(areal);
487: }
488: #else
489: PetscDrawLGAddCommonPoint(lgctx->lg, ptime, a);
490: #endif
491: VecRestoreArrayRead(v, &a);
492: if (((lgctx->howoften > 0) && (!(step % lgctx->howoften))) || ((lgctx->howoften == -1) && ts->reason)) {
493: PetscDrawLGDraw(lgctx->lg);
494: }
495: return(0);
496: }