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