]> AND Private Git Repository - GMRES2stage.git/commitdiff
Logo AND Algorithmique Numérique Distribuée

Private GIT Repository
ajout ex45 LSQR pas rapide...
authorraphael couturier <couturie@extinction>
Wed, 13 Aug 2014 14:34:14 +0000 (16:34 +0200)
committerraphael couturier <couturie@extinction>
Wed, 13 Aug 2014 14:34:14 +0000 (16:34 +0200)
code/ex45.c [new file with mode: 0644]

diff --git a/code/ex45.c b/code/ex45.c
new file mode 100644 (file)
index 0000000..a44a2d0
--- /dev/null
@@ -0,0 +1,637 @@
+
+// /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
+
+
+/*
+Laplacian in 3D. Modeled by the partial differential equation
+
+   - Laplacian u = 1,0 < x,y,z < 1,
+
+with boundary conditions
+
+   u = 1 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1.
+
+   This uses multigrid to solve the linear system
+
+   See src/snes/examples/tutorials/ex50.c
+
+   Can also be run with -pc_type exotic -ksp_type fgmres
+
+*/
+
+static char help[] = "Solves 3D Laplacian using multigrid.\n\n";
+
+#include <petscksp.h>
+#include <petscdm.h>
+#include <petscdmda.h>
+
+extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,void*);
+extern PetscErrorCode ComputeRHS(KSP,Vec,void*);
+extern PetscErrorCode ComputeInitialGuess(KSP,Vec,void*);
+
+#undef __FUNCT__
+#define __FUNCT__ "main"
+
+
+
+
+
+
+
+int KrylovMinimize(Mat A, Vec b, Vec x) {
+  
+
+  //Variables
+
+  PetscScalar  gamma, alpha, oldgamma, beta;
+  PetscReal norm=20, Eprecision=5e-5, cgprec=1e-40;     
+  PetscInt giter=0, ColS=12, col=0, Emaxiter=50000000, iter=0, iterations=15, Iiter=0;
+  PetscErrorCode ierr;
+  PetscScalar T1, T2;
+  KSP ksp;
+  PetscInt total=0;  
+  PetscInt size;
+  PetscInt Istart,Iend;
+  PetscInt i,its;
+  Vec       x_old, residu;
+  Mat S, AS;
+  PetscScalar *array;
+  PetscInt *ind_row;
+  Vec Alpha, p, ss, vect, r, q, Ax;
+
+
+  PetscInt first=1;
+
+  ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
+  ierr = KSPSetOperators(ksp,A,A);CHKERRQ(ierr);
+  ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
+
+
+
+
+  VecGetSize(b,&size);
+
+  ierr = PetscPrintf(PETSC_COMM_WORLD, "Size of vector %D\n", size); CHKERRQ(ierr);  
+
+  PetscInt aa,bb;
+  MatGetOwnershipRange(A,&aa,&bb);
+
+  //  ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%D %D\n", aa,bb); CHKERRQ(ierr);  
+  //PetscSynchronizedFlush(PETSC_COMM_WORLD);
+
+
+  ierr = MatCreate(PETSC_COMM_WORLD, &S);CHKERRQ(ierr);
+  ierr = MatSetSizes(S, bb-aa,  PETSC_DECIDE, size, ColS); CHKERRQ(ierr);
+  ierr = MatSetType(S, MATMPIDENSE); CHKERRQ(ierr);
+  ierr = MatSetUp(S); CHKERRQ(ierr);
+
+  ierr = MatGetOwnershipRange(S, &Istart, &Iend); CHKERRQ(ierr);
+
+  ierr = VecCreate(PETSC_COMM_WORLD, &Alpha); CHKERRQ(ierr);
+  ierr = VecSetSizes(Alpha, PETSC_DECIDE, ColS); CHKERRQ(ierr);
+  ierr = VecSetFromOptions(Alpha); CHKERRQ(ierr);
+  ierr = VecDuplicate(Alpha, &vect); CHKERRQ(ierr);
+  ierr = VecDuplicate(Alpha, &p); CHKERRQ(ierr);
+  ierr = VecDuplicate(Alpha, &ss); CHKERRQ(ierr);
+  ierr = VecDuplicate(b, &r); CHKERRQ(ierr);
+  ierr = VecDuplicate(b, &q); CHKERRQ(ierr);
+  ierr = VecDuplicate(b, &Ax); CHKERRQ(ierr);
+
+  ierr = VecDuplicate(b,&x_old);CHKERRQ(ierr);
+  ierr = VecDuplicate(b,&residu);CHKERRQ(ierr);
+
+
+
+  //indexes of row (these indexes are global)
+  ind_row = (PetscInt*)malloc(sizeof(PetscInt)*(Iend-Istart));
+  for(i=0; i<Iend-Istart; i++) ind_row[i] = i + Istart;
+
+  //Initializations
+  //  ierr = KSPGMRESSetRestart(ksp, 16); CHKERRQ(ierr);
+  ierr = KSPSetTolerances(ksp, 1e-10, 1e-10, PETSC_DEFAULT, 16); CHKERRQ(ierr);
+  ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE); CHKERRQ(ierr);
+
+
+
+  //GMRES WITH MINIMIZATION
+  T1 = MPI_Wtime();
+  while(giter<Emaxiter && norm>Eprecision ){
+    for(col=0; col<ColS  &&  norm>Eprecision; col++){
+
+
+      //Solve
+      ierr = KSPSolve(ksp, b, x); CHKERRQ(ierr);
+
+      ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
+      total += its; Iiter ++;
+
+
+
+      //Build S'
+      ierr = VecGetArray(x, &array);
+      ierr = MatSetValues(S, Iend-Istart, ind_row, 1, &col, array, INSERT_VALUES);
+      VecRestoreArray(x, &array);
+
+
+
+      //Error
+      ierr = VecCopy(x, residu); CHKERRQ(ierr);
+      ierr = VecAXPY(residu, -1, x_old); CHKERRQ(ierr);
+      ierr = VecNorm(residu, NORM_INFINITY, &norm); CHKERRQ(ierr);
+
+
+
+      ierr = PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g, outer iteration %D\n", norm, giter); CHKERRQ(ierr);
+      ierr = VecCopy(x, x_old); CHKERRQ(ierr);
+
+
+    }
+
+
+    //minimization step
+    if( norm>Eprecision) {
+
+      MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY);
+      MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY);
+
+
+      //Build AS
+      if(first) {
+        MatMatMult(A,S, MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AS);
+        first=0;
+      }
+      else
+        MatMatMult(A,S, MAT_REUSE_MATRIX,PETSC_DEFAULT,&AS);
+
+
+
+
+      //Minimization with CGLS 
+      MatMult(AS, Alpha, r);  VecAYPX(r, -1, b); //r_0 = b-AS*x_0
+
+
+      MatMultTranspose(AS, r, p); //p_0 = AS'*r_0 
+
+
+
+
+      ierr = VecCopy(p, ss); CHKERRQ(ierr); //p_0 = s_0
+      VecNorm(ss, NORM_2, &gamma); gamma = gamma * gamma; //gamma = norm2(s)^2
+      while(gamma>cgprec && iter<iterations){
+        MatMult(AS, p, q); //q = AS*p
+        VecNorm(q, NORM_2, &alpha); alpha = alpha * alpha; //alpha = norm2(q)^2
+        alpha = gamma / alpha; 
+        VecAXPY(Alpha, alpha, p); //Alpha += alpha*p
+        VecAXPY(r, -alpha, q); //r -= alpha*q
+        MatMultTranspose(AS, r, ss); // ss = AS'*r
+        oldgamma = gamma;
+        VecNorm(ss, NORM_2, &gamma); gamma = gamma * gamma; //gamma = norm2(s)^2
+        beta = gamma / oldgamma;
+        VecAYPX(p, beta, ss); //p = s+beta*p;
+        iter ++;
+      }        
+
+      iter = 0; giter ++;
+      //Minimizer
+      MatMult(S, Alpha, x); //x = S*Alpha;
+    }
+
+  }
+  T2 = MPI_Wtime();
+  ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Execution time : %g (s)\n", T2-T1); CHKERRQ(ierr);
+  ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Total number of iterations : %D\n", total); CHKERRQ(ierr);
+
+  return 0;
+
+}
+
+
+
+
+
+
+
+
+
+int KrylovMinimizeLSQR(Mat A, Vec b, Vec x) {
+  
+
+  //Variables
+
+  PetscScalar  alpha, beta;
+  PetscReal norm=20, Eprecision=5e-5, tol=1e-40;     
+  PetscInt giter=0, ColS=12, col=0, Emaxiter=50000000, iter=0, iterations=15, Iiter=0;
+  PetscErrorCode ierr;
+  PetscScalar T1, T2;
+  KSP ksp;
+  PetscInt total=0;  
+  PetscInt size;
+  PetscInt Istart,Iend;
+  PetscInt i,its;
+  Vec       x_old, residu;
+  Mat S, AS;
+  PetscScalar *array;
+  PetscInt *ind_row;
+  Vec  Ax;
+  PetscScalar norm_ksp;
+  Vec u,v,d,uu,vt,zero_long,zero_short,x_lsqr;
+
+  PetscInt first=1;
+
+  ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
+  ierr = KSPSetOperators(ksp,A,A);CHKERRQ(ierr);
+  ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
+
+
+
+
+  VecGetSize(b,&size);
+
+  ierr = PetscPrintf(PETSC_COMM_WORLD, "Size of vector %D\n", size); CHKERRQ(ierr);  
+
+  PetscInt aa,bb;
+  MatGetOwnershipRange(A,&aa,&bb);
+
+  //  ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%D %D\n", aa,bb); CHKERRQ(ierr);  
+  //PetscSynchronizedFlush(PETSC_COMM_WORLD);
+
+
+  ierr = MatCreate(PETSC_COMM_WORLD, &S);CHKERRQ(ierr);
+  ierr = MatSetSizes(S, bb-aa,  PETSC_DECIDE, size, ColS); CHKERRQ(ierr);
+  ierr = MatSetType(S, MATMPIDENSE); CHKERRQ(ierr);
+  ierr = MatSetUp(S); CHKERRQ(ierr);
+
+  ierr = MatGetOwnershipRange(S, &Istart, &Iend); CHKERRQ(ierr);
+
+
+
+  ierr = VecDuplicate(b, &Ax); CHKERRQ(ierr);
+
+  ierr = VecDuplicate(b,&x_old);CHKERRQ(ierr);
+  ierr = VecDuplicate(b,&residu);CHKERRQ(ierr);
+
+
+  //long vector
+  ierr = VecDuplicate(b,&u);CHKERRQ(ierr);
+
+
+  ierr = VecDuplicate(b,&uu);CHKERRQ(ierr);
+  ierr = VecDuplicate(b,&zero_long);CHKERRQ(ierr);
+  ierr = VecSet(zero_long,0);CHKERRQ(ierr);
+
+  //small vector
+  ierr = VecCreate(PETSC_COMM_WORLD, &v); CHKERRQ(ierr);
+  ierr = VecSetSizes(v, PETSC_DECIDE, ColS); CHKERRQ(ierr);
+  ierr = VecSetFromOptions(v); CHKERRQ(ierr);
+  ierr = VecDuplicate(v,&zero_short);CHKERRQ(ierr);
+  ierr = VecSet(zero_short,0);CHKERRQ(ierr);
+  ierr = VecDuplicate(v,&d);CHKERRQ(ierr);
+  ierr = VecDuplicate(v,&vt);CHKERRQ(ierr);
+  ierr = VecDuplicate(v,&x_lsqr);CHKERRQ(ierr);
+
+
+  //indexes of row (these indexes are global)
+  ind_row = (PetscInt*)malloc(sizeof(PetscInt)*(Iend-Istart));
+  for(i=0; i<Iend-Istart; i++) ind_row[i] = i + Istart;
+
+  //Initializations
+  //  ierr = KSPGMRESSetRestart(ksp, 16); CHKERRQ(ierr);
+  ierr = KSPSetTolerances(ksp, 1e-10, 1e-10, PETSC_DEFAULT, 16); CHKERRQ(ierr);
+  ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE); CHKERRQ(ierr);
+
+
+
+
+  //GMRES WITH MINIMIZATION
+  T1 = MPI_Wtime();
+  while(giter<Emaxiter && norm>Eprecision ){
+    for(col=0; col<ColS  &&  norm>Eprecision; col++){
+
+
+      //Solve
+      ierr = KSPSolve(ksp, b, x); CHKERRQ(ierr);
+
+      ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
+      total += its; Iiter ++;
+
+
+
+      //Build S'
+      ierr = VecGetArray(x, &array);
+      ierr = MatSetValues(S, Iend-Istart, ind_row, 1, &col, array, INSERT_VALUES);
+      VecRestoreArray(x, &array);
+
+
+
+      //Error
+      ierr = VecCopy(x, residu); CHKERRQ(ierr);
+      ierr = VecAXPY(residu, -1, x_old); CHKERRQ(ierr);
+      ierr = VecNorm(residu, NORM_INFINITY, &norm); CHKERRQ(ierr);
+
+
+
+      ierr = PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g, outer iteration %D\n", norm, giter); CHKERRQ(ierr);
+      ierr = VecCopy(x, x_old); CHKERRQ(ierr);
+
+
+    }
+
+
+    //minimization step
+    if( norm>Eprecision) {
+
+      MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY);
+      MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY);
+
+
+
+
+      //Build AS
+      if(first) {
+        MatMatMult(A,S, MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AS);
+        
+        first=0;
+      }
+      else
+        MatMatMult(A,S, MAT_REUSE_MATRIX,PETSC_DEFAULT,&AS);
+
+
+
+
+
+      //LSQR
+      //LSQR
+      //LSQR
+
+
+
+      PetscScalar n2b,tolb,normr,c,s,phibar,normar,norma,thet,rhot,rho,phi;
+      PetscInt stag;
+      tolb = tol * n2b;
+      VecNorm(b, NORM_2, &n2b); //n2b = norm(b);
+      ierr = VecCopy(b, u); CHKERRQ(ierr);   //u=b
+      VecNorm(u, NORM_2, &beta); // beta=norm(u)
+      normr=beta;
+      if (beta != 0) {
+        VecAYPX(u,1/beta,zero_long); //  u = u / beta;
+      }
+      c=1;
+      s=0;
+      phibar=beta;
+      MatMultTranspose(AS, u, v); //v=A'*u
+      ierr = VecSet(x_lsqr,0);CHKERRQ(ierr);
+      VecNorm(v, NORM_2, &alpha); // alpha=norm(v)
+      if (alpha != 0) {
+        VecAYPX(v,1/alpha,zero_short); //  v = v / alpha;
+      }
+      ierr = VecSet(d,0);CHKERRQ(ierr);
+      normar = alpha * beta;
+      norma=0;
+      //stag=0;
+      for(i=0;i<iterations;i++) {
+        MatMult(AS, v, uu); //uu=A*v
+        VecAYPX(u, -alpha, uu); //u = uu-alpha*u;
+        VecNorm(u, NORM_2, &beta); // beta=norm(u)
+        VecAYPX(u,1/beta,zero_long); //  u = u / beta;
+        norma=sqrt(norma*norma+alpha*alpha+beta*beta); // norma = norm([norma alpha beta]);
+        thet = - s * alpha;
+        rhot = c * alpha;
+        rho = sqrt(rhot*rhot + beta*beta);
+        c = rhot / rho;
+        s = - beta / rho;
+        phi = c * phibar;
+        if (phi == 0) {             // stagnation of the method
+          stag = 1;
+        }
+        phibar = s * phibar;
+        VecAYPX(d,-thet,v);       //d = (v - thet * d);
+        VecAYPX(d,1/rho,zero_short);     //d=d/ rho;
+
+
+        if (normar/(norma*normr) <= tol) { // check for convergence in min{|b-A*x|}
+          break;
+        }
+        if (normr <= tolb) {           // check for convergence in A*x=b
+          break;
+        }
+
+
+        VecAXPY(x_lsqr,phi,d);     // x_lsqr=x_lsqr+phi*d
+        normr = abs(s) * normr;
+        MatMultTranspose(AS, u, vt); //vt=A'*u;
+        VecAYPX(v,-beta,vt);  // v = vt - beta * v;
+        VecNorm(v, NORM_2, &alpha); // alpha=norm(v)
+        VecAYPX(v,1/alpha,zero_short); //  v = v / alpha;
+        normar = alpha * abs( s * phi);
+      }        
+      
+
+
+      iter = 0; giter ++;
+      //Minimizer
+      MatMult(S, x_lsqr, x); //x = S*x_small;
+    }
+
+  }
+  T2 = MPI_Wtime();
+  ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Execution time LSQR : %g (s)\n", T2-T1); CHKERRQ(ierr);
+  ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Total number of iterations LSQR : %D\n", total); CHKERRQ(ierr);
+
+  return 0;
+
+}
+
+
+
+
+
+
+
+int main(int argc,char **argv)
+{
+  PetscErrorCode ierr;
+  KSP            ksp;
+  PetscReal      norm;
+  DM             da;
+  Vec            x,b,r;
+  Mat            A;
+  PetscInt total;
+  PetscScalar T1,T2;
+
+  PetscInitialize(&argc,&argv,(char*)0,help);
+
+  ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
+  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);
+  ierr = KSPSetDM(ksp,da);CHKERRQ(ierr);
+  ierr = KSPSetComputeInitialGuess(ksp,ComputeInitialGuess,NULL);CHKERRQ(ierr);
+  ierr = KSPSetComputeRHS(ksp,ComputeRHS,NULL);CHKERRQ(ierr);
+  ierr = KSPSetComputeOperators(ksp,ComputeMatrix,NULL);CHKERRQ(ierr);
+  ierr = DMDestroy(&da);CHKERRQ(ierr);
+
+  ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
+
+
+  ierr = KSPSetTolerances(ksp, 1e-7, 1e-7, PETSC_DEFAULT, 50000000); CHKERRQ(ierr);
+  T1 = MPI_Wtime();
+  ierr = KSPSolve(ksp,NULL,NULL);CHKERRQ(ierr);
+  T2 = MPI_Wtime();
+
+
+
+  ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
+  ierr = KSPGetRhs(ksp,&b);CHKERRQ(ierr);
+  ierr = VecDuplicate(b,&r);CHKERRQ(ierr);
+  ierr = KSPGetOperators(ksp,&A,NULL);CHKERRQ(ierr);
+
+  ierr = MatMult(A,x,r);CHKERRQ(ierr);
+  ierr = VecAXPY(r,-1.0,b);CHKERRQ(ierr);
+  ierr = VecNorm(r,NORM_2,&norm);CHKERRQ(ierr);
+  ierr = PetscPrintf(PETSC_COMM_WORLD,"Residual norm %g\n",(double)norm);CHKERRQ(ierr);
+  ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Execution time : %g (s)\n", T2-T1); CHKERRQ(ierr);
+
+  ierr = KSPGetIterationNumber(ksp, &total); CHKERRQ(ierr);
+  ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Total number of iterations : %D\n", total); CHKERRQ(ierr);
+
+
+
+  //  ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
+
+
+{
+
+    Vec x;
+    Vec sol;
+    VecDuplicate(b,&x);
+    VecDuplicate(b,&sol);
+    
+    KrylovMinimize(A, b, x);
+
+
+
+    MatMult(A,x,sol);
+
+
+    VecAXPY(sol,-1,b);
+    VecNorm(sol, NORM_2, &norm);
+    ierr = PetscPrintf(PETSC_COMM_WORLD, "Error2 %g\n",norm);
+
+    ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
+
+  }
+
+
+ {
+
+    Vec x;
+    Vec sol;
+    VecDuplicate(b,&x);
+    VecDuplicate(b,&sol);
+    
+    KrylovMinimizeLSQR(A, b, x);
+
+
+
+    MatMult(A,x,sol);
+
+
+    VecAXPY(sol,-1,b);
+    VecNorm(sol, NORM_2, &norm);
+    ierr = PetscPrintf(PETSC_COMM_WORLD, "Error LSQR %g\n",norm);
+
+
+    ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
+  }
+
+  ierr = VecDestroy(&r);CHKERRQ(ierr);
+
+  ierr = PetscFinalize();
+
+  return 0;
+}
+
+#undef __FUNCT__
+#define __FUNCT__ "ComputeRHS"
+PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
+{
+  PetscErrorCode ierr;
+  PetscInt       i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs;
+  DM             dm;
+  PetscScalar    Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy;
+  PetscScalar    ***barray;
+
+  PetscFunctionBeginUser;
+  ierr    = KSPGetDM(ksp,&dm);CHKERRQ(ierr);
+  ierr    = DMDAGetInfo(dm,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
+  Hx      = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1);
+  HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx;
+  ierr    = DMDAGetCorners(dm,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
+  ierr    = DMDAVecGetArray(dm,b,&barray);CHKERRQ(ierr);
+
+  for (k=zs; k<zs+zm; k++) {
+    for (j=ys; j<ys+ym; j++) {
+      for (i=xs; i<xs+xm; i++) {
+        if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1) {
+          barray[k][j][i] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);
+        } else {
+          barray[k][j][i] = Hx*Hy*Hz;
+        }
+      }
+    }
+  }
+  ierr = DMDAVecRestoreArray(dm,b,&barray);CHKERRQ(ierr);
+  PetscFunctionReturn(0);
+}
+
+#undef __FUNCT__
+#define __FUNCT__ "ComputeInitialGuess"
+PetscErrorCode ComputeInitialGuess(KSP ksp,Vec b,void *ctx)
+{
+  PetscErrorCode ierr;
+
+  PetscFunctionBeginUser;
+  ierr = VecSet(b,0);CHKERRQ(ierr);
+  PetscFunctionReturn(0);
+}
+
+#undef __FUNCT__
+#define __FUNCT__ "ComputeMatrix"
+PetscErrorCode ComputeMatrix(KSP ksp,Mat jac,Mat B,void *ctx)
+{
+  DM             da;
+  PetscErrorCode ierr;
+  PetscInt       i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs;
+  PetscScalar    v[7],Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy;
+  MatStencil     row,col[7];
+
+  PetscFunctionBeginUser;
+  ierr    = KSPGetDM(ksp,&da);CHKERRQ(ierr);
+  ierr    = DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
+  Hx      = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1);
+  HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx;
+  ierr    = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
+
+  for (k=zs; k<zs+zm; k++) {
+    for (j=ys; j<ys+ym; j++) {
+      for (i=xs; i<xs+xm; i++) {
+        row.i = i; row.j = j; row.k = k;
+        if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1) {
+          v[0] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);
+          ierr = MatSetValuesStencil(B,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr);
+        } else {
+          v[0] = -HxHydHz;col[0].i = i; col[0].j = j; col[0].k = k-1;
+          v[1] = -HxHzdHy;col[1].i = i; col[1].j = j-1; col[1].k = k;
+          v[2] = -HyHzdHx;col[2].i = i-1; col[2].j = j; col[2].k = k;
+          v[3] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);col[3].i = row.i; col[3].j = row.j; col[3].k = row.k;
+          v[4] = -HyHzdHx;col[4].i = i+1; col[4].j = j; col[4].k = k;
+          v[5] = -HxHzdHy;col[5].i = i; col[5].j = j+1; col[5].k = k;
+          v[6] = -HxHydHz;col[6].i = i; col[6].j = j; col[6].k = k+1;
+          ierr = MatSetValuesStencil(B,1,&row,7,col,v,INSERT_VALUES);CHKERRQ(ierr);
+        }
+      }
+    }
+  }
+  ierr   = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
+  ierr   = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
+  PetscFunctionReturn(0);
+}