Actual source code: ex19.c
petsc-3.12.5 2020-03-29
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" ;
10: /* in HTML, '<' = '<' and '>' = '>' */
12: /*
13: See src/ksp/ksp/examples/tutorials/ex45.c
14: */
16: /*T
17: Concepts: SNES ^solving a system of nonlinear equations (parallel multicomponent example);
18: Concepts: DMDA ^using distributed arrays;
19: Concepts: multicomponent
20: Processors: n
21: T*/
We thank David E. Keyes for contributing the driven cavity discretization within this example code.
This problem is modeled by the partial differential equation system
\begin{eqnarray}
- \triangle U - \nabla_y \Omega & = & 0 \\
- \triangle V + \nabla_x\Omega & = & 0 \\
- \triangle \Omega + \nabla \cdot ([U*\Omega,V*\Omega]) - GR* \nabla_x T & = & 0 \\
- \triangle T + PR* \nabla \cdot ([U*T,V*T]) & = & 0
\end{eqnarray}
in the unit square, which is uniformly discretized in each of x and y in this simple encoding.
No-slip, rigid-wall Dirichlet conditions are used for $ [U,V]$.
Dirichlet conditions are used for Omega, based on the definition of
vorticity: $ \Omega = - \nabla_y U + \nabla_x V$, where along each
constant coordinate boundary, the tangential derivative is zero.
Dirichlet conditions are used for T on the left and right walls,
and insulation homogeneous Neumann conditions are used for T on
the top and bottom walls.
A finite difference approximation with the usual 5-point stencil
is used to discretize the boundary value problem to obtain a
nonlinear system of equations. Upwinding is used for the divergence
(convective) terms and central for the gradient (source) terms.
The Jacobian can be either
* formed via finite differencing using coloring (the default), or
* applied matrix-free via the option -snes_mf
(for larger grid problems this variant may not converge
without a preconditioner due to ill-conditioning).
60: /*
61: Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
62: Include "petscsnes.h" so that we can use SNES solvers. Note that this
63: file automatically includes:
64: petscsys.h - base PETSc routines petscvec.h - vectors
65: petscmat.h - matrices
66: petscis.h - index sets petscksp.h - Krylov subspace methods
67: petscviewer.h - viewers petscpc.h - preconditioners
68: petscksp.h - linear solvers
69: */
70: #if defined(PETSC_APPLE_FRAMEWORK)
71: #import <PETSc/petscsnes.h>
72: #import <PETSc/petscdmda.h>
73: #else
74: #include <petscsnes.h>
75: #include <petscdm.h>
76: #include <petscdmda.h>
77: #endif
79: /*
80: User-defined routines and data structures
81: */
82: typedef struct {
83: PetscScalar u,v,omega,temp;
84: } Field;
86: PetscErrorCode FormFunctionLocal(DMDALocalInfo *,Field**,Field**,void*) ;
88: typedef struct {
89: PetscReal lidvelocity,prandtl,grashof; /* physical parameters */
90: PetscBool draw_contours; /* flag - 1 indicates drawing contours */
91: } AppCtx;
93: extern PetscErrorCode FormInitialGuess(AppCtx*,DM ,Vec ) ;
94: extern PetscErrorCode NonlinearGS(SNES ,Vec ,Vec ,void*) ;
96: int main(int argc,char **argv)
97: {
98: AppCtx user; /* user-defined work context */
99: PetscInt mx,my,its;
101: MPI_Comm comm;
102: SNES snes;
103: DM da;
104: Vec x;
106: PetscInitialize (&argc,&argv,(char*)0,help);if (ierr) return ierr;
109: comm = PETSC_COMM_WORLD ;
110: SNESCreate (comm,&snes);
112: /*
113: Create distributed array object to manage parallel grid and vectors
114: for principal unknowns (x) and governing residuals (f)
115: */
116: DMDACreate2d (PETSC_COMM_WORLD ,DM_BOUNDARY_NONE ,DM_BOUNDARY_NONE ,DMDA_STENCIL_STAR ,4,4,PETSC_DECIDE ,PETSC_DECIDE ,4,1,0,0,&da);
117: DMSetFromOptions (da);
118: DMSetUp (da);
119: SNESSetDM (snes,(DM )da);
120: SNESSetNGS (snes, NonlinearGS, (void*)&user);
122: DMDAGetInfo (da,0,&mx,&my,PETSC_IGNORE ,PETSC_IGNORE ,PETSC_IGNORE ,PETSC_IGNORE ,PETSC_IGNORE ,PETSC_IGNORE ,PETSC_IGNORE ,PETSC_IGNORE ,PETSC_IGNORE ,PETSC_IGNORE );
123: /*
124: Problem parameters (velocity of lid, prandtl, and grashof numbers)
125: */
126: user.lidvelocity = 1.0/(mx*my);
127: user.prandtl = 1.0;
128: user.grashof = 1.0;
130: PetscOptionsGetReal (NULL,NULL,"-lidvelocity" ,&user.lidvelocity,NULL);
131: PetscOptionsGetReal (NULL,NULL,"-prandtl" ,&user.prandtl,NULL);
132: PetscOptionsGetReal (NULL,NULL,"-grashof" ,&user.grashof,NULL);
133: PetscOptionsHasName (NULL,NULL,"-contours" ,&user.draw_contours);
135: DMDASetFieldName (da,0,"x_velocity" );
136: DMDASetFieldName (da,1,"y_velocity" );
137: DMDASetFieldName (da,2,"Omega" );
138: DMDASetFieldName (da,3,"temperature" );
140: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141: Create user context, set problem data, create vector data structures.
142: Also, compute the initial guess.
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146: Create nonlinear solver context
148: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149: DMSetApplicationContext (da,&user);
150: DMDASNESSetFunctionLocal (da,INSERT_VALUES ,(PetscErrorCode (*)(DMDALocalInfo *,void*,void*,void*))FormFunctionLocal,&user);
151: SNESSetFromOptions (snes);
152: PetscPrintf (comm,"lid velocity = %g, prandtl # = %g, grashof # = %g\n" ,(double)user.lidvelocity,(double)user.prandtl,(double)user.grashof);
155: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156: Solve the nonlinear system
157: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
158: DMCreateGlobalVector (da,&x);
159: FormInitialGuess(&user,da,x);
161: SNESSolve (snes,NULL,x);
163: SNESGetIterationNumber (snes,&its);
164: PetscPrintf (comm,"Number of SNES iterations = %D\n" , its);
166: /*
167: Visualize solution
168: */
169: if (user.draw_contours) {
170: VecView (x,PETSC_VIEWER_DRAW_WORLD );
171: }
173: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
174: Free work space. All PETSc objects should be destroyed when they
175: are no longer needed.
176: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
177: VecDestroy (&x);
178: DMDestroy (&da);
179: SNESDestroy (&snes);
180: PetscFinalize ();
181: return ierr;
182: }
184: /* ------------------------------------------------------------------- */
186: /*
187: FormInitialGuess - Forms initial approximation.
189: Input Parameters:
190: user - user-defined application context
191: X - vector
193: Output Parameter:
194: X - vector
195: */
196: PetscErrorCode FormInitialGuess(AppCtx *user,DM da,Vec X)
197: {
198: PetscInt i,j,mx,xs,ys,xm,ym;
200: PetscReal grashof,dx;
201: Field **x;
204: grashof = user->grashof;
206: DMDAGetInfo (da,0,&mx,0,0,0,0,0,0,0,0,0,0,0);
207: dx = 1.0/(mx-1);
209: /*
210: Get local grid boundaries (for 2-dimensional DMDA ):
211: xs, ys - starting grid indices (no ghost points)
212: xm, ym - widths of local grid (no ghost points)
213: */
214: DMDAGetCorners (da,&xs,&ys,NULL,&xm,&ym,NULL);
216: /*
217: Get a pointer to vector data.
218: - For default PETSc vectors, VecGetArray () returns a pointer to
219: the data array. Otherwise, the routine is implementation dependent.
220: - You MUST call VecRestoreArray () when you no longer need access to
221: the array.
222: */
223: DMDAVecGetArrayWrite (da,X,&x);
225: /*
226: Compute initial guess over the locally owned part of the grid
227: Initial condition is motionless fluid and equilibrium temperature
228: */
229: for (j=ys; j<ys+ym; j++) {
230: for (i=xs; i<xs+xm; i++) {
231: x[j][i].u = 0.0;
232: x[j][i].v = 0.0;
233: x[j][i].omega = 0.0;
234: x[j][i].temp = (grashof>0)*i*dx;
235: }
236: }
238: /*
239: Restore vector
240: */
241: DMDAVecRestoreArrayWrite (da,X,&x);
242: return (0);
243: }
245: PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,Field **x,Field **f,void *ptr)
246: {
247: AppCtx *user = (AppCtx*)ptr;
249: PetscInt xints,xinte,yints,yinte,i,j;
250: PetscReal hx,hy,dhx,dhy,hxdhy,hydhx;
251: PetscReal grashof,prandtl,lid;
252: PetscScalar u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;
255: grashof = user->grashof;
256: prandtl = user->prandtl;
257: lid = user->lidvelocity;
259: /*
260: Define mesh intervals ratios for uniform grid.
262: Note: FD formulae below are normalized by multiplying through by
263: local volume element (i.e. hx*hy) to obtain coefficients O(1) in two dimensions.
266: */
267: dhx = (PetscReal )(info->mx-1); dhy = (PetscReal )(info->my-1);
268: hx = 1.0/dhx; hy = 1.0/dhy;
269: hxdhy = hx*dhy; hydhx = hy*dhx;
271: xints = info->xs; xinte = info->xs+info->xm; yints = info->ys; yinte = info->ys+info->ym;
273: /* Test whether we are on the bottom edge of the global array */
274: if (yints == 0) {
275: j = 0;
276: yints = yints + 1;
277: /* bottom edge */
278: for (i=info->xs; i<info->xs+info->xm; i++) {
279: f[j][i].u = x[j][i].u;
280: f[j][i].v = x[j][i].v;
281: f[j][i].omega = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
282: f[j][i].temp = x[j][i].temp-x[j+1][i].temp;
283: }
284: }
286: /* Test whether we are on the top edge of the global array */
287: if (yinte == info->my) {
288: j = info->my - 1;
289: yinte = yinte - 1;
290: /* top edge */
291: for (i=info->xs; i<info->xs+info->xm; i++) {
292: f[j][i].u = x[j][i].u - lid;
293: f[j][i].v = x[j][i].v;
294: f[j][i].omega = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
295: f[j][i].temp = x[j][i].temp-x[j-1][i].temp;
296: }
297: }
299: /* Test whether we are on the left edge of the global array */
300: if (xints == 0) {
301: i = 0;
302: xints = xints + 1;
303: /* left edge */
304: for (j=info->ys; j<info->ys+info->ym; j++) {
305: f[j][i].u = x[j][i].u;
306: f[j][i].v = x[j][i].v;
307: f[j][i].omega = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
308: f[j][i].temp = x[j][i].temp;
309: }
310: }
312: /* Test whether we are on the right edge of the global array */
313: if (xinte == info->mx) {
314: i = info->mx - 1;
315: xinte = xinte - 1;
316: /* right edge */
317: for (j=info->ys; j<info->ys+info->ym; j++) {
318: f[j][i].u = x[j][i].u;
319: f[j][i].v = x[j][i].v;
320: f[j][i].omega = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
321: f[j][i].temp = x[j][i].temp - (PetscReal )(grashof>0);
322: }
323: }
325: /* Compute over the interior points */
326: for (j=yints; j<yinte; j++) {
327: for (i=xints; i<xinte; i++) {
329: /*
330: convective coefficients for upwinding
331: */
332: vx = x[j][i].u; avx = PetscAbsScalar(vx);
333: vxp = .5*(vx+avx); vxm = .5*(vx-avx);
334: vy = x[j][i].v; avy = PetscAbsScalar(vy);
335: vyp = .5*(vy+avy); vym = .5*(vy-avy);
337: /* U velocity */
338: u = x[j][i].u;
339: uxx = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
340: uyy = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
341: f[j][i].u = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;
343: /* V velocity */
344: u = x[j][i].v;
345: uxx = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
346: uyy = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
347: f[j][i].v = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;
349: /* Omega */
350: u = x[j][i].omega;
351: uxx = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
352: uyy = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
353: f[j][i].omega = uxx + uyy + (vxp*(u - x[j][i-1].omega) + vxm*(x[j][i+1].omega - u))*hy +
354: (vyp*(u - x[j-1][i].omega) + vym*(x[j+1][i].omega - u))*hx -
355: .5*grashof*(x[j][i+1].temp - x[j][i-1].temp)*hy;
357: /* Temperature */
358: u = x[j][i].temp;
359: uxx = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
360: uyy = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
361: f[j][i].temp = uxx + uyy + prandtl*((vxp*(u - x[j][i-1].temp) + vxm*(x[j][i+1].temp - u))*hy +
362: (vyp*(u - x[j-1][i].temp) + 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.0*info->ym*info->xm);
370: return (0);
371: }
373: /*
374: Performs sweeps of point block nonlinear Gauss-Seidel on all the local grid points
375: */
376: PetscErrorCode NonlinearGS(SNES snes, Vec X, Vec B, void *ctx)
377: {
378: DMDALocalInfo info;
379: Field **x,**b;
381: Vec localX, localB;
382: DM da;
383: PetscInt xints,xinte,yints,yinte,i,j,k,l;
384: PetscInt max_its,tot_its;
385: PetscInt sweeps;
386: PetscReal rtol,atol,stol;
387: PetscReal hx,hy,dhx,dhy,hxdhy,hydhx;
388: PetscReal grashof,prandtl,lid;
389: PetscScalar u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;
390: PetscScalar fu, fv, fomega, ftemp;
391: PetscScalar dfudu;
392: PetscScalar dfvdv;
393: PetscScalar dfodu, dfodv, dfodo;
394: PetscScalar dftdu, dftdv, dftdt;
395: PetscScalar yu=0, yv=0, yo=0, yt=0;
396: PetscScalar bjiu, bjiv, bjiomega, bjitemp;
397: PetscBool ptconverged;
398: PetscReal pfnorm,pfnorm0,pynorm,pxnorm;
399: AppCtx *user = (AppCtx*)ctx;
402: grashof = user->grashof;
403: prandtl = user->prandtl;
404: lid = user->lidvelocity;
405: tot_its = 0;
406: SNESNGSGetTolerances (snes,&rtol,&atol,&stol,&max_its);
407: SNESNGSGetSweeps (snes,&sweeps);
408: SNESGetDM (snes,(DM *)&da);
409: DMGetLocalVector (da,&localX);
410: if (B) {
411: DMGetLocalVector (da,&localB);
412: }
413: /*
414: Scatter ghost points to local vector, using the 2-step process
415: DMGlobalToLocalBegin (), DMGlobalToLocalEnd ().
416: */
417: DMGlobalToLocalBegin (da,X,INSERT_VALUES ,localX);
418: DMGlobalToLocalEnd (da,X,INSERT_VALUES ,localX);
419: if (B) {
420: DMGlobalToLocalBegin (da,B,INSERT_VALUES ,localB);
421: DMGlobalToLocalEnd (da,B,INSERT_VALUES ,localB);
422: }
423: DMDAGetLocalInfo (da,&info);
424: DMDAVecGetArrayWrite (da,localX,&x);
425: if (B) {
426: DMDAVecGetArrayRead (da,localB,&b);
427: }
428: /* looks like a combination of the formfunction / formjacobian routines */
429: dhx = (PetscReal )(info.mx-1);dhy = (PetscReal )(info.my-1);
430: hx = 1.0/dhx; hy = 1.0/dhy;
431: hxdhy = hx*dhy; hydhx = hy*dhx;
433: xints = info.xs; xinte = info.xs+info.xm; yints = info.ys; yinte = info.ys+info.ym;
435: /* Set the boundary conditions on the momentum equations */
436: /* Test whether we are on the bottom edge of the global array */
437: if (yints == 0) {
438: j = 0;
439: /* bottom edge */
440: for (i=info.xs; i<info.xs+info.xm; i++) {
442: if (B) {
443: bjiu = b[j][i].u;
444: bjiv = b[j][i].v;
445: } else {
446: bjiu = 0.0;
447: bjiv = 0.0;
448: }
449: x[j][i].u = 0.0 + bjiu;
450: x[j][i].v = 0.0 + bjiv;
451: }
452: }
454: /* Test whether we are on the top edge of the global array */
455: if (yinte == info.my) {
456: j = info.my - 1;
457: /* top edge */
458: for (i=info.xs; i<info.xs+info.xm; i++) {
459: if (B) {
460: bjiu = b[j][i].u;
461: bjiv = b[j][i].v;
462: } else {
463: bjiu = 0.0;
464: bjiv = 0.0;
465: }
466: x[j][i].u = lid + bjiu;
467: x[j][i].v = bjiv;
468: }
469: }
471: /* Test whether we are on the left edge of the global array */
472: if (xints == 0) {
473: i = 0;
474: /* left edge */
475: for (j=info.ys; j<info.ys+info.ym; j++) {
476: if (B) {
477: bjiu = b[j][i].u;
478: bjiv = b[j][i].v;
479: } else {
480: bjiu = 0.0;
481: bjiv = 0.0;
482: }
483: x[j][i].u = 0.0 + bjiu;
484: x[j][i].v = 0.0 + bjiv;
485: }
486: }
488: /* Test whether we are on the right edge of the global array */
489: if (xinte == info.mx) {
490: i = info.mx - 1;
491: /* right edge */
492: for (j=info.ys; j<info.ys+info.ym; j++) {
493: if (B) {
494: bjiu = b[j][i].u;
495: bjiv = b[j][i].v;
496: } else {
497: bjiu = 0.0;
498: bjiv = 0.0;
499: }
500: x[j][i].u = 0.0 + bjiu;
501: x[j][i].v = 0.0 + bjiv;
502: }
503: }
505: for (k=0; k < sweeps; k++) {
506: for (j=info.ys; j<info.ys + info.ym; j++) {
507: for (i=info.xs; i<info.xs + info.xm; i++) {
508: ptconverged = PETSC_FALSE ;
509: pfnorm0 = 0.0;
510: fu = 0.0;
511: fv = 0.0;
512: fomega = 0.0;
513: ftemp = 0.0;
514: /* Run Newton's method on a single grid point */
515: for (l = 0; l < max_its && !ptconverged; l++) {
516: if (B) {
517: bjiu = b[j][i].u;
518: bjiv = b[j][i].v;
519: bjiomega = b[j][i].omega;
520: bjitemp = b[j][i].temp;
521: } else {
522: bjiu = 0.0;
523: bjiv = 0.0;
524: bjiomega = 0.0;
525: bjitemp = 0.0;
526: }
528: if (i != 0 && i != info.mx - 1 && j != 0 && j != info.my-1) {
529: /* U velocity */
530: u = x[j][i].u;
531: uxx = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
532: uyy = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
533: fu = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx - bjiu;
534: dfudu = 2.0*(hydhx + hxdhy);
535: /* V velocity */
536: u = x[j][i].v;
537: uxx = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
538: uyy = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
539: fv = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy - bjiv;
540: dfvdv = 2.0*(hydhx + hxdhy);
541: /*
542: convective coefficients for upwinding
543: */
544: vx = x[j][i].u; avx = PetscAbsScalar(vx);
545: vxp = .5*(vx+avx); vxm = .5*(vx-avx);
546: vy = x[j][i].v; avy = PetscAbsScalar(vy);
547: vyp = .5*(vy+avy); vym = .5*(vy-avy);
548: /* Omega */
549: u = x[j][i].omega;
550: uxx = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
551: uyy = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
552: fomega = uxx + uyy + (vxp*(u - x[j][i-1].omega) + vxm*(x[j][i+1].omega - u))*hy +
553: (vyp*(u - x[j-1][i].omega) + vym*(x[j+1][i].omega - u))*hx -
554: .5*grashof*(x[j][i+1].temp - x[j][i-1].temp)*hy - bjiomega;
555: /* convective coefficient derivatives */
556: dfodo = 2.0*(hydhx + hxdhy) + ((vxp - vxm)*hy + (vyp - vym)*hx);
557: if (PetscRealPart (vx) > 0.0) dfodu = (u - x[j][i-1].omega)*hy;
558: else dfodu = (x[j][i+1].omega - u)*hy;
560: if (PetscRealPart (vy) > 0.0) dfodv = (u - x[j-1][i].omega)*hx;
561: else dfodv = (x[j+1][i].omega - u)*hx;
563: /* Temperature */
564: u = x[j][i].temp;
565: uxx = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
566: uyy = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
567: ftemp = uxx + uyy + prandtl*((vxp*(u - x[j][i-1].temp) + vxm*(x[j][i+1].temp - u))*hy + (vyp*(u - x[j-1][i].temp) + vym*(x[j+1][i].temp - u))*hx) - bjitemp;
568: dftdt = 2.0*(hydhx + hxdhy) + prandtl*((vxp - vxm)*hy + (vyp - vym)*hx);
569: if (PetscRealPart (vx) > 0.0) dftdu = prandtl*(u - x[j][i-1].temp)*hy;
570: else dftdu = prandtl*(x[j][i+1].temp - u)*hy;
572: if (PetscRealPart (vy) > 0.0) dftdv = prandtl*(u - x[j-1][i].temp)*hx;
573: else dftdv = prandtl*(x[j+1][i].temp - u)*hx;
575: /* invert the system:
576: [ dfu / du 0 0 0 ][yu] = [fu]
577: [ 0 dfv / dv 0 0 ][yv] [fv]
578: [ dfo / du dfo / dv dfo / do 0 ][yo] [fo]
579: [ dft / du dft / dv 0 dft / dt ][yt] [ft]
580: by simple back-substitution
581: */
582: yu = fu / dfudu;
583: yv = fv / dfvdv;
584: yo = (fomega - (dfodu*yu + dfodv*yv)) / dfodo;
585: yt = (ftemp - (dftdu*yu + dftdv*yv)) / dftdt;
587: x[j][i].u = x[j][i].u - yu;
588: x[j][i].v = x[j][i].v - yv;
589: x[j][i].temp = x[j][i].temp - yt;
590: x[j][i].omega = x[j][i].omega - yo;
591: }
592: if (i == 0) {
593: fomega = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx - bjiomega;
594: ftemp = x[j][i].temp - bjitemp;
595: yo = fomega;
596: yt = ftemp;
597: x[j][i].omega = x[j][i].omega - fomega;
598: x[j][i].temp = x[j][i].temp - ftemp;
599: }
600: if (i == info.mx - 1) {
601: fomega = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx - bjiomega;
602: ftemp = x[j][i].temp - (PetscReal )(grashof>0) - bjitemp;
603: yo = fomega;
604: yt = ftemp;
605: x[j][i].omega = x[j][i].omega - fomega;
606: x[j][i].temp = x[j][i].temp - ftemp;
607: }
608: if (j == 0) {
609: fomega = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy - bjiomega;
610: ftemp = x[j][i].temp-x[j+1][i].temp - bjitemp;
611: yo = fomega;
612: yt = ftemp;
613: x[j][i].omega = x[j][i].omega - fomega;
614: x[j][i].temp = x[j][i].temp - ftemp;
615: }
616: if (j == info.my - 1) {
617: fomega = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy - bjiomega;
618: ftemp = x[j][i].temp-x[j-1][i].temp - bjitemp;
619: yo = fomega;
620: yt = ftemp;
621: x[j][i].omega = x[j][i].omega - fomega;
622: x[j][i].temp = x[j][i].temp - ftemp;
623: }
624: tot_its++;
625: pfnorm = PetscRealPart (fu*fu + fv*fv + fomega*fomega + ftemp*ftemp);
626: pfnorm = PetscSqrtReal(pfnorm);
627: pynorm = PetscRealPart (yu*yu + yv*yv + yo*yo + yt*yt);
628: pynorm = PetscSqrtReal(pynorm);
629: pxnorm = PetscRealPart (x[j][i].u*x[j][i].u + x[j][i].v*x[j][i].v + x[j][i].omega*x[j][i].omega + x[j][i].temp*x[j][i].temp);
630: pxnorm = PetscSqrtReal(pxnorm);
631: if (l == 0) pfnorm0 = pfnorm;
632: if (rtol*pfnorm0 >pfnorm || atol > pfnorm || pxnorm*stol > pynorm) ptconverged = PETSC_TRUE ;
633: }
634: }
635: }
636: }
637: DMDAVecRestoreArrayWrite (da,localX,&x);
638: if (B) {
639: DMDAVecRestoreArrayRead (da,localB,&b);
640: }
641: DMLocalToGlobalBegin (da,localX,INSERT_VALUES ,X);
642: DMLocalToGlobalEnd (da,localX,INSERT_VALUES ,X);
643: PetscLogFlops (tot_its*(84.0 + 41.0 + 26.0));
644: DMRestoreLocalVector (da,&localX);
645: if (B) {
646: DMRestoreLocalVector (da,&localB);
647: }
648: return (0);
649: }
652: /*TEST
654: test:
655: nsize: 2
656: args: -da_refine 3 -snes_monitor_short -pc_type mg -ksp_type fgmres -pc_mg_type full
657: requires: !single
659: test:
660: suffix: 10
661: nsize: 3
662: args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -pc_fieldsplit_type symmetric_multiplicative -snes_view -da_refine 1 -ksp_type fgmres
663: requires: !single
665: test:
666: suffix: 11
667: nsize: 4
668: requires: pastix
669: args: -snes_monitor_short -pc_type redundant -dm_mat_type mpiaij -redundant_pc_factor_mat_solver_type pastix -pc_redundant_number 2 -da_refine 4 -ksp_type fgmres
671: test:
672: suffix: 12
673: nsize: 12
674: requires: pastix
675: args: -snes_monitor_short -pc_type redundant -dm_mat_type mpiaij -redundant_pc_factor_mat_solver_type pastix -pc_redundant_number 5 -da_refine 4 -ksp_type fgmres
677: test:
678: suffix: 13
679: nsize: 3
680: args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -pc_fieldsplit_type multiplicative -snes_view -da_refine 1 -ksp_type fgmres -snes_mf_operator
681: requires: !single
683: test:
684: suffix: 14
685: nsize: 4
686: args: -snes_monitor_short -pc_type mg -dm_mat_type baij -mg_coarse_pc_type bjacobi -da_refine 3 -ksp_type fgmres
687: requires: !single
689: test:
690: suffix: 14_ds
691: nsize: 4
692: args: -snes_converged_reason -pc_type mg -dm_mat_type baij -mg_coarse_pc_type bjacobi -da_refine 3 -ksp_type fgmres -mat_fd_type ds
693: output_file: output/ex19_2.out
694: requires: !single
696: test:
697: suffix: 17
698: args: -snes_monitor_short -ksp_pc_side right
699: requires: !single
701: test:
702: suffix: 18
703: args: -ksp_monitor_snes_lg -ksp_pc_side right
704: requires: x !single
706: test:
707: suffix: 2
708: nsize: 4
709: args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds
710: requires: !single
712: test:
713: suffix: 2_bcols1
714: nsize: 4
715: args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds -mat_fd_coloring_bcols
716: output_file: output/ex19_2.out
717: requires: !single
719: test:
720: suffix: 3
721: nsize: 4
722: requires: mumps
723: args: -da_refine 3 -snes_monitor_short -pc_type redundant -dm_mat_type mpiaij -redundant_ksp_type preonly -redundant_pc_factor_mat_solver_type mumps -pc_redundant_number 2
725: test:
726: suffix: 4
727: nsize: 12
728: requires: mumps
729: args: -da_refine 3 -snes_monitor_short -pc_type redundant -dm_mat_type mpiaij -redundant_ksp_type preonly -redundant_pc_factor_mat_solver_type mumps -pc_redundant_number 5
730: output_file: output/ex19_3.out
732: test:
733: suffix: 6
734: args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -snes_view -ksp_type fgmres -da_refine 1
735: requires: !single
737: test:
738: suffix: 7
739: nsize: 3
740: args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -snes_view -da_refine 1 -ksp_type fgmres
742: requires: !single
743: test:
744: suffix: 8
745: args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -pc_fieldsplit_block_size 2 -pc_fieldsplit_0_fields 0,1 -pc_fieldsplit_1_fields 0,1 -pc_fieldsplit_type multiplicative -snes_view -fieldsplit_pc_type lu -da_refine 1 -ksp_type fgmres
746: requires: !single
748: test:
749: suffix: 9
750: nsize: 3
751: args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -pc_fieldsplit_type multiplicative -snes_view -da_refine 1 -ksp_type fgmres
752: requires: !single
754: test:
755: suffix: aspin
756: nsize: 4
757: args: -da_refine 3 -da_overlap 2 -snes_monitor_short -snes_type aspin -grashof 4e4 -lidvelocity 100 -ksp_monitor_short
758: requires: !single
760: test:
761: suffix: bcgsl
762: nsize: 2
763: args: -ksp_type bcgsl -ksp_monitor_short -da_refine 2 -ksp_bcgsl_ell 3 -snes_view
764: requires: !single
766: test:
767: suffix: bcols1
768: nsize: 2
769: args: -da_refine 3 -snes_monitor_short -pc_type mg -ksp_type fgmres -pc_mg_type full -mat_fd_coloring_bcols 1
770: output_file: output/ex19_1.out
771: requires: !single
773: test:
774: suffix: bjacobi
775: nsize: 4
776: args: -da_refine 4 -ksp_type fgmres -pc_type bjacobi -pc_bjacobi_blocks 2 -sub_ksp_type gmres -sub_ksp_max_it 2 -sub_pc_type bjacobi -sub_sub_ksp_type preonly -sub_sub_pc_type ilu -snes_monitor_short
777: requires: !single
779: test:
780: suffix: cgne
781: args: -da_refine 2 -pc_type lu -ksp_type cgne -ksp_monitor_short -ksp_converged_reason -ksp_view -ksp_norm_type unpreconditioned
782: filter: grep -v HERMITIAN
783: requires: !single
785: test:
786: suffix: cgs
787: args: -da_refine 1 -ksp_monitor_short -ksp_type cgs
788: requires: !single
790: test:
791: suffix: composite_fieldsplit
792: args: -ksp_type fgmres -pc_type composite -pc_composite_type MULTIPLICATIVE -pc_composite_pcs fieldsplit,none -sub_0_pc_fieldsplit_block_size 4 -sub_0_pc_fieldsplit_type additive -sub_0_pc_fieldsplit_0_fields 0,1,2 -sub_0_pc_fieldsplit_1_fields 3 -snes_monitor_short -ksp_monitor_short
793: requires: !single
795: test:
796: suffix: composite_fieldsplit_bjacobi
797: args: -ksp_type fgmres -pc_type composite -pc_composite_type MULTIPLICATIVE -pc_composite_pcs fieldsplit,bjacobi -sub_0_pc_fieldsplit_block_size 4 -sub_0_pc_fieldsplit_type additive -sub_0_pc_fieldsplit_0_fields 0,1,2 -sub_0_pc_fieldsplit_1_fields 3 -sub_1_pc_bjacobi_blocks 16 -sub_1_sub_pc_type lu -snes_monitor_short -ksp_monitor_short
798: requires: !single
800: test:
801: suffix: composite_fieldsplit_bjacobi_2
802: nsize: 4
803: args: -ksp_type fgmres -pc_type composite -pc_composite_type MULTIPLICATIVE -pc_composite_pcs fieldsplit,bjacobi -sub_0_pc_fieldsplit_block_size 4 -sub_0_pc_fieldsplit_type additive -sub_0_pc_fieldsplit_0_fields 0,1,2 -sub_0_pc_fieldsplit_1_fields 3 -sub_1_pc_bjacobi_blocks 16 -sub_1_sub_pc_type lu -snes_monitor_short -ksp_monitor_short
804: requires: !single
806: test:
807: suffix: composite_gs_newton
808: nsize: 2
809: args: -da_refine 3 -grashof 4e4 -lidvelocity 100 -snes_monitor_short -snes_type composite -snes_composite_type additiveoptimal -snes_composite_sneses ngs,newtonls -sub_0_snes_max_it 20 -sub_1_pc_type mg
810: requires: !single
812: test:
813: suffix: cuda
814: requires: cuda !single
815: args: -dm_vec_type cuda -dm_mat_type aijcusparse -pc_type none -ksp_type fgmres -snes_monitor_short -snes_rtol 1.e-5
817: test:
818: suffix: draw
819: args: -pc_type fieldsplit -snes_view draw -fieldsplit_x_velocity_pc_type mg -fieldsplit_x_velocity_pc_mg_galerkin pmat -fieldsplit_x_velocity_pc_mg_levels 2 -da_refine 1 -fieldsplit_x_velocity_mg_coarse_pc_type svd
820: requires: x !single
822: test:
823: suffix: drawports
824: args: -snes_monitor_solution draw::draw_ports -da_refine 1
825: output_file: output/ex19_draw.out
826: requires: x !single
828: test:
829: suffix: fas
830: args: -da_refine 4 -snes_monitor_short -snes_type fas -fas_levels_snes_type ngs -fas_levels_snes_ngs_sweeps 3 -fas_levels_snes_ngs_atol 0.0 -fas_levels_snes_ngs_stol 0.0 -grashof 4e4 -snes_fas_smoothup 6 -snes_fas_smoothdown 6 -lidvelocity 100
831: requires: !single
833: test:
834: suffix: fas_full
835: args: -da_refine 4 -snes_monitor_short -snes_type fas -snes_fas_type full -snes_fas_full_downsweep -fas_levels_snes_type ngs -fas_levels_snes_ngs_sweeps 3 -fas_levels_snes_ngs_atol 0.0 -fas_levels_snes_ngs_stol 0.0 -grashof 4e4 -snes_fas_smoothup 6 -snes_fas_smoothdown 6 -lidvelocity 100
836: requires: !single
838: test:
839: suffix: fdcoloring_ds
840: args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds
841: output_file: output/ex19_2.out
842: requires: !single
844: test:
845: suffix: fdcoloring_ds_baij
846: args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds -dm_mat_type baij
847: output_file: output/ex19_2.out
848: requires: !single
850: test:
851: suffix: fdcoloring_ds_bcols1
852: args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds -mat_fd_coloring_bcols 1
853: output_file: output/ex19_2.out
854: requires: !single
856: test:
857: suffix: fdcoloring_wp
858: args: -da_refine 3 -snes_monitor_short -pc_type mg
859: requires: !single
861: test:
862: suffix: fdcoloring_wp_baij
863: args: -da_refine 3 -snes_monitor_short -pc_type mg -dm_mat_type baij
864: output_file: output/ex19_fdcoloring_wp.out
865: requires: !single
867: test:
868: suffix: fdcoloring_wp_bcols1
869: args: -da_refine 3 -snes_monitor_short -pc_type mg -mat_fd_coloring_bcols 1
870: output_file: output/ex19_fdcoloring_wp.out
871: requires: !single
873: test:
874: suffix: fieldsplit_2
875: args: -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type additive -pc_fieldsplit_0_fields 0,1,2 -pc_fieldsplit_1_fields 3 -snes_monitor_short -ksp_monitor_short
876: requires: !single
878: test:
879: suffix: fieldsplit_3
880: args: -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type additive -pc_fieldsplit_0_fields 0,1,2 -pc_fieldsplit_1_fields 3 -fieldsplit_0_pc_type lu -fieldsplit_1_pc_type lu -snes_monitor_short -ksp_monitor_short
881: requires: !single
883: test:
884: suffix: fieldsplit_4
885: args: -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type SCHUR -pc_fieldsplit_0_fields 0,1,2 -pc_fieldsplit_1_fields 3 -fieldsplit_0_pc_type lu -fieldsplit_1_pc_type lu -snes_monitor_short -ksp_monitor_short
886: requires: !single
888: # HYPRE PtAP broken with complex numbers
889: test:
890: suffix: fieldsplit_hypre
891: nsize: 2
892: requires: hypre mumps !complex
893: args: -pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type SCHUR -pc_fieldsplit_0_fields 0,1,2 -pc_fieldsplit_1_fields 3 -fieldsplit_0_pc_type lu -fieldsplit_0_pc_factor_mat_solver_type mumps -fieldsplit_1_pc_type hypre -fieldsplit_1_pc_hypre_type boomeramg -snes_monitor_short -ksp_monitor_short
895: test:
896: suffix: fieldsplit_mumps
897: nsize: 2
898: requires: mumps
899: args: -pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type SCHUR -pc_fieldsplit_0_fields 0,1,2 -pc_fieldsplit_1_fields 3 -fieldsplit_0_pc_type lu -fieldsplit_1_pc_type lu -snes_monitor_short -ksp_monitor_short -fieldsplit_0_pc_factor_mat_solver_type mumps -fieldsplit_1_pc_factor_mat_solver_type mumps
900: output_file: output/ex19_fieldsplit_5.out
902: test:
903: suffix: greedy_coloring
904: nsize: 2
905: args: -da_refine 3 -snes_monitor_short -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -mat_coloring_weight_type lf -mat_coloring_view> ex19_greedy_coloring.tmp 2>&1
906: requires: !single
908: # HYPRE PtAP broken with complex numbers
909: test:
910: suffix: hypre
911: nsize: 2
912: requires: hypre !complex
913: args: -da_refine 3 -snes_monitor_short -pc_type hypre
915: test:
916: suffix: ibcgs
917: nsize: 2
918: args: -ksp_type ibcgs -ksp_monitor_short -da_refine 2 -snes_view
919: requires: !complex !single
921: test:
922: suffix: kaczmarz
923: nsize: 2
924: args: -pc_type kaczmarz -ksp_monitor_short -snes_monitor_short -snes_view
925: requires: !single
927: test:
928: suffix: klu
929: requires: suitesparse
930: args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type klu
931: output_file: output/ex19_superlu.out
933: test:
934: suffix: klu_2
935: requires: suitesparse
936: args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type klu -mat_klu_ordering PETSC
937: output_file: output/ex19_superlu.out
939: test:
940: suffix: klu_3
941: requires: suitesparse
942: args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type klu -mat_klu_use_btf 0
943: output_file: output/ex19_superlu.out
945: test:
946: suffix: ml
947: nsize: 2
948: requires: ml
949: args: -da_refine 3 -snes_monitor_short -pc_type ml
951: test:
952: suffix: ngmres_fas
953: args: -da_refine 4 -snes_monitor_short -snes_type ngmres -npc_fas_levels_snes_type ngs -npc_fas_levels_snes_ngs_sweeps 3 -npc_fas_levels_snes_ngs_atol 0.0 -npc_fas_levels_snes_ngs_stol 0.0 -npc_snes_type fas -npc_fas_levels_snes_type ngs -npc_snes_max_it 1 -npc_snes_fas_smoothup 6 -npc_snes_fas_smoothdown 6 -lidvelocity 100 -grashof 4e4
954: requires: !single
956: test:
957: suffix: ngmres_fas_gssecant
958: args: -da_refine 3 -snes_monitor_short -snes_type ngmres -npc_snes_type fas -npc_fas_levels_snes_type ngs -npc_fas_levels_snes_max_it 6 -npc_fas_levels_snes_ngs_secant -npc_fas_levels_snes_ngs_max_it 1 -npc_fas_coarse_snes_max_it 1 -lidvelocity 100 -grashof 4e4
959: requires: !single
961: test:
962: suffix: ngmres_fas_ms
963: nsize: 2
964: args: -snes_grid_sequence 2 -lidvelocity 200 -grashof 1e4 -snes_monitor_short -snes_view -snes_converged_reason -snes_type ngmres -npc_snes_type fas -npc_fas_coarse_snes_type newtonls -npc_fas_coarse_ksp_type preonly -npc_snes_max_it 1
965: requires: !single
967: test:
968: suffix: ngmres_nasm
969: nsize: 4
970: args: -da_refine 4 -da_overlap 2 -snes_monitor_short -snes_type ngmres -snes_max_it 10 -npc_snes_type nasm -npc_snes_nasm_type basic -grashof 4e4 -lidvelocity 100
971: requires: !single
973: test:
974: suffix: ngs
975: args: -snes_type ngs -snes_view -snes_monitor -snes_rtol 1e-4
976: requires: !single
978: test:
979: suffix: ngs_fd
980: args: -snes_type ngs -snes_ngs_secant -snes_view -snes_monitor -snes_rtol 1e-4
981: requires: !single
983: test:
984: suffix: parms
985: nsize: 2
986: requires: parms
987: args: -pc_type parms -ksp_monitor_short -snes_view
989: test:
990: suffix: superlu
991: requires: superlu
992: args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu
994: test:
995: suffix: superlu_sell
996: requires: superlu
997: args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu -dm_mat_type sell -pc_factor_mat_ordering_type natural
998: output_file: output/ex19_superlu.out
1000: test:
1001: suffix: superlu_dist
1002: requires: superlu_dist
1003: args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu_dist
1004: output_file: output/ex19_superlu.out
1006: test:
1007: suffix: superlu_dist_2
1008: nsize: 2
1009: requires: superlu_dist
1010: args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu_dist
1011: output_file: output/ex19_superlu.out
1013: test:
1014: suffix: superlu_equil
1015: requires: superlu
1016: args: -da_grid_x 20 -da_grid_y 20 -{snes,ksp}_monitor_short -pc_type lu -pc_factor_mat_solver_type superlu -mat_superlu_equil
1018: test:
1019: suffix: superlu_equil_sell
1020: requires: superlu
1021: args: -da_grid_x 20 -da_grid_y 20 -{snes,ksp}_monitor_short -pc_type lu -pc_factor_mat_solver_type superlu -mat_superlu_equil -dm_mat_type sell -pc_factor_mat_ordering_type natural
1022: output_file: output/ex19_superlu_equil.out
1024: test:
1025: suffix: tcqmr
1026: args: -da_refine 1 -ksp_monitor_short -ksp_type tcqmr
1027: requires: !single
1029: test:
1030: suffix: tfqmr
1031: args: -da_refine 1 -ksp_monitor_short -ksp_type tfqmr
1032: requires: !single
1034: test:
1035: suffix: umfpack
1036: requires: suitesparse
1037: args: -da_refine 2 -pc_type lu -pc_factor_mat_solver_type umfpack -snes_view -snes_monitor_short -ksp_monitor_short
1039: test:
1040: suffix: tut_1
1041: nsize: 4
1042: requires: !single
1043: args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view
1045: test:
1046: suffix: tut_2
1047: nsize: 4
1048: requires: !single
1049: args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view -pc_type mg
1051: # HYPRE PtAP broken with complex numbers
1052: test:
1053: suffix: tut_3
1054: nsize: 4
1055: requires: hypre !single !complex
1056: args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view -pc_type hypre
1058: test:
1059: suffix: tut_8
1060: nsize: 4
1061: requires: ml !single
1062: args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view -pc_type ml
1064: test:
1065: suffix: tut_4
1066: nsize: 1
1067: requires: !single
1068: args: -da_refine 5 -log_view
1069: filter: head -n 2
1070: filter_output: head -n 2
1072: test:
1073: suffix: tut_5
1074: nsize: 1
1075: requires: !single
1076: args: -da_refine 5 -log_view -pc_type mg
1077: filter: head -n 2
1078: filter_output: head -n 2
1080: test:
1081: suffix: tut_6
1082: nsize: 4
1083: requires: !single
1084: args: -da_refine 5 -log_view
1085: filter: head -n 2
1086: filter_output: head -n 2
1088: test:
1089: suffix: tut_7
1090: nsize: 4
1091: requires: !single
1092: args: -da_refine 5 -log_view -pc_type mg
1093: filter: head -n 2
1094: filter_output: head -n 2
1096: test:
1097: suffix: cuda_1
1098: nsize: 1
1099: requires: cuda
1100: args: -snes_monitor -dm_mat_type seqaijcusparse -dm_vec_type seqcuda -pc_type gamg -ksp_monitor -mg_levels_ksp_max_it 3
1103: test:
1104: suffix: cuda_2
1105: nsize: 3
1106: requires: cuda !single
1107: args: -snes_monitor -dm_mat_type mpiaijcusparse -dm_vec_type mpicuda -pc_type gamg -ksp_monitor -mg_levels_ksp_max_it 3
1109: test:
1110: suffix: seqbaijmkl
1111: nsize: 1
1112: requires: define(PETSC_HAVE_MKL_SPARSE_OPTIMIZE)
1113: args: -dm_mat_type baij -snes_monitor -ksp_monitor -snes_view
1115: test:
1116: suffix: mpibaijmkl
1117: nsize: 2
1118: requires: define(PETSC_HAVE_MKL_SPARSE_OPTIMIZE)
1119: args: -dm_mat_type baij -snes_monitor -ksp_monitor -snes_view
1121: test:
1122: suffix: cpardiso
1123: nsize: 4
1124: requires: mkl_cpardiso
1125: args: -pc_type lu -pc_factor_mat_solver_type mkl_cpardiso -ksp_monitor
1127: test:
1128: suffix: logviewmemory
1129: requires: define(PETSC_USE_LOG) !define(PETSC_HAVE_VALGRIND)
1130: args: -log_view -log_view_memory -da_refine 4
1131: filter: grep MatFDColorSetUp | wc -w | xargs -I % sh -c "expr % \> 21"
1133: TEST*/