Actual source code: biharmonic3.c
petsc-dev 2014-02-02
2: static char help[] = "Solves biharmonic equation in 1d.\n";
4: /*
5: Solves the equation biharmonic equation in split form
7: w = -kappa \Delta u
8: u_t = \Delta w
9: -1 <= u <= 1
10: Periodic boundary conditions
12: Evolve the biharmonic heat equation with bounds: (same as biharmonic)
13: ---------------
14: ./biharmonic2 -ts_monitor -snes_monitor -ts_monitor_draw_solution -pc_type lu -draw_pause .1 -snes_converged_reason --wait -ts_type beuler -da_refine 5 -draw_fields 1 -ts_dt 9.53674e-9
16: w = -kappa \Delta u + u^3 - u
17: u_t = \Delta w
18: -1 <= u <= 1
19: Periodic boundary conditions
21: Evolve the Cahn-Hillard equations:
22: ---------------
23: ./biharmonic2 -ts_monitor -snes_monitor -ts_monitor_draw_solution -pc_type lu -draw_pause .1 -snes_converged_reason --wait -ts_type beuler -da_refine 6 -vi -draw_fields 1 -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard
26: */
27: #include <petscdmda.h>
28: #include <petscts.h>
29: #include <petscdraw.h>
31: /*
32: User-defined routines
33: */
34: extern PetscErrorCode FormFunction(TS,PetscReal,Vec,Vec,Vec,void*),FormInitialSolution(DM,Vec,PetscReal);
35: typedef struct {PetscBool cahnhillard;PetscReal kappa;PetscInt energy;PetscReal tol;PetscReal theta;PetscReal theta_c;} UserCtx;
39: int main(int argc,char **argv)
40: {
41: TS ts; /* nonlinear solver */
42: Vec x,r; /* solution, residual vectors */
43: Mat J; /* Jacobian matrix */
44: PetscInt steps,Mx,maxsteps = 10000000;
46: DM da;
47: MatFDColoring matfdcoloring;
48: ISColoring iscoloring;
49: PetscReal dt;
50: PetscReal vbounds[] = {-100000,100000,-1.1,1.1};
51: PetscBool wait;
52: Vec ul,uh;
53: SNES snes;
54: UserCtx ctx;
56: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: Initialize program
58: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59: PetscInitialize(&argc,&argv,(char*)0,help);
60: ctx.kappa = 1.0;
61: PetscOptionsGetReal(NULL,"-kappa",&ctx.kappa,NULL);
62: ctx.cahnhillard = PETSC_FALSE;
63: PetscOptionsGetBool(NULL,"-cahn-hillard",&ctx.cahnhillard,NULL);
64: PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),2,vbounds);
65: PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),600,600);
66: ctx.energy = 1;
67: /* PetscOptionsGetInt(NULL,"-energy",&ctx.energy,NULL); */
68: PetscOptionsInt("-energy","type of energy (1=double well, 2=double obstacle, 3=logarithmic, 4=degenerate mobility and weird splitting)","",ctx.energy,&ctx.energy,NULL);
69: ctx.tol = 1.0e-8;
70: PetscOptionsGetReal(NULL,"-tol",&ctx.tol,NULL);
71: ctx.theta = .001;
72: ctx.theta_c = 1.0;
73: PetscOptionsGetReal(NULL,"-theta",&ctx.theta,NULL);
74: PetscOptionsGetReal(NULL,"-theta_c",&ctx.theta_c,NULL);
76: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77: Create distributed array (DMDA) to manage parallel grid and vectors
78: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79: DMDACreate1d(PETSC_COMM_WORLD, DMDA_BOUNDARY_PERIODIC, -10,2,2,NULL,&da);
80: DMDASetFieldName(da,0,"Biharmonic heat equation: w = -kappa*u_xx");
81: DMDASetFieldName(da,1,"Biharmonic heat equation: u");
82: DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);
83: dt = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx);
85: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86: Extract global vectors from DMDA; then duplicate for remaining
87: vectors that are the same types
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
89: DMCreateGlobalVector(da,&x);
90: VecDuplicate(x,&r);
92: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93: Create timestepping solver context
94: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95: TSCreate(PETSC_COMM_WORLD,&ts);
96: TSSetDM(ts,da);
97: TSSetProblemType(ts,TS_NONLINEAR);
98: TSSetIFunction(ts,NULL,FormFunction,&ctx);
99: TSSetDuration(ts,maxsteps,.02);
100: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
102: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: Create matrix data structure; set Jacobian evaluation routine
105: < Set Jacobian matrix data structure and default Jacobian evaluation
106: routine. User can override with:
107: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
108: (unless user explicitly sets preconditioner)
109: -snes_mf_operator : form preconditioning matrix as set by the user,
110: but use matrix-free approx for Jacobian-vector
111: products within Newton-Krylov method
113: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114: TSGetSNES(ts,&snes);
115: DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring);
116: DMSetMatType(da,MATAIJ);
117: DMCreateMatrix(da,&J);
118: MatFDColoringCreate(J,iscoloring,&matfdcoloring);
119: ISColoringDestroy(&iscoloring);
120: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);
121: MatFDColoringSetFromOptions(matfdcoloring);
122: MatFDColoringSetUp(J,iscoloring,matfdcoloring);
123: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);
125: {
126: VecDuplicate(x,&ul);
127: VecDuplicate(x,&uh);
128: VecStrideSet(ul,0,PETSC_NINFINITY);
129: VecStrideSet(ul,1,-1.0);
130: VecStrideSet(uh,0,PETSC_INFINITY);
131: VecStrideSet(uh,1,1.0);
132: TSVISetVariableBounds(ts,ul,uh);
133: }
135: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136: Customize nonlinear solver
137: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138: TSSetType(ts,TSBEULER);
140: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141: Set initial conditions
142: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143: FormInitialSolution(da,x,ctx.kappa);
144: TSSetInitialTimeStep(ts,0.0,dt);
145: TSSetSolution(ts,x);
147: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148: Set runtime options
149: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150: TSSetFromOptions(ts);
152: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153: Solve nonlinear system
154: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155: TSSolve(ts,x);
156: wait = PETSC_FALSE;
157: PetscOptionsGetBool(NULL,"-wait",&wait,NULL);
158: if (wait) {
159: PetscSleep(-1);
160: }
161: TSGetTimeStepNumber(ts,&steps);
163: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164: Free work space. All PETSc objects should be destroyed when they
165: are no longer needed.
166: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
167: {
168: VecDestroy(&ul);
169: VecDestroy(&uh);
170: }
171: MatDestroy(&J);
172: MatFDColoringDestroy(&matfdcoloring);
173: VecDestroy(&x);
174: VecDestroy(&r);
175: TSDestroy(&ts);
176: DMDestroy(&da);
178: PetscFinalize();
179: return(0);
180: }
182: typedef struct {PetscScalar w,u;} Field;
183: /* ------------------------------------------------------------------- */
186: /*
187: FormFunction - Evaluates nonlinear function, F(x).
189: Input Parameters:
190: . ts - the TS context
191: . X - input vector
192: . ptr - optional user-defined context, as set by SNESSetFunction()
194: Output Parameter:
195: . F - function vector
196: */
197: PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec Xdot,Vec F,void *ptr)
198: {
199: DM da;
201: PetscInt i,Mx,xs,xm;
202: PetscReal hx,sx;
203: PetscScalar r,l;
204: Field *x,*xdot,*f;
205: Vec localX,localXdot;
206: UserCtx *ctx = (UserCtx*)ptr;
209: TSGetDM(ts,&da);
210: DMGetLocalVector(da,&localX);
211: DMGetLocalVector(da,&localXdot);
212: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
213: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
215: hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
217: /*
218: Scatter ghost points to local vector,using the 2-step process
219: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
220: By placing code between these two statements, computations can be
221: done while messages are in transition.
222: */
223: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
224: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
225: DMGlobalToLocalBegin(da,Xdot,INSERT_VALUES,localXdot);
226: DMGlobalToLocalEnd(da,Xdot,INSERT_VALUES,localXdot);
228: /*
229: Get pointers to vector data
230: */
231: DMDAVecGetArray(da,localX,&x);
232: DMDAVecGetArray(da,localXdot,&xdot);
233: DMDAVecGetArray(da,F,&f);
235: /*
236: Get local grid boundaries
237: */
238: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
240: /*
241: Compute function over the locally owned part of the grid
242: */
243: for (i=xs; i<xs+xm; i++) {
244: f[i].w = x[i].w + ctx->kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx;
245: if (ctx->cahnhillard) {
246: switch (ctx->energy) {
247: case 1: /* double well */
248: f[i].w += -x[i].u*x[i].u*x[i].u + x[i].u;
249: break;
250: case 2: /* double obstacle */
251: f[i].w += x[i].u;
252: break;
253: case 3: /* logarithmic */
254: if (x[i].u < -1.0 + 2.0*ctx->tol) f[i].w += .5*ctx->theta*(-log(ctx->tol) + log((1.0-x[i].u)/2.0)) + ctx->theta_c*x[i].u;
255: else if (x[i].u > 1.0 - 2.0*ctx->tol) f[i].w += .5*ctx->theta*(-log((1.0+x[i].u)/2.0) + log(ctx->tol)) + ctx->theta_c*x[i].u;
256: else f[i].w += .5*ctx->theta*(-log((1.0+x[i].u)/2.0) + log((1.0-x[i].u)/2.0)) + ctx->theta_c*x[i].u;
257: break;
258: case 4:
259: break;
260: }
261: }
262: f[i].u = xdot[i].u - (x[i-1].w + x[i+1].w - 2.0*x[i].w)*sx;
263: if (ctx->energy==4) {
264: f[i].u = xdot[i].u;
265: /* approximation of \grad (M(u) \grad w), where M(u) = (1-u^2) */
266: r = (1.0 - x[i+1].u*x[i+1].u)*(x[i+2].w-x[i].w)*.5/hx;
267: l = (1.0 - x[i-1].u*x[i-1].u)*(x[i].w-x[i-2].w)*.5/hx;
268: f[i].u -= (r - l)*.5/hx;
269: f[i].u += 2.0*ctx->theta_c*x[i].u*(x[i+1].u-x[i-1].u)*(x[i+1].u-x[i-1].u)*.25*sx - (ctx->theta - ctx->theta_c*(1-x[i].u*x[i].u))*(x[i+1].u + x[i-1].u - 2.0*x[i].u)*sx;
270: }
271: }
273: /*
274: Restore vectors
275: */
276: DMDAVecRestoreArray(da,localXdot,&xdot);
277: DMDAVecRestoreArray(da,localX,&x);
278: DMDAVecRestoreArray(da,F,&f);
279: DMRestoreLocalVector(da,&localX);
280: DMRestoreLocalVector(da,&localXdot);
281: return(0);
282: }
284: /* ------------------------------------------------------------------- */
287: PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa)
288: {
290: PetscInt i,xs,xm,Mx;
291: Field *x;
292: PetscReal hx,xx,r,sx;
295: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
296: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
298: hx = 1.0/(PetscReal)Mx;
299: sx = 1.0/(hx*hx);
301: /*
302: Get pointers to vector data
303: */
304: DMDAVecGetArray(da,X,&x);
306: /*
307: Get local grid boundaries
308: */
309: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
311: /*
312: Compute function over the locally owned part of the grid
313: */
314: for (i=xs; i<xs+xm; i++) {
315: xx = i*hx;
316: r = PetscSqrtScalar((xx-.5)*(xx-.5));
317: if (r < .125) x[i].u = 1.0;
318: else x[i].u = -.50;
319: /* u[i] = PetscPowScalar(x - .5,4.0); */
320: }
321: for (i=xs; i<xs+xm; i++) x[i].w = -kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx;
323: /*
324: Restore vectors
325: */
326: DMDAVecRestoreArray(da,X,&x);
327: return(0);
328: }