xref: /petsc/src/mat/tests/ex129.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1c4762a1bSJed Brown /*
2c4762a1bSJed Brown   Laplacian in 3D. Use for testing MatSolve routines.
3c4762a1bSJed Brown   Modeled by the partial differential equation
4c4762a1bSJed Brown 
5c4762a1bSJed Brown    - Laplacian u = 1,0 < x,y,z < 1,
6c4762a1bSJed Brown 
7c4762a1bSJed Brown    with boundary conditions
8c4762a1bSJed Brown    u = 1 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1.
9c4762a1bSJed Brown */
10c4762a1bSJed Brown 
11c4762a1bSJed Brown static char help[] = "This example is for testing different MatSolve routines :MatSolve(), MatSolveAdd(), MatSolveTranspose(), MatSolveTransposeAdd(), and MatMatSolve().\n\
12c4762a1bSJed Brown Example usage: ./ex129 -mat_type aij -dof 2\n\n";
13c4762a1bSJed Brown 
14c4762a1bSJed Brown #include <petscdm.h>
15c4762a1bSJed Brown #include <petscdmda.h>
16c4762a1bSJed Brown 
17c4762a1bSJed Brown extern PetscErrorCode ComputeMatrix(DM, Mat);
18c4762a1bSJed Brown extern PetscErrorCode ComputeRHS(DM, Vec);
19c4762a1bSJed Brown extern PetscErrorCode ComputeRHSMatrix(PetscInt, PetscInt, Mat *);
20c4762a1bSJed Brown 
main(int argc,char ** args)21d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
22d71ae5a4SJacob Faibussowitsch {
23c4762a1bSJed Brown   PetscMPIInt   size;
24c4762a1bSJed Brown   Vec           x, b, y, b1;
25c4762a1bSJed Brown   DM            da;
26c4762a1bSJed Brown   Mat           A, F, RHS, X, C1;
27c4762a1bSJed Brown   MatFactorInfo info;
28c4762a1bSJed Brown   IS            perm, iperm;
29c4762a1bSJed Brown   PetscInt      dof = 1, M = 8, m, n, nrhs;
30c4762a1bSJed Brown   PetscScalar   one = 1.0;
31c4762a1bSJed Brown   PetscReal     norm, tol = 1000 * PETSC_MACHINE_EPSILON;
32c4762a1bSJed Brown   PetscBool     InplaceLU = PETSC_FALSE;
33c4762a1bSJed Brown 
34327415f7SBarry Smith   PetscFunctionBeginUser;
35c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
37be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only");
389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dof", &dof, NULL));
399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL));
40c4762a1bSJed Brown 
419566063dSJacob Faibussowitsch   PetscCall(DMDACreate(PETSC_COMM_WORLD, &da));
429566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(da, 3));
439566063dSJacob Faibussowitsch   PetscCall(DMDASetBoundaryType(da, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE));
449566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilType(da, DMDA_STENCIL_STAR));
459566063dSJacob Faibussowitsch   PetscCall(DMDASetSizes(da, M, M, M));
469566063dSJacob Faibussowitsch   PetscCall(DMDASetNumProcs(da, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE));
479566063dSJacob Faibussowitsch   PetscCall(DMDASetDof(da, dof));
489566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilWidth(da, 1));
499566063dSJacob Faibussowitsch   PetscCall(DMDASetOwnershipRanges(da, NULL, NULL, NULL));
509566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(da, MATBAIJ));
519566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
529566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &x));
559566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &b));
569566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(b, &y));
579566063dSJacob Faibussowitsch   PetscCall(ComputeRHS(da, b));
589566063dSJacob Faibussowitsch   PetscCall(VecSet(y, one));
599566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(da, &A));
609566063dSJacob Faibussowitsch   PetscCall(ComputeMatrix(da, A));
619566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &m, &n));
62c4762a1bSJed Brown   nrhs = 2;
639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nrhs", &nrhs, NULL));
649566063dSJacob Faibussowitsch   PetscCall(ComputeRHSMatrix(m, nrhs, &RHS));
659566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(RHS, MAT_DO_NOT_COPY_VALUES, &X));
66c4762a1bSJed Brown 
679566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(A, MATORDERINGND, &perm, &iperm));
68c4762a1bSJed Brown 
699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-inplacelu", &InplaceLU, NULL));
709566063dSJacob Faibussowitsch   PetscCall(MatFactorInfoInitialize(&info));
71c4762a1bSJed Brown   if (!InplaceLU) {
729566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_LU, &F));
73c4762a1bSJed Brown     info.fill = 5.0;
749566063dSJacob Faibussowitsch     PetscCall(MatLUFactorSymbolic(F, A, perm, iperm, &info));
759566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(F, A, &info));
76c4762a1bSJed Brown   } else { /* Test inplace factorization */
779566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &F));
789566063dSJacob Faibussowitsch     PetscCall(MatLUFactor(F, perm, iperm, &info));
79c4762a1bSJed Brown   }
80c4762a1bSJed Brown 
819566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(y, &b1));
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /* MatSolve */
849566063dSJacob Faibussowitsch   PetscCall(MatSolve(F, b, x));
859566063dSJacob Faibussowitsch   PetscCall(MatMult(A, x, b1));
869566063dSJacob Faibussowitsch   PetscCall(VecAXPY(b1, -1.0, b));
879566063dSJacob Faibussowitsch   PetscCall(VecNorm(b1, NORM_2, &norm));
8848a46eb9SPierre Jolivet   if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSolve              : Error of norm %g\n", (double)norm));
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   /* MatSolveTranspose */
919566063dSJacob Faibussowitsch   PetscCall(MatSolveTranspose(F, b, x));
929566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A, x, b1));
939566063dSJacob Faibussowitsch   PetscCall(VecAXPY(b1, -1.0, b));
949566063dSJacob Faibussowitsch   PetscCall(VecNorm(b1, NORM_2, &norm));
9548a46eb9SPierre Jolivet   if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSolveTranspose     : Error of norm %g\n", (double)norm));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* MatSolveAdd */
989566063dSJacob Faibussowitsch   PetscCall(MatSolveAdd(F, b, y, x));
999566063dSJacob Faibussowitsch   PetscCall(MatMult(A, y, b1));
1009566063dSJacob Faibussowitsch   PetscCall(VecScale(b1, -1.0));
1019566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A, x, b1, b1));
1029566063dSJacob Faibussowitsch   PetscCall(VecAXPY(b1, -1.0, b));
1039566063dSJacob Faibussowitsch   PetscCall(VecNorm(b1, NORM_2, &norm));
10448a46eb9SPierre Jolivet   if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSolveAdd           : Error of norm %g\n", (double)norm));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* MatSolveTransposeAdd */
1079566063dSJacob Faibussowitsch   PetscCall(MatSolveTransposeAdd(F, b, y, x));
1089566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A, y, b1));
1099566063dSJacob Faibussowitsch   PetscCall(VecScale(b1, -1.0));
1109566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(A, x, b1, b1));
1119566063dSJacob Faibussowitsch   PetscCall(VecAXPY(b1, -1.0, b));
1129566063dSJacob Faibussowitsch   PetscCall(VecNorm(b1, NORM_2, &norm));
11348a46eb9SPierre Jolivet   if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSolveTransposeAdd  : Error of norm %g\n", (double)norm));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* MatMatSolve */
1169566063dSJacob Faibussowitsch   PetscCall(MatMatSolve(F, RHS, X));
1179566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A, X, MAT_INITIAL_MATRIX, 2.0, &C1));
1189566063dSJacob Faibussowitsch   PetscCall(MatAXPY(C1, -1.0, RHS, SAME_NONZERO_PATTERN));
1199566063dSJacob Faibussowitsch   PetscCall(MatNorm(C1, NORM_FROBENIUS, &norm));
12048a46eb9SPierre Jolivet   if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMatSolve           : Error of norm %g\n", (double)norm));
121c4762a1bSJed Brown 
1229566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
1249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b1));
1259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
1269566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1279566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F));
1289566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&RHS));
1299566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C1));
1309566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&X));
1319566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
1329566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iperm));
1339566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1349566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
135b122ec5aSJacob Faibussowitsch   return 0;
136c4762a1bSJed Brown }
137c4762a1bSJed Brown 
ComputeRHS(DM da,Vec b)138d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeRHS(DM da, Vec b)
139d71ae5a4SJacob Faibussowitsch {
140c4762a1bSJed Brown   PetscInt    mx, my, mz;
141c4762a1bSJed Brown   PetscScalar h;
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   PetscFunctionBegin;
1449566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
145c4762a1bSJed Brown   h = 1.0 / ((mx - 1) * (my - 1) * (mz - 1));
1469566063dSJacob Faibussowitsch   PetscCall(VecSet(b, h));
1473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
148c4762a1bSJed Brown }
149c4762a1bSJed Brown 
ComputeRHSMatrix(PetscInt m,PetscInt nrhs,Mat * C)150d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeRHSMatrix(PetscInt m, PetscInt nrhs, Mat *C)
151d71ae5a4SJacob Faibussowitsch {
152c4762a1bSJed Brown   PetscRandom  rand;
153c4762a1bSJed Brown   Mat          RHS;
154c4762a1bSJed Brown   PetscScalar *array, rval;
155c4762a1bSJed Brown   PetscInt     i, k;
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   PetscFunctionBegin;
1589566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &RHS));
1599566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(RHS, m, PETSC_DECIDE, PETSC_DECIDE, nrhs));
1609566063dSJacob Faibussowitsch   PetscCall(MatSetType(RHS, MATSEQDENSE));
1619566063dSJacob Faibussowitsch   PetscCall(MatSetUp(RHS));
162c4762a1bSJed Brown 
1639566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
1649566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
1659566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(RHS, &array));
166c4762a1bSJed Brown   for (i = 0; i < m; i++) {
1679566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rand, &rval));
168c4762a1bSJed Brown     array[i] = rval;
169c4762a1bSJed Brown   }
170c4762a1bSJed Brown   if (nrhs > 1) {
171c4762a1bSJed Brown     for (k = 1; k < nrhs; k++) {
172ad540459SPierre Jolivet       for (i = 0; i < m; i++) array[m * k + i] = array[i];
173c4762a1bSJed Brown     }
174c4762a1bSJed Brown   }
1759566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(RHS, &array));
1769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(RHS, MAT_FINAL_ASSEMBLY));
1779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(RHS, MAT_FINAL_ASSEMBLY));
178c4762a1bSJed Brown   *C = RHS;
1799566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
1803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
181c4762a1bSJed Brown }
182c4762a1bSJed Brown 
ComputeMatrix(DM da,Mat B)183d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeMatrix(DM da, Mat B)
184d71ae5a4SJacob Faibussowitsch {
185c4762a1bSJed Brown   PetscInt     i, j, k, mx, my, mz, xm, ym, zm, xs, ys, zs, dof, k1, k2, k3;
186c4762a1bSJed Brown   PetscScalar *v, *v_neighbor, Hx, Hy, Hz, HxHydHz, HyHzdHx, HxHzdHy, r1, r2;
187c4762a1bSJed Brown   MatStencil   row, col;
188c4762a1bSJed Brown   PetscRandom  rand;
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   PetscFunctionBegin;
1919566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
1929566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetSeed(rand, 1));
1939566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rand, -.001, .001));
1949566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
195c4762a1bSJed Brown 
1969566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
197c4762a1bSJed Brown   /* For simplicity, this example only works on mx=my=mz */
198e00437b9SBarry Smith   PetscCheck(mx == my && mx == mz, PETSC_COMM_SELF, PETSC_ERR_SUP, "This example only works with mx %" PetscInt_FMT " = my %" PetscInt_FMT " = mz %" PetscInt_FMT, mx, my, mz);
199c4762a1bSJed Brown 
2009371c9d4SSatish Balay   Hx      = 1.0 / (PetscReal)(mx - 1);
2019371c9d4SSatish Balay   Hy      = 1.0 / (PetscReal)(my - 1);
2029371c9d4SSatish Balay   Hz      = 1.0 / (PetscReal)(mz - 1);
2039371c9d4SSatish Balay   HxHydHz = Hx * Hy / Hz;
2049371c9d4SSatish Balay   HxHzdHy = Hx * Hz / Hy;
2059371c9d4SSatish Balay   HyHzdHx = Hy * Hz / Hx;
206c4762a1bSJed Brown 
2079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * dof * dof + 1, &v));
208c4762a1bSJed Brown   v_neighbor = v + dof * dof;
2099566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(v, 2 * dof * dof + 1));
210c4762a1bSJed Brown   k3 = 0;
211c4762a1bSJed Brown   for (k1 = 0; k1 < dof; k1++) {
212c4762a1bSJed Brown     for (k2 = 0; k2 < dof; k2++) {
213c4762a1bSJed Brown       if (k1 == k2) {
214c4762a1bSJed Brown         v[k3]          = 2.0 * (HxHydHz + HxHzdHy + HyHzdHx);
215c4762a1bSJed Brown         v_neighbor[k3] = -HxHydHz;
216c4762a1bSJed Brown       } else {
2179566063dSJacob Faibussowitsch         PetscCall(PetscRandomGetValue(rand, &r1));
2189566063dSJacob Faibussowitsch         PetscCall(PetscRandomGetValue(rand, &r2));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown         v[k3]          = r1;
221c4762a1bSJed Brown         v_neighbor[k3] = r2;
222c4762a1bSJed Brown       }
223c4762a1bSJed Brown       k3++;
224c4762a1bSJed Brown     }
225c4762a1bSJed Brown   }
2269566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
229c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
230c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
2319371c9d4SSatish Balay         row.i = i;
2329371c9d4SSatish Balay         row.j = j;
2339371c9d4SSatish Balay         row.k = k;
234a5b23f4aSJose E. Roman         if (i == 0 || j == 0 || k == 0 || i == mx - 1 || j == my - 1 || k == mz - 1) { /* boundary points */
2359566063dSJacob Faibussowitsch           PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &row, v, INSERT_VALUES));
236c4762a1bSJed Brown         } else { /* interior points */
237c4762a1bSJed Brown           /* center */
2389371c9d4SSatish Balay           col.i = i;
2399371c9d4SSatish Balay           col.j = j;
2409371c9d4SSatish Balay           col.k = k;
2419566063dSJacob Faibussowitsch           PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &col, v, INSERT_VALUES));
242c4762a1bSJed Brown 
243c4762a1bSJed Brown           /* x neighbors */
2449371c9d4SSatish Balay           col.i = i - 1;
2459371c9d4SSatish Balay           col.j = j;
2469371c9d4SSatish Balay           col.k = k;
2479566063dSJacob Faibussowitsch           PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &col, v_neighbor, INSERT_VALUES));
2489371c9d4SSatish Balay           col.i = i + 1;
2499371c9d4SSatish Balay           col.j = j;
2509371c9d4SSatish Balay           col.k = k;
2519566063dSJacob Faibussowitsch           PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &col, v_neighbor, INSERT_VALUES));
252c4762a1bSJed Brown 
253c4762a1bSJed Brown           /* y neighbors */
2549371c9d4SSatish Balay           col.i = i;
2559371c9d4SSatish Balay           col.j = j - 1;
2569371c9d4SSatish Balay           col.k = k;
2579566063dSJacob Faibussowitsch           PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &col, v_neighbor, INSERT_VALUES));
2589371c9d4SSatish Balay           col.i = i;
2599371c9d4SSatish Balay           col.j = j + 1;
2609371c9d4SSatish Balay           col.k = k;
2619566063dSJacob Faibussowitsch           PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &col, v_neighbor, INSERT_VALUES));
262c4762a1bSJed Brown 
263c4762a1bSJed Brown           /* z neighbors */
2649371c9d4SSatish Balay           col.i = i;
2659371c9d4SSatish Balay           col.j = j;
2669371c9d4SSatish Balay           col.k = k - 1;
2679566063dSJacob Faibussowitsch           PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &col, v_neighbor, INSERT_VALUES));
2689371c9d4SSatish Balay           col.i = i;
2699371c9d4SSatish Balay           col.j = j;
2709371c9d4SSatish Balay           col.k = k + 1;
2719566063dSJacob Faibussowitsch           PetscCall(MatSetValuesBlockedStencil(B, 1, &row, 1, &col, v_neighbor, INSERT_VALUES));
272c4762a1bSJed Brown         }
273c4762a1bSJed Brown       }
274c4762a1bSJed Brown     }
275c4762a1bSJed Brown   }
2769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
2789566063dSJacob Faibussowitsch   PetscCall(PetscFree(v));
2799566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
2803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
281c4762a1bSJed Brown }
282c4762a1bSJed Brown 
283c4762a1bSJed Brown /*TEST
284c4762a1bSJed Brown 
285c4762a1bSJed Brown    test:
286c4762a1bSJed Brown       args: -dm_mat_type aij -dof 1
287*3886731fSPierre Jolivet       output_file: output/empty.out
288c4762a1bSJed Brown 
289c4762a1bSJed Brown    test:
290c4762a1bSJed Brown       suffix: 2
291c4762a1bSJed Brown       args: -dm_mat_type aij -dof 1 -inplacelu
292*3886731fSPierre Jolivet       output_file: output/empty.out
293c4762a1bSJed Brown 
294c4762a1bSJed Brown TEST*/
295