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-2011, 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 private/ipimpl.h
26: #include slepcblaslapack.h
27: #include ../src/eps/impls/davidson/common/davidson.h
29: /*
30: IPOrthogonalizeCGS1 - Compute |v'| (estimated), |v| and one step of CGS with only one global synchronization
31: */
34: PetscErrorCode IPBOrthogonalizeCGS1(IP ip,PetscInt nds,Vec *DS,Vec *BDS,PetscInt n,PetscBool *which,Vec *V,Vec *BV,Vec v,Vec Bv,PetscScalar *H,PetscReal *onorm,PetscReal *norm)
35: {
37: PetscInt i,j;
38: PetscScalar alpha;
41: /* h = [DS V]^* Bv ; alpha = (v , v) */
42: VecsMultIa(H,0,nds,DS,0,nds,&Bv,0,1); j = nds;
43: if (!which) {
44: VecsMultIa(H+j,0,n,V,0,n,&Bv,0,1); j+= n;
45: } else {
46: for (i=0; i<n; i++) {
47: if (which[i]) {
48: VecsMultIa(H+j,0,1,V+i,0,1,&Bv,0,1); j++;
49: }
50: }
51: }
52: if (onorm || norm) {
53: VecsMultIa(H+j,0,1,&v,0,1,&Bv,0,1); j++;
54: }
55: VecsMultIb(H,0,j,j,1,PETSC_NULL,v);
56: if (onorm || norm) alpha = H[j-1];
58: /* v = v - V h */
59: SlepcVecMAXPBY(v,1.0,-1.0,nds,H,DS);
60: if (which) {
61: for (j=0; j<n; j++)
62: if (which[j]) { VecAXPBY(v,-H[nds+j],1.0,V[j]); }
63: } else {
64: SlepcVecMAXPBY(v,1.0,-1.0,n,H+nds,V);
65: }
67: /* Bv = Bv - BV h */
68: SlepcVecMAXPBY(Bv,1.0,-1.0,nds,H,BDS);
69: if (which) {
70: for (j=0; j<n; j++)
71: if (which[j]) { VecAXPBY(Bv,-H[nds+j],1.0,BV[j]); }
72: } else {
73: SlepcVecMAXPBY(Bv,1.0,-1.0,n,H+nds,BV);
74: }
75:
76: /* compute |v| */
77: if (onorm) *onorm = PetscSqrtReal(PetscRealPart(alpha));
79: /* compute |v'| */
80: if (norm) {
81: VecDot(Bv, v, &alpha);
82: *norm = PetscSqrtReal(PetscRealPart(alpha));
83: }
84: return(0);
85: }
87: /*
88: IPOrthogonalizeCGS - Orthogonalize with classical Gram-Schmidt
89: */
92: static PetscErrorCode IPBOrthogonalizeCGS(IP ip,PetscInt nds,Vec *DS,Vec *BDS,PetscInt n,PetscBool *which,Vec *V,Vec *BV,Vec v,Vec Bv,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
93: {
95: PetscScalar lh[100],*h,lc[100],*c,alpha;
96: PetscBool allocatedh = PETSC_FALSE,allocatedc = PETSC_FALSE;
97: PetscReal onrm,nrm;
98: PetscInt j,k;
101: /* allocate h and c if needed */
102: if (!H || nds>0) {
103: if (nds+n+1<=100) h = lh;
104: else {
105: PetscMalloc((nds+n+1)*sizeof(PetscScalar),&h);
106: allocatedh = PETSC_TRUE;
107: }
108: } else h = H;
109: if (ip->orthog_ref != IP_ORTHOG_REFINE_NEVER) {
110: if (nds+n+1<=100) c = lc;
111: else {
112: PetscMalloc((nds+n+1)*sizeof(PetscScalar),&c);
113: allocatedc = PETSC_TRUE;
114: }
115: }
117: /* orthogonalize and compute onorm */
118: switch (ip->orthog_ref) {
119:
120: case IP_ORTHOG_REFINE_NEVER:
121: IPBOrthogonalizeCGS1(ip,nds,DS,BDS,n,which,V,BV,v,Bv,h,PETSC_NULL,PETSC_NULL);
122: /* compute |v| */
123: if (norm) {
124: VecDot(Bv,v,&alpha);
125: *norm = PetscSqrtReal(PetscRealPart(alpha));
126: }
127: /* linear dependence check does not work without refinement */
128: if (lindep) *lindep = PETSC_FALSE;
129: break;
130:
131: case IP_ORTHOG_REFINE_ALWAYS:
132: IPBOrthogonalizeCGS1(ip,nds,DS,BDS,n,which,V,BV,v,Bv,h,PETSC_NULL,PETSC_NULL);
133: if (lindep) {
134: IPBOrthogonalizeCGS1(ip,nds,DS,BDS,n,which,V,BV,v,Bv,c,&onrm,&nrm);
135: if (norm) *norm = nrm;
136: if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
137: else *lindep = PETSC_FALSE;
138: } else {
139: IPBOrthogonalizeCGS1(ip,nds,DS,BDS,n,which,V,BV,v,Bv,c,PETSC_NULL,norm);
140: }
141: for (j=0;j<n;j++)
142: if (!which || which[j]) h[nds+j] += c[nds+j];
143: break;
144:
145: case IP_ORTHOG_REFINE_IFNEEDED:
146: IPBOrthogonalizeCGS1(ip,nds,DS,BDS,n,which,V,BV,v,Bv,h,&onrm,&nrm);
147: /* ||q|| < eta ||h|| */
148: k = 1;
149: while (k<3 && nrm < ip->orthog_eta * onrm) {
150: k++;
151: if (!ip->matrix) {
152: IPBOrthogonalizeCGS1(ip,nds,DS,BDS,n,which,V,BV,v,Bv,c,&onrm,&nrm);
153: } else {
154: onrm = nrm;
155: IPBOrthogonalizeCGS1(ip,nds,DS,BDS,n,which,V,BV,v,Bv,c,PETSC_NULL,&nrm);
156: }
157: for (j=0;j<n;j++)
158: if (!which || which[j]) h[nds+j] += c[nds+j];
159: }
160: if (norm) *norm = nrm;
161: if (lindep) {
162: if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
163: else *lindep = PETSC_FALSE;
164: }
165: break;
167: default:
168: SETERRQ(((PetscObject)ip)->comm,PETSC_ERR_ARG_WRONG,"Unknown orthogonalization refinement");
169: }
171: /* recover H from workspace */
172: if (H && nds>0) {
173: for (j=0;j<n;j++)
174: if (!which || which[j]) H[j] = h[nds+j];
175: }
177: /* free work space */
178: if (allocatedc) { PetscFree(c); }
179: if (allocatedh) { PetscFree(h); }
180: return(0);
181: }
185: /*@
186: IPBOrthogonalize - B-Orthogonalize a vector with respect to two set of vectors.
188: Collective on IP
190: Input Parameters:
191: + ip - the inner product (IP) context
192: . nds - number of columns of DS
193: . DS - first set of vectors
194: . BDS - B * DS
195: . n - number of columns of V
196: . which - logical array indicating columns of V to be used
197: . V - second set of vectors
198: - BV - B * V
200: Input/Output Parameter:
201: + v - (input) vector to be orthogonalized and (output) result of
202: orthogonalization
203: - Bv - (input/output) B * v
205: Output Parameter:
206: + H - coefficients computed during orthogonalization with V, of size nds+n
207: if norm == PETSC_NULL, and nds+n+1 otherwise.
208: . norm - norm of the vector after being orthogonalized
209: - lindep - flag indicating that refinement did not improve the quality
210: of orthogonalization
212: Notes:
213: This function applies an orthogonal projector to project vector v onto the
214: orthogonal complement of the span of the columns of DS and V.
216: On exit, v0 = [V v]*H, where v0 is the original vector v.
218: This routine does not normalize the resulting vector.
220: Level: developer
222: .seealso: IPSetOrthogonalization(), IPBiOrthogonalize()
223: @*/
224: PetscErrorCode IPBOrthogonalize(IP ip,PetscInt nds,Vec *DS, Vec *BDS,PetscInt n,PetscBool *which,Vec *V,Vec *BV,Vec v,Vec Bv,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
225: {
227: PetscScalar alpha;
233: PetscLogEventBegin(IP_Orthogonalize,ip,0,0,0);
234:
235: /* Bv <- B * v */
236: PetscLogEventBegin(IP_ApplyMatrix,ip,0,0,0);
237: MatMult(ip->matrix, v, Bv);
238: PetscLogEventEnd(IP_ApplyMatrix,ip,0,0,0);
239:
240: if (nds==0 && n==0) {
241: if (norm) {
242: VecDot(Bv, v, &alpha);
243: *norm = PetscSqrtReal(PetscRealPart(alpha));
244: }
245: if (lindep) *lindep = PETSC_FALSE;
246: } else {
247: switch (ip->orthog_type) {
248: case IP_ORTHOG_CGS:
249: IPBOrthogonalizeCGS(ip,nds,DS,BDS,n,which,V,BV,v,Bv,H,norm,lindep);
250: break;
251: case IP_ORTHOG_MGS:
252: SETERRQ(((PetscObject)ip)->comm,PETSC_ERR_ARG_WRONG,"Unsupported orthogonalization type");
253: break;
254: default:
255: SETERRQ(((PetscObject)ip)->comm,PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
256: }
257: }
258: PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);
259: return(0);
260: }