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