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