Actual source code: ipbiorthog.c

  1: /*
  2:      Routines related to bi-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-2012, Universitat 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 <slepc-private/ipimpl.h>      /*I "slepcip.h" I*/
 26: #include <slepcblaslapack.h>

 28: /*
 29:     Biorthogonalization routine using classical Gram-Schmidt with refinement.
 30:  */
 33: static PetscErrorCode IPCGSBiOrthogonalization(IP ip,PetscInt n_,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *hnorm,PetscReal *norm)
 34: {
 35: #if defined(SLEPC_MISSING_LAPACK_GELQF) || defined(SLEPC_MISSING_LAPACK_ORMLQ)
 37:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GELQF/ORMLQ - Lapack routine is unavailable");
 38: #else
 40:   PetscBLASInt   j,ione=1,lwork,info,n=n_;
 41:   PetscScalar    shh[100],*lhh,*vw,*tau,one=1.0,*work;

 44:   /* Don't allocate small arrays */
 45:   if (n<=100) lhh = shh;
 46:   else { PetscMalloc(n*sizeof(PetscScalar),&lhh); }
 47:   PetscMalloc(n*n*sizeof(PetscScalar),&vw);
 48: 
 49:   for (j=0;j<n;j++) {
 50:     IPMInnerProduct(ip,V[j],n,W,vw+j*n);
 51:   }
 52:   lwork = n;
 53:   PetscMalloc(n*sizeof(PetscScalar),&tau);
 54:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
 55:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
 56:   LAPACKgelqf_(&n,&n,vw,&n,tau,work,&lwork,&info);
 57:   PetscFPTrapPop();
 58:   if (info) SETERRQ1(((PetscObject)ip)->comm,PETSC_ERR_LIB,"Error in Lapack xGELQF %d",info);
 59: 
 60:   /*** First orthogonalization ***/

 62:   /* h = W^* v */
 63:   /* q = v - V h */
 64:   IPMInnerProduct(ip,v,n,W,H);
 65:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
 66:   BLAStrsm_("L","L","N","N",&n,&ione,&one,vw,&n,H,&n);
 67:   LAPACKormlq_("L","N",&n,&ione,&n,vw,&n,tau,H,&n,work,&lwork,&info);
 68:   PetscFPTrapPop();
 69:   if (info) SETERRQ1(((PetscObject)ip)->comm,PETSC_ERR_LIB,"Error in Lapack xORMLQ %d",info);
 70:   SlepcVecMAXPBY(v,1.0,-1.0,n,H,V);

 72:   /* compute norm of v */
 73:   if (norm) { IPNorm(ip,v,norm); }
 74: 
 75:   if (n>100) { PetscFree(lhh); }
 76:   PetscFree(vw);
 77:   PetscFree(tau);
 78:   PetscFree(work);
 79:   return(0);
 80: #endif
 81: }

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

 88:    Collective on IP and Vec

 90:    Input Parameters:
 91: +  ip - the inner product context
 92: .  n - number of columns of V
 93: .  V - set of vectors
 94: -  W - set of vectors

 96:    Input/Output Parameter:
 97: .  v - vector to be orthogonalized

 99:    Output Parameter:
100: +  H  - coefficients computed during orthogonalization
101: -  norm - norm of the vector after being orthogonalized

103:    Notes:
104:    This function applies an oblique projector to project vector v onto the
105:    span of the columns of V along the orthogonal complement of the column
106:    space of W. 

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

110:    This routine does not normalize the resulting vector.

112:    Level: developer

114: .seealso: IPSetOrthogonalization(), IPOrthogonalize()
115: @*/
116: PetscErrorCode IPBiOrthogonalize(IP ip,PetscInt n,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *norm)
117: {
119:   PetscScalar    lh[100],*h;
120:   PetscBool      allocated = PETSC_FALSE;
121:   PetscReal      lhnrm,*hnrm,lnrm,*nrm;

126:   if (n==0) {
127:     if (norm) { IPNorm(ip,v,norm); }
128:   } else {
129:     PetscLogEventBegin(IP_Orthogonalize,ip,0,0,0);
130:     /* allocate H if needed */
131:     if (!H) {
132:       if (n<=100) h = lh;
133:       else {
134:         PetscMalloc(n*sizeof(PetscScalar),&h);
135:         allocated = PETSC_TRUE;
136:       }
137:     } else h = H;
138: 
139:     /* retrieve hnrm and nrm for linear dependence check or conditional refinement */
140:     if (ip->orthog_ref == IP_ORTHOG_REFINE_IFNEEDED) {
141:       hnrm = &lhnrm;
142:       if (norm) nrm = norm;
143:       else nrm = &lnrm;
144:     } else {
145:       hnrm = PETSC_NULL;
146:       nrm = norm;
147:     }
148: 
149:     switch (ip->orthog_type) {
150:       case IP_ORTHOG_CGS:
151:         IPCGSBiOrthogonalization(ip,n,V,W,v,h,hnrm,nrm);
152:         break;
153:       default:
154:         SETERRQ(((PetscObject)ip)->comm,PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
155:     }
156: 
157:     if (allocated) { PetscFree(h); }
158:     PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);
159:   }
160:   return(0);
161: }

163: /* 
164:    IPPseudoOrthogonalizeCGS1 - Compute |v'| (estimated), |v| and one step of CGS with only one global synchronization (indefinite)
165: */
168: PetscErrorCode IPPseudoOrthogonalizeCGS1(IP ip,PetscInt n,Vec *V,PetscReal* omega,Vec v,PetscScalar *H,PetscReal *onorm,PetscReal *norm)
169: {
171:   PetscInt       j;
172:   PetscScalar    alpha;
173:   PetscReal      sum;

176:   /* h = W^* v ; alpha = (v , v) */
177:   if (!onorm && !norm) {
178:     /* use simpler function */
179:     IPMInnerProduct(ip,v,n,V,H);
180:   } else {
181:     /* merge comunications */
182:     IPMInnerProductBegin(ip,v,n,V,H);
183:     if (onorm || (norm && !ip->matrix)) {
184:       IPInnerProductBegin(ip,v,v,&alpha);
185:     }

187:     IPMInnerProductEnd(ip,v,n,V,H);
188:     if (onorm || (norm && !ip->matrix)) {
189:       IPInnerProductEnd(ip,v,v,&alpha);
190:     }
191:   }

193:   /* q = v - V h */
194:   for (j=0;j<n;j++) H[j] /= omega[j];  /* apply inverse of signature */
195:   SlepcVecMAXPBY(v,1.0,-1.0,n,H,V);
196:   for (j=0;j<n;j++) H[j] *= omega[j];  /* revert signature */
197: 
198:   /* compute |v| */
199:   if (onorm) {
200:     if (PetscRealPart(alpha)>0.0) *onorm = PetscSqrtReal(PetscRealPart(alpha));
201:     else *onorm = -PetscSqrtReal(-PetscRealPart(alpha));
202:   }

204:   if (norm) {
205:     if (!ip->matrix) {
206:       /* estimate |v'| from |v| */
207:       sum = 0.0;
208:       for (j=0; j<n; j++)
209:         sum += PetscRealPart(H[j] * PetscConj(H[j]));
210:       *norm = PetscRealPart(alpha)-sum;
211:       if (*norm <= 0.0) {
212:         IPNorm(ip,v,norm);
213:       } else *norm = PetscSqrtReal(*norm);
214:     } else {
215:       /* compute |v'| */
216:       IPNorm(ip,v,norm);
217:     }
218:   }
219:   return(0);
220: }

222: /*
223:   IPPseudoOrthogonalizeCGS - Orthogonalize with classical Gram-Schmidt (indefinite)
224: */
227: static PetscErrorCode IPPseudoOrthogonalizeCGS(IP ip,PetscInt n,Vec *V,PetscReal *omega,Vec v,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
228: {
230:   PetscScalar    lh[100],*h,lc[100],*c;
231:   PetscBool      allocatedh = PETSC_FALSE,allocatedc = PETSC_FALSE;
232:   PetscReal      onrm,nrm;
233:   PetscInt       j,k;

236:   /* allocate h and c if needed */
237:   if (!H) {
238:     if (n<=100) h = lh;
239:     else {
240:       PetscMalloc(n*sizeof(PetscScalar),&h);
241:       allocatedh = PETSC_TRUE;
242:     }
243:   } else h = H;
244:   if (ip->orthog_ref != IP_ORTHOG_REFINE_NEVER) {
245:     if (n<=100) c = lc;
246:     else {
247:       PetscMalloc(n*sizeof(PetscScalar),&c);
248:       allocatedc = PETSC_TRUE;
249:     }
250:   }

252:   /* orthogonalize and compute onorm */
253:   switch (ip->orthog_ref) {
254: 
255:   case IP_ORTHOG_REFINE_NEVER:
256:     IPPseudoOrthogonalizeCGS1(ip,n,V,omega,v,h,PETSC_NULL,PETSC_NULL);
257:     /* compute |v| */
258:     if (norm) { IPNorm(ip,v,norm); }
259:     /* linear dependence check does not work without refinement */
260:     if (lindep) *lindep = PETSC_FALSE;
261:     break;
262: 
263:   case IP_ORTHOG_REFINE_ALWAYS:
264:     IPPseudoOrthogonalizeCGS1(ip,n,V,omega,v,h,PETSC_NULL,PETSC_NULL);
265:     if (lindep) {
266:       IPPseudoOrthogonalizeCGS1(ip,n,V,omega,v,c,&onrm,&nrm);
267:       if (norm) *norm = nrm;
268:       if (PetscAbs(nrm) < ip->orthog_eta * PetscAbs(onrm)) *lindep = PETSC_TRUE;
269:       else *lindep = PETSC_FALSE;
270:     } else {
271:       IPPseudoOrthogonalizeCGS1(ip,n,V,omega,v,c,PETSC_NULL,norm);
272:     }
273:     for (j=0;j<n;j++)
274:       h[j] += c[j];
275:     break;
276: 
277:   case IP_ORTHOG_REFINE_IFNEEDED:
278:     IPPseudoOrthogonalizeCGS1(ip,n,V,omega,v,h,&onrm,&nrm);
279:     /* ||q|| < eta ||h|| */
280:     k = 1;
281:     while (k<3 && PetscAbs(nrm) < ip->orthog_eta * PetscAbs(onrm)) {
282:       k++;
283:       if (!ip->matrix) {
284:         IPPseudoOrthogonalizeCGS1(ip,n,V,omega,v,c,&onrm,&nrm);
285:       } else {
286:         onrm = nrm;
287:         IPPseudoOrthogonalizeCGS1(ip,n,V,omega,v,c,PETSC_NULL,&nrm);
288:       }
289:       for (j=0;j<n;j++)
290:         h[j] += c[j];
291:     }
292:     if (norm) *norm = nrm;
293:     if (lindep) {
294:       if (PetscAbs(nrm) < ip->orthog_eta * PetscAbs(onrm)) *lindep = PETSC_TRUE;
295:       else *lindep = PETSC_FALSE;
296:     }
297:     break;

299:   default:
300:     SETERRQ(((PetscObject)ip)->comm,PETSC_ERR_ARG_WRONG,"Unknown orthogonalization refinement");
301:   }

303:   /* recover H from workspace */
304:   if (H) {
305:     for (j=0;j<n;j++)
306:       H[j] = h[j];
307:   }

309:   /* free work space */
310:   if (allocatedc) { PetscFree(c); }
311:   if (allocatedh) { PetscFree(h); }
312:   return(0);
313: }

317: /*@
318:    IPPseudoOrthogonalize - Orthogonalize a vector with respect to two set of vectors
319:    in the sense of a pseudo-inner product.

321:    Collective on IP and Vec

323:    Input Parameters:
324: +  ip     - the inner product (IP) context
325: .  n      - number of columns of V
326: .  V      - set of vectors
327: -  omega  - set of signs that define a signature matrix

329:    Input/Output Parameter:
330: .  v      - (input) vector to be orthogonalized and (output) result of 
331:             orthogonalization

333:    Output Parameter:
334: +  H      - coefficients computed during orthogonalization
335: .  norm   - norm of the vector after being orthogonalized
336: -  lindep - flag indicating that refinement did not improve the quality
337:             of orthogonalization

339:    Notes:
340:    This function is the analogue of IPOrthogonalize, but for the indefinite
341:    case. When using an indefinite IP the norm is not well defined, so we
342:    take the convention of having negative norms in such cases. The 
343:    orthogonalization is then defined by a set of vectors V satisfying 
344:    V'*B*V=Omega, where Omega is a signature matrix diag([+/-1,...,+/-1]).

346:    On exit, v = v0 - V*Omega*H, where v0 is the original vector v.

348:    This routine does not normalize the resulting vector. The output
349:    argument 'norm' may be negative.

351:    Level: developer

353: .seealso: IPSetOrthogonalization(), IPOrthogonalize()
354: @*/
355: PetscErrorCode IPPseudoOrthogonalize(IP ip,PetscInt n,Vec *V,PetscReal *omega,Vec v,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
356: {

362:   PetscLogEventBegin(IP_Orthogonalize,ip,0,0,0);
363:   if (n==0) {
364:     if (norm) { IPNorm(ip,v,norm); }
365:     if (lindep) *lindep = PETSC_FALSE;
366:   } else {
367:     switch (ip->orthog_type) {
368:     case IP_ORTHOG_CGS:
369:       IPPseudoOrthogonalizeCGS(ip,n,V,omega,v,H,norm,lindep);
370:       break;
371:     case IP_ORTHOG_MGS:
372:       SETERRQ(((PetscObject)ip)->comm,PETSC_ERR_SUP,"Modified Gram-Schmidt not implemented for indefinite case");
373:       break;
374:     default:
375:       SETERRQ(((PetscObject)ip)->comm,PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
376:     }
377:   }
378:   PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);
379:   return(0);
380: }