Actual source code: borthog.c

  1: #define PETSCKSP_DLL

  3: /*
  4:     Routines used for the orthogonalization of the Hessenberg matrix.

  6:     Note that for the complex numbers version, the VecDot() and
  7:     VecMDot() arguments within the code MUST remain in the order
  8:     given for correct computation of inner products.
  9: */
 10:  #include src/ksp/ksp/impls/gmres/gmresp.h

 12: /*@C
 13:      KSPGMRESModifiedGramSchmidtOrthogonalization -  This is the basic orthogonalization routine 
 14:                 using modified Gram-Schmidt.

 16:      Collective on KSP

 18:   Input Parameters:
 19: +   ksp - KSP object, must be associated with GMRES, FGMRES, or LGMRES Krylov method
 20: -   its - one less then the current GMRES restart iteration, i.e. the size of the Krylov space

 22:    Options Database Keys:
 23: .  -ksp_gmres_modifiedgramschmidt - Activates KSPGMRESModifiedGramSchmidtOrthogonalization()

 25:    Level: intermediate

 27: .seealso:  KSPGMRESSetOrthogonalization(), KSPGMRESClassicalGramSchmidtOrthogonalization()

 29: @*/
 32: PetscErrorCode  KSPGMRESModifiedGramSchmidtOrthogonalization(KSP ksp,PetscInt it)
 33: {
 34:   KSP_GMRES      *gmres = (KSP_GMRES *)(ksp->data);
 36:   PetscInt       j;
 37:   PetscScalar    *hh,*hes;

 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:     *hes++ = *hh;
 48:     /* vv(it+1) <- vv(it+1) - hh[it+1][j] vv(j) */
 49:     VecAXPY(VEC_VV(it+1),-(*hh++),VEC_VV(j));
 50:   }
 52:   return(0);
 53: }