Actual source code: borthog.c

petsc-3.3-p7 2013-05-11
  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:    Level: intermediate

 26: .seealso:  KSPGMRESSetOrthogonalization(), KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESGetOrthogonalization()

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

 39:   PetscLogEventBegin(KSP_GMRESOrthogonalization,ksp,0,0,0);
 40:   /* update Hessenberg matrix and do Gram-Schmidt */
 41:   hh  = HH(0,it);
 42:   hes = HES(0,it);
 43:   for (j=0; j<=it; j++) {
 44:     /* (vv(it+1), vv(j)) */
 45:     VecDot(VEC_VV(it+1),VEC_VV(j),hh);
 46:     *hes++ = *hh;
 47:     /* vv(it+1) <- vv(it+1) - hh[it+1][j] vv(j) */
 48:     VecAXPY(VEC_VV(it+1),-(*hh++),VEC_VV(j));
 49:   }
 50:   PetscLogEventEnd(KSP_GMRESOrthogonalization,ksp,0,0,0);
 51:   return(0);
 52: }