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-2013, Universitat Politecnica de Valencia, Spain
9: This file is part of SLEPc.
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);
45: j = nds;
46: if (!which) {
47: VecsMultIa(H+j,0,n,V,0,n,&Bv,0,1);
48: j+= n;
49: } else {
50: for (i=0; i<n; i++) {
51: if (which[i]) {
52: VecsMultIa(H+j,0,1,V+i,0,1,&Bv,0,1);
53: j++;
54: }
55: }
56: }
57: if (onorm || norm) {
58: VecsMultIa(H+j,0,1,&v,0,1,&Bv,0,1);
59: j++;
60: }
61: VecsMultIb(H,0,j,j,1,NULL,v);
62: if (onorm || norm) alpha = H[j-1];
64: /* h = J * h */
65: if (BDSnorms && defl) for (i=0; i<nds; i++) H[i]*= BDSnorms[i];
66: if (BVnorms && V) {
67: if (!which) {
68: for (i=0; i<n; i++) H[i+nds]*= BVnorms[i];
69: } else {
70: for (i=j=0; i<n; i++) {
71: if (which[i]) H[j++]*= BVnorms[i];
72: }
73: }
74: }
76: /* v = v - V h */
77: SlepcVecMAXPBY(v,1.0,-1.0,nds,H,defl);
78: if (which) {
79: for (j=0; j<n; j++)
80: if (which[j]) { VecAXPBY(v,-H[nds+j],1.0,V[j]); }
81: } else {
82: SlepcVecMAXPBY(v,1.0,-1.0,n,H+nds,V);
83: }
85: /* Bv = Bv - BV h */
86: SlepcVecMAXPBY(Bv,1.0,-1.0,nds,H,BDS);
87: if (which) {
88: for (j=0; j<n; j++)
89: if (which[j]) { VecAXPBY(Bv,-H[nds+j],1.0,BV[j]); }
90: } else {
91: SlepcVecMAXPBY(Bv,1.0,-1.0,n,H+nds,BV);
92: }
94: /* compute |v| */
95: if (onorm) *onorm = MyPetscSqrtReal(alpha);
97: /* compute |v'| */
98: if (norm) {
99: VecDot(Bv,v,&alpha);
100: *norm = MyPetscSqrtReal(alpha);
101: }
102: return(0);
103: }
105: /*
106: IPOrthogonalizeCGS - Orthogonalize with classical Gram-Schmidt
107: */
110: 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)
111: {
113: PetscScalar *h,*c,alpha;
114: PetscReal onrm,nrm;
115: PetscInt sz=0,sz1,j,k;
118: /* allocate h and c if needed */
119: if (!H || nds>0) sz = nds+n;
120: sz1 = sz;
121: if (ip->orthog_ref != IP_ORTHOG_REFINE_NEVER) sz += nds+n;
122: if (sz>ip->lwork) {
123: PetscFree(ip->work);
124: PetscMalloc(sz*sizeof(PetscScalar),&ip->work);
125: PetscLogObjectMemory(ip,(sz-ip->lwork)*sizeof(PetscScalar));
126: ip->lwork = sz;
127: }
128: if (!H || nds>0) h = ip->work;
129: else h = H;
130: if (ip->orthog_ref != IP_ORTHOG_REFINE_NEVER) c = ip->work + sz1;
132: /* orthogonalize and compute onorm */
133: switch (ip->orthog_ref) {
135: case IP_ORTHOG_REFINE_NEVER:
136: IPBOrthogonalizeCGS1(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,h,NULL,NULL);
137: /* compute |v| */
138: if (norm) {
139: VecDot(Bv,v,&alpha);
140: *norm = MyPetscSqrtReal(alpha);
141: }
142: /* linear dependence check does not work without refinement */
143: if (lindep) *lindep = PETSC_FALSE;
144: break;
146: case IP_ORTHOG_REFINE_ALWAYS:
147: IPBOrthogonalizeCGS1(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,h,NULL,NULL);
148: if (lindep) {
149: IPBOrthogonalizeCGS1(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,c,&onrm,&nrm);
150: if (norm) *norm = nrm;
151: if (PetscAbs(nrm) < ip->orthog_eta * PetscAbs(onrm)) *lindep = PETSC_TRUE;
152: else *lindep = PETSC_FALSE;
153: } else {
154: IPBOrthogonalizeCGS1(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,c,NULL,norm);
155: }
156: for (j=0;j<n;j++)
157: if (!which || which[j]) h[nds+j] += c[nds+j];
158: break;
160: case IP_ORTHOG_REFINE_IFNEEDED:
161: IPBOrthogonalizeCGS1(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,h,&onrm,&nrm);
162: /* ||q|| < eta ||h|| */
163: k = 1;
164: while (k<3 && PetscAbs(nrm) < ip->orthog_eta * PetscAbs(onrm)) {
165: k++;
166: if (!ip->matrix) {
167: IPBOrthogonalizeCGS1(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,c,&onrm,&nrm);
168: } else {
169: onrm = nrm;
170: IPBOrthogonalizeCGS1(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,c,NULL,&nrm);
171: }
172: for (j=0;j<n;j++)
173: if (!which || which[j]) h[nds+j] += c[nds+j];
174: }
175: if (norm) *norm = nrm;
176: if (lindep) {
177: if (PetscAbs(nrm) < ip->orthog_eta * PetscAbs(onrm)) *lindep = PETSC_TRUE;
178: else *lindep = PETSC_FALSE;
179: }
180: break;
182: default:
183: SETERRQ(PetscObjectComm((PetscObject)ip),PETSC_ERR_ARG_WRONG,"Unknown orthogonalization refinement");
184: }
186: /* recover H from workspace */
187: if (H && nds>0) {
188: for (j=0;j<n;j++)
189: if (!which || which[j]) H[j] = h[nds+j];
190: }
191: return(0);
192: }
196: /*@
197: IPBOrthogonalize - B-Orthogonalize a vector with respect to two set of vectors.
199: Collective on IP
201: Input Parameters:
202: + ip - the inner product (IP) context
203: . nds - number of columns of defl
204: . defl - first set of vectors
205: . BDS - B * defl
206: . BDSnorms - DS_i' * B * DS_i
207: . n - number of columns of V
208: . which - logical array indicating columns of V to be used
209: . V - second set of vectors
210: . BV - B * V
211: - BVnorms - V_i' * B * V_i
213: Input/Output Parameter:
214: + v - (input) vector to be orthogonalized and (output) result of
215: orthogonalization
216: - Bv - (input/output) B * v
218: Output Parameter:
219: + H - coefficients computed during orthogonalization with V, of size nds+n
220: if norm == NULL, and nds+n+1 otherwise.
221: . norm - norm of the vector after being orthogonalized
222: - lindep - flag indicating that refinement did not improve the quality
223: of orthogonalization
225: Notes:
226: This function applies an orthogonal projector to project vector v onto the
227: orthogonal complement of the span of the columns of defl and V.
229: On exit, v0 = [V v]*H, where v0 is the original vector v.
231: This routine does not normalize the resulting vector.
233: Level: developer
235: .seealso: IPSetOrthogonalization(), IPBiOrthogonalize()
236: @*/
237: 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)
238: {
240: PetscScalar alpha;
246: PetscLogEventBegin(IP_Orthogonalize,ip,0,0,0);
248: /* Bv <- B * v */
249: PetscLogEventBegin(IP_ApplyMatrix,ip,0,0,0);
250: MatMult(ip->matrix,v,Bv);
251: PetscLogEventEnd(IP_ApplyMatrix,ip,0,0,0);
253: if (!nds && !n) {
254: if (norm) {
255: VecDot(Bv,v,&alpha);
256: *norm = MyPetscSqrtReal(alpha);
257: }
258: if (lindep) *lindep = PETSC_FALSE;
259: } else {
260: switch (ip->orthog_type) {
261: case IP_ORTHOG_CGS:
262: IPBOrthogonalizeCGS(ip,nds,defl,BDS,BDSnorms,n,which,V,BV,BVnorms,v,Bv,H,norm,lindep);
263: break;
264: case IP_ORTHOG_MGS:
265: SETERRQ(PetscObjectComm((PetscObject)ip),PETSC_ERR_ARG_WRONG,"Unsupported orthogonalization type");
266: break;
267: default:
268: SETERRQ(PetscObjectComm((PetscObject)ip),PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
269: }
270: }
271: PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);
272: return(0);
273: }