Actual source code: hpddm.cxx

  1: #define HPDDM_MIXED_PRECISION 1
  2: #include <petsc/private/petschpddm.h>

  4: const char *const KSPHPDDMTypes[]          = {KSPGMRES, "bgmres", KSPCG, "bcg", "gcrodr", "bgcrodr", "bfbcg", KSPPREONLY};
  5: const char *const KSPHPDDMPrecisionTypes[] = {"HALF", "SINGLE", "DOUBLE", "QUADRUPLE", "KSPHPDDMPrecisionType", "KSP_HPDDM_PRECISION_", NULL};
  6: const char *const HPDDMOrthogonalization[] = {"cgs", "mgs"};
  7: const char *const HPDDMQR[]                = {"cholqr", "cgs", "mgs"};
  8: const char *const HPDDMVariant[]           = {"left", "right", "flexible"};
  9: const char *const HPDDMRecycleTarget[]     = {"SM", "LM", "SR", "LR", "SI", "LI"};
 10: const char *const HPDDMRecycleStrategy[]   = {"A", "B"};

 12: PetscBool  HPDDMCite       = PETSC_FALSE;
 13: const char HPDDMCitation[] = "@article{jolivet2020petsc,\n"
 14:                              "  Author = {Jolivet, Pierre and Roman, Jose E. and Zampini, Stefano},\n"
 15:                              "  Title = {{KSPHPDDM} and {PCHPDDM}: Extending {PETSc} with Robust Overlapping {Schwarz} Preconditioners and Advanced {Krylov} Methods},\n"
 16:                              "  Year = {2021},\n"
 17:                              "  Publisher = {Elsevier},\n"
 18:                              "  Journal = {Computer \\& Mathematics with Applications},\n"
 19:                              "  Volume = {84},\n"
 20:                              "  Pages = {277--295},\n"
 21:                              "  Url = {https://github.com/prj-/jolivet2020petsc}\n"
 22:                              "}\n";

 24: #if PetscDefined(HAVE_SLEPC) && PetscDefined(HAVE_DYNAMIC_LIBRARIES) && PetscDefined(USE_SHARED_LIBRARIES)
 25: static PetscBool loadedDL = PETSC_FALSE;
 26: #endif

 28: static PetscErrorCode KSPSetFromOptions_HPDDM(KSP ksp, PetscOptionItems *PetscOptionsObject)
 29: {
 30:   KSP_HPDDM  *data = (KSP_HPDDM *)ksp->data;
 31:   PetscInt    i, j;
 32:   PetscMPIInt size;

 34:   PetscFunctionBegin;
 35:   PetscOptionsHeadBegin(PetscOptionsObject, "KSPHPDDM options, cf. https://github.com/hpddm/hpddm");
 36:   i = (data->cntl[0] == static_cast<char>(PETSC_DECIDE) ? HPDDM_KRYLOV_METHOD_GMRES : data->cntl[0]);
 37:   PetscCall(PetscOptionsEList("-ksp_hpddm_type", "Type of Krylov method", "KSPHPDDMGetType", KSPHPDDMTypes, PETSC_STATIC_ARRAY_LENGTH(KSPHPDDMTypes), KSPHPDDMTypes[HPDDM_KRYLOV_METHOD_GMRES], &i, NULL));
 38:   if (i == PETSC_STATIC_ARRAY_LENGTH(KSPHPDDMTypes) - 1) i = HPDDM_KRYLOV_METHOD_NONE; /* need to shift the value since HPDDM_KRYLOV_METHOD_RICHARDSON is not registered in PETSc */
 39:   data->cntl[0] = i;
 40:   PetscCall(PetscOptionsEnum("-ksp_hpddm_precision", "Precision in which Krylov bases are stored", "KSPHPDDM", KSPHPDDMPrecisionTypes, (PetscEnum)data->precision, (PetscEnum *)&data->precision, NULL));
 41:   PetscCheck(data->precision != KSP_HPDDM_PRECISION_QUADRUPLE || PetscDefined(HAVE_REAL___FLOAT128), PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP_SYS, "Unsupported %s precision", KSPHPDDMPrecisionTypes[data->precision]);
 42:   PetscCheck(std::abs(data->precision - PETSC_KSPHPDDM_DEFAULT_PRECISION) <= 1, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Unhandled mixed %s and %s precisions", KSPHPDDMPrecisionTypes[data->precision], KSPHPDDMPrecisionTypes[PETSC_KSPHPDDM_DEFAULT_PRECISION]);
 43:   if (data->cntl[0] != HPDDM_KRYLOV_METHOD_NONE) {
 44:     if (data->cntl[0] != HPDDM_KRYLOV_METHOD_BCG && data->cntl[0] != HPDDM_KRYLOV_METHOD_BFBCG) {
 45:       i = (data->cntl[1] == static_cast<char>(PETSC_DECIDE) ? HPDDM_VARIANT_LEFT : data->cntl[1]);
 46:       if (ksp->pc_side_set == PC_SIDE_DEFAULT)
 47:         PetscCall(PetscOptionsEList("-ksp_hpddm_variant", "Left, right, or variable preconditioning", "KSPHPDDM", HPDDMVariant, PETSC_STATIC_ARRAY_LENGTH(HPDDMVariant), HPDDMVariant[HPDDM_VARIANT_LEFT], &i, NULL));
 48:       else if (ksp->pc_side_set == PC_RIGHT) i = HPDDM_VARIANT_RIGHT;
 49:       data->cntl[1] = i;
 50:       if (i > 0) PetscCall(KSPSetPCSide(ksp, PC_RIGHT));
 51:     }
 52:     if (data->cntl[0] == HPDDM_KRYLOV_METHOD_BGMRES || data->cntl[0] == HPDDM_KRYLOV_METHOD_BGCRODR || data->cntl[0] == HPDDM_KRYLOV_METHOD_BFBCG) {
 53:       data->rcntl[0] = (PetscAbsReal(data->rcntl[0] - static_cast<PetscReal>(PETSC_DECIDE)) < PETSC_SMALL ? -1.0 : data->rcntl[0]);
 54:       PetscCall(PetscOptionsReal("-ksp_hpddm_deflation_tol", "Tolerance when deflating right-hand sides inside block methods", "KSPHPDDM", data->rcntl[0], data->rcntl, NULL));
 55:       i = (data->scntl[data->cntl[0] != HPDDM_KRYLOV_METHOD_BFBCG] == static_cast<unsigned short>(PETSC_DECIDE) ? 1 : PetscMax(1, data->scntl[data->cntl[0] != HPDDM_KRYLOV_METHOD_BFBCG]));
 56:       PetscCall(PetscOptionsRangeInt("-ksp_hpddm_enlarge_krylov_subspace", "Split the initial right-hand side into multiple vectors", "KSPHPDDM", i, &i, NULL, 1, std::numeric_limits<unsigned short>::max() - 1));
 57:       data->scntl[data->cntl[0] != HPDDM_KRYLOV_METHOD_BFBCG] = i;
 58:     } else data->scntl[data->cntl[0] != HPDDM_KRYLOV_METHOD_BCG] = 0;
 59:     if (data->cntl[0] == HPDDM_KRYLOV_METHOD_GMRES || data->cntl[0] == HPDDM_KRYLOV_METHOD_BGMRES || data->cntl[0] == HPDDM_KRYLOV_METHOD_GCRODR || data->cntl[0] == HPDDM_KRYLOV_METHOD_BGCRODR) {
 60:       i = (data->cntl[2] == static_cast<char>(PETSC_DECIDE) ? HPDDM_ORTHOGONALIZATION_CGS : data->cntl[2] & 3);
 61:       PetscCall(PetscOptionsEList("-ksp_hpddm_orthogonalization", "Classical (faster) or Modified (more robust) Gram--Schmidt process", "KSPHPDDM", HPDDMOrthogonalization, PETSC_STATIC_ARRAY_LENGTH(HPDDMOrthogonalization), HPDDMOrthogonalization[HPDDM_ORTHOGONALIZATION_CGS], &i, NULL));
 62:       j = (data->cntl[2] == static_cast<char>(PETSC_DECIDE) ? HPDDM_QR_CHOLQR : ((data->cntl[2] >> 2) & 7));
 63:       PetscCall(PetscOptionsEList("-ksp_hpddm_qr", "Distributed QR factorizations computed with Cholesky QR, Classical or Modified Gram--Schmidt process", "KSPHPDDM", HPDDMQR, PETSC_STATIC_ARRAY_LENGTH(HPDDMQR), HPDDMQR[HPDDM_QR_CHOLQR], &j, NULL));
 64:       data->cntl[2] = static_cast<char>(i) + (static_cast<char>(j) << 2);
 65:       i             = (data->scntl[0] == static_cast<unsigned short>(PETSC_DECIDE) ? PetscMin(30, ksp->max_it) : data->scntl[0]);
 66:       PetscCall(PetscOptionsRangeInt("-ksp_gmres_restart", "Maximum number of Arnoldi vectors generated per cycle", "KSPHPDDM", i, &i, NULL, PetscMin(1, ksp->max_it), PetscMin(ksp->max_it, std::numeric_limits<unsigned short>::max() - 1)));
 67:       data->scntl[0] = i;
 68:     }
 69:     if (data->cntl[0] == HPDDM_KRYLOV_METHOD_BCG || data->cntl[0] == HPDDM_KRYLOV_METHOD_BFBCG) {
 70:       j = (data->cntl[1] == static_cast<char>(PETSC_DECIDE) ? HPDDM_QR_CHOLQR : data->cntl[1]);
 71:       PetscCall(PetscOptionsEList("-ksp_hpddm_qr", "Distributed QR factorizations computed with Cholesky QR, Classical or Modified Gram--Schmidt process", "KSPHPDDM", HPDDMQR, PETSC_STATIC_ARRAY_LENGTH(HPDDMQR), HPDDMQR[HPDDM_QR_CHOLQR], &j, NULL));
 72:       data->cntl[1] = j;
 73:     }
 74:     if (data->cntl[0] == HPDDM_KRYLOV_METHOD_GCRODR || data->cntl[0] == HPDDM_KRYLOV_METHOD_BGCRODR) {
 75:       i = (data->icntl[0] == static_cast<int>(PETSC_DECIDE) ? PetscMin(20, data->scntl[0] - 1) : data->icntl[0]);
 76:       PetscCall(PetscOptionsRangeInt("-ksp_hpddm_recycle", "Number of harmonic Ritz vectors to compute", "KSPHPDDM", i, &i, NULL, 1, data->scntl[0] - 1));
 77:       data->icntl[0] = i;
 78:       if (!PetscDefined(HAVE_SLEPC) || !PetscDefined(USE_SHARED_LIBRARIES) || data->cntl[0] == HPDDM_KRYLOV_METHOD_GCRODR) {
 79:         i = (data->cntl[3] == static_cast<char>(PETSC_DECIDE) ? HPDDM_RECYCLE_TARGET_SM : data->cntl[3]);
 80:         PetscCall(PetscOptionsEList("-ksp_hpddm_recycle_target", "Criterion to select harmonic Ritz vectors", "KSPHPDDM", HPDDMRecycleTarget, PETSC_STATIC_ARRAY_LENGTH(HPDDMRecycleTarget), HPDDMRecycleTarget[HPDDM_RECYCLE_TARGET_SM], &i, NULL));
 81:         data->cntl[3] = i;
 82:       } else {
 83:         PetscCheck(data->precision == PETSC_KSPHPDDM_DEFAULT_PRECISION, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_INCOMP, "Cannot use SLEPc with a different precision than PETSc for harmonic Ritz eigensolves");
 84:         PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ksp), &size));
 85:         i = (data->cntl[3] == static_cast<char>(PETSC_DECIDE) ? 1 : data->cntl[3]);
 86:         PetscCall(PetscOptionsRangeInt("-ksp_hpddm_recycle_redistribute", "Number of processes used to solve eigenvalue problems when recycling in BGCRODR", "KSPHPDDM", i, &i, NULL, 1, PetscMin(size, 192)));
 87:         data->cntl[3] = i;
 88:       }
 89:       i = (data->cntl[4] == static_cast<char>(PETSC_DECIDE) ? HPDDM_RECYCLE_STRATEGY_A : data->cntl[4]);
 90:       PetscCall(PetscOptionsEList("-ksp_hpddm_recycle_strategy", "Generalized eigenvalue problem to solve for recycling", "KSPHPDDM", HPDDMRecycleStrategy, PETSC_STATIC_ARRAY_LENGTH(HPDDMRecycleStrategy), HPDDMRecycleStrategy[HPDDM_RECYCLE_STRATEGY_A], &i, NULL));
 91:       data->cntl[4] = i;
 92:     }
 93:   } else {
 94:     data->cntl[0]  = HPDDM_KRYLOV_METHOD_NONE;
 95:     data->scntl[1] = 1;
 96:   }
 97:   PetscCheck(ksp->nmax >= std::numeric_limits<int>::min() && ksp->nmax <= std::numeric_limits<int>::max(), PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "KSPMatSolve() block size %" PetscInt_FMT " not representable by an integer, which is not handled by KSPHPDDM",
 98:              ksp->nmax);
 99:   data->icntl[1] = static_cast<int>(ksp->nmax);
100:   PetscOptionsHeadEnd();
101:   PetscFunctionReturn(PETSC_SUCCESS);
102: }

104: static PetscErrorCode KSPView_HPDDM(KSP ksp, PetscViewer viewer)
105: {
106:   KSP_HPDDM            *data  = (KSP_HPDDM *)ksp->data;
107:   HPDDM::PETScOperator *op    = data->op;
108:   const PetscScalar    *array = op ? op->storage() : NULL;
109:   PetscBool             ascii;

111:   PetscFunctionBegin;
112:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &ascii));
113:   if (op && ascii) {
114:     PetscCall(PetscViewerASCIIPrintf(viewer, "HPDDM type: %s\n", KSPHPDDMTypes[std::min(static_cast<PetscInt>(data->cntl[0]), static_cast<PetscInt>(PETSC_STATIC_ARRAY_LENGTH(KSPHPDDMTypes) - 1))]));
115:     PetscCall(PetscViewerASCIIPrintf(viewer, "precision: %s\n", KSPHPDDMPrecisionTypes[data->precision]));
116:     if (data->cntl[0] == HPDDM_KRYLOV_METHOD_BGMRES || data->cntl[0] == HPDDM_KRYLOV_METHOD_BGCRODR || data->cntl[0] == HPDDM_KRYLOV_METHOD_BFBCG) {
117:       if (PetscAbsReal(data->rcntl[0] - static_cast<PetscReal>(PETSC_DECIDE)) < PETSC_SMALL) PetscCall(PetscViewerASCIIPrintf(viewer, "no deflation at restarts\n"));
118:       else PetscCall(PetscViewerASCIIPrintf(viewer, "deflation tolerance: %g\n", static_cast<double>(data->rcntl[0])));
119:     }
120:     if (data->cntl[0] == HPDDM_KRYLOV_METHOD_GCRODR || data->cntl[0] == HPDDM_KRYLOV_METHOD_BGCRODR) {
121:       PetscCall(PetscViewerASCIIPrintf(viewer, "deflation subspace attached? %s\n", PetscBools[array ? PETSC_TRUE : PETSC_FALSE]));
122:       if (!PetscDefined(HAVE_SLEPC) || !PetscDefined(USE_SHARED_LIBRARIES) || data->cntl[0] == HPDDM_KRYLOV_METHOD_GCRODR) PetscCall(PetscViewerASCIIPrintf(viewer, "deflation target: %s\n", HPDDMRecycleTarget[static_cast<PetscInt>(data->cntl[3])]));
123:       else PetscCall(PetscViewerASCIIPrintf(viewer, "redistribution size: %d\n", static_cast<PetscMPIInt>(data->cntl[3])));
124:     }
125:     if (data->icntl[1] != static_cast<int>(PETSC_DECIDE)) PetscCall(PetscViewerASCIIPrintf(viewer, "  block size is %d\n", data->icntl[1]));
126:   }
127:   PetscFunctionReturn(PETSC_SUCCESS);
128: }

130: static PetscErrorCode KSPSetUp_HPDDM(KSP ksp)
131: {
132:   KSP_HPDDM *data = (KSP_HPDDM *)ksp->data;
133:   Mat        A;
134:   PetscInt   n, bs;
135:   PetscBool  match;

137:   PetscFunctionBegin;
138:   PetscCall(KSPGetOperators(ksp, &A, NULL));
139:   PetscCall(MatGetLocalSize(A, &n, NULL));
140:   PetscCall(MatGetBlockSize(A, &bs));
141:   PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &match, MATSEQKAIJ, MATMPIKAIJ, ""));
142:   if (match) n /= bs;
143:   data->op = new HPDDM::PETScOperator(ksp, n);
144:   if (PetscUnlikely(!ksp->setfromoptionscalled || data->cntl[0] == static_cast<char>(PETSC_DECIDE))) { /* what follows is basically a copy/paste of KSPSetFromOptions_HPDDM, with no call to PetscOptions() */
145:     PetscCall(PetscInfo(ksp, "KSPSetFromOptions() not called or uninitialized internal structure, hardwiring default KSPHPDDM options\n"));
146:     if (data->cntl[0] == static_cast<char>(PETSC_DECIDE)) data->cntl[0] = 0; /* GMRES by default */
147:     if (data->cntl[0] != HPDDM_KRYLOV_METHOD_NONE) {                         /* following options do not matter with PREONLY */
148:       if (data->cntl[0] != HPDDM_KRYLOV_METHOD_BCG && data->cntl[0] != HPDDM_KRYLOV_METHOD_BFBCG) {
149:         data->cntl[1] = HPDDM_VARIANT_LEFT; /* left preconditioning by default */
150:         if (ksp->pc_side_set == PC_RIGHT) data->cntl[1] = HPDDM_VARIANT_RIGHT;
151:         if (data->cntl[1] > 0) PetscCall(KSPSetPCSide(ksp, PC_RIGHT));
152:       }
153:       if (data->cntl[0] == HPDDM_KRYLOV_METHOD_BGMRES || data->cntl[0] == HPDDM_KRYLOV_METHOD_BGCRODR || data->cntl[0] == HPDDM_KRYLOV_METHOD_BFBCG) {
154:         data->rcntl[0]                                          = -1.0; /* no deflation by default */
155:         data->scntl[data->cntl[0] != HPDDM_KRYLOV_METHOD_BFBCG] = 1;    /* Krylov subspace not enlarged by default */
156:       } else data->scntl[data->cntl[0] != HPDDM_KRYLOV_METHOD_BCG] = 0;
157:       if (data->cntl[0] == HPDDM_KRYLOV_METHOD_GMRES || data->cntl[0] == HPDDM_KRYLOV_METHOD_BGMRES || data->cntl[0] == HPDDM_KRYLOV_METHOD_GCRODR || data->cntl[0] == HPDDM_KRYLOV_METHOD_BGCRODR) {
158:         data->cntl[2]  = static_cast<char>(HPDDM_ORTHOGONALIZATION_CGS) + (static_cast<char>(HPDDM_QR_CHOLQR) << 2); /* CGS and CholQR by default */
159:         data->scntl[0] = PetscMin(30, ksp->max_it);                                                                  /* restart parameter of 30 by default */
160:       }
161:       if (data->cntl[0] == HPDDM_KRYLOV_METHOD_BCG || data->cntl[0] == HPDDM_KRYLOV_METHOD_BFBCG) { data->cntl[1] = HPDDM_QR_CHOLQR; /* CholQR by default */ }
162:       if (data->cntl[0] == HPDDM_KRYLOV_METHOD_GCRODR || data->cntl[0] == HPDDM_KRYLOV_METHOD_BGCRODR) {
163:         data->icntl[0] = PetscMin(20, data->scntl[0] - 1); /* recycled subspace of size 20 by default */
164:         if (!PetscDefined(HAVE_SLEPC) || !PetscDefined(USE_SHARED_LIBRARIES) || data->cntl[0] == HPDDM_KRYLOV_METHOD_GCRODR) {
165:           data->cntl[3] = HPDDM_RECYCLE_TARGET_SM; /* default recycling target */
166:         } else {
167:           data->cntl[3] = 1; /* redistribution parameter of 1 by default */
168:         }
169:         data->cntl[4] = HPDDM_RECYCLE_STRATEGY_A; /* default recycling strategy */
170:       }
171:     } else data->scntl[1] = 1;
172:   }
173:   PetscCheck(ksp->nmax >= std::numeric_limits<int>::min() && ksp->nmax <= std::numeric_limits<int>::max(), PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "KSPMatSolve() block size %" PetscInt_FMT " not representable by an integer, which is not handled by KSPHPDDM",
174:              ksp->nmax);
175:   data->icntl[1] = static_cast<int>(ksp->nmax);
176:   PetscFunctionReturn(PETSC_SUCCESS);
177: }

179: static inline PetscErrorCode KSPReset_HPDDM_Private(KSP ksp)
180: {
181:   KSP_HPDDM *data = (KSP_HPDDM *)ksp->data;

183:   PetscFunctionBegin;
184:   /* cast PETSC_DECIDE into the appropriate types to avoid compiler warnings */
185:   std::fill_n(data->rcntl, PETSC_STATIC_ARRAY_LENGTH(data->rcntl), static_cast<PetscReal>(PETSC_DECIDE));
186:   std::fill_n(data->icntl, PETSC_STATIC_ARRAY_LENGTH(data->icntl), static_cast<int>(PETSC_DECIDE));
187:   std::fill_n(data->scntl, PETSC_STATIC_ARRAY_LENGTH(data->scntl), static_cast<unsigned short>(PETSC_DECIDE));
188:   std::fill_n(data->cntl, PETSC_STATIC_ARRAY_LENGTH(data->cntl), static_cast<char>(PETSC_DECIDE));
189:   data->precision = PETSC_KSPHPDDM_DEFAULT_PRECISION;
190:   PetscFunctionReturn(PETSC_SUCCESS);
191: }

193: static PetscErrorCode KSPReset_HPDDM(KSP ksp)
194: {
195:   KSP_HPDDM *data = (KSP_HPDDM *)ksp->data;

197:   PetscFunctionBegin;
198:   delete data->op;
199:   data->op = NULL;
200:   PetscCall(KSPReset_HPDDM_Private(ksp));
201:   PetscFunctionReturn(PETSC_SUCCESS);
202: }

204: static PetscErrorCode KSPDestroy_HPDDM(KSP ksp)
205: {
206:   PetscFunctionBegin;
207:   PetscCall(KSPReset_HPDDM(ksp));
208:   PetscCall(KSPDestroyDefault(ksp));
209:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPHPDDMSetDeflationMat_C", NULL));
210:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPHPDDMGetDeflationMat_C", NULL));
211:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPHPDDMSetType_C", NULL));
212:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPHPDDMGetType_C", NULL));
213:   PetscFunctionReturn(PETSC_SUCCESS);
214: }

216: template <PetscMemType type = PETSC_MEMTYPE_HOST>
217: static inline PetscErrorCode KSPSolve_HPDDM_Private(KSP ksp, const PetscScalar *b, PetscScalar *x, PetscInt n)
218: {
219:   KSP_HPDDM              *data = (KSP_HPDDM *)ksp->data;
220:   KSPConvergedDefaultCtx *ctx  = (KSPConvergedDefaultCtx *)ksp->cnvP;
221:   const PetscInt          N    = data->op->getDof() * n;
222:   PetscBool               flg;
223: #if !PetscDefined(USE_REAL_DOUBLE) || PetscDefined(HAVE_F2CBLASLAPACK___FLOAT128_BINDINGS)
224:   HPDDM::upscaled_type<PetscScalar> *high[2];
225: #endif
226: #if !PetscDefined(USE_REAL_SINGLE) || PetscDefined(HAVE_F2CBLASLAPACK___FP16_BINDINGS)
227:   typedef HPDDM::downscaled_type<PetscScalar> PetscDownscaledScalar PETSC_ATTRIBUTE_MAY_ALIAS;
228:   PetscDownscaledScalar                                            *low[2];
229: #endif
230: #if PetscDefined(HAVE_CUDA)
231:   Mat     A;
232:   VecType vtype;
233: #endif

235:   PetscFunctionBegin;
236: #if PetscDefined(HAVE_CUDA)
237:   PetscCall(KSPGetOperators(ksp, &A, NULL));
238:   PetscCall(MatGetVecType(A, &vtype));
239:   std::initializer_list<std::string>                 list = {VECCUDA, VECSEQCUDA, VECMPICUDA};
240:   std::initializer_list<std::string>::const_iterator it   = std::find(list.begin(), list.end(), std::string(vtype));
241:   PetscCheck(type != PETSC_MEMTYPE_HOST || it == list.end(), PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "MatGetVecType() must return a Vec with the same PetscMemType as the right-hand side and solution, PetscMemType(%s) != %s", vtype, PetscMemTypeToString(type));
242: #endif
243:   PetscCall(PCGetDiagonalScale(ksp->pc, &flg));
244:   PetscCheck(!flg, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);
245:   if (n > 1) {
246:     if (ksp->converged == KSPConvergedDefault) {
247:       PetscCheck(!ctx->mininitialrtol, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Krylov method %s does not support KSPConvergedDefaultSetUMIRNorm()", ((PetscObject)ksp)->type_name);
248:       if (!ctx->initialrtol) {
249:         PetscCall(PetscInfo(ksp, "Forcing KSPConvergedDefaultSetUIRNorm() since KSPConvergedDefault() cannot handle multiple norms\n"));
250:         ctx->initialrtol = PETSC_TRUE;
251:       }
252:     } else PetscCall(PetscInfo(ksp, "Using a special \"converged\" callback, be careful, it is used in KSPHPDDM to track blocks of residuals\n"));
253:   }
254:   /* initial guess is always nonzero with recycling methods if there is a deflation subspace available */
255:   if ((data->cntl[0] == HPDDM_KRYLOV_METHOD_GCRODR || data->cntl[0] == HPDDM_KRYLOV_METHOD_BGCRODR) && data->op->storage()) ksp->guess_zero = PETSC_FALSE;
256:   ksp->its    = 0;
257:   ksp->reason = KSP_CONVERGED_ITERATING;
258:   if (data->precision > PETSC_KSPHPDDM_DEFAULT_PRECISION) { /* Krylov basis stored in higher precision than PetscScalar */
259: #if !PetscDefined(USE_REAL_DOUBLE) || PetscDefined(HAVE_F2CBLASLAPACK___FLOAT128_BINDINGS)
260:     if (type == PETSC_MEMTYPE_HOST) {
261:       PetscCall(PetscMalloc2(N, high, N, high + 1));
262:       HPDDM::copy_n(b, N, high[0]);
263:       HPDDM::copy_n(x, N, high[1]);
264:       PetscCall(HPDDM::IterativeMethod::solve(*data->op, high[0], high[1], n, PetscObjectComm((PetscObject)ksp)));
265:       HPDDM::copy_n(high[1], N, x);
266:       PetscCall(PetscFree2(high[0], high[1]));
267:     } else {
268:       PetscCheck(PetscDefined(HAVE_CUDA) && PetscDefined(USE_REAL_SINGLE), PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "CUDA in PETSc has no support for precisions other than single or double");
269:   #if PetscDefined(HAVE_CUDA)
270:     #if PetscDefined(HAVE_HPDDM)
271:       PetscCall(KSPSolve_HPDDM_CUDA_Private(data, b, x, n, PetscObjectComm((PetscObject)ksp)));
272:     #else
273:       SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "No CUDA support with --download-hpddm from SLEPc");
274:     #endif
275:   #endif
276:     }
277: #else
278:     PetscCheck(data->precision != KSP_HPDDM_PRECISION_QUADRUPLE, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Reconfigure with --download-f2cblaslapack --with-f2cblaslapack-float128-bindings");
279: #endif
280:   } else if (data->precision < PETSC_KSPHPDDM_DEFAULT_PRECISION) { /* Krylov basis stored in lower precision than PetscScalar */
281: #if !PetscDefined(USE_REAL_SINGLE) || PetscDefined(HAVE_F2CBLASLAPACK___FP16_BINDINGS)
282:     if (type == PETSC_MEMTYPE_HOST) {
283:       PetscCall(PetscMalloc1(N, low));
284:       low[1] = reinterpret_cast<PetscDownscaledScalar *>(x);
285:       std::copy_n(b, N, low[0]);
286:       for (PetscInt i = 0; i < N; ++i) low[1][i] = x[i];
287:       PetscCall(HPDDM::IterativeMethod::solve(*data->op, low[0], low[1], n, PetscObjectComm((PetscObject)ksp)));
288:       for (PetscInt i = N; i-- > 0;) x[i] = low[1][i];
289:       PetscCall(PetscFree(low[0]));
290:     } else {
291:       PetscCheck(PetscDefined(HAVE_CUDA) && PetscDefined(USE_REAL_DOUBLE), PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "CUDA in PETSc has no support for precisions other than single or double");
292:   #if PetscDefined(HAVE_CUDA)
293:     #if PetscDefined(HAVE_HPDDM)
294:       PetscCall(KSPSolve_HPDDM_CUDA_Private(data, b, x, n, PetscObjectComm((PetscObject)ksp)));
295:     #else
296:       SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "No CUDA support with --download-hpddm from SLEPc");
297:     #endif
298:   #endif
299:     }
300: #else
301:     PetscCheck(data->precision != KSP_HPDDM_PRECISION_HALF, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Reconfigure with --download-f2cblaslapack --with-f2cblaslapack-fp16-bindings");
302: #endif
303:   } else { /* Krylov basis stored in the same precision as PetscScalar */
304:     if (type == PETSC_MEMTYPE_HOST) PetscCall(HPDDM::IterativeMethod::solve(*data->op, b, x, n, PetscObjectComm((PetscObject)ksp)));
305:     else {
306:       PetscCheck(PetscDefined(USE_REAL_SINGLE) || PetscDefined(USE_REAL_DOUBLE), PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "CUDA in PETSc has no support for precisions other than single or double");
307: #if PetscDefined(HAVE_CUDA)
308:   #if PetscDefined(HAVE_HPDDM)
309:       PetscCall(KSPSolve_HPDDM_CUDA_Private(data, b, x, n, PetscObjectComm((PetscObject)ksp)));
310:   #else
311:       SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "No CUDA support with --download-hpddm from SLEPc");
312:   #endif
313: #endif
314:     }
315:   }
316:   if (!ksp->reason) { /* KSPConvergedDefault() is still returning 0 (= KSP_CONVERGED_ITERATING) */
317:     if (ksp->its >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
318:     else ksp->reason = KSP_CONVERGED_RTOL; /* early exit by HPDDM, which only happens on breakdowns or convergence */
319:   }
320:   ksp->its = PetscMin(ksp->its, ksp->max_it);
321:   PetscFunctionReturn(PETSC_SUCCESS);
322: }

324: static PetscErrorCode KSPSolve_HPDDM(KSP ksp)
325: {
326:   KSP_HPDDM         *data = (KSP_HPDDM *)ksp->data;
327:   Mat                A, B;
328:   PetscScalar       *x, *bt = NULL, **ptr;
329:   const PetscScalar *b;
330:   PetscInt           i, j, n;
331:   PetscBool          flg;
332:   PetscMemType       type[2];

334:   PetscFunctionBegin;
335:   PetscCall(PetscCitationsRegister(HPDDMCitation, &HPDDMCite));
336:   PetscCall(KSPGetOperators(ksp, &A, NULL));
337:   PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQKAIJ, MATMPIKAIJ, ""));
338:   PetscCall(VecGetArrayWriteAndMemType(ksp->vec_sol, &x, type));
339:   PetscCall(VecGetArrayReadAndMemType(ksp->vec_rhs, &b, type + 1));
340:   PetscCheck(type[0] == type[1], PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Right-hand side and solution vectors must have the same PetscMemType, %s != %s", PetscMemTypeToString(type[0]), PetscMemTypeToString(type[1]));
341:   if (!flg) {
342:     if (PetscMemTypeCUDA(type[0])) PetscCall(KSPSolve_HPDDM_Private<PETSC_MEMTYPE_CUDA>(ksp, b, x, 1));
343:     else {
344:       PetscCheck(PetscMemTypeHost(type[0]), PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "PetscMemType (%s) is neither PETSC_MEMTYPE_HOST nor PETSC_MEMTYPE_CUDA", PetscMemTypeToString(type[0]));
345:       PetscCall(KSPSolve_HPDDM_Private(ksp, b, x, 1));
346:     }
347:   } else {
348:     PetscCheck(PetscMemTypeHost(type[0]), PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "PetscMemType (%s) is not PETSC_MEMTYPE_HOST", PetscMemTypeToString(type[0]));
349:     PetscCall(MatKAIJGetScaledIdentity(A, &flg));
350:     PetscCall(MatKAIJGetAIJ(A, &B));
351:     PetscCall(MatGetBlockSize(A, &n));
352:     PetscCall(MatGetLocalSize(B, &i, NULL));
353:     j = data->op->getDof();
354:     if (!flg) i *= n; /* S and T are not scaled identities, cannot use block methods */
355:     if (i != j) {     /* switching between block and standard methods */
356:       delete data->op;
357:       data->op = new HPDDM::PETScOperator(ksp, i);
358:     }
359:     if (flg && n > 1) {
360:       PetscCall(PetscMalloc1(i * n, &bt));
361:       /* from row- to column-major to be consistent with HPDDM */
362:       HPDDM::Wrapper<PetscScalar>::omatcopy<'T'>(i, n, b, n, bt, i);
363:       ptr = const_cast<PetscScalar **>(&b);
364:       std::swap(*ptr, bt);
365:       HPDDM::Wrapper<PetscScalar>::imatcopy<'T'>(i, n, x, n, i);
366:     }
367:     PetscCall(KSPSolve_HPDDM_Private(ksp, b, x, flg ? n : 1));
368:     if (flg && n > 1) {
369:       std::swap(*ptr, bt);
370:       PetscCall(PetscFree(bt));
371:       /* from column- to row-major to be consistent with MatKAIJ format */
372:       HPDDM::Wrapper<PetscScalar>::imatcopy<'T'>(n, i, x, i, n);
373:     }
374:   }
375:   PetscCall(VecRestoreArrayReadAndMemType(ksp->vec_rhs, &b));
376:   PetscCall(VecRestoreArrayWriteAndMemType(ksp->vec_sol, &x));
377:   PetscFunctionReturn(PETSC_SUCCESS);
378: }

380: /*@
381:      KSPHPDDMSetDeflationMat - Sets the deflation space used by Krylov methods in `KSPHPDDM` with recycling. This space is viewed as a set of vectors stored in
382:      a `MATDENSE` (column major).

384:    Input Parameters:
385: +     ksp - iterative context
386: -     U - deflation space to be used during KSPSolve()

388:    Level: intermediate

390: .seealso: [](chapter_ksp), `KSPHPDDM`, `KSPCreate()`, `KSPType`, `KSPHPDDMGetDeflationMat()`
391: @*/
392: PetscErrorCode KSPHPDDMSetDeflationMat(KSP ksp, Mat U)
393: {
394:   PetscFunctionBegin;
397:   PetscCheckSameComm(ksp, 1, U, 2);
398:   PetscUseMethod(ksp, "KSPHPDDMSetDeflationMat_C", (KSP, Mat), (ksp, U));
399:   PetscFunctionReturn(PETSC_SUCCESS);
400: }

402: /*@
403:      KSPHPDDMGetDeflationMat - Gets the deflation space computed by Krylov methods in `KSPHPDDM`  with recycling or NULL if `KSPSolve()` has not been called yet.
404:      This space is viewed as a set of vectors stored in a `MATDENSE` (column major). It is the responsibility of the user to free the returned `Mat`.

406:    Input Parameter:
407: .     ksp - iterative context

409:    Output Parameter:
410: .     U - deflation space generated during `KSPSolve()`

412:    Level: intermediate

414: .seealso: [](chapter_ksp), `KSPHPDDM`, `KSPCreate()`, `KSPType`, `KSPHPDDMSetDeflationMat()`
415: @*/
416: PetscErrorCode KSPHPDDMGetDeflationMat(KSP ksp, Mat *U)
417: {
418:   PetscFunctionBegin;
420:   if (U) {
422:     PetscUseMethod(ksp, "KSPHPDDMGetDeflationMat_C", (KSP, Mat *), (ksp, U));
423:   }
424:   PetscFunctionReturn(PETSC_SUCCESS);
425: }

427: static PetscErrorCode KSPHPDDMSetDeflationMat_HPDDM(KSP ksp, Mat U)
428: {
429:   KSP_HPDDM            *data = (KSP_HPDDM *)ksp->data;
430:   HPDDM::PETScOperator *op   = data->op;
431:   Mat                   A;
432:   const PetscScalar    *array;
433:   PetscScalar          *copy;
434:   PetscInt              m1, M1, m2, M2, n2, N2, ldu;
435:   PetscBool             match;

437:   PetscFunctionBegin;
438:   if (!op) {
439:     PetscCall(KSPSetUp(ksp));
440:     op = data->op;
441:   }
442:   PetscCheck(data->precision == PETSC_KSPHPDDM_DEFAULT_PRECISION, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "%s != %s", KSPHPDDMPrecisionTypes[data->precision], KSPHPDDMPrecisionTypes[PETSC_KSPHPDDM_DEFAULT_PRECISION]);
443:   PetscCall(KSPGetOperators(ksp, &A, NULL));
444:   PetscCall(MatGetLocalSize(A, &m1, NULL));
445:   PetscCall(MatGetLocalSize(U, &m2, &n2));
446:   PetscCall(MatGetSize(A, &M1, NULL));
447:   PetscCall(MatGetSize(U, &M2, &N2));
448:   PetscCheck(m1 == m2 && M1 == M2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot use a deflation space with (m2,M2) = (%" PetscInt_FMT ",%" PetscInt_FMT ") for a linear system with (m1,M1) = (%" PetscInt_FMT ",%" PetscInt_FMT ")", m2, M2, m1, M1);
449:   PetscCall(PetscObjectTypeCompareAny((PetscObject)U, &match, MATSEQDENSE, MATMPIDENSE, ""));
450:   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided deflation space not stored in a dense Mat");
451:   PetscCall(MatDenseGetArrayRead(U, &array));
452:   copy = op->allocate(m2, 1, N2);
453:   PetscCheck(copy, PETSC_COMM_SELF, PETSC_ERR_POINTER, "Memory allocation error");
454:   PetscCall(MatDenseGetLDA(U, &ldu));
455:   HPDDM::Wrapper<PetscScalar>::omatcopy<'N'>(N2, m2, array, ldu, copy, m2);
456:   PetscCall(MatDenseRestoreArrayRead(U, &array));
457:   PetscFunctionReturn(PETSC_SUCCESS);
458: }

460: static PetscErrorCode KSPHPDDMGetDeflationMat_HPDDM(KSP ksp, Mat *U)
461: {
462:   KSP_HPDDM            *data = (KSP_HPDDM *)ksp->data;
463:   HPDDM::PETScOperator *op   = data->op;
464:   Mat                   A;
465:   const PetscScalar    *array;
466:   PetscScalar          *copy;
467:   PetscInt              m1, M1, N2;

469:   PetscFunctionBegin;
470:   if (!op) {
471:     PetscCall(KSPSetUp(ksp));
472:     op = data->op;
473:   }
474:   PetscCheck(data->precision == PETSC_KSPHPDDM_DEFAULT_PRECISION, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "%s != %s", KSPHPDDMPrecisionTypes[data->precision], KSPHPDDMPrecisionTypes[PETSC_KSPHPDDM_DEFAULT_PRECISION]);
475:   array = op->storage();
476:   N2    = op->k().first * op->k().second;
477:   if (!array) *U = NULL;
478:   else {
479:     PetscCall(KSPGetOperators(ksp, &A, NULL));
480:     PetscCall(MatGetLocalSize(A, &m1, NULL));
481:     PetscCall(MatGetSize(A, &M1, NULL));
482:     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)ksp), m1, PETSC_DECIDE, M1, N2, NULL, U));
483:     PetscCall(MatDenseGetArrayWrite(*U, &copy));
484:     PetscCall(PetscArraycpy(copy, array, m1 * N2));
485:     PetscCall(MatDenseRestoreArrayWrite(*U, &copy));
486:   }
487:   PetscFunctionReturn(PETSC_SUCCESS);
488: }

490: static PetscErrorCode KSPMatSolve_HPDDM(KSP ksp, Mat B, Mat X)
491: {
492:   KSP_HPDDM         *data = (KSP_HPDDM *)ksp->data;
493:   Mat                A;
494:   const PetscScalar *b;
495:   PetscScalar       *x;
496:   PetscInt           n, lda;
497:   PetscMemType       type[2];

499:   PetscFunctionBegin;
500:   PetscCall(PetscCitationsRegister(HPDDMCitation, &HPDDMCite));
501:   if (!data->op) PetscCall(KSPSetUp(ksp));
502:   PetscCall(KSPGetOperators(ksp, &A, NULL));
503:   PetscCall(MatGetLocalSize(B, &n, NULL));
504:   PetscCall(MatDenseGetLDA(B, &lda));
505:   PetscCheck(n == lda, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Unhandled leading dimension lda = %" PetscInt_FMT " with n = %" PetscInt_FMT, lda, n);
506:   PetscCall(MatGetLocalSize(A, &n, NULL));
507:   PetscCall(MatDenseGetLDA(X, &lda));
508:   PetscCheck(n == lda, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Unhandled leading dimension lda = %" PetscInt_FMT " with n = %" PetscInt_FMT, lda, n);
509:   PetscCall(MatGetSize(X, NULL, &n));
510:   PetscCall(MatDenseGetArrayWriteAndMemType(X, &x, type));
511:   PetscCall(MatDenseGetArrayReadAndMemType(B, &b, type + 1));
512:   PetscCheck(type[0] == type[1], PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Right-hand side and solution matrices must have the same PetscMemType, %s != %s", PetscMemTypeToString(type[0]), PetscMemTypeToString(type[1]));
513:   if (PetscMemTypeCUDA(type[0])) PetscCall(KSPSolve_HPDDM_Private<PETSC_MEMTYPE_CUDA>(ksp, b, x, n));
514:   else {
515:     PetscCheck(PetscMemTypeHost(type[0]), PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "PetscMemType (%s) is neither PETSC_MEMTYPE_HOST nor PETSC_MEMTYPE_CUDA", PetscMemTypeToString(type[0]));
516:     PetscCall(KSPSolve_HPDDM_Private(ksp, b, x, n));
517:   }
518:   PetscCall(MatDenseRestoreArrayReadAndMemType(B, &b));
519:   PetscCall(MatDenseRestoreArrayWriteAndMemType(X, &x));
520:   PetscFunctionReturn(PETSC_SUCCESS);
521: }

523: /*@
524:      KSPHPDDMSetType - Sets the type of Krylov method used in `KSPHPDDM`.

526:    Collective

528:    Input Parameters:
529: +     ksp - iterative context
530: -     type - any of gmres, bgmres, cg, bcg, gcrodr, bgcrodr, bfbcg, or preonly

532:    Level: intermediate

534:    Notes:
535:      Unlike `KSPReset()`, this function does not destroy any deflation space attached to the `KSP`.

537:      As an example, in the following sequence:
538: .vb
539:      KSPHPDDMSetType(ksp, KSPGCRODR);
540:      KSPSolve(ksp, b, x);
541:      KSPHPDDMSetType(ksp, KSPGMRES);
542:      KSPHPDDMSetType(ksp, KSPGCRODR);
543:      KSPSolve(ksp, b, x);
544: .ve
545:     the recycled space is reused in the second `KSPSolve()`.

547: .seealso: [](chapter_ksp), `KSPCreate()`, `KSPType`, `KSPHPDDMType`, `KSPHPDDMGetType()`
548: @*/
549: PetscErrorCode KSPHPDDMSetType(KSP ksp, KSPHPDDMType type)
550: {
551:   PetscFunctionBegin;
554:   PetscUseMethod(ksp, "KSPHPDDMSetType_C", (KSP, KSPHPDDMType), (ksp, type));
555:   PetscFunctionReturn(PETSC_SUCCESS);
556: }

558: /*@
559:      KSPHPDDMGetType - Gets the type of Krylov method used in `KSPHPDDM`.

561:    Input Parameter:
562: .     ksp - iterative context

564:    Output Parameter:
565: .     type - any of gmres, bgmres, cg, bcg, gcrodr, bgcrodr, bfbcg, or preonly

567:    Level: intermediate

569: .seealso: [](chapter_ksp), `KSPCreate()`, `KSPType`, `KSPHPDDMType`, `KSPHPDDMSetType()`
570: @*/
571: PetscErrorCode KSPHPDDMGetType(KSP ksp, KSPHPDDMType *type)
572: {
573:   PetscFunctionBegin;
575:   if (type) {
577:     PetscUseMethod(ksp, "KSPHPDDMGetType_C", (KSP, KSPHPDDMType *), (ksp, type));
578:   }
579:   PetscFunctionReturn(PETSC_SUCCESS);
580: }

582: static PetscErrorCode KSPHPDDMSetType_HPDDM(KSP ksp, KSPHPDDMType type)
583: {
584:   KSP_HPDDM *data = (KSP_HPDDM *)ksp->data;
585:   PetscInt   i;
586:   PetscBool  flg = PETSC_FALSE;

588:   PetscFunctionBegin;
589:   for (i = 0; i < static_cast<PetscInt>(PETSC_STATIC_ARRAY_LENGTH(KSPHPDDMTypes)); ++i) {
590:     PetscCall(PetscStrcmp(KSPHPDDMTypes[type], KSPHPDDMTypes[i], &flg));
591:     if (flg) break;
592:   }
593:   PetscCheck(i != PETSC_STATIC_ARRAY_LENGTH(KSPHPDDMTypes), PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown KSPHPDDMType %d", type);
594:   if (data->cntl[0] != static_cast<char>(PETSC_DECIDE) && data->cntl[0] != i) PetscCall(KSPReset_HPDDM_Private(ksp));
595:   data->cntl[0] = i;
596:   PetscFunctionReturn(PETSC_SUCCESS);
597: }

599: static PetscErrorCode KSPHPDDMGetType_HPDDM(KSP ksp, KSPHPDDMType *type)
600: {
601:   KSP_HPDDM *data = (KSP_HPDDM *)ksp->data;

603:   PetscFunctionBegin;
604:   PetscCheck(data->cntl[0] != static_cast<char>(PETSC_DECIDE), PETSC_COMM_SELF, PETSC_ERR_ORDER, "KSPHPDDMType not set yet");
605:   /* need to shift by -1 for HPDDM_KRYLOV_METHOD_NONE */
606:   *type = static_cast<KSPHPDDMType>(PetscMin(data->cntl[0], static_cast<char>(PETSC_STATIC_ARRAY_LENGTH(KSPHPDDMTypes) - 1)));
607:   PetscFunctionReturn(PETSC_SUCCESS);
608: }

610: /*MC
611:      KSPHPDDM - Interface with the HPDDM library. This `KSP` may be used to further select methods that are currently not implemented natively in PETSc, e.g.,
612:      GCRODR [2006], a recycled Krylov method which is similar to `KSPLGMRES`, see [2016] for a comparison. ex75.c shows how to reproduce the results
613:      from the aforementioned paper [2006]. A chronological bibliography of relevant publications linked with `KSP` available in HPDDM through `KSPHPDDM`,
614:      and not available directly in PETSc, may be found below. The interface is explained in details in [2021].

616:    Options Database Keys:
617: +   -ksp_gmres_restart <restart, default=30> - see `KSPGMRES`
618: .   -ksp_hpddm_type <type, default=gmres> - any of gmres, bgmres, cg, bcg, gcrodr, bgcrodr, bfbcg, or preonly, see `KSPHPDDMType`
619: .   -ksp_hpddm_precision <value, default=same as PetscScalar> - any of half, single, double or quadruple, see `KSPHPDDMPrecision`
620: .   -ksp_hpddm_deflation_tol <eps, default=\-1.0> - tolerance when deflating right-hand sides inside block methods (no deflation by default, only relevant with block methods)
621: .   -ksp_hpddm_enlarge_krylov_subspace <p, default=1> - split the initial right-hand side into multiple vectors (only relevant with nonblock methods)
622: .   -ksp_hpddm_orthogonalization <type, default=cgs> - any of cgs or mgs, see KSPGMRES
623: .   -ksp_hpddm_qr <type, default=cholqr> - distributed QR factorizations with any of cholqr, cgs, or mgs (only relevant with block methods)
624: .   -ksp_hpddm_variant <type, default=left> - any of left, right, or flexible (this option is superseded by `KSPSetPCSide()`)
625: .   -ksp_hpddm_recycle <n, default=0> - number of harmonic Ritz vectors to compute (only relevant with GCRODR or BGCRODR)
626: .   -ksp_hpddm_recycle_target <type, default=SM> - criterion to select harmonic Ritz vectors using either SM, LM, SR, LR, SI, or LI (only relevant with GCRODR or BGCRODR).
627:      For BGCRODR, if PETSc is compiled with SLEPc, this option is not relevant, since SLEPc is used instead. Options are set with the prefix -ksp_hpddm_recycle_eps_
628: .   -ksp_hpddm_recycle_strategy <type, default=A> - generalized eigenvalue problem A or B to solve for recycling (only relevant with flexible GCRODR or BGCRODR)
629: -   -ksp_hpddm_recycle_symmetric <true, default=false> - symmetric generalized eigenproblems in BGCRODR, useful to switch to distributed solvers like EPSELEMENTAL or EPSSCALAPACK
630:      (only relevant when PETSc is compiled with SLEPc)

632:    Level: intermediate

634:    References:
635: +   1980 - The block conjugate gradient algorithm and related methods. O'Leary. Linear Algebra and its Applications.
636: .   2006 - Recycling Krylov subspaces for sequences of linear systems. Parks, de Sturler, Mackey, Johnson, and Maiti. SIAM Journal on Scientific Computing
637: .   2013 - A modified block flexible GMRES method with deflation at each iteration for the solution of non-Hermitian linear systems with multiple right-hand sides.
638:            Calandra, Gratton, Lago, Vasseur, and Carvalho. SIAM Journal on Scientific Computing.
639: .   2016 - Block iterative methods and recycling for improved scalability of linear solvers. Jolivet and Tournier. SC16.
640: .   2017 - A breakdown-free block conjugate gradient method. Ji and Li. BIT Numerical Mathematics.
641: -   2021 - KSPHPDDM and PCHPDDM: extending PETSc with advanced Krylov methods and robust multilevel overlapping Schwarz preconditioners. Jolivet, Roman, and Zampini.
642:            Computer & Mathematics with Applications.

644: .seealso: [](chapter_ksp), [](sec_flexibleksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPGMRES`, `KSPCG`, `KSPLGMRES`, `KSPDGMRES`
645: M*/

647: PETSC_EXTERN PetscErrorCode KSPCreate_HPDDM(KSP ksp)
648: {
649:   KSP_HPDDM  *data;
650:   PetscInt    i;
651:   const char *common[] = {KSPGMRES, KSPCG, KSPPREONLY};
652:   PetscBool   flg      = PETSC_FALSE;

654:   PetscFunctionBegin;
655:   PetscCall(PetscNew(&data));
656:   ksp->data = (void *)data;
657:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 2));
658:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 1));
659:   ksp->ops->solve          = KSPSolve_HPDDM;
660:   ksp->ops->matsolve       = KSPMatSolve_HPDDM;
661:   ksp->ops->setup          = KSPSetUp_HPDDM;
662:   ksp->ops->setfromoptions = KSPSetFromOptions_HPDDM;
663:   ksp->ops->destroy        = KSPDestroy_HPDDM;
664:   ksp->ops->view           = KSPView_HPDDM;
665:   ksp->ops->reset          = KSPReset_HPDDM;
666:   PetscCall(KSPReset_HPDDM_Private(ksp));
667:   for (i = 0; i < static_cast<PetscInt>(PETSC_STATIC_ARRAY_LENGTH(common)); ++i) {
668:     PetscCall(PetscStrcmp(((PetscObject)ksp)->type_name, common[i], &flg));
669:     if (flg) break;
670:   }
671:   if (!i) data->cntl[0] = HPDDM_KRYLOV_METHOD_GMRES;
672:   else if (i == 1) data->cntl[0] = HPDDM_KRYLOV_METHOD_CG;
673:   else if (i == 2) data->cntl[0] = HPDDM_KRYLOV_METHOD_NONE;
674:   if (data->cntl[0] != static_cast<char>(PETSC_DECIDE)) PetscCall(PetscInfo(ksp, "Using the previously set KSPType %s\n", common[i]));
675:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPHPDDMSetDeflationMat_C", KSPHPDDMSetDeflationMat_HPDDM));
676:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPHPDDMGetDeflationMat_C", KSPHPDDMGetDeflationMat_HPDDM));
677:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPHPDDMSetType_C", KSPHPDDMSetType_HPDDM));
678:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPHPDDMGetType_C", KSPHPDDMGetType_HPDDM));
679: #if PetscDefined(HAVE_SLEPC) && PetscDefined(HAVE_DYNAMIC_LIBRARIES) && PetscDefined(USE_SHARED_LIBRARIES)
680:   if (!loadedDL) PetscCall(HPDDMLoadDL_Private(&loadedDL));
681: #endif
682:   data->precision = PETSC_KSPHPDDM_DEFAULT_PRECISION;
683:   PetscFunctionReturn(PETSC_SUCCESS);
684: }