]> AND Private Git Repository - GMRES2stage.git/blob - code/ex45.c
Logo AND Algorithmique Numérique Distribuée

Private GIT Repository
petites modifs
[GMRES2stage.git] / code / ex45.c
1
2 // /home/couturie/work/petsc-3.5.1/arch-linux2-c-debug/bin/mpirun -np 3    ./ex45  -da_grid_x 160 -da_grid_y 160 -da_grid_z 160  -ksp_type fgmres
3
4
5 /*
6 Laplacian in 3D. Modeled by the partial differential equation
7
8    - Laplacian u = 1,0 < x,y,z < 1,
9
10 with boundary conditions
11
12    u = 1 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1.
13
14    This uses multigrid to solve the linear system
15
16    See src/snes/examples/tutorials/ex50.c
17
18    Can also be run with -pc_type exotic -ksp_type fgmres
19
20 */
21
22 static char help[] = "Solves 3D Laplacian using multigrid.\n\n";
23
24 #include <petscksp.h>
25 #include <petscdm.h>
26 #include <petscdmda.h>
27
28 extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,void*);
29 extern PetscErrorCode ComputeRHS(KSP,Vec,void*);
30 extern PetscErrorCode ComputeInitialGuess(KSP,Vec,void*);
31
32 #undef __FUNCT__
33 #define __FUNCT__ "main"
34
35
36
37
38
39
40
41 int KrylovMinimize(Mat A, Vec b, Vec x) {
42   
43
44   //Variables
45
46   PetscScalar  gamma, alpha, oldgamma, beta;
47   PetscReal norm=20, Eprecision=1e-8, cgprec=1e-40;     
48   PetscInt giter=0, ColS=12, col=0, Emaxiter=50000000, iter=0, iterations=15, Iiter=0;
49   PetscErrorCode ierr;
50   PetscScalar T1, T2;
51   KSP ksp;
52   PetscInt total=0;  
53   PetscInt size;
54   PetscInt Istart,Iend;
55   PetscInt i,its;
56   Vec       x_old, residu;
57   Mat S, AS;
58   PetscScalar *array;
59   PetscInt *ind_row;
60   Vec Alpha, p, ss, vect, r, q, Ax;
61
62
63   PetscInt first=1;
64
65   ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
66   ierr = KSPSetOperators(ksp,A,A);CHKERRQ(ierr);
67   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
68
69
70
71
72   VecGetSize(b,&size);
73
74   ierr = PetscPrintf(PETSC_COMM_WORLD, "Size of vector %D\n", size); CHKERRQ(ierr);  
75
76   PetscInt aa,bb;
77   MatGetOwnershipRange(A,&aa,&bb);
78
79   //  ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%D %D\n", aa,bb); CHKERRQ(ierr);  
80   //PetscSynchronizedFlush(PETSC_COMM_WORLD);
81
82
83   ierr = MatCreate(PETSC_COMM_WORLD, &S);CHKERRQ(ierr);
84   ierr = MatSetSizes(S, bb-aa,  PETSC_DECIDE, size, ColS); CHKERRQ(ierr);
85   ierr = MatSetType(S, MATMPIDENSE); CHKERRQ(ierr);
86   ierr = MatSetUp(S); CHKERRQ(ierr);
87
88   ierr = MatGetOwnershipRange(S, &Istart, &Iend); CHKERRQ(ierr);
89
90   ierr = VecCreate(PETSC_COMM_WORLD, &Alpha); CHKERRQ(ierr);
91   ierr = VecSetSizes(Alpha, PETSC_DECIDE, ColS); CHKERRQ(ierr);
92   ierr = VecSetFromOptions(Alpha); CHKERRQ(ierr);
93   ierr = VecDuplicate(Alpha, &vect); CHKERRQ(ierr);
94   ierr = VecDuplicate(Alpha, &p); CHKERRQ(ierr);
95   ierr = VecDuplicate(Alpha, &ss); CHKERRQ(ierr);
96   ierr = VecDuplicate(b, &r); CHKERRQ(ierr);
97   ierr = VecDuplicate(b, &q); CHKERRQ(ierr);
98   ierr = VecDuplicate(b, &Ax); CHKERRQ(ierr);
99
100   ierr = VecDuplicate(b,&x_old);CHKERRQ(ierr);
101   ierr = VecDuplicate(b,&residu);CHKERRQ(ierr);
102
103
104
105   //indexes of row (these indexes are global)
106   ind_row = (PetscInt*)malloc(sizeof(PetscInt)*(Iend-Istart));
107   for(i=0; i<Iend-Istart; i++) ind_row[i] = i + Istart;
108
109   //Initializations
110   //  ierr = KSPGMRESSetRestart(ksp, 16); CHKERRQ(ierr);
111   ierr = KSPSetTolerances(ksp, 1e-13, 1e-13, PETSC_DEFAULT, 16); CHKERRQ(ierr);
112   ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE); CHKERRQ(ierr);
113
114
115
116   //GMRES WITH MINIMIZATION
117   T1 = MPI_Wtime();
118   while(giter<Emaxiter && norm>Eprecision ){
119     for(col=0; col<ColS  &&  norm>Eprecision; col++){
120
121
122       //Solve
123       ierr = KSPSolve(ksp, b, x); CHKERRQ(ierr);
124
125       ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
126       total += its; Iiter ++;
127
128
129
130       //Build S'
131       ierr = VecGetArray(x, &array);
132       ierr = MatSetValues(S, Iend-Istart, ind_row, 1, &col, array, INSERT_VALUES);
133       VecRestoreArray(x, &array);
134
135
136
137       KSPGetResidualNorm(ksp,&norm);
138       /*      //Error
139       ierr = VecCopy(x, residu); CHKERRQ(ierr);
140       ierr = VecAXPY(residu, -1, x_old); CHKERRQ(ierr);
141       ierr = VecNorm(residu, NORM_INFINITY, &norm); CHKERRQ(ierr);
142        */
143
144
145       ierr = PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g, outer iteration %D\n", norm, giter); CHKERRQ(ierr);
146       ierr = VecCopy(x, x_old); CHKERRQ(ierr);
147
148
149     }
150
151
152     //minimization step
153     if( norm>Eprecision) {
154
155       MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY);
156       MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY);
157
158
159       //Build AS
160       if(first) {
161         MatMatMult(A,S, MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AS);
162         first=0;
163       }
164       else
165         MatMatMult(A,S, MAT_REUSE_MATRIX,PETSC_DEFAULT,&AS);
166
167
168
169
170       //Minimization with CGLS  
171       MatMult(AS, Alpha, r);  VecAYPX(r, -1, b); //r_0 = b-AS*x_0
172
173
174       MatMultTranspose(AS, r, p); //p_0 = AS'*r_0 
175
176
177
178
179       ierr = VecCopy(p, ss); CHKERRQ(ierr); //p_0 = s_0
180       VecNorm(ss, NORM_2, &gamma); gamma = gamma * gamma; //gamma = norm2(s)^2
181       while(gamma>cgprec && iter<iterations){
182         MatMult(AS, p, q); //q = AS*p
183         VecNorm(q, NORM_2, &alpha); alpha = alpha * alpha; //alpha = norm2(q)^2
184         alpha = gamma / alpha; 
185         VecAXPY(Alpha, alpha, p); //Alpha += alpha*p
186         VecAXPY(r, -alpha, q); //r -= alpha*q
187         MatMultTranspose(AS, r, ss); // ss = AS'*r
188         oldgamma = gamma;
189         VecNorm(ss, NORM_2, &gamma); gamma = gamma * gamma; //gamma = norm2(s)^2
190         beta = gamma / oldgamma;
191         VecAYPX(p, beta, ss); //p = s+beta*p;
192         iter ++;
193       } 
194
195       iter = 0; giter ++;
196       //Minimizer
197       MatMult(S, Alpha, x); //x = S*Alpha;
198     }
199
200   }
201   T2 = MPI_Wtime();
202   ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Execution time : %g (s)\n", T2-T1); CHKERRQ(ierr);
203   ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Total number of iterations : %D\n", total); CHKERRQ(ierr);
204
205   return 0;
206
207 }
208
209
210
211
212
213
214
215
216
217 int KrylovMinimizeLSQR(Mat A, Vec b, Vec x) {
218   
219
220   //Variables
221
222   PetscScalar  alpha, beta;
223   PetscReal norm=20, Eprecision=1e-8, tol=1e-40;     
224   PetscInt giter=0, ColS=12, col=0, Emaxiter=50000000, iter=0, iterations=15, Iiter=0;
225   PetscErrorCode ierr;
226   PetscScalar T1, T2;
227   KSP ksp;
228   PetscInt total=0;  
229   PetscInt size;
230   PetscInt Istart,Iend;
231   PetscInt i,its;
232   Vec       x_old, residu;
233   Mat S, AS;
234   PetscScalar *array;
235   PetscInt *ind_row;
236   Vec  Ax;
237   PetscScalar norm_ksp;
238   Vec u,v,d,uu,vt,zero_long,zero_short,x_lsqr;
239
240   PetscInt first=1;
241
242   ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
243   ierr = KSPSetOperators(ksp,A,A);CHKERRQ(ierr);
244   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
245
246
247
248
249   VecGetSize(b,&size);
250
251   ierr = PetscPrintf(PETSC_COMM_WORLD, "Size of vector %D\n", size); CHKERRQ(ierr);  
252
253   PetscInt aa,bb;
254   MatGetOwnershipRange(A,&aa,&bb);
255
256   //  ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%D %D\n", aa,bb); CHKERRQ(ierr);  
257   //PetscSynchronizedFlush(PETSC_COMM_WORLD);
258
259
260   ierr = MatCreate(PETSC_COMM_WORLD, &S);CHKERRQ(ierr);
261   ierr = MatSetSizes(S, bb-aa,  PETSC_DECIDE, size, ColS); CHKERRQ(ierr);
262   ierr = MatSetType(S, MATMPIDENSE); CHKERRQ(ierr);
263   ierr = MatSetUp(S); CHKERRQ(ierr);
264
265   ierr = MatGetOwnershipRange(S, &Istart, &Iend); CHKERRQ(ierr);
266
267
268
269   ierr = VecDuplicate(b, &Ax); CHKERRQ(ierr);
270
271   ierr = VecDuplicate(b,&x_old);CHKERRQ(ierr);
272   ierr = VecDuplicate(b,&residu);CHKERRQ(ierr);
273
274
275   //long vector
276   ierr = VecDuplicate(b,&u);CHKERRQ(ierr);
277
278
279   ierr = VecDuplicate(b,&uu);CHKERRQ(ierr);
280   ierr = VecDuplicate(b,&zero_long);CHKERRQ(ierr);
281   ierr = VecSet(zero_long,0);CHKERRQ(ierr);
282
283   //small vector
284   ierr = VecCreate(PETSC_COMM_WORLD, &v); CHKERRQ(ierr);
285   ierr = VecSetSizes(v, PETSC_DECIDE, ColS); CHKERRQ(ierr);
286   ierr = VecSetFromOptions(v); CHKERRQ(ierr);
287   ierr = VecDuplicate(v,&zero_short);CHKERRQ(ierr);
288   ierr = VecSet(zero_short,0);CHKERRQ(ierr);
289   ierr = VecDuplicate(v,&d);CHKERRQ(ierr);
290   ierr = VecDuplicate(v,&vt);CHKERRQ(ierr);
291   ierr = VecDuplicate(v,&x_lsqr);CHKERRQ(ierr);
292
293
294   //indexes of row (these indexes are global)
295   ind_row = (PetscInt*)malloc(sizeof(PetscInt)*(Iend-Istart));
296   for(i=0; i<Iend-Istart; i++) ind_row[i] = i + Istart;
297
298   //Initializations
299   //  ierr = KSPGMRESSetRestart(ksp, 16); CHKERRQ(ierr);
300   ierr = KSPSetTolerances(ksp, 1e-13, 1e-13, PETSC_DEFAULT, 16); CHKERRQ(ierr);
301   ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE); CHKERRQ(ierr);
302
303
304
305
306   //GMRES WITH MINIMIZATION
307   T1 = MPI_Wtime();
308   while(giter<Emaxiter && norm>Eprecision ){
309     for(col=0; col<ColS  &&  norm>Eprecision; col++){
310
311
312       //Solve
313       ierr = KSPSolve(ksp, b, x); CHKERRQ(ierr);
314
315       ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
316       total += its; Iiter ++;
317
318
319
320       //Build S'
321       ierr = VecGetArray(x, &array);
322       ierr = MatSetValues(S, Iend-Istart, ind_row, 1, &col, array, INSERT_VALUES);
323       VecRestoreArray(x, &array);
324
325
326       KSPGetResidualNorm(ksp,&norm);
327
328         /*
329       //Error
330       ierr = VecCopy(x, residu); CHKERRQ(ierr);
331       ierr = VecAXPY(residu, -1, x_old); CHKERRQ(ierr);
332       ierr = VecNorm(residu, NORM_INFINITY, &norm); CHKERRQ(ierr);
333          */
334
335
336       ierr = PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g, outer iteration %D\n", norm, giter); CHKERRQ(ierr);
337       ierr = VecCopy(x, x_old); CHKERRQ(ierr);
338
339
340     }
341
342
343     //minimization step
344     if( norm>Eprecision) {
345
346       MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY);
347       MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY);
348
349
350
351
352       //Build AS
353       if(first) {
354         MatMatMult(A,S, MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AS);
355         
356         first=0;
357       }
358       else
359         MatMatMult(A,S, MAT_REUSE_MATRIX,PETSC_DEFAULT,&AS);
360
361
362
363
364
365       //LSQR
366       //LSQR
367       //LSQR
368
369
370
371       PetscScalar n2b,tolb,normr,c,s,phibar,normar,norma,thet,rhot,rho,phi;
372       PetscInt stag;
373       tolb = tol * n2b;
374       VecNorm(b, NORM_2, &n2b); //n2b = norm(b);
375       ierr = VecCopy(b, u); CHKERRQ(ierr);   //u=b
376       VecNorm(u, NORM_2, &beta); // beta=norm(u)
377       normr=beta;
378       if (beta != 0) {
379         VecAYPX(u,1/beta,zero_long); //  u = u / beta;
380       }
381       c=1;
382       s=0;
383       phibar=beta;
384       MatMultTranspose(AS, u, v); //v=A'*u
385       ierr = VecSet(x_lsqr,0);CHKERRQ(ierr);
386       VecNorm(v, NORM_2, &alpha); // alpha=norm(v)
387       if (alpha != 0) {
388         VecAYPX(v,1/alpha,zero_short); //  v = v / alpha;
389       }
390       ierr = VecSet(d,0);CHKERRQ(ierr);
391       normar = alpha * beta;
392       norma=0;
393       //stag=0;
394       for(i=0;i<iterations;i++) {
395         MatMult(AS, v, uu); //uu=A*v
396         VecAYPX(u, -alpha, uu); //u = uu-alpha*u;
397         VecNorm(u, NORM_2, &beta); // beta=norm(u)
398         VecAYPX(u,1/beta,zero_long); //  u = u / beta;
399         norma=sqrt(norma*norma+alpha*alpha+beta*beta); // norma = norm([norma alpha beta]);
400         thet = - s * alpha;
401         rhot = c * alpha;
402         rho = sqrt(rhot*rhot + beta*beta);
403         c = rhot / rho;
404         s = - beta / rho;
405         phi = c * phibar;
406         if (phi == 0) {             // stagnation of the method
407           stag = 1;
408         }
409         phibar = s * phibar;
410         VecAYPX(d,-thet,v);       //d = (v - thet * d);
411         VecAYPX(d,1/rho,zero_short);     //d=d/ rho;
412
413
414         VecAXPY(x_lsqr,phi,d);     // x_lsqr=x_lsqr+phi*d
415         normr = abs(s) * normr;
416         MatMultTranspose(AS, u, vt); //vt=A'*u;
417         VecAYPX(v,-beta,vt);  // v = vt - beta * v;
418         VecNorm(v, NORM_2, &alpha); // alpha=norm(v)
419         VecAYPX(v,1/alpha,zero_short); //  v = v / alpha;
420         normar = alpha * abs( s * phi);
421       } 
422       
423
424
425       iter = 0; giter ++;
426       //Minimizer
427       MatMult(S, x_lsqr, x); //x = S*x_small;
428     }
429
430   }
431   T2 = MPI_Wtime();
432   ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Execution time LSQR : %g (s)\n", T2-T1); CHKERRQ(ierr);
433   ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Total number of iterations LSQR : %D\n", total); CHKERRQ(ierr);
434
435   return 0;
436
437 }
438
439
440
441
442
443
444
445 int main(int argc,char **argv)
446 {
447   PetscErrorCode ierr;
448   KSP            ksp;
449   PetscReal      norm;
450   DM             da;
451   Vec            x,b,r;
452   Mat            A;
453   PetscInt total;
454   PetscScalar T1,T2;
455
456   PetscInitialize(&argc,&argv,(char*)0,help);
457
458   ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
459   ierr = DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,-7,-7,-7,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,0,&da);CHKERRQ(ierr);
460   ierr = KSPSetDM(ksp,da);CHKERRQ(ierr);
461   ierr = KSPSetComputeInitialGuess(ksp,ComputeInitialGuess,NULL);CHKERRQ(ierr);
462   ierr = KSPSetComputeRHS(ksp,ComputeRHS,NULL);CHKERRQ(ierr);
463   ierr = KSPSetComputeOperators(ksp,ComputeMatrix,NULL);CHKERRQ(ierr);
464   ierr = DMDestroy(&da);CHKERRQ(ierr);
465
466   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
467
468
469   ierr = KSPSetTolerances(ksp, 1e-10, 1e-10, PETSC_DEFAULT, 50000000); CHKERRQ(ierr);
470   T1 = MPI_Wtime();
471   ierr = KSPSolve(ksp,NULL,NULL);CHKERRQ(ierr);
472   T2 = MPI_Wtime();
473
474
475
476   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
477   ierr = KSPGetRhs(ksp,&b);CHKERRQ(ierr);
478   ierr = VecDuplicate(b,&r);CHKERRQ(ierr);
479   ierr = KSPGetOperators(ksp,&A,NULL);CHKERRQ(ierr);
480
481   ierr = MatMult(A,x,r);CHKERRQ(ierr);
482   ierr = VecAXPY(r,-1.0,b);CHKERRQ(ierr);
483   ierr = VecNorm(r,NORM_2,&norm);CHKERRQ(ierr);
484   ierr = PetscPrintf(PETSC_COMM_WORLD,"Residual norm %g\n",(double)norm);CHKERRQ(ierr);
485   ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Execution time : %g (s)\n", T2-T1); CHKERRQ(ierr);
486
487   ierr = KSPGetIterationNumber(ksp, &total); CHKERRQ(ierr);
488   ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Total number of iterations : %D\n", total); CHKERRQ(ierr);
489
490
491
492   //  ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
493
494
495 {
496
497     Vec x;
498     Vec sol;
499     VecDuplicate(b,&x);
500     VecDuplicate(b,&sol);
501     
502     KrylovMinimize(A, b, x);
503
504
505
506     MatMult(A,x,sol);
507
508
509     VecAXPY(sol,-1,b);
510     VecNorm(sol, NORM_2, &norm);
511     ierr = PetscPrintf(PETSC_COMM_WORLD, "Error2 %g\n",norm);
512
513     ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
514
515   }
516
517
518  {
519
520     Vec x;
521     Vec sol;
522     VecDuplicate(b,&x);
523     VecDuplicate(b,&sol);
524     
525     KrylovMinimizeLSQR(A, b, x);
526
527
528
529     MatMult(A,x,sol);
530
531
532     VecAXPY(sol,-1,b);
533     VecNorm(sol, NORM_2, &norm);
534     ierr = PetscPrintf(PETSC_COMM_WORLD, "Error LSQR %g\n",norm);
535
536
537     ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
538   }
539  
540
541   ierr = VecDestroy(&r);CHKERRQ(ierr);
542
543   ierr = PetscFinalize();
544
545   return 0;
546 }
547
548 #undef __FUNCT__
549 #define __FUNCT__ "ComputeRHS"
550 PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
551 {
552   PetscErrorCode ierr;
553   PetscInt       i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs;
554   DM             dm;
555   PetscScalar    Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy;
556   PetscScalar    ***barray;
557
558   PetscFunctionBeginUser;
559   ierr    = KSPGetDM(ksp,&dm);CHKERRQ(ierr);
560   ierr    = DMDAGetInfo(dm,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
561   Hx      = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1);
562   HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx;
563   ierr    = DMDAGetCorners(dm,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
564   ierr    = DMDAVecGetArray(dm,b,&barray);CHKERRQ(ierr);
565
566   for (k=zs; k<zs+zm; k++) {
567     for (j=ys; j<ys+ym; j++) {
568       for (i=xs; i<xs+xm; i++) {
569         if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1) {
570           barray[k][j][i] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);
571         } else {
572           barray[k][j][i] = Hx*Hy*Hz;
573         }
574       }
575     }
576   }
577   ierr = DMDAVecRestoreArray(dm,b,&barray);CHKERRQ(ierr);
578   PetscFunctionReturn(0);
579 }
580
581 #undef __FUNCT__
582 #define __FUNCT__ "ComputeInitialGuess"
583 PetscErrorCode ComputeInitialGuess(KSP ksp,Vec b,void *ctx)
584 {
585   PetscErrorCode ierr;
586
587   PetscFunctionBeginUser;
588   ierr = VecSet(b,0);CHKERRQ(ierr);
589   PetscFunctionReturn(0);
590 }
591
592 #undef __FUNCT__
593 #define __FUNCT__ "ComputeMatrix"
594 PetscErrorCode ComputeMatrix(KSP ksp,Mat jac,Mat B,void *ctx)
595 {
596   DM             da;
597   PetscErrorCode ierr;
598   PetscInt       i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs;
599   PetscScalar    v[7],Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy;
600   MatStencil     row,col[7];
601
602   PetscFunctionBeginUser;
603   ierr    = KSPGetDM(ksp,&da);CHKERRQ(ierr);
604   ierr    = DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
605   Hx      = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1);
606   HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx;
607   ierr    = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
608
609   for (k=zs; k<zs+zm; k++) {
610     for (j=ys; j<ys+ym; j++) {
611       for (i=xs; i<xs+xm; i++) {
612         row.i = i; row.j = j; row.k = k;
613         if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1) {
614           v[0] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);
615           ierr = MatSetValuesStencil(B,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr);
616         } else {
617           v[0] = -HxHydHz;col[0].i = i; col[0].j = j; col[0].k = k-1;
618           v[1] = -HxHzdHy;col[1].i = i; col[1].j = j-1; col[1].k = k;
619           v[2] = -HyHzdHx;col[2].i = i-1; col[2].j = j; col[2].k = k;
620           v[3] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);col[3].i = row.i; col[3].j = row.j; col[3].k = row.k;
621           v[4] = -HyHzdHx;col[4].i = i+1; col[4].j = j; col[4].k = k;
622           v[5] = -HxHzdHy;col[5].i = i; col[5].j = j+1; col[5].k = k;
623           v[6] = -HxHydHz;col[6].i = i; col[6].j = j; col[6].k = k+1;
624           ierr = MatSetValuesStencil(B,1,&row,7,col,v,INSERT_VALUES);CHKERRQ(ierr);
625         }
626       }
627     }
628   }
629   ierr   = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
630   ierr   = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
631   PetscFunctionReturn(0);
632 }