Actual source code: tsirm.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  1: /*
  2:     This file implements TSIRM, the Two-Stage Iteration with least-squares Residual Minimization method. 
  3:     It is an iterative method to solve large sparse linear systems of the form Ax=b, and it improves the convergence of Krylov based iterative methods.
  4:     The principle is to build an external iteration over a Krylov method (for example GMRES), and to frequently store its current residual in a matrix S. After a given number of outer iterations, a least-squares minimization step (with CGLS or LSQR) is applied on S, in order to compute a better solution and to make new iterations if required.
  5: */
  6: #include <petsc/private/kspimpl.h>        /*I "petscksp.h" I*/
  7: #include <petscksp.h>

  9: typedef struct {
 10:   PetscReal tol_ls;
 11:   PetscInt  size_ls,maxiter_ls,cgls,size,Istart,Iend;
 12:   Mat       A,S;
 13:   Vec       Alpha,r;
 14: } KSP_TSIRM;

 18: static PetscErrorCode KSPSetUp_TSIRM(KSP ksp)
 19: {
 21:   KSP_TSIRM      *tsirm = (KSP_TSIRM*)ksp->data;
 22: 
 24:   /* Initialization */
 25:   tsirm->tol_ls     = 1e-40;
 26:   tsirm->size_ls    = 12;
 27:   tsirm->maxiter_ls = 15;
 28:   tsirm->cgls       = 0;
 29: 
 30:   /* Matrix of the system */
 31:   KSPGetOperators(ksp,&tsirm->A,NULL);    /* Matrix of the system   */
 32:   MatGetSize(tsirm->A,&tsirm->size,NULL); /* Size of the system     */
 33:   MatGetOwnershipRange(tsirm->A,&tsirm->Istart,&tsirm->Iend);
 34: 
 35:   /* Matrix S of residuals */
 36:   MatCreate(PETSC_COMM_WORLD,&tsirm->S);
 37:   MatSetSizes(tsirm->S,tsirm->Iend-tsirm->Istart,PETSC_DECIDE,tsirm->size,tsirm->size_ls);
 38:   MatSetType(tsirm->S,MATMPIDENSE);
 39:   MatSetUp(tsirm->S);
 40: 
 41:   /* Residual and vector Alpha computed in the minimization step */
 42:   MatCreateVecs(tsirm->S,&tsirm->Alpha,&tsirm->r);
 43:   return(0);
 44: }

 48: PetscErrorCode KSPSolve_TSIRM(KSP ksp)
 49: {
 51:   KSP_TSIRM      *tsirm = (KSP_TSIRM*)ksp->data;
 52:   KSP            sub_ksp;
 53:   PC             pc;
 54:   Mat            AS;
 55:   Vec            x,b;
 56:   PetscScalar    *array;
 57:   PetscReal      norm = 20;
 58:   PetscInt       i,*ind_row,first_iteration = 1,its = 0,total = 0,col = 0;
 59:   PetscInt       iter_minimization = 0,restart = 30;
 60:   KSP            ksp_min;  /* KSP for minimization */
 61:   PC             pc_min;    /* PC for minimization */
 62: 
 64:   x = ksp->vec_sol; /* Solution vector        */
 65:   b = ksp->vec_rhs; /* Right-hand side vector */
 66: 
 67:   /* Row indexes (these indexes are global) */
 68:   PetscMalloc(tsirm->Iend-tsirm->Istart,&ind_row);
 69:   for (i=0;i<tsirm->Iend-tsirm->Istart;i++) ind_row[i] = i+tsirm->Istart;
 70: 
 71:   /* Inner solver */
 72:   KSPGetPC(ksp,&pc);
 73:   KSPSetFromOptions(ksp);
 74:   PCKSPGetKSP(pc,&sub_ksp);
 75:   KSPSetTolerances(sub_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,restart);
 76:   KSPSetFromOptions(sub_ksp);
 77: 
 78:   /* previously it seemed good but with SNES it seems not good... */
 79:   KSP_MatMult(sub_ksp,tsirm->A,x,tsirm->r);
 80:   VecAXPY(tsirm->r,-1,b);
 81:   VecNorm(tsirm->r,NORM_2,&norm);
 82:   ksp->its = 0;
 83:   KSPConvergedDefault(ksp,ksp->its,norm,&ksp->reason,ksp->cnvP);
 84:   KSPSetInitialGuessNonzero(sub_ksp,PETSC_TRUE);
 85:   do {
 86:     for (col=0;col<tsirm->size_ls && ksp->reason==0;col++) {
 87:       /* Solve (inner ietration) */
 88:       KSPSolve(sub_ksp,b,x);
 89:       KSPGetIterationNumber(sub_ksp,&its);
 90:       total += its;
 91: 
 92:       /* Build S^T */
 93:       VecGetArray(x,&array);
 94:       MatSetValues(tsirm->S,tsirm->Iend-tsirm->Istart,ind_row,1,&col,array,INSERT_VALUES);
 95:       VecRestoreArray(x,&array);
 96: 
 97:       KSPGetResidualNorm(sub_ksp,&norm);
 98:       ksp->rnorm = norm;
 99:       ksp->its ++;
100:       KSPConvergedDefault(ksp,ksp->its,norm,&ksp->reason,ksp->cnvP);
101:       KSPMonitor(ksp,ksp->its,norm);
102:     }
103: 
104:     /* Minimization step */
105:     if (!ksp->reason) {
106:       MatAssemblyBegin(tsirm->S,MAT_FINAL_ASSEMBLY);
107:       MatAssemblyEnd(tsirm->S,MAT_FINAL_ASSEMBLY);
108:       if (first_iteration) {
109:         MatMatMult(tsirm->A,tsirm->S,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AS);
110:         first_iteration = 0;
111:       } else {
112:         MatMatMult(tsirm->A,tsirm->S,MAT_REUSE_MATRIX,PETSC_DEFAULT,&AS);
113:       }
114: 
115:       /* CGLS or LSQR method to minimize the residuals*/
116: 
117:       KSPCreate(PETSC_COMM_WORLD,&ksp_min);
118:       if (tsirm->cgls) {
119:         KSPSetType(ksp_min,KSPCGLS);
120:       } else {
121:         KSPSetType(ksp_min,KSPLSQR);
122:       }
123:       KSPSetOperators(ksp_min,AS,AS);
124:       KSPSetTolerances(ksp_min,tsirm->tol_ls,PETSC_DEFAULT,PETSC_DEFAULT,tsirm->maxiter_ls);
125:       KSPGetPC(ksp_min,&pc_min);
126:       PCSetType(pc_min,PCNONE);
127:       KSPSolve(ksp_min,b,tsirm->Alpha);    /* Find Alpha such that ||AS Alpha = b|| */
128:       KSPGetIterationNumber(ksp_min,&iter_minimization);
129:       KSPGetResidualNorm(ksp_min,&norm);
130:       KSPMonitor(ksp_min,iter_minimization,norm);

132:       /* Apply minimization */
133:       MatMult(tsirm->S,tsirm->Alpha,x); /* x = S * Alpha */
134:     }
135:   } while (ksp->its<ksp->max_it && !ksp->reason);
136:   PetscFree(ind_row);
137:   ksp->its = total;
138:   return(0);
139: }

143: PetscErrorCode KSPSetFromOptions_TSIRM(PetscOptionItems *PetscOptionsObject,KSP ksp)
144: {
146:   KSP_TSIRM      *tsirm = (KSP_TSIRM*)ksp->data;

149:   PetscOptionsHead(PetscOptionsObject,"KSP TSIRM options");
150:   PetscOptionsInt("-ksp_tsirm_cgls","Method used for the minimization step","",tsirm->cgls,&tsirm->cgls,NULL); /*0:LSQR, 1:CGLS*/
151:   PetscOptionsReal("-ksp_tsirm_tol_ls","Tolerance threshold for the minimization step","",tsirm->tol_ls,&tsirm->tol_ls,NULL);
152:   PetscOptionsInt("-ksp_tsirm_max_it_ls","Maximum number of iterations for the minimization step","",tsirm->maxiter_ls,&tsirm->maxiter_ls,NULL);
153:   PetscOptionsInt("-ksp_tsirm_size_ls","Number of residuals for minimization","",tsirm->size_ls,&tsirm->size_ls,NULL);
154:   PetscOptionsTail();
155:   return(0);
156: }

160: PetscErrorCode KSPDestroy_TSIRM(KSP ksp)
161: {
162:   KSP_TSIRM       *tsirm = (KSP_TSIRM*)ksp->data;

166:   /* Free work vectors */
167:   MatDestroy(&tsirm->S);
168:   VecDestroy(&tsirm->Alpha);
169:   VecDestroy(&tsirm->r);
170:   return(0);
171: }

175: PETSC_EXTERN PetscErrorCode KSPCreate_TSIRM(KSP ksp)
176: {
178:   KSP_TSIRM      *tsirm;
179: 
181:   PetscNewLog(ksp,&tsirm);
182:   ksp->data                = (void*)tsirm;
183:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
184:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,1);
185:   ksp->ops->setup          = KSPSetUp_TSIRM;
186:   ksp->ops->solve          = KSPSolve_TSIRM;
187:   ksp->ops->destroy        = KSPDestroy_TSIRM;
188:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
189:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
190:   ksp->ops->setfromoptions = KSPSetFromOptions_TSIRM;
191:   ksp->ops->view           = 0;
192: #if defined(PETSC_USE_COMPLEX)
193:   SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"This is not supported for complex numbers");
194: #endif
195:   return(0);
196: }