Actual source code: biharmonic3.c
petsc-3.8.4 2018-03-24
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 <petscdm.h>
28: #include <petscdmda.h>
29: #include <petscts.h>
30: #include <petscdraw.h>
32: /*
33: User-defined routines
34: */
35: extern PetscErrorCode FormFunction(TS,PetscReal,Vec,Vec,Vec,void*),FormInitialSolution(DM,Vec,PetscReal);
36: typedef struct {PetscBool cahnhillard;PetscReal kappa;PetscInt energy;PetscReal tol;PetscReal theta;PetscReal theta_c;} UserCtx;
38: int main(int argc,char **argv)
39: {
40: TS ts; /* nonlinear solver */
41: Vec x,r; /* solution, residual vectors */
42: Mat J; /* Jacobian matrix */
43: PetscInt steps,Mx,maxsteps = 10000000;
45: DM da;
46: MatFDColoring matfdcoloring;
47: ISColoring iscoloring;
48: PetscReal dt;
49: PetscReal vbounds[] = {-100000,100000,-1.1,1.1};
50: PetscBool wait;
51: Vec ul,uh;
52: SNES snes;
53: UserCtx ctx;
55: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56: Initialize program
57: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
59: ctx.kappa = 1.0;
60: PetscOptionsGetReal(NULL,"-kappa",&ctx.kappa,NULL);
61: ctx.cahnhillard = PETSC_FALSE;
62: PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL);
63: PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),2,vbounds);
64: PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),600,600);
65: ctx.energy = 1;
66: /* PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL); */
67: PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL);
68: ctx.tol = 1.0e-8;
69: PetscOptionsGetReal(NULL,"-tol",&ctx.tol,NULL);
70: ctx.theta = .001;
71: ctx.theta_c = 1.0;
72: PetscOptionsGetReal(NULL,"-theta",&ctx.theta,NULL);
73: PetscOptionsGetReal(NULL,"-theta_c",&ctx.theta_c,NULL);
75: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76: Create distributed array (DMDA) to manage parallel grid and vectors
77: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
78: DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, -10,2,2,NULL,&da);
79: DMSetFromOptions(da);
80: DMSetUp(da);
81: DMDASetFieldName(da,0,"Biharmonic heat equation: w = -kappa*u_xx");
82: DMDASetFieldName(da,1,"Biharmonic heat equation: u");
83: DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);
84: dt = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx);
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Extract global vectors from DMDA; then duplicate for remaining
88: vectors that are the same types
89: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90: DMCreateGlobalVector(da,&x);
91: VecDuplicate(x,&r);
93: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94: Create timestepping solver context
95: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: TSCreate(PETSC_COMM_WORLD,&ts);
97: TSSetDM(ts,da);
98: TSSetProblemType(ts,TS_NONLINEAR);
99: TSSetIFunction(ts,NULL,FormFunction,&ctx);
100: TSSetMaxTime(ts,.02);
101: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
103: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: Create matrix data structure; set Jacobian evaluation routine
106: < Set Jacobian matrix data structure and default Jacobian evaluation
107: routine. User can override with:
108: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
109: (unless user explicitly sets preconditioner)
110: -snes_mf_operator : form preconditioning matrix as set by the user,
111: but use matrix-free approx for Jacobian-vector
112: products within Newton-Krylov method
114: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115: TSGetSNES(ts,&snes);
116: DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring);
117: DMSetMatType(da,MATAIJ);
118: DMCreateMatrix(da,&J);
119: MatFDColoringCreate(J,iscoloring,&matfdcoloring);
120: ISColoringDestroy(&iscoloring);
121: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);
122: MatFDColoringSetFromOptions(matfdcoloring);
123: MatFDColoringSetUp(J,iscoloring,matfdcoloring);
124: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);
126: {
127: VecDuplicate(x,&ul);
128: VecDuplicate(x,&uh);
129: VecStrideSet(ul,0,PETSC_NINFINITY);
130: VecStrideSet(ul,1,-1.0);
131: VecStrideSet(uh,0,PETSC_INFINITY);
132: VecStrideSet(uh,1,1.0);
133: TSVISetVariableBounds(ts,ul,uh);
134: }
136: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137: Customize nonlinear solver
138: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139: TSSetType(ts,TSBEULER);
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: Set initial conditions
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144: FormInitialSolution(da,x,ctx.kappa);
145: TSSetTimeStep(ts,dt);
146: TSSetSolution(ts,x);
148: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149: Set runtime options
150: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151: TSSetFromOptions(ts);
153: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154: Solve nonlinear system
155: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
156: TSSolve(ts,x);
157: wait = PETSC_FALSE;
158: PetscOptionsGetBool(NULL,NULL,"-wait",&wait,NULL);
159: if (wait) {
160: PetscSleep(-1);
161: }
162: TSGetStepNumber(ts,&steps);
164: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165: Free work space. All PETSc objects should be destroyed when they
166: are no longer needed.
167: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
168: {
169: VecDestroy(&ul);
170: VecDestroy(&uh);
171: }
172: MatDestroy(&J);
173: MatFDColoringDestroy(&matfdcoloring);
174: VecDestroy(&x);
175: VecDestroy(&r);
176: TSDestroy(&ts);
177: DMDestroy(&da);
179: PetscFinalize();
180: return(0);
181: }
183: typedef struct {PetscScalar w,u;} Field;
184: /* ------------------------------------------------------------------- */
185: /*
186: FormFunction - Evaluates nonlinear function, F(x).
188: Input Parameters:
189: . ts - the TS context
190: . X - input vector
191: . ptr - optional user-defined context, as set by SNESSetFunction()
193: Output Parameter:
194: . F - function vector
195: */
196: PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec Xdot,Vec F,void *ptr)
197: {
198: DM da;
200: PetscInt i,Mx,xs,xm;
201: PetscReal hx,sx;
202: PetscScalar r,l;
203: Field *x,*xdot,*f;
204: Vec localX,localXdot;
205: UserCtx *ctx = (UserCtx*)ptr;
208: TSGetDM(ts,&da);
209: DMGetLocalVector(da,&localX);
210: DMGetLocalVector(da,&localXdot);
211: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
212: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
214: hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
216: /*
217: Scatter ghost points to local vector,using the 2-step process
218: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
219: By placing code between these two statements, computations can be
220: done while messages are in transition.
221: */
222: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
223: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
224: DMGlobalToLocalBegin(da,Xdot,INSERT_VALUES,localXdot);
225: DMGlobalToLocalEnd(da,Xdot,INSERT_VALUES,localXdot);
227: /*
228: Get pointers to vector data
229: */
230: DMDAVecGetArrayRead(da,localX,&x);
231: DMDAVecGetArrayRead(da,localXdot,&xdot);
232: DMDAVecGetArray(da,F,&f);
234: /*
235: Get local grid boundaries
236: */
237: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
239: /*
240: Compute function over the locally owned part of the grid
241: */
242: for (i=xs; i<xs+xm; i++) {
243: f[i].w = x[i].w + ctx->kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx;
244: if (ctx->cahnhillard) {
245: switch (ctx->energy) {
246: case 1: /* double well */
247: f[i].w += -x[i].u*x[i].u*x[i].u + x[i].u;
248: break;
249: case 2: /* double obstacle */
250: f[i].w += x[i].u;
251: break;
252: case 3: /* logarithmic */
253: 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;
254: 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;
255: 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;
256: break;
257: case 4:
258: break;
259: }
260: }
261: f[i].u = xdot[i].u - (x[i-1].w + x[i+1].w - 2.0*x[i].w)*sx;
262: if (ctx->energy==4) {
263: f[i].u = xdot[i].u;
264: /* approximation of \grad (M(u) \grad w), where M(u) = (1-u^2) */
265: r = (1.0 - x[i+1].u*x[i+1].u)*(x[i+2].w-x[i].w)*.5/hx;
266: l = (1.0 - x[i-1].u*x[i-1].u)*(x[i].w-x[i-2].w)*.5/hx;
267: f[i].u -= (r - l)*.5/hx;
268: 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;
269: }
270: }
272: /*
273: Restore vectors
274: */
275: DMDAVecRestoreArrayRead(da,localXdot,&xdot);
276: DMDAVecRestoreArrayRead(da,localX,&x);
277: DMDAVecRestoreArray(da,F,&f);
278: DMRestoreLocalVector(da,&localX);
279: DMRestoreLocalVector(da,&localXdot);
280: return(0);
281: }
283: /* ------------------------------------------------------------------- */
284: PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa)
285: {
287: PetscInt i,xs,xm,Mx;
288: Field *x;
289: PetscReal hx,xx,r,sx;
292: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
293: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
295: hx = 1.0/(PetscReal)Mx;
296: sx = 1.0/(hx*hx);
298: /*
299: Get pointers to vector data
300: */
301: DMDAVecGetArray(da,X,&x);
303: /*
304: Get local grid boundaries
305: */
306: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
308: /*
309: Compute function over the locally owned part of the grid
310: */
311: for (i=xs; i<xs+xm; i++) {
312: xx = i*hx;
313: r = PetscSqrtScalar((xx-.5)*(xx-.5));
314: if (r < .125) x[i].u = 1.0;
315: else x[i].u = -.50;
316: /* u[i] = PetscPowScalar(x - .5,4.0); */
317: }
318: 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;
320: /*
321: Restore vectors
322: */
323: DMDAVecRestoreArray(da,X,&x);
324: return(0);
325: }