Actual source code: biharmonic2.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: ./biharmonic2 -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: (this fails after a few timesteps 12/17/2017)
22: ---------------
23: ./biharmonic2 -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;
43: DM da;
44: MatFDColoring matfdcoloring;
45: ISColoring iscoloring;
46: PetscReal dt;
47: PetscReal vbounds[] = {-100000,100000,-1.1,1.1};
48: SNES snes;
49: UserCtx ctx;
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Initialize program
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54: PetscInitialize(&argc,&argv,(char*)0,help);
55: ctx.kappa = 1.0;
56: PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL);
57: 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: PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL);
65: ctx.tol = 1.0e-8;
66: PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL);
67: ctx.theta = .001;
68: ctx.theta_c = 1.0;
69: PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL);
70: PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL);
72: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73: Create distributed array (DMDA) to manage parallel grid and vectors
74: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75: DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,2,2,NULL,&da);
76: DMSetFromOptions(da);
77: DMSetUp(da);
78: DMDASetFieldName(da,0,"Biharmonic heat equation: w = -kappa*u_xx");
79: DMDASetFieldName(da,1,"Biharmonic heat equation: u");
80: DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);
81: dt = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx);
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Extract global vectors from DMDA; then duplicate for remaining
85: vectors that are the same types
86: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: DMCreateGlobalVector(da,&x);
88: VecDuplicate(x,&r);
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Create timestepping solver context
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93: TSCreate(PETSC_COMM_WORLD,&ts);
94: TSSetDM(ts,da);
95: TSSetProblemType(ts,TS_NONLINEAR);
96: TSSetIFunction(ts,NULL,FormFunction,&ctx);
97: TSSetMaxTime(ts,.02);
98: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_INTERPOLATE);
100: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101: Create matrix data structure; set Jacobian evaluation routine
103: < Set Jacobian matrix data structure and default Jacobian evaluation
104: routine. User can override with:
105: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
106: (unless user explicitly sets preconditioner)
107: -snes_mf_operator : form preconditioning matrix as set by the user,
108: but use matrix-free approx for Jacobian-vector
109: products within Newton-Krylov method
111: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
112: TSGetSNES(ts,&snes);
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 0;
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;
177: PetscInt i,Mx,xs,xm;
178: PetscReal hx,sx;
179: Field *x,*xdot,*f;
180: Vec localX,localXdot;
181: UserCtx *ctx = (UserCtx*)ptr;
183: TSGetDM(ts,&da);
184: DMGetLocalVector(da,&localX);
185: DMGetLocalVector(da,&localXdot);
186: 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);
188: hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
190: /*
191: Scatter ghost points to local vector,using the 2-step process
192: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
193: By placing code between these two statements, computations can be
194: done while messages are in transition.
195: */
196: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
197: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
198: DMGlobalToLocalBegin(da,Xdot,INSERT_VALUES,localXdot);
199: DMGlobalToLocalEnd(da,Xdot,INSERT_VALUES,localXdot);
201: /*
202: Get pointers to vector data
203: */
204: DMDAVecGetArrayRead(da,localX,&x);
205: DMDAVecGetArrayRead(da,localXdot,&xdot);
206: DMDAVecGetArray(da,F,&f);
208: /*
209: Get local grid boundaries
210: */
211: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
213: /*
214: Compute function over the locally owned part of the grid
215: */
216: for (i=xs; i<xs+xm; i++) {
217: f[i].w = x[i].w + ctx->kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx;
218: if (ctx->cahnhillard) {
219: switch (ctx->energy) {
220: case 1: /* double well */
221: f[i].w += -x[i].u*x[i].u*x[i].u + x[i].u;
222: break;
223: case 2: /* double obstacle */
224: f[i].w += x[i].u;
225: break;
226: case 3: /* logarithmic */
227: if (PetscRealPart(x[i].u) < -1.0 + 2.0*ctx->tol) f[i].w += .5*ctx->theta*(-PetscLogReal(ctx->tol) + PetscLogScalar((1.0-x[i].u)/2.0)) + ctx->theta_c*x[i].u;
228: else if (PetscRealPart(x[i].u) > 1.0 - 2.0*ctx->tol) f[i].w += .5*ctx->theta*(-PetscLogScalar((1.0+x[i].u)/2.0) + PetscLogReal(ctx->tol)) + ctx->theta_c*x[i].u;
229: 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;
230: break;
231: }
232: }
233: f[i].u = xdot[i].u - (x[i-1].w + x[i+1].w - 2.0*x[i].w)*sx;
234: }
236: /*
237: Restore vectors
238: */
239: DMDAVecRestoreArrayRead(da,localXdot,&xdot);
240: DMDAVecRestoreArrayRead(da,localX,&x);
241: DMDAVecRestoreArray(da,F,&f);
242: DMRestoreLocalVector(da,&localX);
243: DMRestoreLocalVector(da,&localXdot);
244: return 0;
245: }
247: /* ------------------------------------------------------------------- */
248: PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa)
249: {
250: PetscInt i,xs,xm,Mx,xgs,xgm;
251: Field *x;
252: PetscReal hx,xx,r,sx;
253: Vec Xg;
255: 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);
257: hx = 1.0/(PetscReal)Mx;
258: sx = 1.0/(hx*hx);
260: /*
261: Get pointers to vector data
262: */
263: DMCreateLocalVector(da,&Xg);
264: DMDAVecGetArray(da,Xg,&x);
266: /*
267: Get local grid boundaries
268: */
269: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
270: DMDAGetGhostCorners(da,&xgs,NULL,NULL,&xgm,NULL,NULL);
272: /*
273: Compute u function over the locally owned part of the grid including ghost points
274: */
275: for (i=xgs; i<xgs+xgm; i++) {
276: xx = i*hx;
277: r = PetscSqrtReal((xx-.5)*(xx-.5));
278: if (r < .125) x[i].u = 1.0;
279: else x[i].u = -.50;
280: /* fill in x[i].w so that valgrind doesn't detect use of uninitialized memory */
281: x[i].w = 0;
282: }
283: 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;
285: /*
286: Restore vectors
287: */
288: DMDAVecRestoreArray(da,Xg,&x);
290: /* Grab only the global part of the vector */
291: VecSet(X,0);
292: DMLocalToGlobalBegin(da,Xg,ADD_VALUES,X);
293: DMLocalToGlobalEnd(da,Xg,ADD_VALUES,X);
294: VecDestroy(&Xg);
295: return 0;
296: }
298: /*TEST
300: build:
301: requires: !complex !single
303: test:
304: 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
305: requires: x
307: TEST*/