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
6 Laplacian in 3D. Modeled by the partial differential equation
8 - Laplacian u = 1,0 < x,y,z < 1,
10 with boundary conditions
12 u = 1 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1.
14 This uses multigrid to solve the linear system
16 See src/snes/examples/tutorials/ex50.c
18 Can also be run with -pc_type exotic -ksp_type fgmres
22 static char help[] = "Solves 3D Laplacian using multigrid.\n\n";
26 #include <petscdmda.h>
28 extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,void*);
29 extern PetscErrorCode ComputeRHS(KSP,Vec,void*);
30 extern PetscErrorCode ComputeInitialGuess(KSP,Vec,void*);
33 #define __FUNCT__ "main"
41 int KrylovMinimize(Mat A, Vec b, Vec x) {
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;
60 Vec Alpha, p, ss, vect, r, q, Ax;
65 ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
66 ierr = KSPSetOperators(ksp,A,A);CHKERRQ(ierr);
67 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
74 ierr = PetscPrintf(PETSC_COMM_WORLD, "Size of vector %D\n", size); CHKERRQ(ierr);
77 MatGetOwnershipRange(A,&aa,&bb);
79 // ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%D %D\n", aa,bb); CHKERRQ(ierr);
80 //PetscSynchronizedFlush(PETSC_COMM_WORLD);
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);
88 ierr = MatGetOwnershipRange(S, &Istart, &Iend); CHKERRQ(ierr);
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);
100 ierr = VecDuplicate(b,&x_old);CHKERRQ(ierr);
101 ierr = VecDuplicate(b,&residu);CHKERRQ(ierr);
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;
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);
116 //GMRES WITH MINIMIZATION
118 while(giter<Emaxiter && norm>Eprecision ){
119 for(col=0; col<ColS && norm>Eprecision; col++){
123 ierr = KSPSolve(ksp, b, x); CHKERRQ(ierr);
125 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
126 total += its; Iiter ++;
131 ierr = VecGetArray(x, &array);
132 ierr = MatSetValues(S, Iend-Istart, ind_row, 1, &col, array, INSERT_VALUES);
133 VecRestoreArray(x, &array);
137 KSPGetResidualNorm(ksp,&norm);
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);
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);
153 if( norm>Eprecision) {
155 MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY);
156 MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY);
161 MatMatMult(A,S, MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AS);
165 MatMatMult(A,S, MAT_REUSE_MATRIX,PETSC_DEFAULT,&AS);
170 //Minimization with CGLS
171 MatMult(AS, Alpha, r); VecAYPX(r, -1, b); //r_0 = b-AS*x_0
174 MatMultTranspose(AS, r, p); //p_0 = AS'*r_0
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
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;
197 MatMult(S, Alpha, x); //x = S*Alpha;
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);
217 int KrylovMinimizeLSQR(Mat A, Vec b, Vec x) {
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;
230 PetscInt Istart,Iend;
237 PetscScalar norm_ksp;
238 Vec u,v,d,uu,vt,zero_long,zero_short,x_lsqr;
242 ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
243 ierr = KSPSetOperators(ksp,A,A);CHKERRQ(ierr);
244 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
251 ierr = PetscPrintf(PETSC_COMM_WORLD, "Size of vector %D\n", size); CHKERRQ(ierr);
254 MatGetOwnershipRange(A,&aa,&bb);
256 // ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%D %D\n", aa,bb); CHKERRQ(ierr);
257 //PetscSynchronizedFlush(PETSC_COMM_WORLD);
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);
265 ierr = MatGetOwnershipRange(S, &Istart, &Iend); CHKERRQ(ierr);
269 ierr = VecDuplicate(b, &Ax); CHKERRQ(ierr);
271 ierr = VecDuplicate(b,&x_old);CHKERRQ(ierr);
272 ierr = VecDuplicate(b,&residu);CHKERRQ(ierr);
276 ierr = VecDuplicate(b,&u);CHKERRQ(ierr);
279 ierr = VecDuplicate(b,&uu);CHKERRQ(ierr);
280 ierr = VecDuplicate(b,&zero_long);CHKERRQ(ierr);
281 ierr = VecSet(zero_long,0);CHKERRQ(ierr);
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);
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;
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);
306 //GMRES WITH MINIMIZATION
308 while(giter<Emaxiter && norm>Eprecision ){
309 for(col=0; col<ColS && norm>Eprecision; col++){
313 ierr = KSPSolve(ksp, b, x); CHKERRQ(ierr);
315 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
316 total += its; Iiter ++;
321 ierr = VecGetArray(x, &array);
322 ierr = MatSetValues(S, Iend-Istart, ind_row, 1, &col, array, INSERT_VALUES);
323 VecRestoreArray(x, &array);
326 KSPGetResidualNorm(ksp,&norm);
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);
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);
344 if( norm>Eprecision) {
346 MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY);
347 MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY);
354 MatMatMult(A,S, MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AS);
359 MatMatMult(A,S, MAT_REUSE_MATRIX,PETSC_DEFAULT,&AS);
371 PetscScalar n2b,tolb,normr,c,s,phibar,normar,norma,thet,rhot,rho,phi;
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)
379 VecAYPX(u,1/beta,zero_long); // u = u / 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)
388 VecAYPX(v,1/alpha,zero_short); // v = v / alpha;
390 ierr = VecSet(d,0);CHKERRQ(ierr);
391 normar = alpha * beta;
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]);
402 rho = sqrt(rhot*rhot + beta*beta);
406 if (phi == 0) { // stagnation of the method
410 VecAYPX(d,-thet,v); //d = (v - thet * d);
411 VecAYPX(d,1/rho,zero_short); //d=d/ rho;
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);
427 MatMult(S, x_lsqr, x); //x = S*x_small;
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);
445 int main(int argc,char **argv)
456 PetscInitialize(&argc,&argv,(char*)0,help);
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);
466 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
469 ierr = KSPSetTolerances(ksp, 1e-10, 1e-10, PETSC_DEFAULT, 50000000); CHKERRQ(ierr);
471 ierr = KSPSolve(ksp,NULL,NULL);CHKERRQ(ierr);
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);
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);
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);
492 // ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
500 VecDuplicate(b,&sol);
502 KrylovMinimize(A, b, x);
510 VecNorm(sol, NORM_2, &norm);
511 ierr = PetscPrintf(PETSC_COMM_WORLD, "Error2 %g\n",norm);
513 ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
523 VecDuplicate(b,&sol);
525 KrylovMinimizeLSQR(A, b, x);
533 VecNorm(sol, NORM_2, &norm);
534 ierr = PetscPrintf(PETSC_COMM_WORLD, "Error LSQR %g\n",norm);
537 ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
541 ierr = VecDestroy(&r);CHKERRQ(ierr);
543 ierr = PetscFinalize();
549 #define __FUNCT__ "ComputeRHS"
550 PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
553 PetscInt i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs;
555 PetscScalar Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy;
556 PetscScalar ***barray;
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);
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);
572 barray[k][j][i] = Hx*Hy*Hz;
577 ierr = DMDAVecRestoreArray(dm,b,&barray);CHKERRQ(ierr);
578 PetscFunctionReturn(0);
582 #define __FUNCT__ "ComputeInitialGuess"
583 PetscErrorCode ComputeInitialGuess(KSP ksp,Vec b,void *ctx)
587 PetscFunctionBeginUser;
588 ierr = VecSet(b,0);CHKERRQ(ierr);
589 PetscFunctionReturn(0);
593 #define __FUNCT__ "ComputeMatrix"
594 PetscErrorCode ComputeMatrix(KSP ksp,Mat jac,Mat B,void *ctx)
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];
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);
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);
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);
629 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
630 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
631 PetscFunctionReturn(0);