xref: /petsc/src/ksp/ksp/impls/gmres/borthog.c (revision 1ed6e3ff8437baa411029a28c2b08f047df9ad9a)
1 /*
2     Routines used for the orthogonalization of the Hessenberg matrix.
3 
4     Note that for the complex numbers version, the VecDot() and
5     VecMDot() arguments within the code MUST remain in the order
6     given for correct computation of inner products.
7 */
8 #include <../src/ksp/ksp/impls/gmres/gmresimpl.h>
9 
10 /*@C
11   KSPGMRESModifiedGramSchmidtOrthogonalization -  This is the basic orthogonalization routine
12   using modified Gram-Schmidt.
13 
14   Collective, No Fortran Support
15 
16   Input Parameters:
17 + ksp - `KSP` object, must be associated with `KSPGMRES`, `KSPFGMRES`, or `KSPLGMRES` Krylov method
18 - it  - one less than the current GMRES restart iteration, i.e. the size of the Krylov space
19 
20   Options Database Key:
21 . -ksp_gmres_modifiedgramschmidt - Activates `KSPGMRESModifiedGramSchmidtOrthogonalization()`
22 
23   Level: intermediate
24 
25   Note:
26   In general this is much slower than `KSPGMRESClassicalGramSchmidtOrthogonalization()` but has better stability properties.
27 
28 .seealso: [](ch_ksp), `KSPGMRESSetOrthogonalization()`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESGetOrthogonalization()`
29 @*/
KSPGMRESModifiedGramSchmidtOrthogonalization(KSP ksp,PetscInt it)30 PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP ksp, PetscInt it)
31 {
32   KSP_GMRES   *gmres = (KSP_GMRES *)ksp->data;
33   PetscInt     j;
34   PetscScalar *hh, *hes;
35 
36   PetscFunctionBegin;
37   PetscCall(PetscLogEventBegin(KSP_GMRESOrthogonalization, ksp, 0, 0, 0));
38   /* update Hessenberg matrix and do Gram-Schmidt */
39   hh  = HH(0, it);
40   hes = HES(0, it);
41   for (j = 0; j <= it; j++) {
42     /* (vv(it+1), vv(j)) */
43     PetscCall(VecDot(VEC_VV(it + 1), VEC_VV(j), hh));
44     KSPCheckDot(ksp, *hh);
45     if (ksp->reason) break;
46     *hes++ = *hh;
47     /* vv(it+1) <- vv(it+1) - hh[it+1][j] vv(j) */
48     PetscCall(VecAXPY(VEC_VV(it + 1), -(*hh++), VEC_VV(j)));
49   }
50   PetscCall(PetscLogEventEnd(KSP_GMRESOrthogonalization, ksp, 0, 0, 0));
51   PetscFunctionReturn(PETSC_SUCCESS);
52 }
53