1c4762a1bSJed Brown 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown Laplacian in 3D. Use for testing MatSolve routines. 4c4762a1bSJed Brown Modeled by the partial differential equation 5c4762a1bSJed Brown 6c4762a1bSJed Brown - Laplacian u = 1,0 < x,y,z < 1, 7c4762a1bSJed Brown 8c4762a1bSJed Brown with boundary conditions 9c4762a1bSJed Brown u = 1 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1. 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown 12c4762a1bSJed Brown static char help[] = "This example is for testing different MatSolve routines :MatSolve(), MatSolveAdd(), MatSolveTranspose(), MatSolveTransposeAdd(), and MatMatSolve().\n\ 13c4762a1bSJed Brown Example usage: ./ex129 -mat_type aij -dof 2\n\n"; 14c4762a1bSJed Brown 15c4762a1bSJed Brown #include <petscdm.h> 16c4762a1bSJed Brown #include <petscdmda.h> 17c4762a1bSJed Brown 18c4762a1bSJed Brown extern PetscErrorCode ComputeMatrix(DM,Mat); 19c4762a1bSJed Brown extern PetscErrorCode ComputeRHS(DM,Vec); 20c4762a1bSJed Brown extern PetscErrorCode ComputeRHSMatrix(PetscInt,PetscInt,Mat*); 21c4762a1bSJed Brown 22c4762a1bSJed Brown int main(int argc,char **args) 23c4762a1bSJed Brown { 24c4762a1bSJed Brown PetscMPIInt size; 25c4762a1bSJed Brown Vec x,b,y,b1; 26c4762a1bSJed Brown DM da; 27c4762a1bSJed Brown Mat A,F,RHS,X,C1; 28c4762a1bSJed Brown MatFactorInfo info; 29c4762a1bSJed Brown IS perm,iperm; 30c4762a1bSJed Brown PetscInt dof =1,M=8,m,n,nrhs; 31c4762a1bSJed Brown PetscScalar one = 1.0; 32c4762a1bSJed Brown PetscReal norm,tol = 1000*PETSC_MACHINE_EPSILON; 33c4762a1bSJed Brown PetscBool InplaceLU=PETSC_FALSE; 34c4762a1bSJed Brown 35*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 365f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 372c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only"); 385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 40c4762a1bSJed Brown 415f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate(PETSC_COMM_WORLD,&da)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(da,3)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetBoundaryType(da,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetStencilType(da,DMDA_STENCIL_STAR)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetSizes(da,M,M,M)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetNumProcs(da,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetDof(da,dof)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetStencilWidth(da,1)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetOwnershipRanges(da,NULL,NULL,NULL)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(da,MATBAIJ)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 53c4762a1bSJed Brown 545f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&x)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&b)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(b,&y)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeRHS(da,b)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(y,one)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(da,&A)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeMatrix(da,A)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&m,&n)); 62c4762a1bSJed Brown nrhs = 2; 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeRHSMatrix(m,nrhs,&RHS)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(RHS,MAT_DO_NOT_COPY_VALUES,&X)); 66c4762a1bSJed Brown 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering(A,MATORDERINGND,&perm,&iperm)); 68c4762a1bSJed Brown 695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-inplacelu",&InplaceLU,NULL)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatFactorInfoInitialize(&info)); 71c4762a1bSJed Brown if (!InplaceLU) { 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&F)); 73c4762a1bSJed Brown info.fill = 5.0; 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorSymbolic(F,A,perm,iperm,&info)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorNumeric(F,A,&info)); 76c4762a1bSJed Brown } else { /* Test inplace factorization */ 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&F)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactor(F,perm,iperm,&info)); 79c4762a1bSJed Brown } 80c4762a1bSJed Brown 815f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(y,&b1)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* MatSolve */ 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(F,b,x)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,x,b1)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(b1,-1.0,b)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(b1,NORM_2,&norm)); 88c4762a1bSJed Brown if (norm > tol) { 895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"MatSolve : Error of norm %g\n",(double)norm)); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* MatSolveTranspose */ 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolveTranspose(F,b,x)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A,x,b1)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(b1,-1.0,b)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(b1,NORM_2,&norm)); 97c4762a1bSJed Brown if (norm > tol) { 985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"MatSolveTranspose : Error of norm %g\n",(double)norm)); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* MatSolveAdd */ 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolveAdd(F,b,y,x)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,y,b1)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(b1,-1.0)); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(A,x,b1,b1)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(b1,-1.0,b)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(b1,NORM_2,&norm)); 108c4762a1bSJed Brown if (norm > tol) { 1095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"MatSolveAdd : Error of norm %g\n",(double)norm)); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* MatSolveTransposeAdd */ 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolveTransposeAdd(F,b,y,x)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A,y,b1)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(b1,-1.0)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(A,x,b1,b1)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(b1,-1.0,b)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(b1,NORM_2,&norm)); 119c4762a1bSJed Brown if (norm > tol) { 1205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"MatSolveTransposeAdd : Error of norm %g\n",(double)norm)); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* MatMatSolve */ 1245f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatSolve(F,RHS,X)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A,X,MAT_INITIAL_MATRIX,2.0,&C1)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(C1,-1.0,RHS,SAME_NONZERO_PATTERN)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C1,NORM_FROBENIUS,&norm)); 128c4762a1bSJed Brown if (norm > tol) { 1295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"MatMatSolve : Error of norm %g\n",(double)norm)); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 1325f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b1)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&F)); 1385f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&RHS)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C1)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&X)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&perm)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iperm)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 144*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 145*b122ec5aSJacob Faibussowitsch return 0; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown 148c4762a1bSJed Brown PetscErrorCode ComputeRHS(DM da,Vec b) 149c4762a1bSJed Brown { 150c4762a1bSJed Brown PetscInt mx,my,mz; 151c4762a1bSJed Brown PetscScalar h; 152c4762a1bSJed Brown 153c4762a1bSJed Brown PetscFunctionBegin; 1545f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0)); 155c4762a1bSJed Brown h = 1.0/((mx-1)*(my-1)*(mz-1)); 1565f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(b,h)); 157c4762a1bSJed Brown PetscFunctionReturn(0); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown PetscErrorCode ComputeRHSMatrix(PetscInt m,PetscInt nrhs,Mat *C) 161c4762a1bSJed Brown { 162c4762a1bSJed Brown PetscRandom rand; 163c4762a1bSJed Brown Mat RHS; 164c4762a1bSJed Brown PetscScalar *array,rval; 165c4762a1bSJed Brown PetscInt i,k; 166c4762a1bSJed Brown 167c4762a1bSJed Brown PetscFunctionBegin; 1685f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&RHS)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(RHS,m,PETSC_DECIDE,PETSC_DECIDE,nrhs)); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(RHS,MATSEQDENSE)); 1715f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(RHS)); 172c4762a1bSJed Brown 1735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rand)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(RHS,&array)); 176c4762a1bSJed Brown for (i=0; i<m; i++) { 1775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rand,&rval)); 178c4762a1bSJed Brown array[i] = rval; 179c4762a1bSJed Brown } 180c4762a1bSJed Brown if (nrhs > 1) { 181c4762a1bSJed Brown for (k=1; k<nrhs; k++) { 182c4762a1bSJed Brown for (i=0; i<m; i++) { 183c4762a1bSJed Brown array[m*k+i] = array[i]; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown } 186c4762a1bSJed Brown } 1875f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(RHS,&array)); 1885f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(RHS,MAT_FINAL_ASSEMBLY)); 1895f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(RHS,MAT_FINAL_ASSEMBLY)); 190c4762a1bSJed Brown *C = RHS; 1915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 192c4762a1bSJed Brown PetscFunctionReturn(0); 193c4762a1bSJed Brown } 194c4762a1bSJed Brown 195c4762a1bSJed Brown PetscErrorCode ComputeMatrix(DM da,Mat B) 196c4762a1bSJed Brown { 197c4762a1bSJed Brown PetscInt i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs,dof,k1,k2,k3; 198c4762a1bSJed Brown PetscScalar *v,*v_neighbor,Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy,r1,r2; 199c4762a1bSJed Brown MatStencil row,col; 200c4762a1bSJed Brown PetscRandom rand; 201c4762a1bSJed Brown 202c4762a1bSJed Brown PetscFunctionBegin; 2035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 2045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetSeed(rand,1)); 2055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rand,-.001,.001)); 2065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rand)); 207c4762a1bSJed Brown 2085f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0)); 209c4762a1bSJed Brown /* For simplicity, this example only works on mx=my=mz */ 2102c71b3e2SJacob Faibussowitsch PetscCheckFalse(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); 211c4762a1bSJed Brown 212c4762a1bSJed Brown Hx = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1); 213c4762a1bSJed Brown HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx; 214c4762a1bSJed Brown 2155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(2*dof*dof+1,&v)); 216c4762a1bSJed Brown v_neighbor = v + dof*dof; 2175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(v,2*dof*dof+1)); 218c4762a1bSJed Brown k3 = 0; 219c4762a1bSJed Brown for (k1=0; k1<dof; k1++) { 220c4762a1bSJed Brown for (k2=0; k2<dof; k2++) { 221c4762a1bSJed Brown if (k1 == k2) { 222c4762a1bSJed Brown v[k3] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx); 223c4762a1bSJed Brown v_neighbor[k3] = -HxHydHz; 224c4762a1bSJed Brown } else { 2255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rand,&r1)); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rand,&r2)); 227c4762a1bSJed Brown 228c4762a1bSJed Brown v[k3] = r1; 229c4762a1bSJed Brown v_neighbor[k3] = r2; 230c4762a1bSJed Brown } 231c4762a1bSJed Brown k3++; 232c4762a1bSJed Brown } 233c4762a1bSJed Brown } 2345f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 235c4762a1bSJed Brown 236c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 237c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 238c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 239c4762a1bSJed Brown row.i = i; row.j = j; row.k = k; 240a5b23f4aSJose E. Roman if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1) { /* boundary points */ 2415f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&row,v,INSERT_VALUES)); 242c4762a1bSJed Brown } else { /* interior points */ 243c4762a1bSJed Brown /* center */ 244c4762a1bSJed Brown col.i = i; col.j = j; col.k = k; 2455f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&col,v,INSERT_VALUES)); 246c4762a1bSJed Brown 247c4762a1bSJed Brown /* x neighbors */ 248c4762a1bSJed Brown col.i = i-1; col.j = j; col.k = k; 2495f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 250c4762a1bSJed Brown col.i = i+1; col.j = j; col.k = k; 2515f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* y neighbors */ 254c4762a1bSJed Brown col.i = i; col.j = j-1; col.k = k; 2555f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 256c4762a1bSJed Brown col.i = i; col.j = j+1; col.k = k; 2575f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 258c4762a1bSJed Brown 259c4762a1bSJed Brown /* z neighbors */ 260c4762a1bSJed Brown col.i = i; col.j = j; col.k = k-1; 2615f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 262c4762a1bSJed Brown col.i = i; col.j = j; col.k = k+1; 2635f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 264c4762a1bSJed Brown } 265c4762a1bSJed Brown } 266c4762a1bSJed Brown } 267c4762a1bSJed Brown } 2685f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 2705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 2715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 272c4762a1bSJed Brown PetscFunctionReturn(0); 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 275c4762a1bSJed Brown /*TEST 276c4762a1bSJed Brown 277c4762a1bSJed Brown test: 278c4762a1bSJed Brown args: -dm_mat_type aij -dof 1 279c4762a1bSJed Brown output_file: output/ex129.out 280c4762a1bSJed Brown 281c4762a1bSJed Brown test: 282c4762a1bSJed Brown suffix: 2 283c4762a1bSJed Brown args: -dm_mat_type aij -dof 1 -inplacelu 284c4762a1bSJed Brown output_file: output/ex129.out 285c4762a1bSJed Brown 286c4762a1bSJed Brown TEST*/ 287