Actual source code: biharmonic3.c
petsc-3.12.2 2019-11-22
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
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;
45: DM da;
46: MatFDColoring matfdcoloring;
47: ISColoring iscoloring;
48: PetscReal dt;
49: PetscReal vbounds[] = {-100000,100000,-1.1,1.1};
50: SNES snes;
51: UserCtx ctx;
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Initialize program
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
57: ctx.kappa = 1.0;
58: PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL);
59: 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: 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_STEPOVER);
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: SNESSetType(snes,SNESVINEWTONRSLS);
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: PetscScalar r,l;
182: Field *x,*xdot,*f;
183: Vec localX,localXdot;
184: UserCtx *ctx = (UserCtx*)ptr;
187: TSGetDM(ts,&da);
188: DMGetLocalVector(da,&localX);
189: DMGetLocalVector(da,&localXdot);
190: 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);
192: hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
194: /*
195: Scatter ghost points to local vector,using the 2-step process
196: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
197: By placing code between these two statements, computations can be
198: done while messages are in transition.
199: */
200: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
201: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
202: DMGlobalToLocalBegin(da,Xdot,INSERT_VALUES,localXdot);
203: DMGlobalToLocalEnd(da,Xdot,INSERT_VALUES,localXdot);
205: /*
206: Get pointers to vector data
207: */
208: DMDAVecGetArrayRead(da,localX,&x);
209: DMDAVecGetArrayRead(da,localXdot,&xdot);
210: DMDAVecGetArray(da,F,&f);
212: /*
213: Get local grid boundaries
214: */
215: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
217: /*
218: Compute function over the locally owned part of the grid
219: */
220: for (i=xs; i<xs+xm; i++) {
221: f[i].w = x[i].w + ctx->kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx;
222: if (ctx->cahnhillard) {
223: switch (ctx->energy) {
224: case 1: /* double well */
225: f[i].w += -x[i].u*x[i].u*x[i].u + x[i].u;
226: break;
227: case 2: /* double obstacle */
228: f[i].w += x[i].u;
229: break;
230: case 3: /* logarithmic */
231: 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;
232: 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;
233: 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;
234: break;
235: case 4:
236: break;
237: }
238: }
239: f[i].u = xdot[i].u - (x[i-1].w + x[i+1].w - 2.0*x[i].w)*sx;
240: if (ctx->energy==4) {
241: f[i].u = xdot[i].u;
242: /* approximation of \grad (M(u) \grad w), where M(u) = (1-u^2) */
243: r = (1.0 - x[i+1].u*x[i+1].u)*(x[i+2].w-x[i].w)*.5/hx;
244: l = (1.0 - x[i-1].u*x[i-1].u)*(x[i].w-x[i-2].w)*.5/hx;
245: f[i].u -= (r - l)*.5/hx;
246: 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;
247: }
248: }
250: /*
251: Restore vectors
252: */
253: DMDAVecRestoreArrayRead(da,localXdot,&xdot);
254: DMDAVecRestoreArrayRead(da,localX,&x);
255: DMDAVecRestoreArray(da,F,&f);
256: DMRestoreLocalVector(da,&localX);
257: DMRestoreLocalVector(da,&localXdot);
258: return(0);
259: }
261: /* ------------------------------------------------------------------- */
262: PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa)
263: {
265: PetscInt i,xs,xm,Mx;
266: Field *x;
267: PetscReal hx,xx,r,sx;
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: DMDAVecGetArray(da,X,&x);
280: /*
281: Get local grid boundaries
282: */
283: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
285: /*
286: Compute function over the locally owned part of the grid
287: */
288: for (i=xs; i<xs+xm; i++) {
289: xx = i*hx;
290: r = PetscSqrtScalar((xx-.5)*(xx-.5));
291: if (r < .125) x[i].u = 1.0;
292: else x[i].u = -.50;
293: /* u[i] = PetscPowScalar(x - .5,4.0); */
294: }
295: 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;
297: /*
298: Restore vectors
299: */
300: DMDAVecRestoreArray(da,X,&x);
301: return(0);
302: }
304: /*TEST
306: build:
307: requires: !complex !single
309: test:
310: 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
311: requires: x
313: TEST*/