Actual source code: default.c
1: /*
2: This file contains some simple default routines for common operations.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2013, 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/epsimpl.h> /*I "slepceps.h" I*/
25: #include <slepcblaslapack.h>
29: PetscErrorCode EPSReset_Default(EPS eps)
30: {
34: VecDestroyVecs(eps->nwork,&eps->work);
35: eps->nwork = 0;
36: EPSFreeSolution(eps);
37: return(0);
38: }
42: PetscErrorCode EPSBackTransform_Default(EPS eps)
43: {
47: STBackTransform(eps->st,eps->nconv,eps->eigr,eps->eigi);
48: return(0);
49: }
53: /*
54: EPSComputeVectors_Default - Compute eigenvectors from the vectors
55: provided by the eigensolver. This version just copies the vectors
56: and is intended for solvers such as power that provide the eigenvector.
57: */
58: PetscErrorCode EPSComputeVectors_Default(EPS eps)
59: {
61: eps->evecsavailable = PETSC_TRUE;
62: return(0);
63: }
67: /*
68: EPSComputeVectors_Hermitian - Copies the Lanczos vectors as eigenvectors
69: using purification for generalized eigenproblems.
70: */
71: PetscErrorCode EPSComputeVectors_Hermitian(EPS eps)
72: {
74: PetscInt i;
75: PetscReal norm;
76: Vec w;
79: if (eps->isgeneralized) {
80: /* Purify eigenvectors */
81: VecDuplicate(eps->V[0],&w);
82: for (i=0;i<eps->nconv;i++) {
83: VecCopy(eps->V[i],w);
84: STApply(eps->st,w,eps->V[i]);
85: IPNorm(eps->ip,eps->V[i],&norm);
86: VecScale(eps->V[i],1.0/norm);
87: }
88: VecDestroy(&w);
89: }
90: eps->evecsavailable = PETSC_TRUE;
91: return(0);
92: }
96: /*
97: EPSComputeVectors_Indefinite - similar to the Schur version but
98: for indefinite problems
99: */
100: PetscErrorCode EPSComputeVectors_Indefinite(EPS eps)
101: {
103: PetscInt n,ld,i;
104: PetscScalar *Z;
105: Vec v;
106: #if !defined(PETSC_USE_COMPLEX)
107: PetscScalar tmp;
108: PetscReal norm,normi;
109: #endif
112: DSGetLeadingDimension(eps->ds,&ld);
113: DSGetDimensions(eps->ds,&n,NULL,NULL,NULL,NULL);
114: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
115: DSGetArray(eps->ds,DS_MAT_X,&Z);
116: SlepcUpdateVectors(n,eps->V,0,n,Z,ld,PETSC_FALSE);
117: DSRestoreArray(eps->ds,DS_MAT_X,&Z);
118: /* purification */
119: VecDuplicate(eps->V[0],&v);
120: for (i=0;i<eps->nconv;i++) {
121: VecCopy(eps->V[i],v);
122: STApply(eps->st,v,eps->V[i]);
123: }
124: VecDestroy(&v);
125: /* normalization */
126: for (i=0;i<n;i++) {
127: #if !defined(PETSC_USE_COMPLEX)
128: if (eps->eigi[i] != 0.0) {
129: VecNorm(eps->V[i],NORM_2,&norm);
130: VecNorm(eps->V[i+1],NORM_2,&normi);
131: tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
132: VecScale(eps->V[i],tmp);
133: VecScale(eps->V[i+1],tmp);
134: i++;
135: } else
136: #endif
137: {
138: VecNormalize(eps->V[i],NULL);
139: }
140: }
141: eps->evecsavailable = PETSC_TRUE;
142: return(0);
143: }
147: /*
148: EPSComputeVectors_Schur - Compute eigenvectors from the vectors
149: provided by the eigensolver. This version is intended for solvers
150: that provide Schur vectors. Given the partial Schur decomposition
151: OP*V=V*T, the following steps are performed:
152: 1) compute eigenvectors of T: T*Z=Z*D
153: 2) compute eigenvectors of OP: X=V*Z
154: If left eigenvectors are required then also do Z'*T=D*Z', Y=W*Z
155: */
156: PetscErrorCode EPSComputeVectors_Schur(EPS eps)
157: {
159: PetscInt n,i,ld;
160: PetscScalar *Z;
161: #if !defined(PETSC_USE_COMPLEX)
162: PetscScalar tmp;
163: PetscReal norm,normi;
164: #endif
165: Vec w;
168: if (eps->ishermitian) {
169: if (eps->isgeneralized && !eps->ispositive) {
170: EPSComputeVectors_Indefinite(eps);
171: } else {
172: EPSComputeVectors_Hermitian(eps);
173: }
174: return(0);
175: }
176: DSGetLeadingDimension(eps->ds,&ld);
177: DSGetDimensions(eps->ds,&n,NULL,NULL,NULL,NULL);
179: /* right eigenvectors */
180: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
182: /* V = V * Z */
183: DSGetArray(eps->ds,DS_MAT_X,&Z);
184: SlepcUpdateVectors(n,eps->V,0,n,Z,ld,PETSC_FALSE);
185: DSRestoreArray(eps->ds,DS_MAT_X,&Z);
187: /* Purify eigenvectors */
188: if (eps->ispositive) {
189: VecDuplicate(eps->V[0],&w);
190: for (i=0;i<n;i++) {
191: VecCopy(eps->V[i],w);
192: STApply(eps->st,w,eps->V[i]);
193: }
194: VecDestroy(&w);
195: }
197: /* Fix eigenvectors if balancing was used */
198: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
199: for (i=0;i<n;i++) {
200: VecPointwiseDivide(eps->V[i],eps->V[i],eps->D);
201: }
202: }
204: /* normalize eigenvectors (when using purification or balancing) */
205: if (eps->ispositive || (eps->balance!=EPS_BALANCE_NONE && eps->D)) {
206: for (i=0;i<n;i++) {
207: #if !defined(PETSC_USE_COMPLEX)
208: if (eps->eigi[i] != 0.0) {
209: VecNorm(eps->V[i],NORM_2,&norm);
210: VecNorm(eps->V[i+1],NORM_2,&normi);
211: tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
212: VecScale(eps->V[i],tmp);
213: VecScale(eps->V[i+1],tmp);
214: i++;
215: } else
216: #endif
217: {
218: VecNormalize(eps->V[i],NULL);
219: }
220: }
221: }
223: /* left eigenvectors */
224: if (eps->leftvecs) {
225: DSVectors(eps->ds,DS_MAT_Y,NULL,NULL);
226: /* W = W * Z */
227: DSGetArray(eps->ds,DS_MAT_Y,&Z);
228: SlepcUpdateVectors(n,eps->W,0,n,Z,ld,PETSC_FALSE);
229: DSRestoreArray(eps->ds,DS_MAT_Y,&Z);
230: }
231: eps->evecsavailable = PETSC_TRUE;
232: return(0);
233: }
237: /*@
238: EPSSetWorkVecs - Sets a number of work vectors into a EPS object
240: Collective on EPS
242: Input Parameters:
243: + eps - eigensolver context
244: - nw - number of work vectors to allocate
246: Developers Note:
247: This is PETSC_EXTERN because it may be required by user plugin EPS
248: implementations.
250: Level: developer
251: @*/
252: PetscErrorCode EPSSetWorkVecs(EPS eps,PetscInt nw)
253: {
257: if (eps->nwork != nw) {
258: VecDestroyVecs(eps->nwork,&eps->work);
259: eps->nwork = nw;
260: VecDuplicateVecs(eps->t,nw,&eps->work);
261: PetscLogObjectParents(eps,nw,eps->work);
262: }
263: return(0);
264: }
268: /*
269: EPSSetWhichEigenpairs_Default - Sets the default value for which,
270: depending on the ST.
271: */
272: PetscErrorCode EPSSetWhichEigenpairs_Default(EPS eps)
273: {
274: PetscBool target;
278: PetscObjectTypeCompareAny((PetscObject)eps->st,&target,STSINVERT,STCAYLEY,STFOLD,"");
279: if (target) eps->which = EPS_TARGET_MAGNITUDE;
280: else eps->which = EPS_LARGEST_MAGNITUDE;
281: return(0);
282: }
286: /*
287: EPSConvergedEigRelative - Checks convergence relative to the eigenvalue.
288: */
289: PetscErrorCode EPSConvergedEigRelative(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
290: {
291: PetscReal w;
294: w = SlepcAbsEigenvalue(eigr,eigi);
295: *errest = res/w;
296: return(0);
297: }
301: /*
302: EPSConvergedAbsolute - Checks convergence absolutely.
303: */
304: PetscErrorCode EPSConvergedAbsolute(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
305: {
307: *errest = res;
308: return(0);
309: }
313: /*
314: EPSConvergedNormRelative - Checks convergence relative to the eigenvalue and
315: the matrix norms.
316: */
317: PetscErrorCode EPSConvergedNormRelative(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
318: {
319: PetscReal w;
322: w = SlepcAbsEigenvalue(eigr,eigi);
323: *errest = res / (eps->nrma + w*eps->nrmb);
324: return(0);
325: }
329: /*
330: EPSComputeRitzVector - Computes the current Ritz vector.
332: Simple case (complex scalars or real scalars with Zi=NULL):
333: x = V*Zr (V is an array of nv vectors, Zr has length nv)
335: Split case:
336: x = V*Zr y = V*Zi (Zr and Zi have length nv)
337: */
338: PetscErrorCode EPSComputeRitzVector(EPS eps,PetscScalar *Zr,PetscScalar *Zi,Vec *V,PetscInt nv,Vec x,Vec y)
339: {
341: PetscReal norm;
342: #if !defined(PETSC_USE_COMPLEX)
343: Vec z;
344: #endif
347: /* compute eigenvector */
348: SlepcVecMAXPBY(x,0.0,1.0,nv,Zr,V);
350: /* purify eigenvector in positive generalized problems */
351: if (eps->ispositive) {
352: STApply(eps->st,x,y);
353: if (eps->ishermitian) {
354: IPNorm(eps->ip,y,&norm);
355: } else {
356: VecNorm(y,NORM_2,&norm);
357: }
358: VecScale(y,1.0/norm);
359: VecCopy(y,x);
360: }
361: /* fix eigenvector if balancing is used */
362: if (!eps->ishermitian && eps->balance!=EPS_BALANCE_NONE && eps->D) {
363: VecPointwiseDivide(x,x,eps->D);
364: VecNormalize(x,&norm);
365: }
366: #if !defined(PETSC_USE_COMPLEX)
367: /* compute imaginary part of eigenvector */
368: if (Zi) {
369: SlepcVecMAXPBY(y,0.0,1.0,nv,Zi,V);
370: if (eps->ispositive) {
371: VecDuplicate(V[0],&z);
372: STApply(eps->st,y,z);
373: VecNorm(z,NORM_2,&norm);
374: VecScale(z,1.0/norm);
375: VecCopy(z,y);
376: VecDestroy(&z);
377: }
378: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
379: VecPointwiseDivide(y,y,eps->D);
380: VecNormalize(y,&norm);
381: }
382: } else
383: #endif
384: { VecSet(y,0.0); }
385: return(0);
386: }
390: /*
391: EPSBuildBalance_Krylov - uses a Krylov subspace method to compute the
392: diagonal matrix to be applied for balancing in non-Hermitian problems.
393: */
394: PetscErrorCode EPSBuildBalance_Krylov(EPS eps)
395: {
396: Vec z,p,r;
397: PetscInt i,j;
398: PetscReal norma;
399: PetscScalar *pz,*pD;
400: const PetscScalar *pr,*pp;
401: PetscErrorCode ierr;
404: VecDuplicate(eps->V[0],&r);
405: VecDuplicate(eps->V[0],&p);
406: VecDuplicate(eps->V[0],&z);
407: VecSet(eps->D,1.0);
409: for (j=0;j<eps->balance_its;j++) {
411: /* Build a random vector of +-1's */
412: SlepcVecSetRandom(z,eps->rand);
413: VecGetArray(z,&pz);
414: for (i=0;i<eps->nloc;i++) {
415: if (PetscRealPart(pz[i])<0.5) pz[i]=-1.0;
416: else pz[i]=1.0;
417: }
418: VecRestoreArray(z,&pz);
420: /* Compute p=DA(D\z) */
421: VecPointwiseDivide(r,z,eps->D);
422: STApply(eps->st,r,p);
423: VecPointwiseMult(p,p,eps->D);
424: if (j==0) {
425: /* Estimate the matrix inf-norm */
426: VecAbs(p);
427: VecMax(p,NULL,&norma);
428: }
429: if (eps->balance == EPS_BALANCE_TWOSIDE) {
430: /* Compute r=D\(A'Dz) */
431: VecPointwiseMult(z,z,eps->D);
432: STApplyTranspose(eps->st,z,r);
433: VecPointwiseDivide(r,r,eps->D);
434: }
436: /* Adjust values of D */
437: VecGetArrayRead(r,&pr);
438: VecGetArrayRead(p,&pp);
439: VecGetArray(eps->D,&pD);
440: for (i=0;i<eps->nloc;i++) {
441: if (eps->balance == EPS_BALANCE_TWOSIDE) {
442: if (PetscAbsScalar(pp[i])>eps->balance_cutoff*norma && pr[i]!=0.0)
443: pD[i] *= PetscSqrtReal(PetscAbsScalar(pr[i]/pp[i]));
444: } else {
445: if (pp[i]!=0.0) pD[i] *= 1.0/PetscAbsScalar(pp[i]);
446: }
447: }
448: VecRestoreArrayRead(r,&pr);
449: VecRestoreArrayRead(p,&pp);
450: VecRestoreArray(eps->D,&pD);
451: }
453: VecDestroy(&r);
454: VecDestroy(&p);
455: VecDestroy(&z);
456: return(0);
457: }