Actual source code: borthog3.c

  1: /*$Id: borthog3.c,v 1.25 2001/08/07 03:03:51 balay Exp $*/
  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/sles/ksp/impls/gmres/gmresp.h

 11: /*
 12:   This version uses 1 iteration of iterative refinement of UNMODIFIED Gram-Schmidt.  
 13:   It can give better performance when running in a parallel 
 14:   environment and in some cases even in a sequential environment (because
 15:   MAXPY has more data reuse).

 17:   Care is taken to accumulate the updated HH/HES values.
 18:  */
 19: int KSPGMRESIROrthogonalization(KSP  ksp,int it)
 20: {
 21:   KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
 22:   int       j,ncnt,ierr;
 23:   PetscScalar    *hh,*hes,shh[100],*lhh;

 26:   PetscLogEventBegin(KSP_GMRESOrthogonalization,ksp,0,0,0);
 27:   /* Don't allocate small arrays */
 28:   if (it < 100) lhh = shh;
 29:   else {
 30:     PetscMalloc((it+1) * sizeof(PetscScalar),&lhh);
 31:   }
 32: 
 33:   /* update Hessenberg matrix and do unmodified Gram-Schmidt */
 34:   hh  = HH(0,it);
 35:   hes = HES(0,it);

 37:   /* Clear hh and hes since we will accumulate values into them */
 38:   for (j=0; j<=it; j++) {
 39:     hh[j]  = 0.0;
 40:     hes[j] = 0.0;
 41:   }

 43:   ncnt = 0;
 44:   do {
 45:     /* 
 46:          This is really a matrix-vector product, with the matrix stored
 47:          as pointer to rows 
 48:     */
 49:     VecMDot(it+1,VEC_VV(it+1),&(VEC_VV(0)),lhh); /* <v,vnew> */

 51:     /*
 52:          This is really a matrix vector product: 
 53:          [h[0],h[1],...]*[ v[0]; v[1]; ...] subtracted from v[it+1].
 54:     */
 55:     for (j=0; j<=it; j++) lhh[j] = - lhh[j];
 56:     VecMAXPY(it+1,lhh,VEC_VV(it+1),&VEC_VV(0));
 57:     for (j=0; j<=it; j++) {
 58:       hh[j]  -= lhh[j];     /* hh += <v,vnew> */
 59:       hes[j] += lhh[j];     /* hes += - <v,vnew> */
 60:     }
 61:   } while (ncnt++ < 2);

 63:   if (it >= 100) {PetscFree(lhh);}
 64:   PetscLogEventEnd(KSP_GMRESOrthogonalization,ksp,0,0,0);
 65:   return(0);
 66: }