Actual source code: iporthog.c

  1: /*
  2:      Routines related to orthogonalization.
  3:      See the SLEPc Technical Report STR-1 for a detailed explanation.

  5:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  6:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  7:    Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain

  9:    This file is part of SLEPc.
 10:       
 11:    SLEPc is free software: you can redistribute it and/or modify it under  the
 12:    terms of version 3 of the GNU Lesser General Public License as published by
 13:    the Free Software Foundation.

 15:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 16:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 17:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 18:    more details.

 20:    You  should have received a copy of the GNU Lesser General  Public  License
 21:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 22:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 23: */

 25:  #include private/ipimpl.h
 26:  #include slepcblaslapack.h

 28: /* 
 29:    IPOrthogonalizeMGS - Compute one step of Modified Gram-Schmidt 
 30: */
 33: static PetscErrorCode IPOrthogonalizeMGS(IP ip,PetscInt n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H)
 34: {
 36:   PetscInt       j;
 37: 
 39:   for (j=0; j<n; j++)
 40:     if (!which || which[j]) {
 41:       /* h_j = ( v, v_j ) */
 42:       IPInnerProduct(ip,v,V[j],&H[j]);
 43:       /* v <- v - h_j v_j */
 44:       VecAXPY(v,-H[j],V[j]);
 45:     }
 46:   return(0);
 47: }

 49: /* 
 50:    IPOrthogonalizeCGS - Compute |v'| (estimated), |v| and one step of CGS with only one global synchronization
 51: */
 54: PetscErrorCode IPOrthogonalizeCGS(IP ip,PetscInt n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H,PetscReal *onorm,PetscReal *norm,Vec work)
 55: {
 57:   PetscInt       j;
 58:   PetscScalar    alpha;
 59:   PetscReal      sum;
 60: 
 62:   /* h = W^* v ; alpha = (v , v) */
 63:   if (which) {
 64:   /* use select array */
 65:     for (j=0; j<n; j++)
 66:       if (which[j]) { IPInnerProductBegin(ip,v,V[j],&H[j]); }
 67:     if (onorm || norm) { IPInnerProductBegin(ip,v,v,&alpha); }
 68:     for (j=0; j<n; j++)
 69:       if (which[j]) { IPInnerProductEnd(ip,v,V[j],&H[j]); }
 70:     if (onorm || norm) { IPInnerProductEnd(ip,v,v,&alpha); }
 71:   } else { /* merge comunications */
 72:     if (onorm || norm) {
 73:       IPMInnerProductBegin(ip,v,n,V,H);
 74:       IPInnerProductBegin(ip,v,v,&alpha);
 75:       IPMInnerProductEnd(ip,v,n,V,H);
 76:       IPInnerProductEnd(ip,v,v,&alpha);
 77:     } else { /* use simpler function */
 78:       IPMInnerProduct(ip,v,n,V,H);
 79:     }
 80:   }

 82:   /* q = v - V h */
 83:   VecSet(work,0.0);
 84:   if (which) {
 85:     for (j=0; j<n; j++)
 86:       if (which[j]) { VecAXPY(work,H[j],V[j]); }
 87:   } else {
 88:     VecMAXPY(work,n,H,V);
 89:   }
 90:   VecAXPY(v,-1.0,work);
 91: 
 92:   /* compute |v| */
 93:   if (onorm) *onorm = sqrt(PetscRealPart(alpha));

 95:   /* compute |v'| */
 96:   if (norm) {
 97:     sum = 0.0;
 98:     for (j=0; j<n; j++)
 99:       if (!which || which[j])
100:         sum += PetscRealPart(H[j] * PetscConj(H[j]));
101:     *norm = PetscRealPart(alpha)-sum;
102:     if (*norm <= 0.0) {
103:       IPNorm(ip,v,norm);
104:     } else *norm = sqrt(*norm);
105:   }
106:   return(0);
107: }

109: /* 
110:    IPOrthogonalizeGS - Compute |v'|, |v| and one step of CGS or MGS 
111: */
114: static PetscErrorCode IPOrthogonalizeGS(IP ip,PetscInt n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H,PetscReal *onorm,PetscReal *norm,Vec work)
115: {
117: 
119:   switch (ip->orthog_type) {
120:   case IP_CGS_ORTH:
121:     IPOrthogonalizeCGS(ip,n,which,V,v,H,onorm,norm,work);
122:     break;
123:   case IP_MGS_ORTH:
124:     /* compute |v| */
125:     if (onorm) { IPNorm(ip,v,onorm); }
126:     IPOrthogonalizeMGS(ip,n,which,V,v,H);
127:     /* compute |v'| */
128:     if (norm) { IPNorm(ip,v,norm); }
129:     break;
130:   default:
131:     SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
132:   }
133:   return(0);
134: }

138: /*@
139:    IPOrthogonalize - Orthogonalize a vector with respect to a set of vectors.

141:    Collective on IP

143:    Input Parameters:
144: +  ip    - the inner product (IP) context
145: .  n      - number of columns of V
146: .  which  - logical array indicating columns of V to be used
147: .  V      - set of vectors
148: .  work   - workspace vector 
149: -  swork  - workspace array 

151:    Input/Output Parameter:
152: .  v      - (input) vector to be orthogonalized and (output) result of 
153:             orthogonalization

155:    Output Parameter:
156: +  H      - coefficients computed during orthogonalization
157: .  norm   - norm of the vector after being orthogonalized
158: -  lindep - flag indicating that refinement did not improve the quality
159:             of orthogonalization

161:    Notes:
162:    This function applies an orthogonal projector to project vector v onto the
163:    orthogonal complement of the span of the columns of V.

165:    On exit, v0 = [V v]*H, where v0 is the original vector v.

167:    This routine does not normalize the resulting vector.

169:    Level: developer

171: .seealso: IPSetOrthogonalization(), IPBiOrthogonalize()
172: @*/
173: PetscErrorCode IPOrthogonalize(IP ip,PetscInt n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscTruth *lindep,Vec work,PetscScalar* swork)
174: {
176:   PetscScalar    lh[100],*h,lc[100],*c;
177:   PetscTruth     allocatedh = PETSC_FALSE,allocatedc = PETSC_FALSE,allocatedw = PETSC_FALSE;
178:   PetscReal      onrm,nrm;
179:   PetscInt       j,k;
181:   if (n==0) {
182:     if (norm) { IPNorm(ip,v,norm); }
183:     if (lindep) *lindep = PETSC_FALSE;
184:     return(0);
185:   }
186:   PetscLogEventBegin(IP_Orthogonalize,ip,0,0,0);

188:   /* allocate H, c and work if needed */
189:   if (!H) {
190:     if (n<=100) h = lh;
191:     else {
192:       PetscMalloc(n*sizeof(PetscScalar),&h);
193:       allocatedh = PETSC_TRUE;
194:     }
195:   } else h = H;
196:   if (!swork) {
197:     if (ip->orthog_ref != IP_ORTH_REFINE_NEVER) {
198:       if (n<=100) c = lc;
199:       else {
200:         PetscMalloc(n*sizeof(PetscScalar),&c);
201:         allocatedc = PETSC_TRUE;
202:       }
203:     }
204:   } else c = swork;
205:   if (!work && ip->orthog_type == IP_CGS_ORTH) {
206:     VecDuplicate(v,&work);
207:     allocatedw = PETSC_TRUE;
208:   }

210:   /* orthogonalize and compute onorm */
211:   switch (ip->orthog_ref) {
212: 
213:   case IP_ORTH_REFINE_NEVER:
214:     IPOrthogonalizeGS(ip,n,which,V,v,h,PETSC_NULL,PETSC_NULL,work);
215:     /* compute |v| */
216:     if (norm) { IPNorm(ip,v,norm); }
217:     /* linear dependence check does not work without refinement */
218:     if (lindep) *lindep = PETSC_FALSE;
219:     break;
220: 
221:   case IP_ORTH_REFINE_ALWAYS:
222:     IPOrthogonalizeGS(ip,n,which,V,v,h,PETSC_NULL,PETSC_NULL,work);
223:     if (lindep) {
224:       IPOrthogonalizeGS(ip,n,which,V,v,c,&onrm,&nrm,work);
225:       if (norm) *norm = nrm;
226:       if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
227:       else *lindep = PETSC_FALSE;
228:     } else {
229:       IPOrthogonalizeGS(ip,n,which,V,v,c,PETSC_NULL,norm,work);
230:     }
231:     for (j=0;j<n;j++)
232:       if (!which || which[j]) h[j] += c[j];
233:     break;
234: 
235:   case IP_ORTH_REFINE_IFNEEDED:
236:     IPOrthogonalizeGS(ip,n,which,V,v,h,&onrm,&nrm,work);
237:     /* ||q|| < eta ||h|| */
238:     k = 1;
239:     while (k<3 && nrm < ip->orthog_eta * onrm) {
240:       k++;
241:       switch (ip->orthog_type) {
242:       case IP_CGS_ORTH:
243:         IPOrthogonalizeGS(ip,n,which,V,v,c,&onrm,&nrm,work);
244:         break;
245:       case IP_MGS_ORTH:
246:         onrm = nrm;
247:         IPOrthogonalizeGS(ip,n,which,V,v,c,PETSC_NULL,&nrm,work);
248:         break;
249:       default:
250:         SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
251:       }
252:       for (j=0;j<n;j++)
253:         if (!which || which[j]) h[j] += c[j];
254:     }
255:     if (norm) *norm = nrm;
256:     if (lindep) {
257:       if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
258:       else *lindep = PETSC_FALSE;
259:     }
260: 
261:     break;
262:   default:
263:     SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization refinement");
264:   }

266:   /* free work space */
267:   if (allocatedc) { PetscFree(c); }
268:   if (allocatedh) { PetscFree(h); }
269:   if (allocatedw) { VecDestroy(work); }
270: 
271:   PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);
272:   return(0);
273: }

277: /*@
278:    IPQRDecomposition - Compute the QR factorization of a set of vectors.

280:    Collective on IP

282:    Input Parameters:
283: +  ip - the eigenproblem solver context
284: .  V - set of vectors
285: .  m - starting column
286: .  n - ending column
287: .  ldr - leading dimension of R
288: -  work - workspace vector 

290:    Output Parameter:
291: .  R  - triangular matrix of QR factorization

293:    Notes:
294:    This routine orthonormalizes the columns of V so that V'*V=I up to 
295:    working precision. It assumes that the first m columns (from 0 to m-1) 
296:    are already orthonormal. The coefficients of the orthogonalization are
297:    stored in R so that V*R is equal to the original V.

299:    In some cases, this routine makes V B-orthonormal, that is, V'*B*V=I.

301:    If one of the vectors is linearly dependent wrt the rest (at working
302:    precision) then it is replaced by a random vector.

304:    Level: developer

306: .seealso: IPOrthogonalize(), IPNorm(), IPInnerProduct().
307: @*/
308: PetscErrorCode IPQRDecomposition(IP ip,Vec *V,PetscInt m,PetscInt n,PetscScalar *R,PetscInt ldr,Vec work)
309: {
311:   PetscInt       k;
312:   PetscReal      norm;
313:   PetscTruth     lindep;
314: 

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

319:     /* orthogonalize v_k with respect to v_0, ..., v_{k-1} */
320:     if (R) { IPOrthogonalize(ip,k,PETSC_NULL,V,V[k],&R[0+ldr*k],&norm,&lindep,work,PETSC_NULL); }
321:     else   { IPOrthogonalize(ip,k,PETSC_NULL,V,V[k],PETSC_NULL,&norm,&lindep,work,PETSC_NULL); }

323:     /* normalize v_k: r_{k,k} = ||v_k||_2; v_k = v_k/r_{k,k} */
324:     if (norm==0.0 || lindep) {
325:       PetscInfo(ip,"Linearly dependent vector found, generating a new random vector\n");
326:       SlepcVecSetRandom(V[k]);
327:       IPNorm(ip,V[k],&norm);
328:     }
329:     VecScale(V[k],1.0/norm);
330:     if (R) R[k+ldr*k] = norm;

332:   }

334:   return(0);
335: }

337: /*
338:     Biorthogonalization routine using classical Gram-Schmidt with refinement.
339:  */
342: static PetscErrorCode IPCGSBiOrthogonalization(IP ip,PetscInt n_,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *hnorm,PetscReal *norm)
343: {
344: #if defined(SLEPC_MISSING_LAPACK_GELQF) || defined(SLEPC_MISSING_LAPACK_ORMLQ)
346:   SETERRQ(PETSC_ERR_SUP,"xGELQF - Lapack routine is unavailable.");
347: #else
349:   PetscBLASInt   j,ione=1,lwork,info,n=n_;
350:   PetscScalar    shh[100],*lhh,*vw,*tau,one=1.0,*work;
351:   Vec            w;


355:   /* Don't allocate small arrays */
356:   if (n<=100) lhh = shh;
357:   else { PetscMalloc(n*sizeof(PetscScalar),&lhh); }
358:   PetscMalloc(n*n*sizeof(PetscScalar),&vw);
359:   VecDuplicate(v,&w);
360: 
361:   for (j=0;j<n;j++) {
362:     IPMInnerProduct(ip,V[j],n,W,vw+j*n);
363:   }
364:   lwork = n;
365:   PetscMalloc(n*sizeof(PetscScalar),&tau);
366:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
367:   LAPACKgelqf_(&n,&n,vw,&n,tau,work,&lwork,&info);
368:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGELQF %i",info);
369: 
370:   /*** First orthogonalization ***/

372:   /* h = W^* v */
373:   /* q = v - V h */
374:   IPMInnerProduct(ip,v,n,W,H);
375:   BLAStrsm_("L","L","N","N",&n,&ione,&one,vw,&n,H,&n);
376:   LAPACKormlq_("L","N",&n,&ione,&n,vw,&n,tau,H,&n,work,&lwork,&info);
377:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xORMLQ %i",info);
378:   VecSet(w,0.0);
379:   VecMAXPY(w,n,H,V);
380:   VecAXPY(v,-1.0,w);

382:   /* compute norm of v */
383:   if (norm) { IPNorm(ip,v,norm); }
384: 
385:   if (n>100) { PetscFree(lhh); }
386:   PetscFree(vw);
387:   PetscFree(tau);
388:   PetscFree(work);
389:   VecDestroy(w);
390:   return(0);
391: #endif
392: }

396: /*@
397:    IPBiOrthogonalize - Bi-orthogonalize a vector with respect to a set of vectors.

399:    Collective on IP

401:    Input Parameters:
402: +  ip - the inner product context
403: .  n - number of columns of V
404: .  V - set of vectors
405: -  W - set of vectors

407:    Input/Output Parameter:
408: .  v - vector to be orthogonalized

410:    Output Parameter:
411: +  H  - coefficients computed during orthogonalization
412: -  norm - norm of the vector after being orthogonalized

414:    Notes:
415:    This function applies an oblique projector to project vector v onto the
416:    span of the columns of V along the orthogonal complement of the column
417:    space of W. 

419:    On exit, v0 = [V v]*H, where v0 is the original vector v.

421:    This routine does not normalize the resulting vector.

423:    Level: developer

425: .seealso: IPSetOrthogonalization(), IPOrthogonalize()
426: @*/
427: PetscErrorCode IPBiOrthogonalize(IP ip,PetscInt n,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *norm)
428: {
430:   PetscScalar    lh[100],*h;
431:   PetscTruth     allocated = PETSC_FALSE;
432:   PetscReal      lhnrm,*hnrm,lnrm,*nrm;
434:   if (n==0) {
435:     if (norm) { IPNorm(ip,v,norm); }
436:   } else {
437:     PetscLogEventBegin(IP_Orthogonalize,ip,0,0,0);
438: 
439:     /* allocate H if needed */
440:     if (!H) {
441:       if (n<=100) h = lh;
442:       else {
443:         PetscMalloc(n*sizeof(PetscScalar),&h);
444:         allocated = PETSC_TRUE;
445:       }
446:     } else h = H;
447: 
448:     /* retrieve hnrm and nrm for linear dependence check or conditional refinement */
449:     if (ip->orthog_ref == IP_ORTH_REFINE_IFNEEDED) {
450:       hnrm = &lhnrm;
451:       if (norm) nrm = norm;
452:       else nrm = &lnrm;
453:     } else {
454:       hnrm = PETSC_NULL;
455:       nrm = norm;
456:     }
457: 
458:     switch (ip->orthog_type) {
459:       case IP_CGS_ORTH:
460:         IPCGSBiOrthogonalization(ip,n,V,W,v,h,hnrm,nrm);
461:         break;
462:       default:
463:         SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
464:     }
465: 
466:     if (allocated) { PetscFree(h); }
467: 
468:     PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);
469:   }
470:   return(0);
471: }