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;
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;
60: PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL);
61: PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),2,vbounds);
62: PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),600,600);
63: ctx.energy = 1;
64: /*PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL);*/
65: PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL);
66: ctx.tol = 1.0e-8;
67: PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL);
68: ctx.theta = .001;
69: ctx.theta_c = 1.0;
70: PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL);
71: PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL);
73: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74: Create distributed array (DMDA) to manage parallel grid and vectors
75: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
76: DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,2,2,NULL,&da);
77: DMSetFromOptions(da);
78: DMSetUp(da);
79: DMDASetFieldName(da,0,"Biharmonic heat equation: w = -kappa*u_xx");
80: DMDASetFieldName(da,1,"Biharmonic heat equation: u");
81: DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);
82: dt = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx);
84: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: Extract global vectors from DMDA; then duplicate for remaining
86: vectors that are the same types
87: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88: DMCreateGlobalVector(da,&x);
89: VecDuplicate(x,&r);
91: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92: Create timestepping solver context
93: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: TSCreate(PETSC_COMM_WORLD,&ts);
95: TSSetDM(ts,da);
96: TSSetProblemType(ts,TS_NONLINEAR);
97: TSSetIFunction(ts,NULL,FormFunction,&ctx);
98: TSSetMaxTime(ts,.02);
99: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_INTERPOLATE);
101: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: Create matrix data structure; set Jacobian evaluation routine
104: < Set Jacobian matrix data structure and default Jacobian evaluation
105: routine. User can override with:
106: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
107: (unless user explicitly sets preconditioner)
108: -snes_mf_operator : form preconditioning matrix as set by the user,
109: but use matrix-free approx for Jacobian-vector
110: products within Newton-Krylov method
112: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113: TSGetSNES(ts,&snes);
114: DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring);
115: DMSetMatType(da,MATAIJ);
116: DMCreateMatrix(da,&J);
117: MatFDColoringCreate(J,iscoloring,&matfdcoloring);
118: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);
119: MatFDColoringSetFromOptions(matfdcoloring);
120: MatFDColoringSetUp(J,iscoloring,matfdcoloring);
121: ISColoringDestroy(&iscoloring);
122: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);
124: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125: Customize nonlinear solver
126: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127: TSSetType(ts,TSBEULER);
129: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130: Set initial conditions
131: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132: FormInitialSolution(da,x,ctx.kappa);
133: TSSetTimeStep(ts,dt);
134: TSSetSolution(ts,x);
136: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137: Set runtime options
138: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139: TSSetFromOptions(ts);
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: Solve nonlinear system
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144: TSSolve(ts,x);
145: TSGetStepNumber(ts,&steps);
147: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148: Free work space. All PETSc objects should be destroyed when they
149: are no longer needed.
150: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151: MatDestroy(&J);
152: MatFDColoringDestroy(&matfdcoloring);
153: VecDestroy(&x);
154: VecDestroy(&r);
155: TSDestroy(&ts);
156: DMDestroy(&da);
158: PetscFinalize();
159: return ierr;
160: }
162: typedef struct {PetscScalar w,u;} Field;
163: /* ------------------------------------------------------------------- */
164: /*
165: FormFunction - Evaluates nonlinear function, F(x).
167: Input Parameters:
168: . ts - the TS context
169: . X - input vector
170: . ptr - optional user-defined context, as set by SNESSetFunction()
172: Output Parameter:
173: . F - function vector
174: */
175: PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec Xdot,Vec F,void *ptr)
176: {
177: DM da;
179: PetscInt i,Mx,xs,xm;
180: PetscReal hx,sx;
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 (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;
231: 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;
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: }
235: }
236: f[i].u = xdot[i].u - (x[i-1].w + x[i+1].w - 2.0*x[i].w)*sx;
237: }
239: /*
240: Restore vectors
241: */
242: DMDAVecRestoreArrayRead(da,localXdot,&xdot);
243: DMDAVecRestoreArrayRead(da,localX,&x);
244: DMDAVecRestoreArray(da,F,&f);
245: DMRestoreLocalVector(da,&localX);
246: DMRestoreLocalVector(da,&localXdot);
247: return(0);
248: }
250: /* ------------------------------------------------------------------- */
251: PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa)
252: {
254: PetscInt i,xs,xm,Mx,xgs,xgm;
255: Field *x;
256: PetscReal hx,xx,r,sx;
257: Vec Xg;
260: 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);
262: hx = 1.0/(PetscReal)Mx;
263: sx = 1.0/(hx*hx);
265: /*
266: Get pointers to vector data
267: */
268: DMCreateLocalVector(da,&Xg);
269: DMDAVecGetArray(da,Xg,&x);
271: /*
272: Get local grid boundaries
273: */
274: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
275: DMDAGetGhostCorners(da,&xgs,NULL,NULL,&xgm,NULL,NULL);
277: /*
278: Compute u function over the locally owned part of the grid including ghost points
279: */
280: for (i=xgs; i<xgs+xgm; i++) {
281: xx = i*hx;
282: r = PetscSqrtReal((xx-.5)*(xx-.5));
283: if (r < .125) x[i].u = 1.0;
284: else x[i].u = -.50;
285: /* fill in x[i].w so that valgrind doesn't detect use of uninitialized memory */
286: x[i].w = 0;
287: }
288: 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;
290: /*
291: Restore vectors
292: */
293: DMDAVecRestoreArray(da,Xg,&x);
295: /* Grab only the global part of the vector */
296: VecSet(X,0);
297: DMLocalToGlobalBegin(da,Xg,ADD_VALUES,X);
298: DMLocalToGlobalEnd(da,Xg,ADD_VALUES,X);
299: VecDestroy(&Xg);
300: return(0);
301: }
303: /*TEST
305: build:
306: requires: !complex !single
308: test:
309: 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
310: requires: x
312: TEST*/