Actual source code: bvorthog.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    BV orthogonalization routines.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

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

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

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

 24: #include <slepc-private/bvimpl.h>          /*I   "slepcbv.h"   I*/

 28: /*
 29:    BVOrthogonalizeMGS1 - Compute one step of Modified Gram-Schmidt
 30: */
 31: static PetscErrorCode BVOrthogonalizeMGS1(BV bv,PetscInt k,Vec v,PetscBool *which,PetscScalar *H)
 32: {
 34:   PetscInt       i;
 35:   PetscScalar    dot;
 36:   Vec            vi,z;

 39:   z = v;
 40:   for (i=-bv->nc;i<k;i++) {
 41:     if (which && i>=0 && !which[i]) continue;
 42:     BVGetColumn(bv,i,&vi);
 43:     /* h_i = ( v, v_i ) */
 44:     if (bv->matrix) {
 45:       BV_IPMatMult(bv,v);
 46:       z = bv->Bx;
 47:     }
 48:     VecDot(z,vi,&dot);
 49:     /* v <- v - h_i v_i */
 50:     if (bv->indef) dot /= bv->omega[bv->nc+i];
 51:     VecAXPY(v,-dot,vi);
 52:     if (bv->indef) dot *= bv->omega[bv->nc+i];
 53:     if (H) H[bv->nc+i] += dot;
 54:     BVRestoreColumn(bv,i,&vi);
 55:   }
 56:   return(0);
 57: }

 61: /*
 62:    BVOrthogonalizeCGS1 - Compute |v'| (estimated), |v| and one step of CGS with
 63:    only one global synchronization
 64: */
 65: PetscErrorCode BVOrthogonalizeCGS1(BV bv,PetscInt j,Vec v,PetscScalar *H,PetscReal *onorm,PetscReal *norm)
 66: {
 68:   PetscInt       i;
 69:   PetscReal      sum,nrm,beta;
 70:   Vec            w=v;

 73:   /* h = W^* v ; alpha = (v, v) */
 74:   bv->k = j;
 75:   if (onorm || norm) {
 76:     if (!v) {
 77:       bv->k++;
 78:       BVGetColumn(bv,j,&w);
 79:     }
 80:     BVDotVec(bv,w,H);
 81:     if (!v) {
 82:       BVRestoreColumn(bv,j,&w);
 83:       bv->k--;
 84:       beta = PetscSqrtReal(PetscRealPart(H[bv->nc+j]));
 85:     } else {
 86:       BVNormVec(bv,w,NORM_2,&beta);
 87:     }
 88:   } else {
 89:     if (!v) { BVDotColumn(bv,j,H); }
 90:     else { BVDotVec(bv,w,H); }
 91:   }

 93:   /* q = v - V h */
 94:   if (bv->indef) {
 95:     for (i=0;i<bv->nc+j;i++) H[i] /= bv->omega[i];  /* apply inverse of signature */
 96:   }
 97:   if (!v) { BVMultColumn(bv,-1.0,1.0,j,H); }
 98:   else { BVMultVec(bv,-1.0,1.0,w,H); }
 99:   if (bv->indef) {
100:     for (i=0;i<bv->nc+j;i++) H[i] *= bv->omega[i];  /* revert signature */
101:   }

103:   /* compute |v| */
104:   if (onorm) *onorm = beta;

106:   if (bv->indef) {
107:     if (!v) { BVNormColumn(bv,j,NORM_2,&nrm); }
108:     else { BVNormVec(bv,w,NORM_2,&nrm); }
109:     if (norm) *norm = nrm;
110:     bv->omega[bv->nc+j] = (nrm<0.0)? -1.0: 1.0;
111:   } else if (norm) {
112:     /* estimate |v'| from |v| */
113:     sum = 0.0;
114:     for (i=0;i<bv->nc+j;i++) sum += PetscRealPart(H[i]*PetscConj(H[i]));
115:     *norm = beta*beta-sum;
116:     if (*norm <= 0.0) {
117:       if (!v) { BVNormColumn(bv,j,NORM_2,norm); }
118:       else { BVNormVec(bv,w,NORM_2,norm); }
119:     } else *norm = PetscSqrtReal(*norm);
120:   }
121:   return(0);
122: }

126: /*
127:   BVOrthogonalizeMGS - Orthogonalize with modified Gram-Schmidt
128: */
129: static PetscErrorCode BVOrthogonalizeMGS(BV bv,PetscInt j,Vec v,PetscBool *which,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
130: {
132:   PetscReal      onrm,nrm;
133:   PetscInt       k,l;
134:   Vec            w;

137:   if (v) {
138:     w = v;
139:     k = bv->k;
140:   } else {
141:     BVGetColumn(bv,j,&w);
142:     k = j;
143:   }
144:   PetscMemzero(bv->h,(bv->nc+k)*sizeof(PetscScalar));
145:   switch (bv->orthog_ref) {

147:   case BV_ORTHOG_REFINE_IFNEEDED:
148:     /* first step */
149:     BVNormVec(bv,w,NORM_2,&onrm);
150:     BVOrthogonalizeMGS1(bv,k,w,which,bv->h);
151:     BVNormVec(bv,w,NORM_2,&nrm);
152:     /* ||q|| < eta ||h|| */
153:     l = 1;
154:     while (l<3 && nrm && nrm < bv->orthog_eta*onrm) {
155:       l++;
156:       onrm = nrm;
157:       BVOrthogonalizeMGS1(bv,k,w,which,bv->c);
158:       BVNormVec(bv,w,NORM_2,&nrm);
159:     }
160:     if (lindep) {
161:       if (nrm < bv->orthog_eta*onrm) *lindep = PETSC_TRUE;
162:       else *lindep = PETSC_FALSE;
163:     }
164:     break;

166:   case BV_ORTHOG_REFINE_NEVER:
167:     BVOrthogonalizeMGS1(bv,k,w,which,bv->h);
168:     /* compute |v| */
169:     if (norm || lindep) {
170:       BVNormVec(bv,w,NORM_2,&nrm);
171:     }
172:     /* linear dependence check: just test for exactly zero norm */
173:     if (lindep) *lindep = nrm? PETSC_FALSE: PETSC_TRUE;
174:     break;

176:   case BV_ORTHOG_REFINE_ALWAYS:
177:     /* first step */
178:     BVOrthogonalizeMGS1(bv,k,w,which,bv->h);
179:     if (lindep) {
180:       BVNormVec(bv,w,NORM_2,&onrm);
181:     }
182:     /* second step */
183:     BVOrthogonalizeMGS1(bv,k,w,which,bv->h);
184:     if (norm || lindep) {
185:       BVNormVec(bv,w,NORM_2,&nrm);
186:     }
187:     if (lindep) {
188:       if (nrm==0.0 || nrm < bv->orthog_eta*onrm) *lindep = PETSC_TRUE;
189:       else *lindep = PETSC_FALSE;
190:     }
191:     break;
192:   }
193:   if (bv->indef) {
194:     BVNormVec(bv,w,NORM_2,&nrm);
195:     bv->omega[bv->nc+j] = (nrm<0.0)? -1.0: 1.0;
196:   }
197:   if (!v) { BVRestoreColumn(bv,j,&w); }
198:   if (norm) *norm = nrm;
199:   return(0);
200: }

204: /*
205:   BVOrthogonalizeCGS - Orthogonalize with classical Gram-Schmidt
206: */
207: static PetscErrorCode BVOrthogonalizeCGS(BV bv,PetscInt j,Vec v,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
208: {
210:   PetscReal      onrm,nrm;
211:   PetscInt       i,k,l;

214:   if (v) k = bv->k;
215:   else k = j;
216:   switch (bv->orthog_ref) {

218:   case BV_ORTHOG_REFINE_IFNEEDED:
219:     BVOrthogonalizeCGS1(bv,k,v,bv->h,&onrm,&nrm);
220:     /* ||q|| < eta ||h|| */
221:     l = 1;
222:     while (l<3 && nrm && nrm < bv->orthog_eta*onrm) {
223:       l++;
224:       BVOrthogonalizeCGS1(bv,k,v,bv->c,&onrm,&nrm);
225:       for (i=0;i<bv->nc+k;i++) bv->h[i] += bv->c[i];
226:     }
227:     if (norm) *norm = nrm;
228:     if (lindep) {
229:       if (nrm < bv->orthog_eta*onrm) *lindep = PETSC_TRUE;
230:       else *lindep = PETSC_FALSE;
231:     }
232:     break;

234:   case BV_ORTHOG_REFINE_NEVER:
235:     BVOrthogonalizeCGS1(bv,k,v,bv->h,NULL,NULL);
236:     /* compute |v| */
237:     if (norm || lindep) {
238:       if (v) { BVNormVec(bv,v,NORM_2,&nrm); }
239:       else { BVNormColumn(bv,k,NORM_2,&nrm); }
240:     }
241:     if (norm) *norm = nrm;
242:     /* linear dependence check: just test for exactly zero norm */
243:     if (lindep) *lindep = nrm? PETSC_FALSE: PETSC_TRUE;
244:     break;

246:   case BV_ORTHOG_REFINE_ALWAYS:
247:     BVOrthogonalizeCGS1(bv,k,v,bv->h,NULL,NULL);
248:     if (lindep) {
249:       BVOrthogonalizeCGS1(bv,k,v,bv->c,&onrm,&nrm);
250:       if (norm) *norm = nrm;
251:       if (nrm==0.0 || nrm < bv->orthog_eta*onrm) *lindep = PETSC_TRUE;
252:       else *lindep = PETSC_FALSE;
253:     } else {
254:       BVOrthogonalizeCGS1(bv,k,v,bv->c,NULL,norm);
255:     }
256:     for (i=0;i<bv->nc+k;i++) bv->h[i] += bv->c[i];
257:     break;
258:   }
259:   return(0);
260: }

264: /*@
265:    BVOrthogonalizeVec - Orthogonalize a given vector with respect to all
266:    active columns.

268:    Collective on BV

270:    Input Parameters:
271: +  bv     - the basis vectors context
272: -  v      - the vector

274:    Output Parameters:
275: +  H      - (optional) coefficients computed during orthogonalization
276: .  norm   - (optional) norm of the vector after being orthogonalized
277: -  lindep - (optional) flag indicating that refinement did not improve the quality
278:             of orthogonalization

280:    Notes:
281:    This function is equivalent to BVOrthogonalizeColumn() but orthogonalizes
282:    a vector as an argument rather than taking one of the BV columns. The
283:    vector is orthogonalized against all active columns.

285:    Level: advanced

287: .seealso: BVOrthogonalizeColumn(), BVSetOrthogonalization(), BVSetActiveColumns()
288: @*/
289: PetscErrorCode BVOrthogonalizeVec(BV bv,Vec v,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
290: {
292:   PetscInt       i,ksave,lsave;

298:   BVCheckSizes(bv,1);

302:   PetscLogEventBegin(BV_Orthogonalize,bv,0,0,0);
303:   ksave = bv->k;
304:   lsave = bv->l;
305:   bv->l = -bv->nc;  /* must also orthogonalize against constraints and leading columns */
306:   BV_AllocateCoeffs(bv);
307:   BV_AllocateSignature(bv);
308:   switch (bv->orthog_type) {
309:   case BV_ORTHOG_CGS:
310:     BVOrthogonalizeCGS(bv,0,v,H,norm,lindep);
311:     break;
312:   case BV_ORTHOG_MGS:
313:     BVOrthogonalizeMGS(bv,0,v,NULL,H,norm,lindep);
314:     break;
315:   }
316:   bv->k = ksave;
317:   bv->l = lsave;
318:   if (H) for (i=bv->l;i<bv->k;i++) H[i-bv->l] = bv->h[bv->nc+i];
319:   PetscLogEventEnd(BV_Orthogonalize,bv,0,0,0);
320:   return(0);
321: }

325: /*@
326:    BVOrthogonalizeColumn - Orthogonalize one of the column vectors with respect to
327:    the previous ones.

329:    Collective on BV

331:    Input Parameters:
332: +  bv     - the basis vectors context
333: -  j      - index of column to be orthogonalized

335:    Output Parameters:
336: +  H      - (optional) coefficients computed during orthogonalization
337: .  norm   - (optional) norm of the vector after being orthogonalized
338: -  lindep - (optional) flag indicating that refinement did not improve the quality
339:             of orthogonalization

341:    Notes:
342:    This function applies an orthogonal projector to project vector V[j] onto
343:    the orthogonal complement of the span of the columns of V[0..j-1],
344:    where V[.] are the vectors of BV. The columns V[0..j-1] are assumed to be
345:    mutually orthonormal.

347:    Leading columns V[0..l-1] also participate in the orthogonalization.

349:    If a non-standard inner product has been specified with BVSetMatrix(),
350:    then the vector is B-orthogonalized, using the non-standard inner product
351:    defined by matrix B. The output vector satisfies V[j]'*B*V[0..j-1] = 0.

353:    This routine does not normalize the resulting vector.

355:    Level: advanced

357: .seealso: BVSetOrthogonalization(), BVSetMatrix(), BVSetActiveColumns(), BVOrthogonalize(), BVOrthogonalizeVec()
358: @*/
359: PetscErrorCode BVOrthogonalizeColumn(BV bv,PetscInt j,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
360: {
362:   PetscInt       i,ksave,lsave;

368:   BVCheckSizes(bv,1);
369:   if (j<0) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
370:   if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Index j=%D but BV only has %D columns",j,bv->m);

372:   PetscLogEventBegin(BV_Orthogonalize,bv,0,0,0);
373:   ksave = bv->k;
374:   lsave = bv->l;
375:   bv->l = -bv->nc;  /* must also orthogonalize against constraints and leading columns */
376:   BV_AllocateCoeffs(bv);
377:   BV_AllocateSignature(bv);
378:   switch (bv->orthog_type) {
379:   case BV_ORTHOG_CGS:
380:     BVOrthogonalizeCGS(bv,j,NULL,H,norm,lindep);
381:     break;
382:   case BV_ORTHOG_MGS:
383:     BVOrthogonalizeMGS(bv,j,NULL,NULL,H,norm,lindep);
384:     break;
385:   }
386:   bv->k = ksave;
387:   bv->l = lsave;
388:   if (H) for (i=bv->l;i<j;i++) H[i-bv->l] = bv->h[bv->nc+i];
389:   PetscLogEventEnd(BV_Orthogonalize,bv,0,0,0);
390:   return(0);
391: }

395: /*@
396:    BVOrthogonalizeSomeColumn - Orthogonalize one of the column vectors with
397:    respect to some of the previous ones.

399:    Collective on BV

401:    Input Parameters:
402: +  bv     - the basis vectors context
403: .  j      - index of column to be orthogonalized
404: -  which  - logical array indicating selected columns

406:    Output Parameters:
407: +  H      - (optional) coefficients computed during orthogonalization
408: .  norm   - (optional) norm of the vector after being orthogonalized
409: -  lindep - (optional) flag indicating that refinement did not improve the quality
410:             of orthogonalization

412:    Notes:
413:    This function is similar to BVOrthogonalizeColumn(), but V[j] is
414:    orthogonalized only against columns V[i] having which[i]=PETSC_TRUE.
415:    The length of array which must be j at least.

417:    The use of this operation is restricted to MGS orthogonalization type.

419:    Level: advanced

421: .seealso: BVOrthogonalizeColumn(), BVSetOrthogonalization()
422: @*/
423: PetscErrorCode BVOrthogonalizeSomeColumn(BV bv,PetscInt j,PetscBool *which,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
424: {
426:   PetscInt       i,ksave,lsave;

433:   BVCheckSizes(bv,1);
434:   if (j<0) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
435:   if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Index j=%D but BV only has %D columns",j,bv->m);
436:   if (bv->orthog_type!=BV_ORTHOG_MGS) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation only available for MGS orthogonalization");

438:   PetscLogEventBegin(BV_Orthogonalize,bv,0,0,0);
439:   ksave = bv->k;
440:   lsave = bv->l;
441:   bv->l = -bv->nc;  /* must also orthogonalize against constraints and leading columns */
442:   BV_AllocateCoeffs(bv);
443:   BV_AllocateSignature(bv);
444:   BVOrthogonalizeMGS(bv,j,NULL,which,H,norm,lindep);
445:   bv->k = ksave;
446:   bv->l = lsave;
447:   if (H) for (i=bv->l;i<j;i++) H[i-bv->l] = bv->h[bv->nc+i];
448:   PetscLogEventEnd(BV_Orthogonalize,bv,0,0,0);
449:   return(0);
450: }

454: static PetscErrorCode BVOrthogonalize_GS(BV V,Mat R)
455: {
457:   PetscScalar    *r=NULL;
458:   PetscReal      norm;
459:   PetscInt       j,ldr;

462:   ldr = V->k;
463:   if (R) {
464:     MatDenseGetArray(R,&r);
465:     PetscMemzero(r+V->l*ldr,ldr*(ldr-V->l)*sizeof(PetscScalar));
466:   }
467:   for (j=V->l;j<V->k;j++) {
468:     if (R) {
469:       BVOrthogonalizeColumn(V,j,r+j*ldr,&norm,NULL);
470:       r[j+j*ldr] = norm;
471:     } else {
472:       BVOrthogonalizeColumn(V,j,NULL,&norm,NULL);
473:     }
474:     BVScaleColumn(V,j,1.0/norm);
475:   }
476:   if (R) { MatDenseRestoreArray(R,&r); }
477:   return(0);
478: }

482: /*@
483:    BVOrthogonalize - Orthogonalize all columns (except leading ones), that is,
484:    compute the QR decomposition.

486:    Collective on BV

488:    Input Parameter:
489: .  V - basis vectors

491:    Output Parameters:
492: +  V - the modified basis vectors
493: -  R - a sequential dense matrix (or NULL)

495:    Notes:
496:    On input, matrix R must be a sequential dense Mat, with at least as many rows
497:    and columns as the number of active columns of V. The output satisfies
498:    V0 = V*R (where V0 represent the input V) and V'*V = I.

500:    If V has leading columns, then they are not modified (are assumed to be already
501:    orthonormal) and the corresponding part of R is not referenced.

503:    Can pass NULL if R is not required.

505:    Level: intermediate

507: .seealso: BVOrthogonalizeColumn(), BVOrthogonalizeVec(), BVSetActiveColumns()
508: @*/
509: PetscErrorCode BVOrthogonalize(BV V,Mat R)
510: {
512:   PetscBool      match;
513:   PetscInt       m,n;

518:   BVCheckSizes(V,1);
519:   if (R) {
522:     PetscObjectTypeCompare((PetscObject)R,MATSEQDENSE,&match);
523:     if (!match) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Mat argument must be of type seqdense");
524:     MatGetSize(R,&m,&n);
525:     if (m!=n) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Mat argument is not square, it has %D rows and %D columns",m,n);
526:     if (n<V->k) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Mat size %D is smaller than the number of BV active columns %D",n,V->k);
527:   }
528:   if (V->matrix) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Not implemented for non-standard inner product, use BVOrthogonalizeColumn() instead");
529:   if (V->nc) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Not implemented for BV with constraints, use BVOrthogonalizeColumn() instead");

531:   PetscLogEventBegin(BV_Orthogonalize,V,R,0,0);
532:   if (*V->ops->orthogonalize) {
533:     (*V->ops->orthogonalize)(V,R);
534:   } else { /* no specific QR function available, so proceed column by column with Gram-Schmidt */
535:     BVOrthogonalize_GS(V,R);
536:   }
537:   PetscLogEventEnd(BV_Orthogonalize,V,R,0,0);
538:   PetscObjectStateIncrease((PetscObject)V);
539:   return(0);
540: }