Actual source code: ex19.c
2: static char help[] = "Nonlinear driven cavity with multigrid in 2d.\n\
3: \n\
4: The 2D driven cavity problem is solved in a velocity-vorticity formulation.\n\
5: The flow can be driven with the lid or with bouyancy or both:\n\
6: -lidvelocity <lid>, where <lid> = dimensionless velocity of lid\n\
7: -grashof <gr>, where <gr> = dimensionless temperature gradent\n\
8: -prandtl <pr>, where <pr> = dimensionless thermal/momentum diffusity ratio\n\
9: -contours : draw contour plots of solution\n\n";
11: /*T
12: Concepts: SNES^solving a system of nonlinear equations (parallel multicomponent example);
13: Concepts: DA^using distributed arrays;
14: Concepts: multicomponent
15: Processors: n
16: T*/
18: /* ------------------------------------------------------------------------
20: We thank David E. Keyes for contributing the driven cavity discretization
21: within this example code.
23: This problem is modeled by the partial differential equation system
24:
25: - Lap(U) - Grad_y(Omega) = 0
26: - Lap(V) + Grad_x(Omega) = 0
27: - Lap(Omega) + Div([U*Omega,V*Omega]) - GR*Grad_x(T) = 0
28: - Lap(T) + PR*Div([U*T,V*T]) = 0
30: in the unit square, which is uniformly discretized in each of x and
31: y in this simple encoding.
33: No-slip, rigid-wall Dirichlet conditions are used for [U,V].
34: Dirichlet conditions are used for Omega, based on the definition of
35: vorticity: Omega = - Grad_y(U) + Grad_x(V), where along each
36: constant coordinate boundary, the tangential derivative is zero.
37: Dirichlet conditions are used for T on the left and right walls,
38: and insulation homogeneous Neumann conditions are used for T on
39: the top and bottom walls.
41: A finite difference approximation with the usual 5-point stencil
42: is used to discretize the boundary value problem to obtain a
43: nonlinear system of equations. Upwinding is used for the divergence
44: (convective) terms and central for the gradient (source) terms.
45:
46: The Jacobian can be either
47: * formed via finite differencing using coloring (the default), or
48: * applied matrix-free via the option -snes_mf
49: (for larger grid problems this variant may not converge
50: without a preconditioner due to ill-conditioning).
52: ------------------------------------------------------------------------- */
54: /*
55: Include "petscda.h" so that we can use distributed arrays (DAs).
56: Include "petscsnes.h" so that we can use SNES solvers. Note that this
57: file automatically includes:
58: petsc.h - base PETSc routines petscvec.h - vectors
59: petscsys.h - system routines petscmat.h - matrices
60: petscis.h - index sets petscksp.h - Krylov subspace methods
61: petscviewer.h - viewers petscpc.h - preconditioners
62: petscksp.h - linear solvers
63: */
64: #include petscsnes.h
65: #include petscda.h
67: /*
68: User-defined routines and data structures
69: */
70: typedef struct {
71: PetscScalar u,v,omega,temp;
72: } Field;
78: typedef struct {
79: PassiveReal lidvelocity,prandtl,grashof; /* physical parameters */
80: PetscTruth draw_contours; /* flag - 1 indicates drawing contours */
81: } AppCtx;
85: int main(int argc,char **argv)
86: {
87: DMMG *dmmg; /* multilevel grid structure */
88: AppCtx user; /* user-defined work context */
89: PetscInt mx,my,its;
91: MPI_Comm comm;
92: SNES snes;
93: DA da;
95: PetscInitialize(&argc,&argv,(char *)0,help);
96: comm = PETSC_COMM_WORLD;
99: PreLoadBegin(PETSC_TRUE,"SetUp");
100: DMMGCreate(comm,2,&user,&dmmg);
103: /*
104: Create distributed array multigrid object (DMMG) to manage parallel grid and vectors
105: for principal unknowns (x) and governing residuals (f)
106: */
107: DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,4,1,0,0,&da);
108: DMMGSetDM(dmmg,(DM)da);
109: DADestroy(da);
111: DAGetInfo(DMMGGetDA(dmmg),0,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
112: PETSC_IGNORE,PETSC_IGNORE);
113: /*
114: Problem parameters (velocity of lid, prandtl, and grashof numbers)
115: */
116: user.lidvelocity = 1.0/(mx*my);
117: user.prandtl = 1.0;
118: user.grashof = 1.0;
119: PetscOptionsGetReal(PETSC_NULL,"-lidvelocity",&user.lidvelocity,PETSC_NULL);
120: PetscOptionsGetReal(PETSC_NULL,"-prandtl",&user.prandtl,PETSC_NULL);
121: PetscOptionsGetReal(PETSC_NULL,"-grashof",&user.grashof,PETSC_NULL);
122: PetscOptionsHasName(PETSC_NULL,"-contours",&user.draw_contours);
124: DASetFieldName(DMMGGetDA(dmmg),0,"x-velocity");
125: DASetFieldName(DMMGGetDA(dmmg),1,"y-velocity");
126: DASetFieldName(DMMGGetDA(dmmg),2,"Omega");
127: DASetFieldName(DMMGGetDA(dmmg),3,"temperature");
129: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130: Create user context, set problem data, create vector data structures.
131: Also, compute the initial guess.
132: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
134: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135: Create nonlinear solver context
137: Process adiC(36): FormFunctionLocal FormFunctionLocali
138: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139: DMMGSetSNESLocal(dmmg,FormFunctionLocal,0,ad_FormFunctionLocal,admf_FormFunctionLocal);
140: DMMGSetSNESLocali(dmmg,FormFunctionLocali,ad_FormFunctionLocali,admf_FormFunctionLocali);
142: PetscPrintf(comm,"lid velocity = %g, prandtl # = %g, grashof # = %g\n",
143: user.lidvelocity,user.prandtl,user.grashof);
146: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147: Solve the nonlinear system
148: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149: DMMGSetInitialGuess(dmmg,FormInitialGuess);
151: PreLoadStage("Solve");
152: DMMGSolve(dmmg);
154: snes = DMMGGetSNES(dmmg);
155: SNESGetIterationNumber(snes,&its);
156: PetscPrintf(comm,"Number of Newton iterations = %D\n", its);
158: /*
159: Visualize solution
160: */
162: if (user.draw_contours) {
163: VecView(DMMGGetx(dmmg),PETSC_VIEWER_DRAW_WORLD);
164: }
166: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167: Free work space. All PETSc objects should be destroyed when they
168: are no longer needed.
169: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
171: DMMGDestroy(dmmg);
172: PreLoadEnd();
174: PetscFinalize();
175: return 0;
176: }
178: /* ------------------------------------------------------------------- */
183: /*
184: FormInitialGuess - Forms initial approximation.
186: Input Parameters:
187: user - user-defined application context
188: X - vector
190: Output Parameter:
191: X - vector
192: */
193: PetscErrorCode FormInitialGuess(DMMG dmmg,Vec X)
194: {
195: AppCtx *user = (AppCtx*)dmmg->user;
196: DA da = (DA)dmmg->dm;
197: PetscInt i,j,mx,xs,ys,xm,ym;
199: PetscReal grashof,dx;
200: Field **x;
202: grashof = user->grashof;
204: DAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0);
205: dx = 1.0/(mx-1);
207: /*
208: Get local grid boundaries (for 2-dimensional DA):
209: xs, ys - starting grid indices (no ghost points)
210: xm, ym - widths of local grid (no ghost points)
211: */
212: DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
214: /*
215: Get a pointer to vector data.
216: - For default PETSc vectors, VecGetArray() returns a pointer to
217: the data array. Otherwise, the routine is implementation dependent.
218: - You MUST call VecRestoreArray() when you no longer need access to
219: the array.
220: */
221: DAVecGetArray(da,X,&x);
223: /*
224: Compute initial guess over the locally owned part of the grid
225: Initial condition is motionless fluid and equilibrium temperature
226: */
227: for (j=ys; j<ys+ym; j++) {
228: for (i=xs; i<xs+xm; i++) {
229: x[j][i].u = 0.0;
230: x[j][i].v = 0.0;
231: x[j][i].omega = 0.0;
232: x[j][i].temp = (grashof>0)*i*dx;
233: }
234: }
236: /*
237: Restore vector
238: */
239: DAVecRestoreArray(da,X,&x);
240: return 0;
241: }
242: PetscErrorCode FormFunctionLocal(DALocalInfo *info,Field **x,Field **f,void *ptr)
243: {
244: AppCtx *user = (AppCtx*)ptr;
246: PetscInt xints,xinte,yints,yinte,i,j;
247: PetscReal hx,hy,dhx,dhy,hxdhy,hydhx;
248: PetscReal grashof,prandtl,lid;
249: PetscScalar u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;
252: grashof = user->grashof;
253: prandtl = user->prandtl;
254: lid = user->lidvelocity;
256: /*
257: Define mesh intervals ratios for uniform grid.
258: [Note: FD formulae below are normalized by multiplying through by
259: local volume element to obtain coefficients O(1) in two dimensions.]
260: */
261: dhx = (PetscReal)(info->mx-1); dhy = (PetscReal)(info->my-1);
262: hx = 1.0/dhx; hy = 1.0/dhy;
263: hxdhy = hx*dhy; hydhx = hy*dhx;
265: xints = info->xs; xinte = info->xs+info->xm; yints = info->ys; yinte = info->ys+info->ym;
267: /* Test whether we are on the bottom edge of the global array */
268: if (yints == 0) {
269: j = 0;
270: yints = yints + 1;
271: /* bottom edge */
272: for (i=info->xs; i<info->xs+info->xm; i++) {
273: f[j][i].u = x[j][i].u;
274: f[j][i].v = x[j][i].v;
275: f[j][i].omega = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
276: f[j][i].temp = x[j][i].temp-x[j+1][i].temp;
277: }
278: }
280: /* Test whether we are on the top edge of the global array */
281: if (yinte == info->my) {
282: j = info->my - 1;
283: yinte = yinte - 1;
284: /* top edge */
285: for (i=info->xs; i<info->xs+info->xm; i++) {
286: f[j][i].u = x[j][i].u - lid;
287: f[j][i].v = x[j][i].v;
288: f[j][i].omega = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
289: f[j][i].temp = x[j][i].temp-x[j-1][i].temp;
290: }
291: }
293: /* Test whether we are on the left edge of the global array */
294: if (xints == 0) {
295: i = 0;
296: xints = xints + 1;
297: /* left edge */
298: for (j=info->ys; j<info->ys+info->ym; j++) {
299: f[j][i].u = x[j][i].u;
300: f[j][i].v = x[j][i].v;
301: f[j][i].omega = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
302: f[j][i].temp = x[j][i].temp;
303: }
304: }
306: /* Test whether we are on the right edge of the global array */
307: if (xinte == info->mx) {
308: i = info->mx - 1;
309: xinte = xinte - 1;
310: /* right edge */
311: for (j=info->ys; j<info->ys+info->ym; j++) {
312: f[j][i].u = x[j][i].u;
313: f[j][i].v = x[j][i].v;
314: f[j][i].omega = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
315: f[j][i].temp = x[j][i].temp - (PetscReal)(grashof>0);
316: }
317: }
319: /* Compute over the interior points */
320: for (j=yints; j<yinte; j++) {
321: for (i=xints; i<xinte; i++) {
323: /*
324: convective coefficients for upwinding
325: */
326: vx = x[j][i].u; avx = PetscAbsScalar(vx);
327: vxp = .5*(vx+avx); vxm = .5*(vx-avx);
328: vy = x[j][i].v; avy = PetscAbsScalar(vy);
329: vyp = .5*(vy+avy); vym = .5*(vy-avy);
331: /* U velocity */
332: u = x[j][i].u;
333: uxx = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
334: uyy = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
335: f[j][i].u = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;
337: /* V velocity */
338: u = x[j][i].v;
339: uxx = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
340: uyy = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
341: f[j][i].v = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;
343: /* Omega */
344: u = x[j][i].omega;
345: uxx = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
346: uyy = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
347: f[j][i].omega = uxx + uyy +
348: (vxp*(u - x[j][i-1].omega) +
349: vxm*(x[j][i+1].omega - u)) * hy +
350: (vyp*(u - x[j-1][i].omega) +
351: vym*(x[j+1][i].omega - u)) * hx -
352: .5 * grashof * (x[j][i+1].temp - x[j][i-1].temp) * hy;
354: /* Temperature */
355: u = x[j][i].temp;
356: uxx = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
357: uyy = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
358: f[j][i].temp = uxx + uyy + prandtl * (
359: (vxp*(u - x[j][i-1].temp) +
360: vxm*(x[j][i+1].temp - u)) * hy +
361: (vyp*(u - x[j-1][i].temp) +
362: vym*(x[j+1][i].temp - u)) * hx);
363: }
364: }
366: /*
367: Flop count (multiply-adds are counted as 2 operations)
368: */
369: PetscLogFlops(84*info->ym*info->xm);
370: return(0);
371: }
373: /*
374: This is an experimental function and can be safely ignored.
375: */
376: PetscErrorCode FormFunctionLocali(DALocalInfo *info,MatStencil *st,Field **x,PetscScalar *f,void *ptr)
377: {
378: AppCtx *user = (AppCtx*)ptr;
379: PetscInt i,j,c;
380: PassiveReal hx,hy,dhx,dhy,hxdhy,hydhx;
381: PassiveReal grashof,prandtl,lid;
382: PetscScalar u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;
385: grashof = user->grashof;
386: prandtl = user->prandtl;
387: lid = user->lidvelocity;
389: /*
390: Define mesh intervals ratios for uniform grid.
391: [Note: FD formulae below are normalized by multiplying through by
392: local volume element to obtain coefficients O(1) in two dimensions.]
393: */
394: dhx = (PetscReal)(info->mx-1); dhy = (PetscReal)(info->my-1);
395: hx = 1.0/dhx; hy = 1.0/dhy;
396: hxdhy = hx*dhy; hydhx = hy*dhx;
398: i = st->i; j = st->j; c = st->c;
400: /* Test whether we are on the right edge of the global array */
401: if (i == info->mx-1) {
402: if (c == 0) *f = x[j][i].u;
403: else if (c == 1) *f = x[j][i].v;
404: else if (c == 2) *f = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
405: else *f = x[j][i].temp - (PetscReal)(grashof>0);
407: /* Test whether we are on the left edge of the global array */
408: } else if (i == 0) {
409: if (c == 0) *f = x[j][i].u;
410: else if (c == 1) *f = x[j][i].v;
411: else if (c == 2) *f = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
412: else *f = x[j][i].temp;
414: /* Test whether we are on the top edge of the global array */
415: } else if (j == info->my-1) {
416: if (c == 0) *f = x[j][i].u - lid;
417: else if (c == 1) *f = x[j][i].v;
418: else if (c == 2) *f = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
419: else *f = x[j][i].temp-x[j-1][i].temp;
421: /* Test whether we are on the bottom edge of the global array */
422: } else if (j == 0) {
423: if (c == 0) *f = x[j][i].u;
424: else if (c == 1) *f = x[j][i].v;
425: else if (c == 2) *f = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
426: else *f = x[j][i].temp - x[j+1][i].temp;
428: /* Compute over the interior points */
429: } else {
430: /*
431: convective coefficients for upwinding
432: */
433: vx = x[j][i].u; avx = PetscAbsScalar(vx);
434: vxp = .5*(vx+avx); vxm = .5*(vx-avx);
435: vy = x[j][i].v; avy = PetscAbsScalar(vy);
436: vyp = .5*(vy+avy); vym = .5*(vy-avy);
438: /* U velocity */
439: if (c == 0) {
440: u = x[j][i].u;
441: uxx = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
442: uyy = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
443: *f = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;
445: /* V velocity */
446: } else if (c == 1) {
447: u = x[j][i].v;
448: uxx = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
449: uyy = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
450: *f = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;
451:
452: /* Omega */
453: } else if (c == 2) {
454: u = x[j][i].omega;
455: uxx = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
456: uyy = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
457: *f = uxx + uyy +
458: (vxp*(u - x[j][i-1].omega) +
459: vxm*(x[j][i+1].omega - u)) * hy +
460: (vyp*(u - x[j-1][i].omega) +
461: vym*(x[j+1][i].omega - u)) * hx -
462: .5 * grashof * (x[j][i+1].temp - x[j][i-1].temp) * hy;
463:
464: /* Temperature */
465: } else {
466: u = x[j][i].temp;
467: uxx = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
468: uyy = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
469: *f = uxx + uyy + prandtl * (
470: (vxp*(u - x[j][i-1].temp) +
471: vxm*(x[j][i+1].temp - u)) * hy +
472: (vyp*(u - x[j-1][i].temp) +
473: vym*(x[j+1][i].temp - u)) * hx);
474: }
475: }
477: return(0);
478: }