Actual source code: ipborthog.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-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
 27:  #include ../src/eps/impls/davidson/common/davidson.h

 29: #define MyPetscSqrtReal(alpha) (PetscSign(PetscRealPart(alpha))*PetscSqrtReal(PetscAbs(PetscRealPart(alpha))))

 31: /* 
 32:    IPOrthogonalizeCGS1 - Compute |v'| (estimated), |v| and one step of CGS with only one global synchronization
 33: */
 36: PetscErrorCode IPBOrthogonalizeCGS1(IP ip,PetscInt nds,Vec *defl,Vec *BDS,PetscReal *BDSnorms,PetscInt n,PetscBool *which,Vec *V,Vec *BV,PetscReal *BVnorms,Vec v,Vec Bv,PetscScalar *H,PetscReal *onorm,PetscReal *norm)
 37: {
 39:   PetscInt       i,j;
 40:   PetscScalar    alpha;

 43:   /* h = [defl V]^* Bv ; alpha = (Bv , v) */
 44:   VecsMultIa(H,0,nds,defl,0,nds,&Bv,0,1); j = nds;
 45:   if (!which) {
 46:     VecsMultIa(H+j,0,n,V,0,n,&Bv,0,1); j+= n;
 47:   } else {
 48:     for (i=0; i<n; i++) {
 49:       if (which[i]) {
 50:         VecsMultIa(H+j,0,1,V+i,0,1,&Bv,0,1); j++;
 51:       }
 52:     }
 53:   }
 54:   if (onorm || norm) {
 55:     VecsMultIa(H+j,0,1,&v,0,1,&Bv,0,1); j++;
 56:   }
 57:   VecsMultIb(H,0,j,j,1,PETSC_NULL,v);
 58:   if (onorm || norm) alpha = H[j-1];

 60:   /* h = J * h */
 61:   if (BDSnorms && defl) for (i=0; i<nds; i++) H[i]*= BDSnorms[i];
 62:   if (BVnorms && V) {
 63:     if (!which) {
 64:       for (i=0; i<n; i++) H[i+nds]*= BVnorms[i];
 65:     } else {
 66:       for (i=j=0; i<n; i++) {
 67:         if (which[i]) H[j++]*= BVnorms[i];
 68:       }
 69:     }
 70:   }

 72:   /* v = v - V h */
 73:   SlepcVecMAXPBY(v,1.0,-1.0,nds,H,defl);
 74:   if (which) {
 75:     for (j=0; j<n; j++)
 76:       if (which[j]) { VecAXPBY(v,-H[nds+j],1.0,V[j]); }
 77:   } else {
 78:     SlepcVecMAXPBY(v,1.0,-1.0,n,H+nds,V);
 79:   }

 81:   /* Bv = Bv - BV h */
 82:   SlepcVecMAXPBY(Bv,1.0,-1.0,nds,H,BDS);
 83:   if (which) {
 84:     for (j=0; j<n; j++)
 85:       if (which[j]) { VecAXPBY(Bv,-H[nds+j],1.0,BV[j]); }
 86:   } else {
 87:     SlepcVecMAXPBY(Bv,1.0,-1.0,n,H+nds,BV);
 88:   }

 90:   /* compute |v| */
 91:   if (onorm) *onorm = MyPetscSqrtReal(alpha);

 93:   /* compute |v'| */
 94:   if (norm) {
 95:     VecDot(Bv, v, &alpha);
 96:     *norm = MyPetscSqrtReal(alpha);
 97:   }
 98:   return(0);
 99: }

101: /*
102:   IPOrthogonalizeCGS - Orthogonalize with classical Gram-Schmidt
103: */
106: static PetscErrorCode IPBOrthogonalizeCGS(IP ip,PetscInt nds,Vec *defl,Vec *BDS,PetscReal *BDSnorms,PetscInt n,PetscBool *which,Vec *V,Vec *BV,PetscReal *BVnorms,Vec v,Vec Bv,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
107: {
109:   PetscScalar    lh[100],*h,lc[100],*c,alpha;
110:   PetscBool      allocatedh = PETSC_FALSE,allocatedc = PETSC_FALSE;
111:   PetscReal      onrm,nrm;
112:   PetscInt       j,k;

115:   /* allocate h and c if needed */
116:   if (!H || nds>0) {
117:     if (nds+n+1<=100) h = lh;
118:     else {
119:       PetscMalloc((nds+n+1)*sizeof(PetscScalar),&h);
120:       allocatedh = PETSC_TRUE;
121:     }
122:   } else h = H;
123:   if (ip->orthog_ref != IP_ORTHOG_REFINE_NEVER) {
124:     if (nds+n+1<=100) c = lc;
125:     else {
126:       PetscMalloc((nds+n+1)*sizeof(PetscScalar),&c);
127:       allocatedc = PETSC_TRUE;
128:     }
129:   }

131:   /* orthogonalize and compute onorm */
132:   switch (ip->orthog_ref) {
133: 
134:   case IP_ORTHOG_REFINE_NEVER:
135:     IPBOrthogonalizeCGS1(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,h,PETSC_NULL,PETSC_NULL);
136:     /* compute |v| */
137:     if (norm) {
138:       VecDot(Bv,v,&alpha);
139:       *norm = MyPetscSqrtReal(alpha);
140:     }
141:     /* linear dependence check does not work without refinement */
142:     if (lindep) *lindep = PETSC_FALSE;
143:     break;
144: 
145:   case IP_ORTHOG_REFINE_ALWAYS:
146:     IPBOrthogonalizeCGS1(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,h,PETSC_NULL,PETSC_NULL);
147:     if (lindep) {
148:       IPBOrthogonalizeCGS1(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,c,&onrm,&nrm);
149:       if (norm) *norm = nrm;
150:       if (PetscAbs(nrm) < ip->orthog_eta * PetscAbs(onrm)) *lindep = PETSC_TRUE;
151:       else *lindep = PETSC_FALSE;
152:     } else {
153:       IPBOrthogonalizeCGS1(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,c,PETSC_NULL,norm);
154:     }
155:     for (j=0;j<n;j++)
156:       if (!which || which[j]) h[nds+j] += c[nds+j];
157:     break;
158: 
159:   case IP_ORTHOG_REFINE_IFNEEDED:
160:     IPBOrthogonalizeCGS1(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,h,&onrm,&nrm);
161:     /* ||q|| < eta ||h|| */
162:     k = 1;
163:     while (k<3 && PetscAbs(nrm) < ip->orthog_eta * PetscAbs(onrm)) {
164:       k++;
165:       if (!ip->matrix) {
166:         IPBOrthogonalizeCGS1(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,c,&onrm,&nrm);
167:       } else {
168:         onrm = nrm;
169:         IPBOrthogonalizeCGS1(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,c,PETSC_NULL,&nrm);
170:       }
171:       for (j=0;j<n;j++)
172:         if (!which || which[j]) h[nds+j] += c[nds+j];
173:     }
174:     if (norm) *norm = nrm;
175:     if (lindep) {
176:       if (PetscAbs(nrm) < ip->orthog_eta * PetscAbs(onrm)) *lindep = PETSC_TRUE;
177:       else *lindep = PETSC_FALSE;
178:     }
179:     break;

181:   default:
182:     SETERRQ(((PetscObject)ip)->comm,PETSC_ERR_ARG_WRONG,"Unknown orthogonalization refinement");
183:   }

185:   /* recover H from workspace */
186:   if (H && nds>0) {
187:     for (j=0;j<n;j++)
188:       if (!which || which[j]) H[j] = h[nds+j];
189:   }

191:   /* free work space */
192:   if (allocatedc) { PetscFree(c); }
193:   if (allocatedh) { PetscFree(h); }
194:   return(0);
195: }

199: /*@
200:    IPBOrthogonalize - B-Orthogonalize a vector with respect to two set of vectors.

202:    Collective on IP

204:    Input Parameters:
205: +  ip     - the inner product (IP) context
206: .  nds    - number of columns of defl
207: .  defl   - first set of vectors
208: .  BDS    - B * defl
209: .  BDSnorms - DS_i' * B * DS_i
210: .  n      - number of columns of V
211: .  which  - logical array indicating columns of V to be used
212: .  V      - second set of vectors
213: .  BV     - B * V
214: -  BVnorms - V_i' * B * V_i

216:    Input/Output Parameter:
217: +  v      - (input) vector to be orthogonalized and (output) result of 
218:             orthogonalization
219: -  Bv     - (input/output) B * v

221:    Output Parameter:
222: +  H      - coefficients computed during orthogonalization with V, of size nds+n
223:             if norm == PETSC_NULL, and nds+n+1 otherwise.
224: .  norm   - norm of the vector after being orthogonalized
225: -  lindep - flag indicating that refinement did not improve the quality
226:             of orthogonalization

228:    Notes:
229:    This function applies an orthogonal projector to project vector v onto the
230:    orthogonal complement of the span of the columns of defl and V.

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

234:    This routine does not normalize the resulting vector.

236:    Level: developer

238: .seealso: IPSetOrthogonalization(), IPBiOrthogonalize()
239: @*/
240: PetscErrorCode IPBOrthogonalize(IP ip,PetscInt nds,Vec *defl, Vec *BDS,PetscReal *BDSnorms,PetscInt n,PetscBool *which,Vec *V,Vec *BV,PetscReal *BVnorms,Vec v,Vec Bv,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
241: {
243:   PetscScalar    alpha;

249:   PetscLogEventBegin(IP_Orthogonalize,ip,0,0,0);
250: 
251:   /* Bv <- B * v */
252:   PetscLogEventBegin(IP_ApplyMatrix,ip,0,0,0);
253:   MatMult(ip->matrix, v, Bv);
254:   PetscLogEventEnd(IP_ApplyMatrix,ip,0,0,0);
255: 
256:   if (nds==0 && n==0) {
257:     if (norm) {
258:       VecDot(Bv, v, &alpha);
259:       *norm = MyPetscSqrtReal(alpha);
260:     }
261:     if (lindep) *lindep = PETSC_FALSE;
262:   } else {
263:     switch (ip->orthog_type) {
264:     case IP_ORTHOG_CGS:
265:       IPBOrthogonalizeCGS(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,H,norm,lindep);
266:       break;
267:     case IP_ORTHOG_MGS:
268:       SETERRQ(((PetscObject)ip)->comm,PETSC_ERR_ARG_WRONG,"Unsupported orthogonalization type");
269:       break;
270:     default:
271:       SETERRQ(((PetscObject)ip)->comm,PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
272:     }
273:   }
274:   PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);
275:   return(0);
276: }