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

Private GIT Repository
ajout ex45 LSQR pas rapide...
[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=5e-5, 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-10, 1e-10, 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       //Error
138       ierr = VecCopy(x, residu); CHKERRQ(ierr);
139       ierr = VecAXPY(residu, -1, x_old); CHKERRQ(ierr);
140       ierr = VecNorm(residu, NORM_INFINITY, &norm); CHKERRQ(ierr);
141
142
143
144       ierr = PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g, outer iteration %D\n", norm, giter); CHKERRQ(ierr);
145       ierr = VecCopy(x, x_old); CHKERRQ(ierr);
146
147
148     }
149
150
151     //minimization step
152     if( norm>Eprecision) {
153
154       MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY);
155       MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY);
156
157
158       //Build AS
159       if(first) {
160         MatMatMult(A,S, MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AS);
161         first=0;
162       }
163       else
164         MatMatMult(A,S, MAT_REUSE_MATRIX,PETSC_DEFAULT,&AS);
165
166
167
168
169       //Minimization with CGLS  
170       MatMult(AS, Alpha, r);  VecAYPX(r, -1, b); //r_0 = b-AS*x_0
171
172
173       MatMultTranspose(AS, r, p); //p_0 = AS'*r_0 
174
175
176
177
178       ierr = VecCopy(p, ss); CHKERRQ(ierr); //p_0 = s_0
179       VecNorm(ss, NORM_2, &gamma); gamma = gamma * gamma; //gamma = norm2(s)^2
180       while(gamma>cgprec && iter<iterations){
181         MatMult(AS, p, q); //q = AS*p
182         VecNorm(q, NORM_2, &alpha); alpha = alpha * alpha; //alpha = norm2(q)^2
183         alpha = gamma / alpha; 
184         VecAXPY(Alpha, alpha, p); //Alpha += alpha*p
185         VecAXPY(r, -alpha, q); //r -= alpha*q
186         MatMultTranspose(AS, r, ss); // ss = AS'*r
187         oldgamma = gamma;
188         VecNorm(ss, NORM_2, &gamma); gamma = gamma * gamma; //gamma = norm2(s)^2
189         beta = gamma / oldgamma;
190         VecAYPX(p, beta, ss); //p = s+beta*p;
191         iter ++;
192       } 
193
194       iter = 0; giter ++;
195       //Minimizer
196       MatMult(S, Alpha, x); //x = S*Alpha;
197     }
198
199   }
200   T2 = MPI_Wtime();
201   ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Execution time : %g (s)\n", T2-T1); CHKERRQ(ierr);
202   ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Total number of iterations : %D\n", total); CHKERRQ(ierr);
203
204   return 0;
205
206 }
207
208
209
210
211
212
213
214
215
216 int KrylovMinimizeLSQR(Mat A, Vec b, Vec x) {
217   
218
219   //Variables
220
221   PetscScalar  alpha, beta;
222   PetscReal norm=20, Eprecision=5e-5, tol=1e-40;     
223   PetscInt giter=0, ColS=12, col=0, Emaxiter=50000000, iter=0, iterations=15, Iiter=0;
224   PetscErrorCode ierr;
225   PetscScalar T1, T2;
226   KSP ksp;
227   PetscInt total=0;  
228   PetscInt size;
229   PetscInt Istart,Iend;
230   PetscInt i,its;
231   Vec       x_old, residu;
232   Mat S, AS;
233   PetscScalar *array;
234   PetscInt *ind_row;
235   Vec  Ax;
236   PetscScalar norm_ksp;
237   Vec u,v,d,uu,vt,zero_long,zero_short,x_lsqr;
238
239   PetscInt first=1;
240
241   ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
242   ierr = KSPSetOperators(ksp,A,A);CHKERRQ(ierr);
243   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
244
245
246
247
248   VecGetSize(b,&size);
249
250   ierr = PetscPrintf(PETSC_COMM_WORLD, "Size of vector %D\n", size); CHKERRQ(ierr);  
251
252   PetscInt aa,bb;
253   MatGetOwnershipRange(A,&aa,&bb);
254
255   //  ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%D %D\n", aa,bb); CHKERRQ(ierr);  
256   //PetscSynchronizedFlush(PETSC_COMM_WORLD);
257
258
259   ierr = MatCreate(PETSC_COMM_WORLD, &S);CHKERRQ(ierr);
260   ierr = MatSetSizes(S, bb-aa,  PETSC_DECIDE, size, ColS); CHKERRQ(ierr);
261   ierr = MatSetType(S, MATMPIDENSE); CHKERRQ(ierr);
262   ierr = MatSetUp(S); CHKERRQ(ierr);
263
264   ierr = MatGetOwnershipRange(S, &Istart, &Iend); CHKERRQ(ierr);
265
266
267
268   ierr = VecDuplicate(b, &Ax); CHKERRQ(ierr);
269
270   ierr = VecDuplicate(b,&x_old);CHKERRQ(ierr);
271   ierr = VecDuplicate(b,&residu);CHKERRQ(ierr);
272
273
274   //long vector
275   ierr = VecDuplicate(b,&u);CHKERRQ(ierr);
276
277
278   ierr = VecDuplicate(b,&uu);CHKERRQ(ierr);
279   ierr = VecDuplicate(b,&zero_long);CHKERRQ(ierr);
280   ierr = VecSet(zero_long,0);CHKERRQ(ierr);
281
282   //small vector
283   ierr = VecCreate(PETSC_COMM_WORLD, &v); CHKERRQ(ierr);
284   ierr = VecSetSizes(v, PETSC_DECIDE, ColS); CHKERRQ(ierr);
285   ierr = VecSetFromOptions(v); CHKERRQ(ierr);
286   ierr = VecDuplicate(v,&zero_short);CHKERRQ(ierr);
287   ierr = VecSet(zero_short,0);CHKERRQ(ierr);
288   ierr = VecDuplicate(v,&d);CHKERRQ(ierr);
289   ierr = VecDuplicate(v,&vt);CHKERRQ(ierr);
290   ierr = VecDuplicate(v,&x_lsqr);CHKERRQ(ierr);
291
292
293   //indexes of row (these indexes are global)
294   ind_row = (PetscInt*)malloc(sizeof(PetscInt)*(Iend-Istart));
295   for(i=0; i<Iend-Istart; i++) ind_row[i] = i + Istart;
296
297   //Initializations
298   //  ierr = KSPGMRESSetRestart(ksp, 16); CHKERRQ(ierr);
299   ierr = KSPSetTolerances(ksp, 1e-10, 1e-10, PETSC_DEFAULT, 16); CHKERRQ(ierr);
300   ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE); CHKERRQ(ierr);
301
302
303
304
305   //GMRES WITH MINIMIZATION
306   T1 = MPI_Wtime();
307   while(giter<Emaxiter && norm>Eprecision ){
308     for(col=0; col<ColS  &&  norm>Eprecision; col++){
309
310
311       //Solve
312       ierr = KSPSolve(ksp, b, x); CHKERRQ(ierr);
313
314       ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
315       total += its; Iiter ++;
316
317
318
319       //Build S'
320       ierr = VecGetArray(x, &array);
321       ierr = MatSetValues(S, Iend-Istart, ind_row, 1, &col, array, INSERT_VALUES);
322       VecRestoreArray(x, &array);
323
324
325
326       //Error
327       ierr = VecCopy(x, residu); CHKERRQ(ierr);
328       ierr = VecAXPY(residu, -1, x_old); CHKERRQ(ierr);
329       ierr = VecNorm(residu, NORM_INFINITY, &norm); CHKERRQ(ierr);
330
331
332
333       ierr = PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g, outer iteration %D\n", norm, giter); CHKERRQ(ierr);
334       ierr = VecCopy(x, x_old); CHKERRQ(ierr);
335
336
337     }
338
339
340     //minimization step
341     if( norm>Eprecision) {
342
343       MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY);
344       MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY);
345
346
347
348
349       //Build AS
350       if(first) {
351         MatMatMult(A,S, MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AS);
352         
353         first=0;
354       }
355       else
356         MatMatMult(A,S, MAT_REUSE_MATRIX,PETSC_DEFAULT,&AS);
357
358
359
360
361
362       //LSQR
363       //LSQR
364       //LSQR
365
366
367
368       PetscScalar n2b,tolb,normr,c,s,phibar,normar,norma,thet,rhot,rho,phi;
369       PetscInt stag;
370       tolb = tol * n2b;
371       VecNorm(b, NORM_2, &n2b); //n2b = norm(b);
372       ierr = VecCopy(b, u); CHKERRQ(ierr);   //u=b
373       VecNorm(u, NORM_2, &beta); // beta=norm(u)
374       normr=beta;
375       if (beta != 0) {
376         VecAYPX(u,1/beta,zero_long); //  u = u / beta;
377       }
378       c=1;
379       s=0;
380       phibar=beta;
381       MatMultTranspose(AS, u, v); //v=A'*u
382       ierr = VecSet(x_lsqr,0);CHKERRQ(ierr);
383       VecNorm(v, NORM_2, &alpha); // alpha=norm(v)
384       if (alpha != 0) {
385         VecAYPX(v,1/alpha,zero_short); //  v = v / alpha;
386       }
387       ierr = VecSet(d,0);CHKERRQ(ierr);
388       normar = alpha * beta;
389       norma=0;
390       //stag=0;
391       for(i=0;i<iterations;i++) {
392         MatMult(AS, v, uu); //uu=A*v
393         VecAYPX(u, -alpha, uu); //u = uu-alpha*u;
394         VecNorm(u, NORM_2, &beta); // beta=norm(u)
395         VecAYPX(u,1/beta,zero_long); //  u = u / beta;
396         norma=sqrt(norma*norma+alpha*alpha+beta*beta); // norma = norm([norma alpha beta]);
397         thet = - s * alpha;
398         rhot = c * alpha;
399         rho = sqrt(rhot*rhot + beta*beta);
400         c = rhot / rho;
401         s = - beta / rho;
402         phi = c * phibar;
403         if (phi == 0) {             // stagnation of the method
404           stag = 1;
405         }
406         phibar = s * phibar;
407         VecAYPX(d,-thet,v);       //d = (v - thet * d);
408         VecAYPX(d,1/rho,zero_short);     //d=d/ rho;
409
410
411         if (normar/(norma*normr) <= tol) { // check for convergence in min{|b-A*x|}
412           break;
413         }
414         if (normr <= tolb) {           // check for convergence in A*x=b
415           break;
416         }
417
418
419         VecAXPY(x_lsqr,phi,d);     // x_lsqr=x_lsqr+phi*d
420         normr = abs(s) * normr;
421         MatMultTranspose(AS, u, vt); //vt=A'*u;
422         VecAYPX(v,-beta,vt);  // v = vt - beta * v;
423         VecNorm(v, NORM_2, &alpha); // alpha=norm(v)
424         VecAYPX(v,1/alpha,zero_short); //  v = v / alpha;
425         normar = alpha * abs( s * phi);
426       } 
427       
428
429
430       iter = 0; giter ++;
431       //Minimizer
432       MatMult(S, x_lsqr, x); //x = S*x_small;
433     }
434
435   }
436   T2 = MPI_Wtime();
437   ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Execution time LSQR : %g (s)\n", T2-T1); CHKERRQ(ierr);
438   ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Total number of iterations LSQR : %D\n", total); CHKERRQ(ierr);
439
440   return 0;
441
442 }
443
444
445
446
447
448
449
450 int main(int argc,char **argv)
451 {
452   PetscErrorCode ierr;
453   KSP            ksp;
454   PetscReal      norm;
455   DM             da;
456   Vec            x,b,r;
457   Mat            A;
458   PetscInt total;
459   PetscScalar T1,T2;
460
461   PetscInitialize(&argc,&argv,(char*)0,help);
462
463   ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
464   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);
465   ierr = KSPSetDM(ksp,da);CHKERRQ(ierr);
466   ierr = KSPSetComputeInitialGuess(ksp,ComputeInitialGuess,NULL);CHKERRQ(ierr);
467   ierr = KSPSetComputeRHS(ksp,ComputeRHS,NULL);CHKERRQ(ierr);
468   ierr = KSPSetComputeOperators(ksp,ComputeMatrix,NULL);CHKERRQ(ierr);
469   ierr = DMDestroy(&da);CHKERRQ(ierr);
470
471   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
472
473
474   ierr = KSPSetTolerances(ksp, 1e-7, 1e-7, PETSC_DEFAULT, 50000000); CHKERRQ(ierr);
475   T1 = MPI_Wtime();
476   ierr = KSPSolve(ksp,NULL,NULL);CHKERRQ(ierr);
477   T2 = MPI_Wtime();
478
479
480
481   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
482   ierr = KSPGetRhs(ksp,&b);CHKERRQ(ierr);
483   ierr = VecDuplicate(b,&r);CHKERRQ(ierr);
484   ierr = KSPGetOperators(ksp,&A,NULL);CHKERRQ(ierr);
485
486   ierr = MatMult(A,x,r);CHKERRQ(ierr);
487   ierr = VecAXPY(r,-1.0,b);CHKERRQ(ierr);
488   ierr = VecNorm(r,NORM_2,&norm);CHKERRQ(ierr);
489   ierr = PetscPrintf(PETSC_COMM_WORLD,"Residual norm %g\n",(double)norm);CHKERRQ(ierr);
490   ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Execution time : %g (s)\n", T2-T1); CHKERRQ(ierr);
491
492   ierr = KSPGetIterationNumber(ksp, &total); CHKERRQ(ierr);
493   ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Total number of iterations : %D\n", total); CHKERRQ(ierr);
494
495
496
497   //  ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
498
499
500 {
501
502     Vec x;
503     Vec sol;
504     VecDuplicate(b,&x);
505     VecDuplicate(b,&sol);
506     
507     KrylovMinimize(A, b, x);
508
509
510
511     MatMult(A,x,sol);
512
513
514     VecAXPY(sol,-1,b);
515     VecNorm(sol, NORM_2, &norm);
516     ierr = PetscPrintf(PETSC_COMM_WORLD, "Error2 %g\n",norm);
517
518     ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
519
520   }
521
522
523  {
524
525     Vec x;
526     Vec sol;
527     VecDuplicate(b,&x);
528     VecDuplicate(b,&sol);
529     
530     KrylovMinimizeLSQR(A, b, x);
531
532
533
534     MatMult(A,x,sol);
535
536
537     VecAXPY(sol,-1,b);
538     VecNorm(sol, NORM_2, &norm);
539     ierr = PetscPrintf(PETSC_COMM_WORLD, "Error LSQR %g\n",norm);
540
541
542     ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
543   }
544  
545
546   ierr = VecDestroy(&r);CHKERRQ(ierr);
547
548   ierr = PetscFinalize();
549
550   return 0;
551 }
552
553 #undef __FUNCT__
554 #define __FUNCT__ "ComputeRHS"
555 PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
556 {
557   PetscErrorCode ierr;
558   PetscInt       i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs;
559   DM             dm;
560   PetscScalar    Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy;
561   PetscScalar    ***barray;
562
563   PetscFunctionBeginUser;
564   ierr    = KSPGetDM(ksp,&dm);CHKERRQ(ierr);
565   ierr    = DMDAGetInfo(dm,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
566   Hx      = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1);
567   HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx;
568   ierr    = DMDAGetCorners(dm,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
569   ierr    = DMDAVecGetArray(dm,b,&barray);CHKERRQ(ierr);
570
571   for (k=zs; k<zs+zm; k++) {
572     for (j=ys; j<ys+ym; j++) {
573       for (i=xs; i<xs+xm; i++) {
574         if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1) {
575           barray[k][j][i] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);
576         } else {
577           barray[k][j][i] = Hx*Hy*Hz;
578         }
579       }
580     }
581   }
582   ierr = DMDAVecRestoreArray(dm,b,&barray);CHKERRQ(ierr);
583   PetscFunctionReturn(0);
584 }
585
586 #undef __FUNCT__
587 #define __FUNCT__ "ComputeInitialGuess"
588 PetscErrorCode ComputeInitialGuess(KSP ksp,Vec b,void *ctx)
589 {
590   PetscErrorCode ierr;
591
592   PetscFunctionBeginUser;
593   ierr = VecSet(b,0);CHKERRQ(ierr);
594   PetscFunctionReturn(0);
595 }
596
597 #undef __FUNCT__
598 #define __FUNCT__ "ComputeMatrix"
599 PetscErrorCode ComputeMatrix(KSP ksp,Mat jac,Mat B,void *ctx)
600 {
601   DM             da;
602   PetscErrorCode ierr;
603   PetscInt       i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs;
604   PetscScalar    v[7],Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy;
605   MatStencil     row,col[7];
606
607   PetscFunctionBeginUser;
608   ierr    = KSPGetDM(ksp,&da);CHKERRQ(ierr);
609   ierr    = DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
610   Hx      = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1);
611   HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx;
612   ierr    = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
613
614   for (k=zs; k<zs+zm; k++) {
615     for (j=ys; j<ys+ym; j++) {
616       for (i=xs; i<xs+xm; i++) {
617         row.i = i; row.j = j; row.k = k;
618         if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1) {
619           v[0] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);
620           ierr = MatSetValuesStencil(B,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr);
621         } else {
622           v[0] = -HxHydHz;col[0].i = i; col[0].j = j; col[0].k = k-1;
623           v[1] = -HxHzdHy;col[1].i = i; col[1].j = j-1; col[1].k = k;
624           v[2] = -HyHzdHx;col[2].i = i-1; col[2].j = j; col[2].k = k;
625           v[3] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);col[3].i = row.i; col[3].j = row.j; col[3].k = row.k;
626           v[4] = -HyHzdHx;col[4].i = i+1; col[4].j = j; col[4].k = k;
627           v[5] = -HxHzdHy;col[5].i = i; col[5].j = j+1; col[5].k = k;
628           v[6] = -HxHydHz;col[6].i = i; col[6].j = j; col[6].k = k+1;
629           ierr = MatSetValuesStencil(B,1,&row,7,col,v,INSERT_VALUES);CHKERRQ(ierr);
630         }
631       }
632     }
633   }
634   ierr   = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
635   ierr   = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
636   PetscFunctionReturn(0);
637 }