3 // /home/couturie/work/petsc-3.5.1/arch-linux2-c-debug/bin/mpirun -np 3 ex29 -da_grid_x 600 -da_grid_y 600 -ksp_type fgmres
8 Concepts: KSP^solving a system of linear equations
9 Concepts: KSP^Laplacian, 2d
14 Added at the request of Marc Garbey.
16 Inhomogeneous Laplacian in 2D. Modeled by the partial differential equation
18 -div \rho grad u = f, 0 < x,y < 1,
22 f = e^{-x^2/\nu} e^{-y^2/\nu}
24 with Dirichlet boundary conditions
26 u = f(x,y) for x = 0, x = 1, y = 0, y = 1
28 or pure Neumman boundary conditions
30 This uses multigrid to solve the linear system
33 static char help[] = "Solves 2D inhomogeneous Laplacian using multigrid.\n\n";
36 #include <petscdmda.h>
39 extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,void*);
40 extern PetscErrorCode ComputeRHS(KSP,Vec,void*);
42 typedef enum {DIRICHLET, NEUMANN} BCType;
51 #define __FUNCT__ "main"
61 int KrylovMinimize(Mat A, Vec b, Vec x) {
66 PetscScalar gamma, alpha, oldgamma, beta;
67 PetscReal norm=20, Eprecision=1e-6, cgprec=1e-40;
68 PetscInt giter=0, ColS=8, col=0, Emaxiter=50000000, iter=0, iterations=15, Iiter=0;
80 Vec Alpha, p, ss, vect, r, q, Ax;
85 ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
86 ierr = KSPSetOperators(ksp,A,A);CHKERRQ(ierr);
87 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
94 ierr = PetscPrintf(PETSC_COMM_WORLD, "Size of vector %D\n", size); CHKERRQ(ierr);
97 MatGetOwnershipRange(A,&aa,&bb);
99 // ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%D %D\n", aa,bb); CHKERRQ(ierr);
100 //PetscSynchronizedFlush(PETSC_COMM_WORLD);
103 ierr = MatCreate(PETSC_COMM_WORLD, &S);CHKERRQ(ierr);
104 ierr = MatSetSizes(S, bb-aa, PETSC_DECIDE, size, ColS); CHKERRQ(ierr);
105 ierr = MatSetType(S, MATMPIDENSE); CHKERRQ(ierr);
106 ierr = MatSetUp(S); CHKERRQ(ierr);
108 ierr = MatGetOwnershipRange(S, &Istart, &Iend); CHKERRQ(ierr);
110 ierr = VecCreate(PETSC_COMM_WORLD, &Alpha); CHKERRQ(ierr);
111 ierr = VecSetSizes(Alpha, PETSC_DECIDE, ColS); CHKERRQ(ierr);
112 ierr = VecSetFromOptions(Alpha); CHKERRQ(ierr);
113 ierr = VecDuplicate(Alpha, &vect); CHKERRQ(ierr);
114 ierr = VecDuplicate(Alpha, &p); CHKERRQ(ierr);
115 ierr = VecDuplicate(Alpha, &ss); CHKERRQ(ierr);
116 ierr = VecDuplicate(b, &r); CHKERRQ(ierr);
117 ierr = VecDuplicate(b, &q); CHKERRQ(ierr);
118 ierr = VecDuplicate(b, &Ax); CHKERRQ(ierr);
120 ierr = VecDuplicate(b,&x_old);CHKERRQ(ierr);
121 ierr = VecDuplicate(b,&residu);CHKERRQ(ierr);
125 //indexes of row (these indexes are global)
126 ind_row = (PetscInt*)malloc(sizeof(PetscInt)*(Iend-Istart));
127 for(i=0; i<Iend-Istart; i++) ind_row[i] = i + Istart;
130 // ierr = KSPGMRESSetRestart(ksp, 16); CHKERRQ(ierr);
131 ierr = KSPSetTolerances(ksp, 1e-13, 1e-13, PETSC_DEFAULT, 16); CHKERRQ(ierr);
132 ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE); CHKERRQ(ierr);
136 //GMRES WITH MINIMIZATION
138 while(giter<Emaxiter && norm>Eprecision ){
139 for(col=0; col<ColS && norm>Eprecision; col++){
143 ierr = KSPSolve(ksp, b, x); CHKERRQ(ierr);
145 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
146 total += its; Iiter ++;
151 ierr = VecGetArray(x, &array);
152 ierr = MatSetValues(S, Iend-Istart, ind_row, 1, &col, array, INSERT_VALUES);
153 VecRestoreArray(x, &array);
157 KSPGetResidualNorm(ksp,&norm);
160 ierr = VecCopy(x, residu); CHKERRQ(ierr);
161 ierr = VecAXPY(residu, -1, x_old); CHKERRQ(ierr);
162 ierr = VecNorm(residu, NORM_INFINITY, &norm); CHKERRQ(ierr);
166 ierr = PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g, outer iteration %D\n", norm, giter); CHKERRQ(ierr);
167 ierr = VecCopy(x, x_old); CHKERRQ(ierr);
174 if( norm>Eprecision) {
176 MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY);
177 MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY);
182 MatMatMult(A,S, MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AS);
186 MatMatMult(A,S, MAT_REUSE_MATRIX,PETSC_DEFAULT,&AS);
191 //Minimization with CGLS
192 MatMult(AS, Alpha, r); VecAYPX(r, -1, b); //r_0 = b-AS*x_0
195 MatMultTranspose(AS, r, p); //p_0 = AS'*r_0
200 ierr = VecCopy(p, ss); CHKERRQ(ierr); //p_0 = s_0
201 VecNorm(ss, NORM_2, &gamma); gamma = gamma * gamma; //gamma = norm2(s)^2
202 while(gamma>cgprec && iter<iterations){
203 MatMult(AS, p, q); //q = AS*p
204 VecNorm(q, NORM_2, &alpha); alpha = alpha * alpha; //alpha = norm2(q)^2
205 alpha = gamma / alpha;
206 VecAXPY(Alpha, alpha, p); //Alpha += alpha*p
207 VecAXPY(r, -alpha, q); //r -= alpha*q
208 MatMultTranspose(AS, r, ss); // ss = AS'*r
210 VecNorm(ss, NORM_2, &gamma); gamma = gamma * gamma; //gamma = norm2(s)^2
211 beta = gamma / oldgamma;
212 VecAYPX(p, beta, ss); //p = s+beta*p;
218 MatMult(S, Alpha, x); //x = S*Alpha;
223 ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Execution time : %g (s)\n", T2-T1); CHKERRQ(ierr);
224 ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Total number of iterations : %D\n", total); CHKERRQ(ierr);
238 int KrylovMinimizeLSQR(Mat A, Vec b, Vec x) {
243 PetscScalar alpha, beta;
244 PetscReal norm=20, Eprecision=1e-6, tol=1e-40;
245 PetscInt giter=0, ColS=8, col=0, Emaxiter=50000000, iter=0, iterations=15, Iiter=0;
251 PetscInt Istart,Iend;
258 PetscScalar norm_ksp;
259 Vec u,v,d,uu,vt,zero_long,zero_short,x_lsqr;
263 ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
264 ierr = KSPSetOperators(ksp,A,A);CHKERRQ(ierr);
265 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
272 ierr = PetscPrintf(PETSC_COMM_WORLD, "Size of vector %D\n", size); CHKERRQ(ierr);
275 MatGetOwnershipRange(A,&aa,&bb);
277 // ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%D %D\n", aa,bb); CHKERRQ(ierr);
278 //PetscSynchronizedFlush(PETSC_COMM_WORLD);
281 ierr = MatCreate(PETSC_COMM_WORLD, &S);CHKERRQ(ierr);
282 ierr = MatSetSizes(S, bb-aa, PETSC_DECIDE, size, ColS); CHKERRQ(ierr);
283 ierr = MatSetType(S, MATMPIDENSE); CHKERRQ(ierr);
284 ierr = MatSetUp(S); CHKERRQ(ierr);
286 ierr = MatGetOwnershipRange(S, &Istart, &Iend); CHKERRQ(ierr);
290 ierr = VecDuplicate(b, &Ax); CHKERRQ(ierr);
292 ierr = VecDuplicate(b,&x_old);CHKERRQ(ierr);
293 ierr = VecDuplicate(b,&residu);CHKERRQ(ierr);
297 ierr = VecDuplicate(b,&u);CHKERRQ(ierr);
300 ierr = VecDuplicate(b,&uu);CHKERRQ(ierr);
301 ierr = VecDuplicate(b,&zero_long);CHKERRQ(ierr);
302 ierr = VecSet(zero_long,0);CHKERRQ(ierr);
305 ierr = VecCreate(PETSC_COMM_WORLD, &v); CHKERRQ(ierr);
306 ierr = VecSetSizes(v, PETSC_DECIDE, ColS); CHKERRQ(ierr);
307 ierr = VecSetFromOptions(v); CHKERRQ(ierr);
308 ierr = VecDuplicate(v,&zero_short);CHKERRQ(ierr);
309 ierr = VecSet(zero_short,0);CHKERRQ(ierr);
310 ierr = VecDuplicate(v,&d);CHKERRQ(ierr);
311 ierr = VecDuplicate(v,&vt);CHKERRQ(ierr);
312 ierr = VecDuplicate(v,&x_lsqr);CHKERRQ(ierr);
315 //indexes of row (these indexes are global)
316 ind_row = (PetscInt*)malloc(sizeof(PetscInt)*(Iend-Istart));
317 for(i=0; i<Iend-Istart; i++) ind_row[i] = i + Istart;
320 // ierr = KSPGMRESSetRestart(ksp, 16); CHKERRQ(ierr);
321 ierr = KSPSetTolerances(ksp, 1e-13, 1e-13, PETSC_DEFAULT, 16); CHKERRQ(ierr);
322 ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE); CHKERRQ(ierr);
327 //GMRES WITH MINIMIZATION
329 while(giter<Emaxiter && norm>Eprecision ){
330 for(col=0; col<ColS && norm>Eprecision; col++){
334 ierr = KSPSolve(ksp, b, x); CHKERRQ(ierr);
336 ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
337 total += its; Iiter ++;
342 ierr = VecGetArray(x, &array);
343 ierr = MatSetValues(S, Iend-Istart, ind_row, 1, &col, array, INSERT_VALUES);
344 VecRestoreArray(x, &array);
348 KSPGetResidualNorm(ksp,&norm);
351 ierr = VecCopy(x, residu); CHKERRQ(ierr);
352 ierr = VecAXPY(residu, -1, x_old); CHKERRQ(ierr);
353 ierr = VecNorm(residu, NORM_INFINITY, &norm); CHKERRQ(ierr);
357 ierr = PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g, outer iteration %D\n", norm, giter); CHKERRQ(ierr);
358 ierr = VecCopy(x, x_old); CHKERRQ(ierr);
365 if( norm>Eprecision) {
367 MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY);
368 MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY);
375 MatMatMult(A,S, MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AS);
380 MatMatMult(A,S, MAT_REUSE_MATRIX,PETSC_DEFAULT,&AS);
392 PetscScalar n2b,tolb,normr,c,s,phibar,normar,norma,thet,rhot,rho,phi;
395 VecNorm(b, NORM_2, &n2b); //n2b = norm(b);
396 ierr = VecCopy(b, u); CHKERRQ(ierr); //u=b
397 VecNorm(u, NORM_2, &beta); // beta=norm(u)
400 VecAYPX(u,1/beta,zero_long); // u = u / beta;
405 MatMultTranspose(AS, u, v); //v=A'*u
406 ierr = VecSet(x_lsqr,0);CHKERRQ(ierr);
407 VecNorm(v, NORM_2, &alpha); // alpha=norm(v)
409 VecAYPX(v,1/alpha,zero_short); // v = v / alpha;
411 ierr = VecSet(d,0);CHKERRQ(ierr);
412 normar = alpha * beta;
415 for(i=0;i<iterations;i++) {
416 MatMult(AS, v, uu); //uu=A*v
417 VecAYPX(u, -alpha, uu); //u = uu-alpha*u;
418 VecNorm(u, NORM_2, &beta); // beta=norm(u)
419 VecAYPX(u,1/beta,zero_long); // u = u / beta;
420 norma=sqrt(norma*norma+alpha*alpha+beta*beta); // norma = norm([norma alpha beta]);
423 rho = sqrt(rhot*rhot + beta*beta);
427 if (phi == 0) { // stagnation of the method
431 VecAYPX(d,-thet,v); //d = (v - thet * d);
432 VecAYPX(d,1/rho,zero_short); //d=d/ rho;
436 VecAXPY(x_lsqr,phi,d); // x_lsqr=x_lsqr+phi*d
437 normr = abs(s) * normr;
438 MatMultTranspose(AS, u, vt); //vt=A'*u;
439 VecAYPX(v,-beta,vt); // v = vt - beta * v;
440 VecNorm(v, NORM_2, &alpha); // alpha=norm(v)
441 VecAYPX(v,1/alpha,zero_short); // v = v / alpha;
442 normar = alpha * abs( s * phi);
449 MatMult(S, x_lsqr, x); //x = S*x_small;
454 ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Execution time LSQR : %g (s)\n", T2-T1); CHKERRQ(ierr);
455 ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Total number of iterations LSQR : %D\n", total); CHKERRQ(ierr);
467 int main(int argc,char **argv)
472 const char *bcTypes[2] = {"dirichlet","neumann"};
477 PetscScalar T1,T2,norm;
481 PetscInitialize(&argc,&argv,(char*)0,help);
485 MPI_Comm_size(PETSC_COMM_WORLD,&size);
486 PetscPrintf(PETSC_COMM_WORLD,"Number of processors = %d\n",size);
488 ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
489 ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,-3,-3,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da);CHKERRQ(ierr);
490 ierr = DMDASetUniformCoordinates(da,0,1,0,1,0,0);CHKERRQ(ierr);
491 ierr = DMDASetFieldName(da,0,"Pressure");CHKERRQ(ierr);
493 ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for the inhomogeneous Poisson equation", "DMqq");
495 ierr = PetscOptionsReal("-rho", "The conductivity", "ex29.c", user.rho, &user.rho, NULL);CHKERRQ(ierr);
497 ierr = PetscOptionsReal("-nu", "The width of the Gaussian source", "ex29.c", user.nu, &user.nu, NULL);CHKERRQ(ierr);
498 bc = (PetscInt)DIRICHLET;
499 ierr = PetscOptionsEList("-bc_type","Type of boundary condition","ex29.c",bcTypes,2,bcTypes[0],&bc,NULL);CHKERRQ(ierr);
500 user.bcType = (BCType)bc;
501 ierr = PetscOptionsEnd();
503 ierr = KSPSetComputeRHS(ksp,ComputeRHS,&user);CHKERRQ(ierr);
504 ierr = KSPSetComputeOperators(ksp,ComputeMatrix,&user);CHKERRQ(ierr);
505 ierr = KSPSetDM(ksp,da);CHKERRQ(ierr);
506 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
507 ierr = KSPSetUp(ksp);CHKERRQ(ierr);
509 ierr = KSPSetTolerances(ksp, 1e-6, 1e-6, PETSC_DEFAULT, 50000000); CHKERRQ(ierr);
511 ierr = KSPSolve(ksp,NULL,NULL);CHKERRQ(ierr);
516 ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
517 ierr = KSPGetRhs(ksp,&b);CHKERRQ(ierr);
518 KSPGetOperators(ksp,&A,NULL);
520 VecDuplicate(x,&sol);
523 VecNorm(sol, NORM_2, &norm);
524 ierr = PetscPrintf(PETSC_COMM_WORLD, "Error %g\n",norm);
525 ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Execution time : %g (s)\n", T2-T1); CHKERRQ(ierr);
527 ierr = KSPGetIterationNumber(ksp, &total); CHKERRQ(ierr);
528 ierr = PetscPrintf(PETSC_COMM_WORLD, "\t\t\t -- Total number of iterations : %D\n", total); CHKERRQ(ierr);
535 KrylovMinimize(A, b, x);
543 VecNorm(sol, NORM_2, &norm);
544 ierr = PetscPrintf(PETSC_COMM_WORLD, "Error2 %g\n",norm);
554 KrylovMinimizeLSQR(A, b, x);
562 VecNorm(sol, NORM_2, &norm);
563 ierr = PetscPrintf(PETSC_COMM_WORLD, "Error LSQR %g\n",norm);
570 ierr = DMDestroy(&da);CHKERRQ(ierr);
571 ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
572 ierr = PetscFinalize();
578 #define __FUNCT__ "ComputeRHS"
579 PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
581 UserContext *user = (UserContext*)ctx;
583 PetscInt i,j,mx,my,xm,ym,xs,ys;
588 PetscFunctionBeginUser;
589 ierr = KSPGetDM(ksp,&da);CHKERRQ(ierr);
590 ierr = DMDAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
591 Hx = 1.0 / (PetscReal)(mx-1);
592 Hy = 1.0 / (PetscReal)(my-1);
593 ierr = DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr);
594 ierr = DMDAVecGetArray(da, b, &array);CHKERRQ(ierr);
595 for (j=ys; j<ys+ym; j++) {
596 for (i=xs; i<xs+xm; i++) {
597 array[j][i] = PetscExpScalar(-((PetscReal)i*Hx)*((PetscReal)i*Hx)/user->nu)*PetscExpScalar(-((PetscReal)j*Hy)*((PetscReal)j*Hy)/user->nu)*Hx*Hy;
600 ierr = DMDAVecRestoreArray(da, b, &array);CHKERRQ(ierr);
601 ierr = VecAssemblyBegin(b);CHKERRQ(ierr);
602 ierr = VecAssemblyEnd(b);CHKERRQ(ierr);
604 /* force right hand side to be consistent for singular matrix */
605 /* note this is really a hack, normally the model would provide you with a consistent right handside */
606 if (user->bcType == NEUMANN) {
607 MatNullSpace nullspace;
609 ierr = MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);CHKERRQ(ierr);
610 ierr = MatNullSpaceRemove(nullspace,b);CHKERRQ(ierr);
611 ierr = MatNullSpaceDestroy(&nullspace);CHKERRQ(ierr);
613 PetscFunctionReturn(0);
618 #define __FUNCT__ "ComputeRho"
619 PetscErrorCode ComputeRho(PetscInt i, PetscInt j, PetscInt mx, PetscInt my, PetscReal centerRho, PetscReal *rho)
621 PetscFunctionBeginUser;
622 if ((i > mx/3.0) && (i < 2.0*mx/3.0) && (j > my/3.0) && (j < 2.0*my/3.0)) {
627 PetscFunctionReturn(0);
631 #define __FUNCT__ "ComputeMatrix"
632 PetscErrorCode ComputeMatrix(KSP ksp,Mat J,Mat jac,void *ctx)
634 UserContext *user = (UserContext*)ctx;
637 PetscInt i,j,mx,my,xm,ym,xs,ys;
639 PetscReal Hx,Hy,HydHx,HxdHy,rho;
640 MatStencil row, col[5];
643 PetscFunctionBeginUser;
644 ierr = KSPGetDM(ksp,&da);CHKERRQ(ierr);
645 centerRho = user->rho;
646 ierr = DMDAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
647 Hx = 1.0 / (PetscReal)(mx-1);
648 Hy = 1.0 / (PetscReal)(my-1);
651 ierr = DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr);
652 for (j=ys; j<ys+ym; j++) {
653 for (i=xs; i<xs+xm; i++) {
654 row.i = i; row.j = j;
655 ierr = ComputeRho(i, j, mx, my, centerRho, &rho);CHKERRQ(ierr);
656 if (i==0 || j==0 || i==mx-1 || j==my-1) {
657 if (user->bcType == DIRICHLET) {
658 v[0] = 2.0*rho*(HxdHy + HydHx);
659 ierr = MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr);
660 } else if (user->bcType == NEUMANN) {
661 PetscInt numx = 0, numy = 0, num = 0;
663 v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j-1;
667 v[num] = -rho*HydHx; col[num].i = i-1; col[num].j = j;
671 v[num] = -rho*HydHx; col[num].i = i+1; col[num].j = j;
675 v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j+1;
678 v[num] = numx*rho*HydHx + numy*rho*HxdHy; col[num].i = i; col[num].j = j;
680 ierr = MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);CHKERRQ(ierr);
683 v[0] = -rho*HxdHy; col[0].i = i; col[0].j = j-1;
684 v[1] = -rho*HydHx; col[1].i = i-1; col[1].j = j;
685 v[2] = 2.0*rho*(HxdHy + HydHx); col[2].i = i; col[2].j = j;
686 v[3] = -rho*HydHx; col[3].i = i+1; col[3].j = j;
687 v[4] = -rho*HxdHy; col[4].i = i; col[4].j = j+1;
688 ierr = MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);CHKERRQ(ierr);
692 ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
693 ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
694 if (user->bcType == NEUMANN) {
695 MatNullSpace nullspace;
697 ierr = MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);CHKERRQ(ierr);
698 ierr = MatSetNullSpace(jac,nullspace);CHKERRQ(ierr);
699 ierr = MatNullSpaceDestroy(&nullspace);CHKERRQ(ierr);
701 PetscFunctionReturn(0);