Actual source code: tsirm.c

  1: #include <petsc/private/kspimpl.h>

  3: typedef struct {
  4:   PetscReal tol_ls;
  5:   PetscInt  size_ls, maxiter_ls, cgls, size, Istart, Iend;
  6:   Mat       A, S;
  7:   Vec       Alpha, r;
  8: } KSP_TSIRM;

 10: static PetscErrorCode KSPSetUp_TSIRM(KSP ksp)
 11: {
 12:   KSP_TSIRM *tsirm = (KSP_TSIRM *)ksp->data;

 14:   PetscFunctionBegin;
 15:   /* Initialization */
 16: #if defined(PETSC_USE_REAL_SINGLE)
 17:   tsirm->tol_ls = 1e-25;
 18: #else
 19:   tsirm->tol_ls = 1e-50;
 20: #endif
 21:   tsirm->size_ls    = 12;
 22:   tsirm->maxiter_ls = 15;
 23:   tsirm->cgls       = 0;

 25:   /* Matrix of the system */
 26:   PetscCall(KSPGetOperators(ksp, &tsirm->A, NULL));    /* Matrix of the system   */
 27:   PetscCall(MatGetSize(tsirm->A, &tsirm->size, NULL)); /* Size of the system     */
 28:   PetscCall(MatGetOwnershipRange(tsirm->A, &tsirm->Istart, &tsirm->Iend));

 30:   /* Matrix S of residuals */
 31:   PetscCall(MatCreate(PETSC_COMM_WORLD, &tsirm->S));
 32:   PetscCall(MatSetSizes(tsirm->S, tsirm->Iend - tsirm->Istart, PETSC_DECIDE, tsirm->size, tsirm->size_ls));
 33:   PetscCall(MatSetType(tsirm->S, MATDENSE));
 34:   PetscCall(MatSetUp(tsirm->S));

 36:   /* Residual and vector Alpha computed in the minimization step */
 37:   PetscCall(MatCreateVecs(tsirm->S, &tsirm->Alpha, &tsirm->r));
 38:   PetscFunctionReturn(PETSC_SUCCESS);
 39: }

 41: static PetscErrorCode KSPSolve_TSIRM(KSP ksp)
 42: {
 43:   KSP_TSIRM   *tsirm = (KSP_TSIRM *)ksp->data;
 44:   KSP          sub_ksp;
 45:   PC           pc;
 46:   Mat          AS = NULL;
 47:   Vec          x, b;
 48:   PetscScalar *array;
 49:   PetscReal    norm = 20;
 50:   PetscInt     i, *ind_row, first_iteration = 1, its = 0, total = 0, col = 0;
 51:   PetscInt     restart = 30;
 52:   KSP          ksp_min; /* KSP for minimization */
 53:   PC           pc_min;  /* PC for minimization */
 54:   PetscBool    isksp;

 56:   PetscFunctionBegin;
 57:   x = ksp->vec_sol; /* Solution vector        */
 58:   b = ksp->vec_rhs; /* Right-hand side vector */

 60:   /* Row indexes (these indexes are global) */
 61:   PetscCall(PetscMalloc1(tsirm->Iend - tsirm->Istart, &ind_row));
 62:   for (i = 0; i < tsirm->Iend - tsirm->Istart; i++) ind_row[i] = i + tsirm->Istart;

 64:   /* Inner solver */
 65:   PetscCall(KSPGetPC(ksp, &pc));
 66:   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCKSP, &isksp));
 67:   PetscCheck(isksp, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "PC must be of type PCKSP");
 68:   PetscCall(PCKSPGetKSP(pc, &sub_ksp));
 69:   PetscCall(KSPSetTolerances(sub_ksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, restart));

 71:   /* previously it seemed good but with SNES it seems not good... */
 72:   PetscCall(KSP_MatMult(sub_ksp, tsirm->A, x, tsirm->r));
 73:   PetscCall(VecAXPY(tsirm->r, -1, b));
 74:   PetscCall(VecNorm(tsirm->r, NORM_2, &norm));
 75:   KSPCheckNorm(ksp, norm);
 76:   ksp->its = 0;
 77:   PetscCall(KSPConvergedDefault(ksp, ksp->its, norm, &ksp->reason, ksp->cnvP));
 78:   PetscCall(KSPSetInitialGuessNonzero(sub_ksp, PETSC_TRUE));
 79:   do {
 80:     for (col = 0; col < tsirm->size_ls && ksp->reason == 0; col++) {
 81:       /* Solve (inner iteration) */
 82:       PetscCall(KSPSolve(sub_ksp, b, x));
 83:       PetscCall(KSPGetIterationNumber(sub_ksp, &its));
 84:       total += its;

 86:       /* Build S^T */
 87:       PetscCall(VecGetArray(x, &array));
 88:       PetscCall(MatSetValues(tsirm->S, tsirm->Iend - tsirm->Istart, ind_row, 1, &col, array, INSERT_VALUES));
 89:       PetscCall(VecRestoreArray(x, &array));

 91:       PetscCall(KSPGetResidualNorm(sub_ksp, &norm));
 92:       ksp->rnorm = norm;
 93:       ksp->its++;
 94:       PetscCall(KSPConvergedDefault(ksp, ksp->its, norm, &ksp->reason, ksp->cnvP));
 95:       PetscCall(KSPMonitor(ksp, ksp->its, norm));
 96:     }

 98:     /* Minimization step */
 99:     if (!ksp->reason) {
100:       PetscCall(MatAssemblyBegin(tsirm->S, MAT_FINAL_ASSEMBLY));
101:       PetscCall(MatAssemblyEnd(tsirm->S, MAT_FINAL_ASSEMBLY));
102:       if (first_iteration) {
103:         PetscCall(MatMatMult(tsirm->A, tsirm->S, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &AS));
104:         first_iteration = 0;
105:       } else {
106:         PetscCall(MatMatMult(tsirm->A, tsirm->S, MAT_REUSE_MATRIX, PETSC_DEFAULT, &AS));
107:       }

109:       /* CGLS or LSQR method to minimize the residuals*/
110:       PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp_min));
111:       if (tsirm->cgls) {
112:         PetscCall(KSPSetType(ksp_min, KSPCGLS));
113:       } else {
114:         PetscCall(KSPSetType(ksp_min, KSPLSQR));
115:       }
116:       PetscCall(KSPSetOperators(ksp_min, AS, AS));
117:       PetscCall(KSPSetTolerances(ksp_min, tsirm->tol_ls, PETSC_DEFAULT, PETSC_DEFAULT, tsirm->maxiter_ls));
118:       PetscCall(KSPGetPC(ksp_min, &pc_min));
119:       PetscCall(PCSetType(pc_min, PCNONE));
120:       PetscCall(KSPSolve(ksp_min, b, tsirm->Alpha)); /* Find Alpha such that ||AS Alpha = b|| */
121:       PetscCall(KSPDestroy(&ksp_min));
122:       /* Apply minimization */
123:       PetscCall(MatMult(tsirm->S, tsirm->Alpha, x)); /* x = S * Alpha */
124:     }
125:   } while (ksp->its < ksp->max_it && !ksp->reason);
126:   PetscCall(MatDestroy(&AS));
127:   PetscCall(PetscFree(ind_row));
128:   ksp->its = total;
129:   PetscFunctionReturn(PETSC_SUCCESS);
130: }

132: static PetscErrorCode KSPSetFromOptions_TSIRM(KSP ksp, PetscOptionItems *PetscOptionsObject)
133: {
134:   KSP_TSIRM *tsirm = (KSP_TSIRM *)ksp->data;

136:   PetscFunctionBegin;
137:   PetscOptionsHeadBegin(PetscOptionsObject, "KSP TSIRM options");
138:   PetscCall(PetscOptionsInt("-ksp_tsirm_cgls", "Method used for the minimization step", "", tsirm->cgls, &tsirm->cgls, NULL)); /*0:LSQR, 1:CGLS*/
139:   PetscCall(PetscOptionsReal("-ksp_tsirm_tol_ls", "Tolerance threshold for the minimization step", "", tsirm->tol_ls, &tsirm->tol_ls, NULL));
140:   PetscCall(PetscOptionsInt("-ksp_tsirm_max_it_ls", "Maximum number of iterations for the minimization step", "", tsirm->maxiter_ls, &tsirm->maxiter_ls, NULL));
141:   PetscCall(PetscOptionsInt("-ksp_tsirm_size_ls", "Number of residuals for minimization", "", tsirm->size_ls, &tsirm->size_ls, NULL));
142:   PetscOptionsHeadEnd();
143:   PetscFunctionReturn(PETSC_SUCCESS);
144: }

146: static PetscErrorCode KSPDestroy_TSIRM(KSP ksp)
147: {
148:   KSP_TSIRM *tsirm = (KSP_TSIRM *)ksp->data;

150:   PetscFunctionBegin;
151:   PetscCall(MatDestroy(&tsirm->S));
152:   PetscCall(VecDestroy(&tsirm->Alpha));
153:   PetscCall(VecDestroy(&tsirm->r));
154:   PetscCall(PetscFree(ksp->data));
155:   PetscFunctionReturn(PETSC_SUCCESS);
156: }

158: /*MC
159:    KSPTSIRM - Implements the two-stage iteration with least-squares residual minimization method {cite}`couturier2016tsirm`

161:    Options Database Keys:
162: +  -ksp_ksp_type <solver>        - the type of the inner solver (GMRES or any of its variants for instance)
163: .  -ksp_pc_type <preconditioner> - the type of the preconditioner applied to the inner solver
164: .  -ksp_ksp_max_it <maxits>      - the maximum number of inner iterations (iterations of the inner solver)
165: .  -ksp_ksp_rtol <tol>           - sets the relative convergence tolerance of the inner solver
166: .  -ksp_tsirm_cgls <number>      - if 1 use CGLS solver in the minimization step, otherwise use LSQR solver
167: .  -ksp_tsirm_max_it_ls <maxits> - the maximum number of iterations for the least-squares minimization solver
168: .  -ksp_tsirm_tol_ls <tol>       - sets the convergence tolerance of the least-squares minimization solver
169: -  -ksp_tsirm_size_ls <size>     - the number of residuals for the least-squares minimization step

171:    Level: advanced

173:    Notes:
174:    `KSPTSIRM` is a two-stage iteration method for solving large sparse linear systems of the form $Ax=b$. The main idea behind this new
175:    method is the use a least-squares residual minimization to improve the convergence of Krylov based iterative methods, typically those of GMRES variants.
176:    The principle of TSIRM algorithm  is to build an outer iteration over a Krylov method, called the inner solver, and to frequently store the current residual
177:    computed by the given Krylov method in a matrix of residuals S. After a few outer iterations, a least-squares minimization step is applied on the matrix
178:    composed by the saved residuals, in order to compute a better solution and to make new iterations if required.
179:    The minimization step consists in solving the least-squares problem $\min||b-ASa||$ to find 'a' which minimizes the
180:    residuals $(b-AS)$. The minimization step is performed using two solvers of linear least-squares problems: `KSPCGLS` or `KSPLSQR`. A new solution x with
181:    a minimal residual is computed with $x=Sa$.

183:    Contributed by:
184:    Lilia Ziane Khodja

186: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPFGMRES`, `KSPLGMRES`,
187:           `KSPGMRESSetRestart()`, `KSPGMRESSetHapTol()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
188:           `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`,
189:           `KSPGMRESCGSRefinementType`, `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESMonitorKrylov()`, `KSPSetPCSide()`
190: M*/
191: PETSC_EXTERN PetscErrorCode KSPCreate_TSIRM(KSP ksp)
192: {
193:   KSP_TSIRM *tsirm;

195:   PetscFunctionBegin;
196:   PetscCall(PetscNew(&tsirm));
197:   ksp->data = (void *)tsirm;
198:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 2));
199:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 1));
200:   ksp->ops->setup          = KSPSetUp_TSIRM;
201:   ksp->ops->solve          = KSPSolve_TSIRM;
202:   ksp->ops->destroy        = KSPDestroy_TSIRM;
203:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
204:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
205:   ksp->ops->setfromoptions = KSPSetFromOptions_TSIRM;
206:   ksp->ops->view           = NULL;
207: #if defined(PETSC_USE_COMPLEX)
208:   SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "This is not supported for complex numbers");
209: #else
210:   PetscFunctionReturn(PETSC_SUCCESS);
211: #endif
212: }