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 PetscErrorCode ierr; 25c4762a1bSJed Brown PetscMPIInt size; 26c4762a1bSJed Brown Vec x,b,y,b1; 27c4762a1bSJed Brown DM da; 28c4762a1bSJed Brown Mat A,F,RHS,X,C1; 29c4762a1bSJed Brown MatFactorInfo info; 30c4762a1bSJed Brown IS perm,iperm; 31c4762a1bSJed Brown PetscInt dof =1,M=8,m,n,nrhs; 32c4762a1bSJed Brown PetscScalar one = 1.0; 33c4762a1bSJed Brown PetscReal norm,tol = 1000*PETSC_MACHINE_EPSILON; 34c4762a1bSJed Brown PetscBool InplaceLU=PETSC_FALSE; 35c4762a1bSJed Brown 36c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 37*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 382c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only"); 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL)); 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 41c4762a1bSJed Brown 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate(PETSC_COMM_WORLD,&da)); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(da,3)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetBoundaryType(da,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetStencilType(da,DMDA_STENCIL_STAR)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetSizes(da,M,M,M)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetNumProcs(da,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetDof(da,dof)); 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetStencilWidth(da,1)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetOwnershipRanges(da,NULL,NULL,NULL)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(da,MATBAIJ)); 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 54c4762a1bSJed Brown 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&x)); 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&b)); 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(b,&y)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeRHS(da,b)); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(y,one)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(da,&A)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeMatrix(da,A)); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&m,&n)); 63c4762a1bSJed Brown nrhs = 2; 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL)); 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeRHSMatrix(m,nrhs,&RHS)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(RHS,MAT_DO_NOT_COPY_VALUES,&X)); 67c4762a1bSJed Brown 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering(A,MATORDERINGND,&perm,&iperm)); 69c4762a1bSJed Brown 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-inplacelu",&InplaceLU,NULL)); 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatFactorInfoInitialize(&info)); 72c4762a1bSJed Brown if (!InplaceLU) { 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&F)); 74c4762a1bSJed Brown info.fill = 5.0; 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorSymbolic(F,A,perm,iperm,&info)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorNumeric(F,A,&info)); 77c4762a1bSJed Brown } else { /* Test inplace factorization */ 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&F)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactor(F,perm,iperm,&info)); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(y,&b1)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* MatSolve */ 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(F,b,x)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,x,b1)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(b1,-1.0,b)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(b1,NORM_2,&norm)); 89c4762a1bSJed Brown if (norm > tol) { 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"MatSolve : Error of norm %g\n",(double)norm)); 91c4762a1bSJed Brown } 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* MatSolveTranspose */ 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolveTranspose(F,b,x)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A,x,b1)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(b1,-1.0,b)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(b1,NORM_2,&norm)); 98c4762a1bSJed Brown if (norm > tol) { 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"MatSolveTranspose : Error of norm %g\n",(double)norm)); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* MatSolveAdd */ 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolveAdd(F,b,y,x)); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,y,b1)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(b1,-1.0)); 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(A,x,b1,b1)); 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(b1,-1.0,b)); 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(b1,NORM_2,&norm)); 109c4762a1bSJed Brown if (norm > tol) { 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"MatSolveAdd : Error of norm %g\n",(double)norm)); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* MatSolveTransposeAdd */ 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolveTransposeAdd(F,b,y,x)); 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A,y,b1)); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(b1,-1.0)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(A,x,b1,b1)); 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(b1,-1.0,b)); 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(b1,NORM_2,&norm)); 120c4762a1bSJed Brown if (norm > tol) { 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"MatSolveTransposeAdd : Error of norm %g\n",(double)norm)); 122c4762a1bSJed Brown } 123c4762a1bSJed Brown 124c4762a1bSJed Brown /* MatMatSolve */ 125*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatSolve(F,RHS,X)); 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A,X,MAT_INITIAL_MATRIX,2.0,&C1)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(C1,-1.0,RHS,SAME_NONZERO_PATTERN)); 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C1,NORM_FROBENIUS,&norm)); 129c4762a1bSJed Brown if (norm > tol) { 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"MatMatSolve : Error of norm %g\n",(double)norm)); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b1)); 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&F)); 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&RHS)); 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C1)); 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&X)); 142*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&perm)); 143*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iperm)); 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 145c4762a1bSJed Brown ierr = PetscFinalize(); 146c4762a1bSJed Brown return ierr; 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 149c4762a1bSJed Brown PetscErrorCode ComputeRHS(DM da,Vec b) 150c4762a1bSJed Brown { 151c4762a1bSJed Brown PetscInt mx,my,mz; 152c4762a1bSJed Brown PetscScalar h; 153c4762a1bSJed Brown 154c4762a1bSJed Brown PetscFunctionBegin; 155*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0)); 156c4762a1bSJed Brown h = 1.0/((mx-1)*(my-1)*(mz-1)); 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(b,h)); 158c4762a1bSJed Brown PetscFunctionReturn(0); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161c4762a1bSJed Brown PetscErrorCode ComputeRHSMatrix(PetscInt m,PetscInt nrhs,Mat *C) 162c4762a1bSJed Brown { 163c4762a1bSJed Brown PetscRandom rand; 164c4762a1bSJed Brown Mat RHS; 165c4762a1bSJed Brown PetscScalar *array,rval; 166c4762a1bSJed Brown PetscInt i,k; 167c4762a1bSJed Brown 168c4762a1bSJed Brown PetscFunctionBegin; 169*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&RHS)); 170*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(RHS,m,PETSC_DECIDE,PETSC_DECIDE,nrhs)); 171*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(RHS,MATSEQDENSE)); 172*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(RHS)); 173c4762a1bSJed Brown 174*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 175*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rand)); 176*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(RHS,&array)); 177c4762a1bSJed Brown for (i=0; i<m; i++) { 178*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rand,&rval)); 179c4762a1bSJed Brown array[i] = rval; 180c4762a1bSJed Brown } 181c4762a1bSJed Brown if (nrhs > 1) { 182c4762a1bSJed Brown for (k=1; k<nrhs; k++) { 183c4762a1bSJed Brown for (i=0; i<m; i++) { 184c4762a1bSJed Brown array[m*k+i] = array[i]; 185c4762a1bSJed Brown } 186c4762a1bSJed Brown } 187c4762a1bSJed Brown } 188*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(RHS,&array)); 189*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(RHS,MAT_FINAL_ASSEMBLY)); 190*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(RHS,MAT_FINAL_ASSEMBLY)); 191c4762a1bSJed Brown *C = RHS; 192*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 193c4762a1bSJed Brown PetscFunctionReturn(0); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196c4762a1bSJed Brown PetscErrorCode ComputeMatrix(DM da,Mat B) 197c4762a1bSJed Brown { 198c4762a1bSJed Brown PetscInt i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs,dof,k1,k2,k3; 199c4762a1bSJed Brown PetscScalar *v,*v_neighbor,Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy,r1,r2; 200c4762a1bSJed Brown MatStencil row,col; 201c4762a1bSJed Brown PetscRandom rand; 202c4762a1bSJed Brown 203c4762a1bSJed Brown PetscFunctionBegin; 204*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 205*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetSeed(rand,1)); 206*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rand,-.001,.001)); 207*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rand)); 208c4762a1bSJed Brown 209*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0)); 210c4762a1bSJed Brown /* For simplicity, this example only works on mx=my=mz */ 2112c71b3e2SJacob 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); 212c4762a1bSJed Brown 213c4762a1bSJed Brown Hx = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1); 214c4762a1bSJed Brown HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx; 215c4762a1bSJed Brown 216*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(2*dof*dof+1,&v)); 217c4762a1bSJed Brown v_neighbor = v + dof*dof; 218*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(v,2*dof*dof+1)); 219c4762a1bSJed Brown k3 = 0; 220c4762a1bSJed Brown for (k1=0; k1<dof; k1++) { 221c4762a1bSJed Brown for (k2=0; k2<dof; k2++) { 222c4762a1bSJed Brown if (k1 == k2) { 223c4762a1bSJed Brown v[k3] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx); 224c4762a1bSJed Brown v_neighbor[k3] = -HxHydHz; 225c4762a1bSJed Brown } else { 226*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rand,&r1)); 227*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rand,&r2)); 228c4762a1bSJed Brown 229c4762a1bSJed Brown v[k3] = r1; 230c4762a1bSJed Brown v_neighbor[k3] = r2; 231c4762a1bSJed Brown } 232c4762a1bSJed Brown k3++; 233c4762a1bSJed Brown } 234c4762a1bSJed Brown } 235*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 236c4762a1bSJed Brown 237c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 238c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 239c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 240c4762a1bSJed Brown row.i = i; row.j = j; row.k = k; 241a5b23f4aSJose E. Roman if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1) { /* boundary points */ 242*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&row,v,INSERT_VALUES)); 243c4762a1bSJed Brown } else { /* interior points */ 244c4762a1bSJed Brown /* center */ 245c4762a1bSJed Brown col.i = i; col.j = j; col.k = k; 246*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&col,v,INSERT_VALUES)); 247c4762a1bSJed Brown 248c4762a1bSJed Brown /* x neighbors */ 249c4762a1bSJed Brown col.i = i-1; col.j = j; col.k = k; 250*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 251c4762a1bSJed Brown col.i = i+1; col.j = j; col.k = k; 252*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* y neighbors */ 255c4762a1bSJed Brown col.i = i; col.j = j-1; col.k = k; 256*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 257c4762a1bSJed Brown col.i = i; col.j = j+1; col.k = k; 258*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 259c4762a1bSJed Brown 260c4762a1bSJed Brown /* z neighbors */ 261c4762a1bSJed Brown col.i = i; col.j = j; col.k = k-1; 262*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 263c4762a1bSJed Brown col.i = i; col.j = j; col.k = k+1; 264*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES)); 265c4762a1bSJed Brown } 266c4762a1bSJed Brown } 267c4762a1bSJed Brown } 268c4762a1bSJed Brown } 269*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 270*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 271*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 272*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 273c4762a1bSJed Brown PetscFunctionReturn(0); 274c4762a1bSJed Brown } 275c4762a1bSJed Brown 276c4762a1bSJed Brown /*TEST 277c4762a1bSJed Brown 278c4762a1bSJed Brown test: 279c4762a1bSJed Brown args: -dm_mat_type aij -dof 1 280c4762a1bSJed Brown output_file: output/ex129.out 281c4762a1bSJed Brown 282c4762a1bSJed Brown test: 283c4762a1bSJed Brown suffix: 2 284c4762a1bSJed Brown args: -dm_mat_type aij -dof 1 -inplacelu 285c4762a1bSJed Brown output_file: output/ex129.out 286c4762a1bSJed Brown 287c4762a1bSJed Brown TEST*/ 288