From: raphael couturier Date: Wed, 13 Aug 2014 12:35:10 +0000 (+0200) Subject: ajout ex29, lsqr ne semble pas bon sur cet exemple... X-Git-Url: https://bilbo.iut-bm.univ-fcomte.fr/and/gitweb/GMRES2stage.git/commitdiff_plain/bc9c5a57e4b0b66714c6ac200b4aff7c26b38608 ajout ex29, lsqr ne semble pas bon sur cet exemple... --- diff --git a/code/ex15.c b/code/ex15.c index 20144dd..0e0a8c4 100644 --- a/code/ex15.c +++ b/code/ex15.c @@ -1,4 +1,8 @@ +// /home/couturie/work/petsc-3.5.1/arch-linux2-c-debug/bin/mpirun -np 4 ./ex15 -m 400 -n 400 + + + static char help[] = "Solves a linear system in parallel with KSP. Also\n\ illustrates setting a user-defined shell preconditioner and using the\n\ macro __FUNCT__ to define routine names for use in error handling.\n\ @@ -51,11 +55,11 @@ int KrylovMinimize(Mat A, Vec b, Vec x) { //Variables - PetscScalar gamma, alpha, oldgamma, beta, t2; + PetscScalar gamma, alpha, oldgamma, beta; PetscReal norm=20, Eprecision=1e-8, cgprec=1e-40; PetscInt giter=0, ColS=12, col=0, Emaxiter=50000000, iter=0, iterations=15, Iiter=0; PetscErrorCode ierr; - PetscScalar T1, T2, t1; + PetscScalar T1, T2; KSP ksp; PetscInt total=0; PetscInt size; @@ -66,7 +70,6 @@ int KrylovMinimize(Mat A, Vec b, Vec x) { PetscScalar *array; PetscInt *ind_row; Vec Alpha, p, ss, vect, r, q, Ax; - PetscScalar norm_ksp; PetscInt first=1; @@ -207,7 +210,6 @@ int KrylovMinimize(Mat A, Vec b, Vec x) { } T2 = MPI_Wtime(); - t1 = T2 - T1; 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); @@ -222,6 +224,242 @@ int KrylovMinimize(Mat A, Vec b, Vec x) { + +int KrylovMinimizeLSQR(Mat A, Vec b, Vec x) { + + + //Variables + + PetscScalar alpha, beta; + PetscReal norm=20, Eprecision=1e-8, 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; iEprecision ){ + for(col=0; colEprecision; 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 +#include +#include + +extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,void*); +extern PetscErrorCode ComputeRHS(KSP,Vec,void*); + +typedef enum {DIRICHLET, NEUMANN} BCType; + +typedef struct { + PetscReal rho; + PetscReal nu; + BCType bcType; +} UserContext; + +#undef __FUNCT__ +#define __FUNCT__ "main" + + + + + + + + + +int KrylovMinimize(Mat A, Vec b, Vec x) { + + + //Variables + + PetscScalar gamma, alpha, oldgamma, beta; + PetscReal norm=20, Eprecision=1e-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; iEprecision ){ + for(col=0; colEprecision; 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 && iterEprecision ){ + for(col=0; colEprecision; 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;inu)*PetscExpScalar(-((PetscReal)j*Hy)*((PetscReal)j*Hy)/user->nu)*Hx*Hy; + } + } + ierr = DMDAVecRestoreArray(da, b, &array);CHKERRQ(ierr); + ierr = VecAssemblyBegin(b);CHKERRQ(ierr); + ierr = VecAssemblyEnd(b);CHKERRQ(ierr); + + /* force right hand side to be consistent for singular matrix */ + /* note this is really a hack, normally the model would provide you with a consistent right handside */ + if (user->bcType == NEUMANN) { + MatNullSpace nullspace; + + ierr = MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);CHKERRQ(ierr); + ierr = MatNullSpaceRemove(nullspace,b);CHKERRQ(ierr); + ierr = MatNullSpaceDestroy(&nullspace);CHKERRQ(ierr); + } + PetscFunctionReturn(0); +} + + +#undef __FUNCT__ +#define __FUNCT__ "ComputeRho" +PetscErrorCode ComputeRho(PetscInt i, PetscInt j, PetscInt mx, PetscInt my, PetscReal centerRho, PetscReal *rho) +{ + PetscFunctionBeginUser; + if ((i > mx/3.0) && (i < 2.0*mx/3.0) && (j > my/3.0) && (j < 2.0*my/3.0)) { + *rho = centerRho; + } else { + *rho = 1.0; + } + PetscFunctionReturn(0); +} + +#undef __FUNCT__ +#define __FUNCT__ "ComputeMatrix" +PetscErrorCode ComputeMatrix(KSP ksp,Mat J,Mat jac,void *ctx) +{ + UserContext *user = (UserContext*)ctx; + PetscReal centerRho; + PetscErrorCode ierr; + PetscInt i,j,mx,my,xm,ym,xs,ys; + PetscScalar v[5]; + PetscReal Hx,Hy,HydHx,HxdHy,rho; + MatStencil row, col[5]; + DM da; + + PetscFunctionBeginUser; + ierr = KSPGetDM(ksp,&da);CHKERRQ(ierr); + centerRho = user->rho; + ierr = DMDAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); + Hx = 1.0 / (PetscReal)(mx-1); + Hy = 1.0 / (PetscReal)(my-1); + HxdHy = Hx/Hy; + HydHx = Hy/Hx; + ierr = DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr); + for (j=ys; jbcType == DIRICHLET) { + v[0] = 2.0*rho*(HxdHy + HydHx); + ierr = MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr); + } else if (user->bcType == NEUMANN) { + PetscInt numx = 0, numy = 0, num = 0; + if (j!=0) { + v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j-1; + numy++; num++; + } + if (i!=0) { + v[num] = -rho*HydHx; col[num].i = i-1; col[num].j = j; + numx++; num++; + } + if (i!=mx-1) { + v[num] = -rho*HydHx; col[num].i = i+1; col[num].j = j; + numx++; num++; + } + if (j!=my-1) { + v[num] = -rho*HxdHy; col[num].i = i; col[num].j = j+1; + numy++; num++; + } + v[num] = numx*rho*HydHx + numy*rho*HxdHy; col[num].i = i; col[num].j = j; + num++; + ierr = MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);CHKERRQ(ierr); + } + } else { + v[0] = -rho*HxdHy; col[0].i = i; col[0].j = j-1; + v[1] = -rho*HydHx; col[1].i = i-1; col[1].j = j; + v[2] = 2.0*rho*(HxdHy + HydHx); col[2].i = i; col[2].j = j; + v[3] = -rho*HydHx; col[3].i = i+1; col[3].j = j; + v[4] = -rho*HxdHy; col[4].i = i; col[4].j = j+1; + ierr = MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);CHKERRQ(ierr); + } + } + } + ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); + ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); + if (user->bcType == NEUMANN) { + MatNullSpace nullspace; + + ierr = MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);CHKERRQ(ierr); + ierr = MatSetNullSpace(jac,nullspace);CHKERRQ(ierr); + ierr = MatNullSpaceDestroy(&nullspace);CHKERRQ(ierr); + } + PetscFunctionReturn(0); +}