Actual source code: heat.c
petsc-3.8.3 2017-12-09
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;
39: int main(int argc,char **argv)
40: {
41: TS ts; /* time integrator */
42: Vec x,r; /* solution, residual vectors */
43: Mat J; /* Jacobian matrix */
44: PetscInt steps,Mx;
46: DM da;
47: PetscReal dt;
48: PetscReal vbounds[] = {-1.1,1.1};
49: PetscBool wait;
50: UserCtx ctx;
51: PetscBool mymonitor;
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Initialize program
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
57: ctx.kappa = 1.0;
58: PetscOptionsGetReal(NULL,"-kappa",&ctx.kappa,NULL);
59: ctx.allencahn = PETSC_FALSE;
60: PetscOptionsHasName(NULL,"-allen-cahn",&ctx.allencahn);
61: PetscOptionsHasName(NULL,"-mymonitor",&mymonitor);
63: PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,vbounds);
64: PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1200,800);
66: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67: Create distributed array (DMDA) to manage parallel grid and vectors
68: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69: DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, -10,1,2,NULL,&da);
70: DMSetFromOptions(da);
71: DMSetUp(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: TSSetTimeStep(ts,dt);
101: TSSetMaxTime(ts,.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: TSGetStepNumber(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: /* ------------------------------------------------------------------- */
142: /*
143: FormFunction - Evaluates nonlinear function, F(x).
145: Input Parameters:
146: . ts - the TS context
147: . X - input vector
148: . ptr - optional user-defined context, as set by SNESSetFunction()
150: Output Parameter:
151: . F - function vector
152: */
153: PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr)
154: {
155: DM da;
157: PetscInt i,Mx,xs,xm;
158: PetscReal hx,sx;
159: PetscScalar *x,*f;
160: Vec localX;
161: UserCtx *ctx = (UserCtx*)ptr;
164: TSGetDM(ts,&da);
165: DMGetLocalVector(da,&localX);
166: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
167: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
169: hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
171: /*
172: Scatter ghost points to local vector,using the 2-step process
173: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
174: By placing code between these two statements, computations can be
175: done while messages are in transition.
176: */
177: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
178: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
180: /*
181: Get pointers to vector data
182: */
183: DMDAVecGetArrayRead(da,localX,&x);
184: DMDAVecGetArray(da,F,&f);
186: /*
187: Get local grid boundaries
188: */
189: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
191: /*
192: Compute function over the locally owned part of the grid
193: */
194: for (i=xs; i<xs+xm; i++) {
195: f[i] = ctx->kappa*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
196: if (ctx->allencahn) f[i] += (x[i] - x[i]*x[i]*x[i]);
197: }
199: /*
200: Restore vectors
201: */
202: DMDAVecRestoreArrayRead(da,localX,&x);
203: DMDAVecRestoreArray(da,F,&f);
204: DMRestoreLocalVector(da,&localX);
205: return(0);
206: }
208: /* ------------------------------------------------------------------- */
209: PetscErrorCode FormInitialSolution(DM da,Vec U)
210: {
211: PetscErrorCode ierr;
212: PetscInt i,xs,xm,Mx,scale,N;
213: PetscScalar *u;
214: const PetscScalar *f;
215: PetscReal hx,x,r;
216: Vec finesolution;
217: PetscViewer viewer;
218: PetscBool flg;
221: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
222: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
224: hx = 1.0/(PetscReal)Mx;
226: /*
227: Get pointers to vector data
228: */
229: DMDAVecGetArray(da,U,&u);
231: /*
232: Get local grid boundaries
233: */
234: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
236: /* InitialSolution is obtained with
237: ./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
238: */
240: /* PetscViewerBinaryOpen(PETSC_COMM_WORLD,"InitialSolution",FILE_MODE_READ,&viewer); */
241: /* VecCreate(PETSC_COMM_WORLD,&finesolution); */
242: /* VecLoad(finesolution,viewer); */
243: /* PetscViewerDestroy(&viewer); */
244: /* VecGetSize(finesolution,&N); */
245: /* scale = N/Mx; */
246: /* VecGetArrayRead(finesolution,&f); */
248: PetscOptionsHasName(NULL,"-square_initial",&flg);
249: if (!flg) {
250: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"InitialSolution.heat",FILE_MODE_READ,&viewer);
251: VecCreate(PETSC_COMM_WORLD,&finesolution);
252: VecLoad(finesolution,viewer);
253: PetscViewerDestroy(&viewer);
254: VecGetSize(finesolution,&N);
255: scale = N/Mx;
256: VecGetArrayRead(finesolution,&f);
257: }
259: /*
260: Compute function over the locally owned part of the grid
261: */
262: for (i=xs; i<xs+xm; i++) {
263: x = i*hx;
264: r = PetscSqrtScalar((x-.5)*(x-.5));
265: if (r < .125) u[i] = 1.0;
266: else u[i] = -.5;
268: /* With the initial condition above the method is first order in space */
269: /* this is a smooth initial condition so the method becomes second order in space */
270: /*u[i] = PetscSinScalar(2*PETSC_PI*x); */
271: /* u[i] = f[scale*i];*/
272: if (!flg) u[i] = f[scale*i];
273: }
274: if (!flg) {
275: VecRestoreArrayRead(finesolution,&f);
276: VecDestroy(&finesolution);
277: }
278: /* VecRestoreArrayRead(finesolution,&f);
279: VecDestroy(&finesolution);*/
281: /*
282: Restore vectors
283: */
284: DMDAVecRestoreArray(da,U,&u);
285: return(0);
286: }
288: /*
289: This routine is not parallel
290: */
291: PetscErrorCode MyMonitor(TS ts,PetscInt step,PetscReal time,Vec U,void *ptr)
292: {
293: UserCtx *ctx = (UserCtx*)ptr;
294: PetscDrawLG lg;
295: PetscErrorCode ierr;
296: PetscScalar *u;
297: PetscInt Mx,i,xs,xm,cnt;
298: PetscReal x,y,hx,pause,sx,len,max,xx[2],yy[2];
299: PetscDraw draw;
300: Vec localU;
301: DM da;
302: int colors[] = {PETSC_DRAW_YELLOW,PETSC_DRAW_RED,PETSC_DRAW_BLUE};
303: const char*const legend[] = {"-kappa (\\grad u,\\grad u)","(1 - u^2)^2"};
304: PetscDrawAxis axis;
305: PetscDrawViewPorts *ports;
309: TSGetDM(ts,&da);
310: DMGetLocalVector(da,&localU);
311: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
312: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
313: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
314: hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
315: DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);
316: DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);
317: DMDAVecGetArrayRead(da,localU,&u);
319: PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,&lg);
320: PetscDrawLGGetDraw(lg,&draw);
321: PetscDrawCheckResizedWindow(draw);
322: if (!ctx->ports) {
323: PetscDrawViewPortsCreateRect(draw,1,3,&ctx->ports);
324: }
325: ports = ctx->ports;
326: PetscDrawLGGetAxis(lg,&axis);
327: PetscDrawLGReset(lg);
329: xx[0] = 0.0; xx[1] = 1.0; cnt = 2;
330: PetscOptionsGetRealArray(NULL,NULL,"-zoom",xx,&cnt,NULL);
331: xs = xx[0]/hx; xm = (xx[1] - xx[0])/hx;
333: /*
334: Plot the energies
335: */
336: PetscDrawLGSetDimension(lg,1 + (ctx->allencahn ? 1 : 0));
337: PetscDrawLGSetColors(lg,colors+1);
338: PetscDrawViewPortsSet(ports,2);
339: x = hx*xs;
340: for (i=xs; i<xs+xm; i++) {
341: xx[0] = xx[1] = x;
342: yy[0] = PetscRealPart(.25*ctx->kappa*(u[i-1] - u[i+1])*(u[i-1] - u[i+1])*sx);
343: if (ctx->allencahn) yy[1] = .25*PetscRealPart((1. - u[i]*u[i])*(1. - u[i]*u[i]));
344: PetscDrawLGAddPoint(lg,xx,yy);
345: x += hx;
346: }
347: PetscDrawGetPause(draw,&pause);
348: PetscDrawSetPause(draw,0.0);
349: PetscDrawAxisSetLabels(axis,"Energy","","");
350: PetscDrawLGSetLegend(lg,legend);
351: PetscDrawLGDraw(lg);
353: /*
354: Plot the forces
355: */
356: PetscDrawViewPortsSet(ports,1);
357: PetscDrawLGReset(lg);
358: x = xs*hx;;
359: max = 0.;
360: for (i=xs; i<xs+xm; i++) {
361: xx[0] = xx[1] = x;
362: yy[0] = PetscRealPart(ctx->kappa*(u[i-1] + u[i+1] - 2.0*u[i])*sx);
363: max = PetscMax(max,PetscAbs(yy[0]));
364: if (ctx->allencahn) {
365: yy[1] = PetscRealPart(u[i] - u[i]*u[i]*u[i]);
366: max = PetscMax(max,PetscAbs(yy[1]));
367: }
368: PetscDrawLGAddPoint(lg,xx,yy);
369: x += hx;
370: }
371: PetscDrawAxisSetLabels(axis,"Right hand side","","");
372: PetscDrawLGSetLegend(lg,NULL);
373: PetscDrawLGDraw(lg);
375: /*
376: Plot the solution
377: */
378: PetscDrawLGSetDimension(lg,1);
379: PetscDrawViewPortsSet(ports,0);
380: PetscDrawLGReset(lg);
381: x = hx*xs;
382: PetscDrawLGSetLimits(lg,x,x+(xm-1)*hx,-1.1,1.1);
383: PetscDrawLGSetColors(lg,colors);
384: for (i=xs; i<xs+xm; i++) {
385: xx[0] = x;
386: yy[0] = PetscRealPart(u[i]);
387: PetscDrawLGAddPoint(lg,xx,yy);
388: x += hx;
389: }
390: PetscDrawAxisSetLabels(axis,"Solution","","");
391: PetscDrawLGDraw(lg);
393: /*
394: Print the forces as arrows on the solution
395: */
396: x = hx*xs;
397: cnt = xm/60;
398: cnt = (!cnt) ? 1 : cnt;
400: for (i=xs; i<xs+xm; i += cnt) {
401: y = PetscRealPart(u[i]);
402: len = .5*PetscRealPart(ctx->kappa*(u[i-1] + u[i+1] - 2.0*u[i])*sx)/max;
403: PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_RED);
404: if (ctx->allencahn) {
405: len = .5*PetscRealPart(u[i] - u[i]*u[i]*u[i])/max;
406: PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_BLUE);
407: }
408: x += cnt*hx;
409: }
410: DMDAVecRestoreArrayRead(da,localU,&x);
411: DMRestoreLocalVector(da,&localU);
412: PetscDrawStringSetSize(draw,.2,.2);
413: PetscDrawFlush(draw);
414: PetscDrawSetPause(draw,pause);
415: PetscDrawPause(draw);
416: return(0);
417: }
419: PetscErrorCode MyDestroy(void **ptr)
420: {
421: UserCtx *ctx = *(UserCtx**)ptr;
425: PetscDrawViewPortsDestroy(ctx->ports);
426: return(0);
427: }