Actual source code: heat.c
petsc-3.7.2 2016-06-05
2: static char help[] = "Solves heat equation in 1d.\n";
4: /*
5: Solves the equation
7: u_t = kappa \Delta u
8: Periodic boundary conditions
10: Evolve the heat equation:
11: ---------------
12: ./heat -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -wait -ts_type cn -da_refine 5 -mymonitor
14: Evolve the Allen-Cahn equation:
15: ---------------
16: ./heat -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -wait -ts_type cn -da_refine 5 -allen-cahn -kappa .001 -ts_final_time 5 -mymonitor
18: Evolve the Allen-Cahn equation: zoom in on part of the domain
19: ---------------
20: ./heat -ts_monitor -snes_monitor -pc_type lu -snes_converged_reason -ts_type cn -da_refine 5 -allen-cahn -kappa .001 -ts_final_time 5 -zoom .25,.45 -wait -mymonitor
23: The option -square_initial indicates it should use a square wave initial condition otherwise it loads the file InitialSolution.heat as the initial solution. You should run with
24: ./heat -square_initial -ts_monitor -snes_monitor -pc_type lu -snes_converged_reason -ts_type cn -da_refine 9 -ts_final_time 1.e-4 -ts_dt .125e-6 -snes_atol 1.e-25 -snes_rtol 1.e-25 -ts_max_steps 15
25: to generate binaryoutput then do mv binaryoutput InitialSolution.heat to obtain the initial solution file
27: */
28: #include <petscdm.h>
29: #include <petscdmda.h>
30: #include <petscts.h>
31: #include <petscdraw.h>
33: /*
34: User-defined routines
35: */
36: extern PetscErrorCode FormFunction(TS,PetscReal,Vec,Vec,void*),FormInitialSolution(DM,Vec),MyMonitor(TS,PetscInt,PetscReal,Vec,void*),MyDestroy(void**);
37: typedef struct {PetscReal kappa;PetscBool allencahn;PetscDrawViewPorts *ports;} UserCtx;
41: int main(int argc,char **argv)
42: {
43: TS ts; /* time integrator */
44: Vec x,r; /* solution, residual vectors */
45: Mat J; /* Jacobian matrix */
46: PetscInt steps,Mx, maxsteps = 1000000;
48: DM da;
49: PetscReal dt;
50: PetscReal vbounds[] = {-1.1,1.1};
51: PetscBool wait;
52: UserCtx ctx;
53: PetscBool mymonitor;
55: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56: Initialize program
57: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58: PetscInitialize(&argc,&argv,(char*)0,help);
59: ctx.kappa = 1.0;
60: PetscOptionsGetReal(NULL,"-kappa",&ctx.kappa,NULL);
61: ctx.allencahn = PETSC_FALSE;
62: PetscOptionsHasName(NULL,"-allen-cahn",&ctx.allencahn);
63: PetscOptionsHasName(NULL,"-mymonitor",&mymonitor);
65: PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,vbounds);
66: PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1200,800);
68: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69: Create distributed array (DMDA) to manage parallel grid and vectors
70: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71: DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, -10,1,2,NULL,&da);
72: DMDASetFieldName(da,0,"Heat equation: u");
73: DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);
74: dt = 1.0/(ctx.kappa*Mx*Mx);
76: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77: Extract global vectors from DMDA; then duplicate for remaining
78: vectors that are the same types
79: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
80: DMCreateGlobalVector(da,&x);
81: VecDuplicate(x,&r);
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Create timestepping solver context
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: TSCreate(PETSC_COMM_WORLD,&ts);
87: TSSetDM(ts,da);
88: TSSetProblemType(ts,TS_NONLINEAR);
89: TSSetRHSFunction(ts,NULL,FormFunction,&ctx);
91: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92: Customize nonlinear solver
93: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: TSSetType(ts,TSCN);
96: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: Set initial conditions
98: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: FormInitialSolution(da,x);
100: TSSetInitialTimeStep(ts,0.0,dt);
101: TSSetDuration(ts,maxsteps,.02);
102: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_INTERPOLATE);
103: TSSetSolution(ts,x);
106: if (mymonitor) {
107: ctx.ports = NULL;
108: TSMonitorSet(ts,MyMonitor,&ctx,MyDestroy);
109: }
111: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112: Set runtime options
113: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114: TSSetFromOptions(ts);
116: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117: Solve nonlinear system
118: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119: TSSolve(ts,x);
120: wait = PETSC_FALSE;
121: PetscOptionsGetBool(NULL,NULL,"-wait",&wait,NULL);
122: if (wait) {
123: PetscSleep(-1);
124: }
125: TSGetTimeStepNumber(ts,&steps);
126: VecView(x,PETSC_VIEWER_BINARY_WORLD);
128: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129: Free work space. All PETSc objects should be destroyed when they
130: are no longer needed.
131: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132: MatDestroy(&J);
133: VecDestroy(&x);
134: VecDestroy(&r);
135: TSDestroy(&ts);
136: DMDestroy(&da);
138: PetscFinalize();
139: return(0);
140: }
141: /* ------------------------------------------------------------------- */
144: /*
145: FormFunction - Evaluates nonlinear function, F(x).
147: Input Parameters:
148: . ts - the TS context
149: . X - input vector
150: . ptr - optional user-defined context, as set by SNESSetFunction()
152: Output Parameter:
153: . F - function vector
154: */
155: PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr)
156: {
157: DM da;
159: PetscInt i,Mx,xs,xm;
160: PetscReal hx,sx;
161: PetscScalar *x,*f;
162: Vec localX;
163: UserCtx *ctx = (UserCtx*)ptr;
166: TSGetDM(ts,&da);
167: DMGetLocalVector(da,&localX);
168: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
169: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
171: hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
173: /*
174: Scatter ghost points to local vector,using the 2-step process
175: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
176: By placing code between these two statements, computations can be
177: done while messages are in transition.
178: */
179: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
180: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
182: /*
183: Get pointers to vector data
184: */
185: DMDAVecGetArrayRead(da,localX,&x);
186: DMDAVecGetArray(da,F,&f);
188: /*
189: Get local grid boundaries
190: */
191: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
193: /*
194: Compute function over the locally owned part of the grid
195: */
196: for (i=xs; i<xs+xm; i++) {
197: f[i] = ctx->kappa*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
198: if (ctx->allencahn) f[i] += (x[i] - x[i]*x[i]*x[i]);
199: }
201: /*
202: Restore vectors
203: */
204: DMDAVecRestoreArrayRead(da,localX,&x);
205: DMDAVecRestoreArray(da,F,&f);
206: DMRestoreLocalVector(da,&localX);
207: return(0);
208: }
210: /* ------------------------------------------------------------------- */
213: PetscErrorCode FormInitialSolution(DM da,Vec U)
214: {
215: PetscErrorCode ierr;
216: PetscInt i,xs,xm,Mx,scale,N;
217: PetscScalar *u;
218: const PetscScalar *f;
219: PetscReal hx,x,r;
220: Vec finesolution;
221: PetscViewer viewer;
222: PetscBool flg;
225: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
226: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
228: hx = 1.0/(PetscReal)Mx;
230: /*
231: Get pointers to vector data
232: */
233: DMDAVecGetArray(da,U,&u);
235: /*
236: Get local grid boundaries
237: */
238: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
240: /* InitialSolution is obtained with
241: ./heat -ts_monitor -snes_monitor -pc_type lu -snes_converged_reason -ts_type cn -da_refine 9 -ts_final_time 1.e-4 -ts_dt .125e-6 -snes_atol 1.e-25 -snes_rtol 1.e-25 -ts_max_steps 15
242: */
244: /* PetscViewerBinaryOpen(PETSC_COMM_WORLD,"InitialSolution",FILE_MODE_READ,&viewer); */
245: /* VecCreate(PETSC_COMM_WORLD,&finesolution); */
246: /* VecLoad(finesolution,viewer); */
247: /* PetscViewerDestroy(&viewer); */
248: /* VecGetSize(finesolution,&N); */
249: /* scale = N/Mx; */
250: /* VecGetArrayRead(finesolution,&f); */
252: PetscOptionsHasName(NULL,"-square_initial",&flg);
253: if (!flg) {
254: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"InitialSolution.heat",FILE_MODE_READ,&viewer);
255: VecCreate(PETSC_COMM_WORLD,&finesolution);
256: VecLoad(finesolution,viewer);
257: PetscViewerDestroy(&viewer);
258: VecGetSize(finesolution,&N);
259: scale = N/Mx;
260: VecGetArrayRead(finesolution,&f);
261: }
263: /*
264: Compute function over the locally owned part of the grid
265: */
266: for (i=xs; i<xs+xm; i++) {
267: x = i*hx;
268: r = PetscSqrtScalar((x-.5)*(x-.5));
269: if (r < .125) u[i] = 1.0;
270: else u[i] = -.5;
272: /* With the initial condition above the method is first order in space */
273: /* this is a smooth initial condition so the method becomes second order in space */
274: /*u[i] = PetscSinScalar(2*PETSC_PI*x); */
275: /* u[i] = f[scale*i];*/
276: if (!flg) u[i] = f[scale*i];
277: }
278: if (!flg) {
279: VecRestoreArrayRead(finesolution,&f);
280: VecDestroy(&finesolution);
281: }
282: /* VecRestoreArrayRead(finesolution,&f);
283: VecDestroy(&finesolution);*/
285: /*
286: Restore vectors
287: */
288: DMDAVecRestoreArray(da,U,&u);
289: return(0);
290: }
294: /*
295: This routine is not parallel
296: */
297: PetscErrorCode MyMonitor(TS ts,PetscInt step,PetscReal time,Vec U,void *ptr)
298: {
299: UserCtx *ctx = (UserCtx*)ptr;
300: PetscDrawLG lg;
301: PetscErrorCode ierr;
302: PetscScalar *u;
303: PetscInt Mx,i,xs,xm,cnt;
304: PetscReal x,y,hx,pause,sx,len,max,xx[2],yy[2];
305: PetscDraw draw;
306: Vec localU;
307: DM da;
308: int colors[] = {PETSC_DRAW_YELLOW,PETSC_DRAW_RED,PETSC_DRAW_BLUE};
309: const char*const legend[] = {"-kappa (\\grad u,\\grad u)","(1 - u^2)^2"};
310: PetscDrawAxis axis;
311: PetscDrawViewPorts *ports;
315: TSGetDM(ts,&da);
316: DMGetLocalVector(da,&localU);
317: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
318: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
319: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
320: hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
321: DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);
322: DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);
323: DMDAVecGetArrayRead(da,localU,&u);
325: PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,&lg);
326: PetscDrawLGGetDraw(lg,&draw);
327: PetscDrawCheckResizedWindow(draw);
328: if (!ctx->ports) {
329: PetscDrawViewPortsCreateRect(draw,1,3,&ctx->ports);
330: }
331: ports = ctx->ports;
332: PetscDrawLGGetAxis(lg,&axis);
333: PetscDrawLGReset(lg);
335: xx[0] = 0.0; xx[1] = 1.0; cnt = 2;
336: PetscOptionsGetRealArray(NULL,NULL,"-zoom",xx,&cnt,NULL);
337: xs = xx[0]/hx; xm = (xx[1] - xx[0])/hx;
339: /*
340: Plot the energies
341: */
342: PetscDrawLGSetDimension(lg,1 + (ctx->allencahn ? 1 : 0));
343: PetscDrawLGSetColors(lg,colors+1);
344: PetscDrawViewPortsSet(ports,2);
345: x = hx*xs;
346: for (i=xs; i<xs+xm; i++) {
347: xx[0] = xx[1] = x;
348: yy[0] = PetscRealPart(.25*ctx->kappa*(u[i-1] - u[i+1])*(u[i-1] - u[i+1])*sx);
349: if (ctx->allencahn) yy[1] = .25*PetscRealPart((1. - u[i]*u[i])*(1. - u[i]*u[i]));
350: PetscDrawLGAddPoint(lg,xx,yy);
351: x += hx;
352: }
353: PetscDrawGetPause(draw,&pause);
354: PetscDrawSetPause(draw,0.0);
355: PetscDrawAxisSetLabels(axis,"Energy","","");
356: PetscDrawLGSetLegend(lg,legend);
357: PetscDrawLGDraw(lg);
359: /*
360: Plot the forces
361: */
362: PetscDrawViewPortsSet(ports,1);
363: PetscDrawLGReset(lg);
364: x = xs*hx;;
365: max = 0.;
366: for (i=xs; i<xs+xm; i++) {
367: xx[0] = xx[1] = x;
368: yy[0] = PetscRealPart(ctx->kappa*(u[i-1] + u[i+1] - 2.0*u[i])*sx);
369: max = PetscMax(max,PetscAbs(yy[0]));
370: if (ctx->allencahn) {
371: yy[1] = PetscRealPart(u[i] - u[i]*u[i]*u[i]);
372: max = PetscMax(max,PetscAbs(yy[1]));
373: }
374: PetscDrawLGAddPoint(lg,xx,yy);
375: x += hx;
376: }
377: PetscDrawAxisSetLabels(axis,"Right hand side","","");
378: PetscDrawLGSetLegend(lg,NULL);
379: PetscDrawLGDraw(lg);
381: /*
382: Plot the solution
383: */
384: PetscDrawLGSetDimension(lg,1);
385: PetscDrawViewPortsSet(ports,0);
386: PetscDrawLGReset(lg);
387: x = hx*xs;
388: PetscDrawLGSetLimits(lg,x,x+(xm-1)*hx,-1.1,1.1);
389: PetscDrawLGSetColors(lg,colors);
390: for (i=xs; i<xs+xm; i++) {
391: xx[0] = x;
392: yy[0] = PetscRealPart(u[i]);
393: PetscDrawLGAddPoint(lg,xx,yy);
394: x += hx;
395: }
396: PetscDrawAxisSetLabels(axis,"Solution","","");
397: PetscDrawLGDraw(lg);
399: /*
400: Print the forces as arrows on the solution
401: */
402: x = hx*xs;
403: cnt = xm/60;
404: cnt = (!cnt) ? 1 : cnt;
406: for (i=xs; i<xs+xm; i += cnt) {
407: y = PetscRealPart(u[i]);
408: len = .5*PetscRealPart(ctx->kappa*(u[i-1] + u[i+1] - 2.0*u[i])*sx)/max;
409: PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_RED);
410: if (ctx->allencahn) {
411: len = .5*PetscRealPart(u[i] - u[i]*u[i]*u[i])/max;
412: PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_BLUE);
413: }
414: x += cnt*hx;
415: }
416: DMDAVecRestoreArrayRead(da,localU,&x);
417: DMRestoreLocalVector(da,&localU);
418: PetscDrawStringSetSize(draw,.2,.2);
419: PetscDrawFlush(draw);
420: PetscDrawSetPause(draw,pause);
421: PetscDrawPause(draw);
422: return(0);
423: }
427: PetscErrorCode MyDestroy(void **ptr)
428: {
429: UserCtx *ctx = *(UserCtx**)ptr;
433: PetscDrawViewPortsDestroy(ctx->ports);
434: return(0);
435: }