Actual source code: ex28.c
petsc-3.14.6 2021-03-30
1: static const char help[] = "1D multiphysics prototype with analytic Jacobians to solve individual problems and a coupled problem.\n\n";
3: /* Solve a PDE coupled to an algebraic system in 1D
4: *
5: * PDE (U):
6: * -(k u_x)_x = 1 on (0,1), subject to u(0) = 0, u(1) = 1
7: * Algebraic (K):
8: * exp(k-1) + k = 1/(1/(1+u) + 1/(1+u_x^2))
9: *
10: * The discretization places k at staggered points, and a separate DMDA is used for each "physics".
11: *
12: * This example is a prototype for coupling in multi-physics problems, therefore residual evaluation and assembly for
13: * each problem (referred to as U and K) are written separately. This permits the same "physics" code to be used for
14: * solving each uncoupled problem as well as the coupled system. In particular, run with -problem_type 0 to solve only
15: * problem U (with K fixed), -problem_type 1 to solve only K (with U fixed), and -problem_type 2 to solve both at once.
16: *
17: * In all cases, a fully-assembled analytic Jacobian is available, so the systems can be solved with a direct solve or
18: * any other standard method. Additionally, by running with
19: *
20: * -pack_dm_mat_type nest
21: *
22: * The same code assembles a coupled matrix where each block is stored separately, which allows the use of PCFieldSplit
23: * without copying values to extract submatrices.
24: */
26: #include <petscsnes.h>
27: #include <petscdm.h>
28: #include <petscdmda.h>
29: #include <petscdmcomposite.h>
31: typedef struct _UserCtx *User;
32: struct _UserCtx {
33: PetscInt ptype;
34: DM pack;
35: Vec Uloc,Kloc;
36: };
38: static PetscErrorCode FormFunctionLocal_U(User user,DMDALocalInfo *info,const PetscScalar u[],const PetscScalar k[],PetscScalar f[])
39: {
40: PetscReal hx = 1./info->mx;
41: PetscInt i;
44: for (i=info->xs; i<info->xs+info->xm; i++) {
45: if (i == 0) f[i] = 1./hx*u[i];
46: else if (i == info->mx-1) f[i] = 1./hx*(u[i] - 1.0);
47: else f[i] = hx*((k[i-1]*(u[i]-u[i-1]) - k[i]*(u[i+1]-u[i]))/(hx*hx) - 1.0);
48: }
49: return(0);
50: }
52: static PetscErrorCode FormFunctionLocal_K(User user,DMDALocalInfo *info,const PetscScalar u[],const PetscScalar k[],PetscScalar f[])
53: {
54: PetscReal hx = 1./info->mx;
55: PetscInt i;
58: for (i=info->xs; i<info->xs+info->xm; i++) {
59: const PetscScalar
60: ubar = 0.5*(u[i+1]+u[i]),
61: gradu = (u[i+1]-u[i])/hx,
62: g = 1. + gradu*gradu,
63: w = 1./(1.+ubar) + 1./g;
64: f[i] = hx*(PetscExpScalar(k[i]-1.0) + k[i] - 1./w);
65: }
66: return(0);
67: }
69: static PetscErrorCode FormFunction_All(SNES snes,Vec X,Vec F,void *ctx)
70: {
71: User user = (User)ctx;
72: DM dau,dak;
73: DMDALocalInfo infou,infok;
74: PetscScalar *u,*k;
75: PetscScalar *fu,*fk;
77: Vec Uloc,Kloc,Fu,Fk;
80: DMCompositeGetEntries(user->pack,&dau,&dak);
81: DMDAGetLocalInfo(dau,&infou);
82: DMDAGetLocalInfo(dak,&infok);
83: DMCompositeGetLocalVectors(user->pack,&Uloc,&Kloc);
84: switch (user->ptype) {
85: case 0:
86: DMGlobalToLocalBegin(dau,X,INSERT_VALUES,Uloc);
87: DMGlobalToLocalEnd (dau,X,INSERT_VALUES,Uloc);
88: DMDAVecGetArray(dau,Uloc,&u);
89: DMDAVecGetArray(dak,user->Kloc,&k);
90: DMDAVecGetArray(dau,F,&fu);
91: FormFunctionLocal_U(user,&infou,u,k,fu);
92: DMDAVecRestoreArray(dau,F,&fu);
93: DMDAVecRestoreArray(dau,Uloc,&u);
94: DMDAVecRestoreArray(dak,user->Kloc,&k);
95: break;
96: case 1:
97: DMGlobalToLocalBegin(dak,X,INSERT_VALUES,Kloc);
98: DMGlobalToLocalEnd (dak,X,INSERT_VALUES,Kloc);
99: DMDAVecGetArray(dau,user->Uloc,&u);
100: DMDAVecGetArray(dak,Kloc,&k);
101: DMDAVecGetArray(dak,F,&fk);
102: FormFunctionLocal_K(user,&infok,u,k,fk);
103: DMDAVecRestoreArray(dak,F,&fk);
104: DMDAVecRestoreArray(dau,user->Uloc,&u);
105: DMDAVecRestoreArray(dak,Kloc,&k);
106: break;
107: case 2:
108: DMCompositeScatter(user->pack,X,Uloc,Kloc);
109: DMDAVecGetArray(dau,Uloc,&u);
110: DMDAVecGetArray(dak,Kloc,&k);
111: DMCompositeGetAccess(user->pack,F,&Fu,&Fk);
112: DMDAVecGetArray(dau,Fu,&fu);
113: DMDAVecGetArray(dak,Fk,&fk);
114: FormFunctionLocal_U(user,&infou,u,k,fu);
115: FormFunctionLocal_K(user,&infok,u,k,fk);
116: DMDAVecRestoreArray(dau,Fu,&fu);
117: DMDAVecRestoreArray(dak,Fk,&fk);
118: DMCompositeRestoreAccess(user->pack,F,&Fu,&Fk);
119: DMDAVecRestoreArray(dau,Uloc,&u);
120: DMDAVecRestoreArray(dak,Kloc,&k);
121: break;
122: }
123: DMCompositeRestoreLocalVectors(user->pack,&Uloc,&Kloc);
124: return(0);
125: }
127: static PetscErrorCode FormJacobianLocal_U(User user,DMDALocalInfo *info,const PetscScalar u[],const PetscScalar k[],Mat Buu)
128: {
129: PetscReal hx = 1./info->mx;
131: PetscInt i;
134: for (i=info->xs; i<info->xs+info->xm; i++) {
135: PetscInt row = i-info->gxs,cols[] = {row-1,row,row+1};
136: PetscScalar val = 1./hx;
137: if (i == 0) {
138: MatSetValuesLocal(Buu,1,&row,1,&row,&val,INSERT_VALUES);
139: } else if (i == info->mx-1) {
140: MatSetValuesLocal(Buu,1,&row,1,&row,&val,INSERT_VALUES);
141: } else {
142: PetscScalar vals[] = {-k[i-1]/hx,(k[i-1]+k[i])/hx,-k[i]/hx};
143: MatSetValuesLocal(Buu,1,&row,3,cols,vals,INSERT_VALUES);
144: }
145: }
146: return(0);
147: }
149: static PetscErrorCode FormJacobianLocal_K(User user,DMDALocalInfo *info,const PetscScalar u[],const PetscScalar k[],Mat Bkk)
150: {
151: PetscReal hx = 1./info->mx;
153: PetscInt i;
156: for (i=info->xs; i<info->xs+info->xm; i++) {
157: PetscInt row = i-info->gxs;
158: PetscScalar vals[] = {hx*(PetscExpScalar(k[i]-1.)+ (PetscScalar)1.)};
159: MatSetValuesLocal(Bkk,1,&row,1,&row,vals,INSERT_VALUES);
160: }
161: return(0);
162: }
164: static PetscErrorCode FormJacobianLocal_UK(User user,DMDALocalInfo *info,DMDALocalInfo *infok,const PetscScalar u[],const PetscScalar k[],Mat Buk)
165: {
166: PetscReal hx = 1./info->mx;
168: PetscInt i;
169: PetscInt row,cols[2];
170: PetscScalar vals[2];
173: if (!Buk) return(0); /* Not assembling this block */
174: for (i=info->xs; i<info->xs+info->xm; i++) {
175: if (i == 0 || i == info->mx-1) continue;
176: row = i-info->gxs;
177: cols[0] = i-1-infok->gxs; vals[0] = (u[i]-u[i-1])/hx;
178: cols[1] = i-infok->gxs; vals[1] = (u[i]-u[i+1])/hx;
179: MatSetValuesLocal(Buk,1,&row,2,cols,vals,INSERT_VALUES);
180: }
181: return(0);
182: }
184: static PetscErrorCode FormJacobianLocal_KU(User user,DMDALocalInfo *info,DMDALocalInfo *infok,const PetscScalar u[],const PetscScalar k[],Mat Bku)
185: {
187: PetscInt i;
188: PetscReal hx = 1./(info->mx-1);
191: if (!Bku) return(0); /* Not assembling this block */
192: for (i=infok->xs; i<infok->xs+infok->xm; i++) {
193: PetscInt row = i-infok->gxs,cols[2];
194: PetscScalar vals[2];
195: const PetscScalar
196: ubar = 0.5*(u[i]+u[i+1]),
197: ubar_L = 0.5,
198: ubar_R = 0.5,
199: gradu = (u[i+1]-u[i])/hx,
200: gradu_L = -1./hx,
201: gradu_R = 1./hx,
202: g = 1. + PetscSqr(gradu),
203: g_gradu = 2.*gradu,
204: w = 1./(1.+ubar) + 1./g,
205: w_ubar = -1./PetscSqr(1.+ubar),
206: w_gradu = -g_gradu/PetscSqr(g),
207: iw = 1./w,
208: iw_ubar = -w_ubar * PetscSqr(iw),
209: iw_gradu = -w_gradu * PetscSqr(iw);
210: cols[0] = i-info->gxs; vals[0] = -hx*(iw_ubar*ubar_L + iw_gradu*gradu_L);
211: cols[1] = i+1-info->gxs; vals[1] = -hx*(iw_ubar*ubar_R + iw_gradu*gradu_R);
212: MatSetValuesLocal(Bku,1,&row,2,cols,vals,INSERT_VALUES);
213: }
214: return(0);
215: }
217: static PetscErrorCode FormJacobian_All(SNES snes,Vec X,Mat J,Mat B,void *ctx)
218: {
219: User user = (User)ctx;
220: DM dau,dak;
221: DMDALocalInfo infou,infok;
222: PetscScalar *u,*k;
224: Vec Uloc,Kloc;
227: DMCompositeGetEntries(user->pack,&dau,&dak);
228: DMDAGetLocalInfo(dau,&infou);
229: DMDAGetLocalInfo(dak,&infok);
230: DMCompositeGetLocalVectors(user->pack,&Uloc,&Kloc);
231: switch (user->ptype) {
232: case 0:
233: DMGlobalToLocalBegin(dau,X,INSERT_VALUES,Uloc);
234: DMGlobalToLocalEnd (dau,X,INSERT_VALUES,Uloc);
235: DMDAVecGetArray(dau,Uloc,&u);
236: DMDAVecGetArray(dak,user->Kloc,&k);
237: FormJacobianLocal_U(user,&infou,u,k,B);
238: DMDAVecRestoreArray(dau,Uloc,&u);
239: DMDAVecRestoreArray(dak,user->Kloc,&k);
240: break;
241: case 1:
242: DMGlobalToLocalBegin(dak,X,INSERT_VALUES,Kloc);
243: DMGlobalToLocalEnd (dak,X,INSERT_VALUES,Kloc);
244: DMDAVecGetArray(dau,user->Uloc,&u);
245: DMDAVecGetArray(dak,Kloc,&k);
246: FormJacobianLocal_K(user,&infok,u,k,B);
247: DMDAVecRestoreArray(dau,user->Uloc,&u);
248: DMDAVecRestoreArray(dak,Kloc,&k);
249: break;
250: case 2: {
251: Mat Buu,Buk,Bku,Bkk;
252: PetscBool nest;
253: IS *is;
254: DMCompositeScatter(user->pack,X,Uloc,Kloc);
255: DMDAVecGetArray(dau,Uloc,&u);
256: DMDAVecGetArray(dak,Kloc,&k);
257: DMCompositeGetLocalISs(user->pack,&is);
258: MatGetLocalSubMatrix(B,is[0],is[0],&Buu);
259: MatGetLocalSubMatrix(B,is[0],is[1],&Buk);
260: MatGetLocalSubMatrix(B,is[1],is[0],&Bku);
261: MatGetLocalSubMatrix(B,is[1],is[1],&Bkk);
262: FormJacobianLocal_U(user,&infou,u,k,Buu);
263: PetscObjectTypeCompare((PetscObject)B,MATNEST,&nest);
264: if (!nest) {
265: /*
266: DMCreateMatrix_Composite() for a nested matrix does not generate off-block matrices that one can call MatSetValuesLocal() on, it just creates dummy
267: matrices with no entries; there cannot insert values into them. Commit b6480e041dd2293a65f96222772d68cdb4ed6306
268: changed Mat_Nest() from returning NULL pointers for these submatrices to dummy matrices because PCFIELDSPLIT could not
269: handle the returned null matrices.
270: */
271: FormJacobianLocal_UK(user,&infou,&infok,u,k,Buk);
272: FormJacobianLocal_KU(user,&infou,&infok,u,k,Bku);
273: }
274: FormJacobianLocal_K(user,&infok,u,k,Bkk);
275: MatRestoreLocalSubMatrix(B,is[0],is[0],&Buu);
276: MatRestoreLocalSubMatrix(B,is[0],is[1],&Buk);
277: MatRestoreLocalSubMatrix(B,is[1],is[0],&Bku);
278: MatRestoreLocalSubMatrix(B,is[1],is[1],&Bkk);
279: DMDAVecRestoreArray(dau,Uloc,&u);
280: DMDAVecRestoreArray(dak,Kloc,&k);
282: ISDestroy(&is[0]);
283: ISDestroy(&is[1]);
284: PetscFree(is);
285: } break;
286: }
287: DMCompositeRestoreLocalVectors(user->pack,&Uloc,&Kloc);
288: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
289: MatAssemblyEnd (B,MAT_FINAL_ASSEMBLY);
290: if (J != B) {
291: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
292: MatAssemblyEnd (J,MAT_FINAL_ASSEMBLY);
293: }
294: return(0);
295: }
297: static PetscErrorCode FormInitial_Coupled(User user,Vec X)
298: {
300: DM dau,dak;
301: DMDALocalInfo infou,infok;
302: Vec Xu,Xk;
303: PetscScalar *u,*k,hx;
304: PetscInt i;
307: DMCompositeGetEntries(user->pack,&dau,&dak);
308: DMCompositeGetAccess(user->pack,X,&Xu,&Xk);
309: DMDAVecGetArray(dau,Xu,&u);
310: DMDAVecGetArray(dak,Xk,&k);
311: DMDAGetLocalInfo(dau,&infou);
312: DMDAGetLocalInfo(dak,&infok);
313: hx = 1./(infok.mx);
314: for (i=infou.xs; i<infou.xs+infou.xm; i++) u[i] = (PetscScalar)i*hx * (1.-(PetscScalar)i*hx);
315: for (i=infok.xs; i<infok.xs+infok.xm; i++) k[i] = 1.0 + 0.5*PetscSinScalar(2*PETSC_PI*i*hx);
316: DMDAVecRestoreArray(dau,Xu,&u);
317: DMDAVecRestoreArray(dak,Xk,&k);
318: DMCompositeRestoreAccess(user->pack,X,&Xu,&Xk);
319: DMCompositeScatter(user->pack,X,user->Uloc,user->Kloc);
320: return(0);
321: }
323: int main(int argc, char *argv[])
324: {
326: DM dau,dak,pack;
327: const PetscInt *lxu;
328: PetscInt *lxk,m,sizes;
329: User user;
330: SNES snes;
331: Vec X,F,Xu,Xk,Fu,Fk;
332: Mat B;
333: IS *isg;
334: PetscBool pass_dm;
336: PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
337: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,10,1,1,NULL,&dau);
338: DMSetOptionsPrefix(dau,"u_");
339: DMSetFromOptions(dau);
340: DMSetUp(dau);
341: DMDAGetOwnershipRanges(dau,&lxu,0,0);
342: DMDAGetInfo(dau,0, &m,0,0, &sizes,0,0, 0,0,0,0,0,0);
343: PetscMalloc1(sizes,&lxk);
344: PetscArraycpy(lxk,lxu,sizes);
345: lxk[0]--;
346: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,m-1,1,1,lxk,&dak);
347: DMSetOptionsPrefix(dak,"k_");
348: DMSetFromOptions(dak);
349: DMSetUp(dak);
350: PetscFree(lxk);
352: DMCompositeCreate(PETSC_COMM_WORLD,&pack);
353: DMSetOptionsPrefix(pack,"pack_");
354: DMCompositeAddDM(pack,dau);
355: DMCompositeAddDM(pack,dak);
356: DMDASetFieldName(dau,0,"u");
357: DMDASetFieldName(dak,0,"k");
358: DMSetFromOptions(pack);
360: DMCreateGlobalVector(pack,&X);
361: VecDuplicate(X,&F);
363: PetscNew(&user);
365: user->pack = pack;
367: DMCompositeGetGlobalISs(pack,&isg);
368: DMCompositeGetLocalVectors(pack,&user->Uloc,&user->Kloc);
369: DMCompositeScatter(pack,X,user->Uloc,user->Kloc);
371: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Coupled problem options","SNES");
372: {
373: user->ptype = 0; pass_dm = PETSC_TRUE;
375: PetscOptionsInt("-problem_type","0: solve for u only, 1: solve for k only, 2: solve for both",0,user->ptype,&user->ptype,NULL);
376: PetscOptionsBool("-pass_dm","Pass the packed DM to SNES to use when determining splits and forward into splits",0,pass_dm,&pass_dm,NULL);
377: }
378: PetscOptionsEnd();
380: FormInitial_Coupled(user,X);
382: SNESCreate(PETSC_COMM_WORLD,&snes);
383: switch (user->ptype) {
384: case 0:
385: DMCompositeGetAccess(pack,X,&Xu,0);
386: DMCompositeGetAccess(pack,F,&Fu,0);
387: DMCreateMatrix(dau,&B);
388: SNESSetFunction(snes,Fu,FormFunction_All,user);
389: SNESSetJacobian(snes,B,B,FormJacobian_All,user);
390: SNESSetFromOptions(snes);
391: SNESSetDM(snes,dau);
392: SNESSolve(snes,NULL,Xu);
393: DMCompositeRestoreAccess(pack,X,&Xu,0);
394: DMCompositeRestoreAccess(pack,F,&Fu,0);
395: break;
396: case 1:
397: DMCompositeGetAccess(pack,X,0,&Xk);
398: DMCompositeGetAccess(pack,F,0,&Fk);
399: DMCreateMatrix(dak,&B);
400: SNESSetFunction(snes,Fk,FormFunction_All,user);
401: SNESSetJacobian(snes,B,B,FormJacobian_All,user);
402: SNESSetFromOptions(snes);
403: SNESSetDM(snes,dak);
404: SNESSolve(snes,NULL,Xk);
405: DMCompositeRestoreAccess(pack,X,0,&Xk);
406: DMCompositeRestoreAccess(pack,F,0,&Fk);
407: break;
408: case 2:
409: DMCreateMatrix(pack,&B);
410: /* This example does not correctly allocate off-diagonal blocks. These options allows new nonzeros (slow). */
411: MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
412: MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
413: SNESSetFunction(snes,F,FormFunction_All,user);
414: SNESSetJacobian(snes,B,B,FormJacobian_All,user);
415: SNESSetFromOptions(snes);
416: if (!pass_dm) { /* Manually provide index sets and names for the splits */
417: KSP ksp;
418: PC pc;
419: SNESGetKSP(snes,&ksp);
420: KSPGetPC(ksp,&pc);
421: PCFieldSplitSetIS(pc,"u",isg[0]);
422: PCFieldSplitSetIS(pc,"k",isg[1]);
423: } else {
424: /* The same names come from the options prefix for dau and dak. This option can support geometric multigrid inside
425: * of splits, but it requires using a DM (perhaps your own implementation). */
426: SNESSetDM(snes,pack);
427: }
428: SNESSolve(snes,NULL,X);
429: break;
430: }
431: VecViewFromOptions(X,NULL,"-view_sol");
433: if (0) {
434: PetscInt col = 0;
435: PetscBool mult_dup = PETSC_FALSE,view_dup = PETSC_FALSE;
436: Mat D;
437: Vec Y;
439: PetscOptionsGetInt(NULL,0,"-col",&col,0);
440: PetscOptionsGetBool(NULL,0,"-mult_dup",&mult_dup,0);
441: PetscOptionsGetBool(NULL,0,"-view_dup",&view_dup,0);
443: VecDuplicate(X,&Y);
444: /* MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); */
445: /* MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); */
446: MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&D);
447: VecZeroEntries(X);
448: VecSetValue(X,col,1.0,INSERT_VALUES);
449: VecAssemblyBegin(X);
450: VecAssemblyEnd(X);
451: MatMult(mult_dup ? D : B,X,Y);
452: MatView(view_dup ? D : B,PETSC_VIEWER_STDOUT_WORLD);
453: /* VecView(X,PETSC_VIEWER_STDOUT_WORLD); */
454: VecView(Y,PETSC_VIEWER_STDOUT_WORLD);
455: MatDestroy(&D);
456: VecDestroy(&Y);
457: }
459: DMCompositeRestoreLocalVectors(pack,&user->Uloc,&user->Kloc);
460: PetscFree(user);
462: ISDestroy(&isg[0]);
463: ISDestroy(&isg[1]);
464: PetscFree(isg);
465: VecDestroy(&X);
466: VecDestroy(&F);
467: MatDestroy(&B);
468: DMDestroy(&dau);
469: DMDestroy(&dak);
470: DMDestroy(&pack);
471: SNESDestroy(&snes);
472: PetscFinalize();
473: return ierr;
474: }
476: /*TEST
478: test:
479: suffix: 0
480: nsize: 3
481: args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -problem_type 0
483: test:
484: suffix: 1
485: nsize: 3
486: args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -problem_type 1
488: test:
489: suffix: 2
490: nsize: 3
491: args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -problem_type 2
493: test:
494: suffix: 3
495: nsize: 3
496: args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -ksp_monitor_short -problem_type 2 -snes_mf_operator -pack_dm_mat_type {{aij nest}} -pc_type fieldsplit -pc_fieldsplit_dm_splits -pc_fieldsplit_type additive -fieldsplit_u_ksp_type gmres -fieldsplit_k_pc_type jacobi
498: test:
499: suffix: 4
500: nsize: 6
501: args: -u_da_grid_x 257 -snes_converged_reason -snes_monitor_short -ksp_monitor_short -problem_type 2 -snes_mf_operator -pack_dm_mat_type aij -pc_type fieldsplit -pc_fieldsplit_type multiplicative -fieldsplit_u_ksp_type gmres -fieldsplit_u_ksp_pc_side right -fieldsplit_u_pc_type mg -fieldsplit_u_pc_mg_levels 4 -fieldsplit_u_mg_levels_ksp_type richardson -fieldsplit_u_mg_levels_ksp_max_it 1 -fieldsplit_u_mg_levels_pc_type sor -fieldsplit_u_pc_mg_galerkin pmat -fieldsplit_u_ksp_converged_reason -fieldsplit_k_pc_type jacobi
502: requires: !single
504: test:
505: suffix: glvis_composite_da_1d
506: args: -u_da_grid_x 20 -problem_type 0 -snes_monitor_solution glvis:
508: test:
509: suffix: cuda
510: nsize: 1
511: requires: cuda
512: args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -problem_type 2 -pack_dm_mat_type aijcusparse -pack_dm_vec_type cuda
514: test:
515: suffix: viennacl
516: nsize: 1
517: requires: viennacl
518: args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -problem_type 2 -pack_dm_mat_type aijviennacl -pack_dm_vec_type viennacl
520: TEST*/