Actual source code: biharmonic3.c
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: ./biharmonic3 -ts_monitor -snes_monitor -ts_monitor_draw_solution -pc_type lu -draw_pause .1 -snes_converged_reason -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: ./biharmonic3 -ts_monitor -snes_monitor -ts_monitor_draw_solution -pc_type lu -draw_pause .1 -snes_converged_reason -ts_type beuler -da_refine 6 -draw_fields 1 -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard
25: */
26: #include <petscdm.h>
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;
37: int main(int argc,char **argv)
38: {
39: TS ts; /* nonlinear solver */
40: Vec x,r; /* solution, residual vectors */
41: Mat J; /* Jacobian matrix */
42: PetscInt steps,Mx;
44: DM da;
45: MatFDColoring matfdcoloring;
46: ISColoring iscoloring;
47: PetscReal dt;
48: PetscReal vbounds[] = {-100000,100000,-1.1,1.1};
49: SNES snes;
50: UserCtx ctx;
52: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53: Initialize program
54: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
56: ctx.kappa = 1.0;
57: PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL);
58: ctx.cahnhillard = PETSC_FALSE;
59: PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL);
60: PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),2,vbounds);
61: PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),600,600);
62: ctx.energy = 1;
63: PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL);
64: ctx.tol = 1.0e-8;
65: PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL);
66: ctx.theta = .001;
67: ctx.theta_c = 1.0;
68: PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL);
69: PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL);
71: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72: Create distributed array (DMDA) to manage parallel grid and vectors
73: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74: DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,2,2,NULL,&da);
75: DMSetFromOptions(da);
76: DMSetUp(da);
77: DMDASetFieldName(da,0,"Biharmonic heat equation: w = -kappa*u_xx");
78: DMDASetFieldName(da,1,"Biharmonic heat equation: u");
79: DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);
80: dt = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Extract global vectors from DMDA; then duplicate for remaining
84: vectors that are the same types
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: DMCreateGlobalVector(da,&x);
87: VecDuplicate(x,&r);
89: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
90: Create timestepping solver context
91: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92: TSCreate(PETSC_COMM_WORLD,&ts);
93: TSSetDM(ts,da);
94: TSSetProblemType(ts,TS_NONLINEAR);
95: TSSetIFunction(ts,NULL,FormFunction,&ctx);
96: TSSetMaxTime(ts,.02);
97: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
99: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100: Create matrix data structure; set Jacobian evaluation routine
102: < Set Jacobian matrix data structure and default Jacobian evaluation
103: routine. User can override with:
104: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
105: (unless user explicitly sets preconditioner)
106: -snes_mf_operator : form preconditioning matrix as set by the user,
107: but use matrix-free approx for Jacobian-vector
108: products within Newton-Krylov method
110: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111: TSGetSNES(ts,&snes);
112: SNESSetType(snes,SNESVINEWTONRSLS);
113: DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring);
114: DMSetMatType(da,MATAIJ);
115: DMCreateMatrix(da,&J);
116: MatFDColoringCreate(J,iscoloring,&matfdcoloring);
117: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);
118: MatFDColoringSetFromOptions(matfdcoloring);
119: MatFDColoringSetUp(J,iscoloring,matfdcoloring);
120: ISColoringDestroy(&iscoloring);
121: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);
123: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124: Customize nonlinear solver
125: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
126: TSSetType(ts,TSBEULER);
128: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129: Set initial conditions
130: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131: FormInitialSolution(da,x,ctx.kappa);
132: TSSetTimeStep(ts,dt);
133: TSSetSolution(ts,x);
135: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136: Set runtime options
137: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138: TSSetFromOptions(ts);
140: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141: Solve nonlinear system
142: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143: TSSolve(ts,x);
144: TSGetStepNumber(ts,&steps);
146: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147: Free work space. All PETSc objects should be destroyed when they
148: are no longer needed.
149: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150: MatDestroy(&J);
151: MatFDColoringDestroy(&matfdcoloring);
152: VecDestroy(&x);
153: VecDestroy(&r);
154: TSDestroy(&ts);
155: DMDestroy(&da);
157: PetscFinalize();
158: return ierr;
159: }
161: typedef struct {PetscScalar w,u;} Field;
162: /* ------------------------------------------------------------------- */
163: /*
164: FormFunction - Evaluates nonlinear function, F(x).
166: Input Parameters:
167: . ts - the TS context
168: . X - input vector
169: . ptr - optional user-defined context, as set by SNESSetFunction()
171: Output Parameter:
172: . F - function vector
173: */
174: PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec Xdot,Vec F,void *ptr)
175: {
176: DM da;
178: PetscInt i,Mx,xs,xm;
179: PetscReal hx,sx;
180: PetscScalar r,l;
181: Field *x,*xdot,*f;
182: Vec localX,localXdot;
183: UserCtx *ctx = (UserCtx*)ptr;
186: TSGetDM(ts,&da);
187: DMGetLocalVector(da,&localX);
188: DMGetLocalVector(da,&localXdot);
189: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
191: hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
193: /*
194: Scatter ghost points to local vector,using the 2-step process
195: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
196: By placing code between these two statements, computations can be
197: done while messages are in transition.
198: */
199: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
200: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
201: DMGlobalToLocalBegin(da,Xdot,INSERT_VALUES,localXdot);
202: DMGlobalToLocalEnd(da,Xdot,INSERT_VALUES,localXdot);
204: /*
205: Get pointers to vector data
206: */
207: DMDAVecGetArrayRead(da,localX,&x);
208: DMDAVecGetArrayRead(da,localXdot,&xdot);
209: DMDAVecGetArray(da,F,&f);
211: /*
212: Get local grid boundaries
213: */
214: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
216: /*
217: Compute function over the locally owned part of the grid
218: */
219: for (i=xs; i<xs+xm; i++) {
220: f[i].w = x[i].w + ctx->kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx;
221: if (ctx->cahnhillard) {
222: switch (ctx->energy) {
223: case 1: /* double well */
224: f[i].w += -x[i].u*x[i].u*x[i].u + x[i].u;
225: break;
226: case 2: /* double obstacle */
227: f[i].w += x[i].u;
228: break;
229: case 3: /* logarithmic */
230: if (x[i].u < -1.0 + 2.0*ctx->tol) f[i].w += .5*ctx->theta*(-PetscLogScalar(ctx->tol) + PetscLogScalar((1.0-x[i].u)/2.0)) + ctx->theta_c*x[i].u;
231: else if (x[i].u > 1.0 - 2.0*ctx->tol) f[i].w += .5*ctx->theta*(-PetscLogScalar((1.0+x[i].u)/2.0) + PetscLogScalar(ctx->tol)) + ctx->theta_c*x[i].u;
232: else f[i].w += .5*ctx->theta*(-PetscLogScalar((1.0+x[i].u)/2.0) + PetscLogScalar((1.0-x[i].u)/2.0)) + ctx->theta_c*x[i].u;
233: break;
234: case 4:
235: break;
236: }
237: }
238: f[i].u = xdot[i].u - (x[i-1].w + x[i+1].w - 2.0*x[i].w)*sx;
239: if (ctx->energy==4) {
240: f[i].u = xdot[i].u;
241: /* approximation of \grad (M(u) \grad w), where M(u) = (1-u^2) */
242: r = (1.0 - x[i+1].u*x[i+1].u)*(x[i+2].w-x[i].w)*.5/hx;
243: l = (1.0 - x[i-1].u*x[i-1].u)*(x[i].w-x[i-2].w)*.5/hx;
244: f[i].u -= (r - l)*.5/hx;
245: 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;
246: }
247: }
249: /*
250: Restore vectors
251: */
252: DMDAVecRestoreArrayRead(da,localXdot,&xdot);
253: DMDAVecRestoreArrayRead(da,localX,&x);
254: DMDAVecRestoreArray(da,F,&f);
255: DMRestoreLocalVector(da,&localX);
256: DMRestoreLocalVector(da,&localXdot);
257: return(0);
258: }
260: /* ------------------------------------------------------------------- */
261: PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa)
262: {
264: PetscInt i,xs,xm,Mx,xgs,xgm;
265: Field *x;
266: PetscReal hx,xx,r,sx;
267: Vec Xg;
270: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
272: hx = 1.0/(PetscReal)Mx;
273: sx = 1.0/(hx*hx);
275: /*
276: Get pointers to vector data
277: */
278: DMCreateLocalVector(da,&Xg);
279: DMDAVecGetArray(da,Xg,&x);
281: /*
282: Get local grid boundaries
283: */
284: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
285: DMDAGetGhostCorners(da,&xgs,NULL,NULL,&xgm,NULL,NULL);
287: /*
288: Compute u function over the locally owned part of the grid including ghost points
289: */
290: for (i=xgs; i<xgs+xgm; i++) {
291: xx = i*hx;
292: r = PetscSqrtReal((xx-.5)*(xx-.5));
293: if (r < .125) x[i].u = 1.0;
294: else x[i].u = -.50;
295: /* fill in x[i].w so that valgrind doesn't detect use of uninitialized memory */
296: x[i].w = 0;
297: }
298: 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;
300: /*
301: Restore vectors
302: */
303: DMDAVecRestoreArray(da,Xg,&x);
305: /* Grab only the global part of the vector */
306: VecSet(X,0);
307: DMLocalToGlobalBegin(da,Xg,ADD_VALUES,X);
308: DMLocalToGlobalEnd(da,Xg,ADD_VALUES,X);
309: VecDestroy(&Xg);
310: return(0);
311: }
313: /*TEST
315: build:
316: requires: !complex !single
318: test:
319: args: -ts_monitor -snes_monitor -pc_type lu -snes_converged_reason -ts_type beuler -da_refine 5 -ts_dt 9.53674e-9 -ts_max_steps 50
320: requires: x
322: TEST*/