Actual source code: borthog2.c
1: /*$Id: borthog2.c,v 1.20 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/ksp/ksp/impls/gmres/gmresp.h
11: /*
12: This version uses classical UNMODIFIED Gram-Schmidt. It has options for using
13: iterative refinement to improve stability.
15: */
18: int KSPGMRESClassicalGramSchmidtOrthogonalization(KSP ksp,int it)
19: {
20: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
21: int j,ierr;
22: PetscScalar *hh,*hes,shh[500],*lhh;
23: PetscReal hnrm, wnrm;
24: PetscTruth refine = (PetscTruth)(gmres->cgstype == KSP_GMRES_CGS_REFINE_ALWAYS);
27: PetscLogEventBegin(KSP_GMRESOrthogonalization,ksp,0,0,0);
28: /* Don't allocate small arrays */
29: if (it < 501) lhh = shh;
30: else {
31: PetscMalloc((it+1) * sizeof(PetscScalar),&lhh);
32: }
33:
34: /* update Hessenberg matrix and do unmodified Gram-Schmidt */
35: hh = HH(0,it);
36: hes = HES(0,it);
38: /* Clear hh and hes since we will accumulate values into them */
39: for (j=0; j<=it; j++) {
40: hh[j] = 0.0;
41: hes[j] = 0.0;
42: }
44: /*
45: This is really a matrix-vector product, with the matrix stored
46: as pointer to rows
47: */
48: VecMDot(it+1,VEC_VV(it+1),&(VEC_VV(0)),lhh); /* <v,vnew> */
49: for (j=0; j<=it; j++) {
50: lhh[j] = - lhh[j];
51: }
53: /*
54: This is really a matrix vector product:
55: [h[0],h[1],...]*[ v[0]; v[1]; ...] subtracted from v[it+1].
56: */
57: VecMAXPY(it+1,lhh,VEC_VV(it+1),&VEC_VV(0));
58: /* note lhh[j] is -<v,vnew> , hence the subtraction */
59: for (j=0; j<=it; j++) {
60: hh[j] -= lhh[j]; /* hh += <v,vnew> */
61: hes[j] -= lhh[j]; /* hes += <v,vnew> */
62: }
64: /*
65: * the second step classical Gram-Schmidt is only necessary
66: * when a simple test criteria is not passed
67: */
68: if (gmres->cgstype == KSP_GMRES_CGS_REFINE_IFNEEDED) {
69: hnrm = 0.0;
70: for (j=0; j<=it; j++) {
71: hnrm += PetscRealPart(lhh[j] * PetscConj(lhh[j]));
72: }
73: hnrm = sqrt(hnrm);
74: VecNorm(VEC_VV(it+1),NORM_2, &wnrm);
75: if (wnrm < 1.0286 * hnrm) {
76: refine = PETSC_TRUE;
77: PetscLogInfo(ksp,"KSPGMRESClassicalGramSchmidtOrthogonalization:Performing iterative refinement wnorm %g hnorm %g\n",wnrm,hnrm);
78: }
79: }
81: if (refine) {
82: VecMDot(it+1,VEC_VV(it+1),&(VEC_VV(0)),lhh); /* <v,vnew> */
83: for (j=0; j<=it; j++) lhh[j] = - lhh[j];
84: VecMAXPY(it+1,lhh,VEC_VV(it+1),&VEC_VV(0));
85: /* note lhh[j] is -<v,vnew> , hence the subtraction */
86: for (j=0; j<=it; j++) {
87: hh[j] -= lhh[j]; /* hh += <v,vnew> */
88: hes[j] -= lhh[j]; /* hes += <v,vnew> */
89: }
90: }
92: if (it >= 501) {PetscFree(lhh);}
93: PetscLogEventEnd(KSP_GMRESOrthogonalization,ksp,0,0,0);
94: return(0);
95: }