Actual source code: borthog.c

  1: /*
  2:     Routines related to orthogonalization.

  4:     See the SLEPc users manual for a detailed explanation.
  5: */
 6:  #include src/eps/epsimpl.h

 10: /*@
 11:    EPSQRDecomposition - Compute the QR factorization of the basis vectors.

 13:    Collective on EPS

 15:    Input Parameters:
 16: +  eps - the eigenproblem solver context
 17: .  m - starting column
 18: .  n - ending column
 19: -  ldr - leading dimension of R

 21:    Output Parameter:
 22: .  R  - triangular matrix of QR factorization

 24:    Level: developer

 26: @*/
 27: int EPSQRDecomposition(EPS eps,int m,int n,PetscScalar *R,int ldr)
 28: {
 29:   int         ierr,k;
 30:   PetscScalar alpha;
 31:   PetscReal   norm;
 32: 

 35:   /* normalize v_m: r_{m,m} = ||v_m||_2; v_m = v_m/r_{m,m} */
 36:   VecNorm(eps->V[m],NORM_2,&norm);
 37:   if (R) { R[m+ldr*m] = norm; }
 38:   if (norm==0.0) SETERRQ( 1,"Zero vector in QR decomposition" );
 39:   alpha = 1.0/norm;
 40:   VecScale(&alpha,eps->V[m]);

 42:   for (k=m+1; k<n; k++) {

 44:     /* orthogonalize v_k with respect to v_0, ..., v_{k-1} */
 45:     PetscLogEventBegin(EPS_Orthogonalization,eps,0,0,0);
 46:     if (R) { (*eps->orthog)(eps,k-1,&R[0+ldr*k],&norm); }
 47:     else   { (*eps->orthog)(eps,k-1,PETSC_NULL,&norm); }
 48:     PetscLogEventEnd(EPS_Orthogonalization,eps,0,0,0);

 50:     /* normalize v_k: r_{k,k} = ||v_k||_2; v_k = v_k/r_{k,k} */
 51:     if (R) { R[k+ldr*k] = norm; }
 52:     if (norm==0.0) SETERRQ( 1,"Zero vector in QR decomposition" );
 53:     alpha = 1.0/norm;
 54:     VecScale(&alpha,eps->V[k]);

 56:   }

 58:   return(0);
 59: }

 63: /*@
 64:    EPSSetOrthogonalization - Specifies the type of orthogonalization technique
 65:    to be used inside the eigensolver.

 67:    Collective on EPS

 69:    Input Parameters:
 70: +  eps        - the eigensolver context 
 71: .  type       - a known type of orthogonalization
 72: .  refinement - type of refinement
 73: -  eta        - parameter for dynamic refinement

 75:    Options Database Keys:
 76: +  -eps_orthog_type <type> -  Where <type> is cgs for Classical Gram-Schmidt
 77:                               orthogonalization (default)
 78:                               or mgs for Modified Gram-Schmidt orthogonalization
 79: .  -eps_orthog_refinement <type> -  Where <type> is one of never, ifneeded
 80:                               (default) or always 
 81: -  -eps_orthog_eta <eta> -  For setting the value of eta (or PETSC_DEFAULT)
 82:     
 83:    Notes:  
 84:    The value of eta is used only when refinement type is "ifneeded". 

 86:    The default orthogonalization technique 
 87:    works well for most problems. MGS is numerically more robust than CGS,
 88:    but CGS may give better scalability.

 90:    Level: intermediate

 92: .seealso: EPSGetOrthogonalization()
 93: @*/
 94: int EPSSetOrthogonalization(EPS eps,EPSOrthogonalizationType type, EPSOrthogonalizationRefinementType refinement, PetscReal eta)
 95: {

 99:   switch (type) {
100:     case EPS_CGS_ORTH:
101:       eps->orthog = EPSClassicalGramSchmidtOrthogonalization;
102:       break;
103:     case EPS_MGS_ORTH:
104:       eps->orthog = EPSModifiedGramSchmidtOrthogonalization;
105:       break;
106:     default:
107:       SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
108:   }
109:   eps->orth_type = type;
110:   switch (refinement) {
111:   case EPS_ORTH_REFINE_NEVER:
112:     eps->orth_eta = 0;
113:     break;
114:   case EPS_ORTH_REFINE_IFNEEDED:
115:     if (eta == PETSC_DEFAULT) eps->orth_eta = 0.7071;
116:     else eps->orth_eta = eta;
117:     break;
118:   case EPS_ORTH_REFINE_ALWAYS:
119:     eps->orth_eta = -1;
120:     break;
121:   default:
122:     SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown refinement type");
123:   }
124: 
125:   return(0);
126: }

130: /*@
131:    EPSGetOrthogonalization - Gets the orthogonalization type from the 
132:    EPS object.

134:    Not Collective

136:    Input Parameter:
137: .  eps - Eigensolver context 

139:    Output Parameter:
140: +  type       - type of orthogonalization technique
141: .  refinement - type of refinement
142: -  eta        - parameter for dynamic refinement

144:    Level: intermediate

146: .seealso: EPSSetOrthogonalization()
147: @*/
148: int EPSGetOrthogonalization(EPS eps,EPSOrthogonalizationType *type,EPSOrthogonalizationRefinementType *refinement, PetscReal *eta)
149: {
152:   if (type) *type = eps->orth_type;
153:   if (eps->orth_eta == 0) {
154:     if (refinement) *refinement = EPS_ORTH_REFINE_NEVER;
155:   } else if (eps->orth_eta < 0) {
156:     if (refinement) *refinement = EPS_ORTH_REFINE_ALWAYS;
157:   } else {
158:     if (refinement) *refinement = EPS_ORTH_REFINE_IFNEEDED;
159:     if (eta) *eta = eps->orth_eta;
160:   }
161:   return(0);
162: }

164: /*
165:     Orthogonalization routine using classical Gram-Schmidt with refinement.
166:  */
169: int EPSClassicalGramSchmidtOrthogonalization(EPS eps,int it,PetscScalar *H,PetscReal *norm)
170: {
171:   int         ierr,j;
172:   PetscScalar shh[100],shh2[100],*lhh;
173:   PetscTruth  alloc = PETSC_FALSE;
174:   PetscReal   hnorm;

177:   if (!H) {
178:     if (it<100) H = shh2;   /* Don't allocate small arrays */
179:     else {
180:       PetscMalloc((it+1)*sizeof(PetscScalar),&H);
181:       alloc = PETSC_TRUE;
182:     }
183:   }
184:   /* Don't allocate small arrays */
185:   if (it<100) lhh = shh;
186:   else { PetscMalloc((it+1)*sizeof(PetscScalar),&lhh); }
187: 
188:   /*** First orthogonalization ***/

190:   /* h = V^* v */
191:   VecMDot(it+1,eps->V[it+1],eps->V,lhh);

193:   /* q = v - V h */
194:   for (j=0;j<=it;j++) {
195:     H[j] = lhh[j];
196:     lhh[j] = -lhh[j];
197:   }
198:   VecMAXPY(it+1,lhh,eps->V[it+1],eps->V);
199: 
200:   VecNorm(eps->V[it+1], NORM_2, norm);

202:   /*** Second orthogonalization if necessary ***/
203:   if (eps->orth_eta != 0.0) {
204:     hnorm = 0.0;
205:     for (j=0; j<=it; j++) {
206:       hnorm  +=  PetscRealPart(H[j] * PetscConj(H[j]));
207:     }
208:     hnorm = sqrt(hnorm);

210:     /* if ||q|| < eta ||h|| */
211:     if (eps->orth_eta < 0 || *norm < eps->orth_eta * hnorm) {
212:       PetscLogInfo(eps,"EPSClassicalGramSchmidtOrthogonalization:Performing iterative refinement wnorm %g hnorm %g\n",*norm,hnorm);
213: 
214:       /* s = V^* q */
215:       VecMDot(it+1,eps->V[it+1],eps->V,lhh); /* <v,vnew> */

217:       /* q = q - V s  ;  h = h + s */
218:       for (j=0;j<=it;j++) {
219:         H[j] += lhh[j];
220:         lhh[j] = -lhh[j];
221:       }
222:       VecMAXPY(it+1,lhh,eps->V[it+1],eps->V);

224:       VecNorm(eps->V[it+1], NORM_2, norm);
225:     }
226:   }

228:   if (it>=100) { PetscFree(lhh); }
229:   if (alloc) { PetscFree(H); }
230:   return(0);
231: }

233: /*
234:     Orthogonalization routine using modified Gram-Schmidt with refinement.
235:  */
238: int EPSModifiedGramSchmidtOrthogonalization(EPS eps,int it,PetscScalar *H,PetscReal *norm)
239: {
240:   int         ierr,j;
241:   PetscScalar alpha;
242:   PetscTruth  refinement,allocated;
243:   PetscScalar lh[100],*h;
244:   PetscReal   hnorm;
245: 
247: 
248:   allocated = PETSC_FALSE;
249:   if (!H) {
250:     if (it<100) h = lh;
251:     else {
252:       PetscMalloc((it+1)*sizeof(PetscScalar),&h);
253:       allocated = PETSC_TRUE;
254:     }
255:   } else h = H;
256: 
257:   for (j=0; j<=it; j++) {
258:     h[j] = 0;
259:   }
260:   refinement = PETSC_FALSE;
261:   do {
262:     for (j=0; j<=it; j++) {
263:       /* alpha = ( v_{it+1}, v_j ) */
264:       VecDot(eps->V[it+1],eps->V[j],&alpha);
265:       /* store coefficients if requested */
266:       h[j] += alpha;
267:       /* v_{it+1} <- v_{it+1} - alpha v_j */
268:       alpha = -alpha;
269:       VecAXPY(&alpha,eps->V[j],eps->V[it+1]);
270:     }
271:     VecNorm(eps->V[it+1], NORM_2, norm);
272:     if (refinement) refinement = PETSC_FALSE;
273:     else if (eps->orth_eta != 0.0) {
274:       hnorm = 0.0;
275:       for (j=0; j<=it; j++) {
276:         hnorm += PetscRealPart(h[j] * PetscConj(h[j]));
277:       }
278:       hnorm = sqrt(hnorm);
279:       if (eps->orth_eta < 0.0 || *norm < eps->orth_eta * hnorm) {
280:         PetscLogInfo(eps,"EPSModifiedGramSchmidtOrthogonalization:Performing iterative refinement wnorm %g hnorm %g\n",*norm,hnorm);
281:         refinement = PETSC_TRUE;
282:       }
283:     }
284:   } while (refinement);
285: 
286:   if (allocated) {
287:     PetscFree(h);
288:   }
289:   return(0);
290: }