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-2009, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
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 private/epsimpl.h
25: #include slepcblaslapack.h
29: PetscErrorCode EPSDestroy_Default(EPS eps)
30: {
35: PetscFree(eps->data);
37: /* free work vectors */
38: EPSDefaultFreeWork(eps);
39: EPSFreeSolution(eps);
40: return(0);
41: }
45: PetscErrorCode EPSBackTransform_Default(EPS eps)
46: {
48: PetscInt i;
52: for (i=0;i<eps->nconv;i++) {
53: STBackTransform(eps->OP,&eps->eigr[i],&eps->eigi[i]);
54: }
55: return(0);
56: }
60: /*
61: EPSComputeVectors_Default - Compute eigenvectors from the vectors
62: provided by the eigensolver. This version just copies the vectors
63: and is intended for solvers such as power that provide the eigenvector.
64: */
65: PetscErrorCode EPSComputeVectors_Default(EPS eps)
66: {
68: eps->evecsavailable = PETSC_TRUE;
69: return(0);
70: }
74: /*
75: EPSComputeVectors_Hermitian - Copies the Lanczos vectors as eigenvectors
76: using purification for generalized eigenproblems.
77: */
78: PetscErrorCode EPSComputeVectors_Hermitian(EPS eps)
79: {
81: PetscInt i;
82: PetscReal norm;
83: Vec w;
86: if (eps->isgeneralized) {
87: /* Purify eigenvectors */
88: VecDuplicate(eps->V[0],&w);
89: for (i=0;i<eps->nconv;i++) {
90: VecCopy(eps->V[i],w);
91: STApply(eps->OP,w,eps->V[i]);
92: VecNormalize(eps->V[i],&norm);
93: }
94: VecDestroy(w);
95: }
96: eps->evecsavailable = PETSC_TRUE;
97: return(0);
98: }
102: /*
103: EPSComputeVectors_Schur - Compute eigenvectors from the vectors
104: provided by the eigensolver. This version is intended for solvers
105: that provide Schur vectors. Given the partial Schur decomposition
106: OP*V=V*T, the following steps are performed:
107: 1) compute eigenvectors of T: T*Z=Z*D
108: 2) compute eigenvectors of OP: X=V*Z
109: If left eigenvectors are required then also do Z'*Tl=D*Z', Y=W*Z
110: */
111: PetscErrorCode EPSComputeVectors_Schur(EPS eps)
112: {
113: #if defined(SLEPC_MISSING_LAPACK_TREVC)
114: SETERRQ(PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable.");
115: #else
117: PetscInt i;
118: PetscBLASInt ncv,nconv,mout,info;
119: PetscScalar *Z,*work;
120: #if defined(PETSC_USE_COMPLEX)
121: PetscReal *rwork;
122: #endif
123: PetscReal norm;
124: Vec w;
125:
127: ncv = PetscBLASIntCast(eps->ncv);
128: nconv = PetscBLASIntCast(eps->nconv);
129: if (eps->ishermitian) {
130: EPSComputeVectors_Hermitian(eps);
131: return(0);
132: }
133: if (eps->ispositive) {
134: VecDuplicate(eps->V[0],&w);
135: }
137: PetscMalloc(nconv*nconv*sizeof(PetscScalar),&Z);
138: PetscMalloc(3*nconv*sizeof(PetscScalar),&work);
139: #if defined(PETSC_USE_COMPLEX)
140: PetscMalloc(nconv*sizeof(PetscReal),&rwork);
141: #endif
143: /* right eigenvectors */
144: #if !defined(PETSC_USE_COMPLEX)
145: LAPACKtrevc_("R","A",PETSC_NULL,&nconv,eps->T,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,&info);
146: #else
147: LAPACKtrevc_("R","A",PETSC_NULL,&nconv,eps->T,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,rwork,&info);
148: #endif
149: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);
151: /* AV = V * Z */
152: SlepcUpdateVectors(eps->nconv,eps->V,0,eps->nconv,Z,eps->nconv,PETSC_FALSE);
153: if (eps->ispositive) {
154: /* Purify eigenvectors */
155: for (i=0;i<eps->nconv;i++) {
156: VecCopy(eps->V[i],w);
157: STApply(eps->OP,w,eps->V[i]);
158: VecNormalize(eps->V[i],&norm);
159: }
160: }
161:
162: /* left eigenvectors */
163: if (eps->solverclass == EPS_TWO_SIDE) {
164: #if !defined(PETSC_USE_COMPLEX)
165: LAPACKtrevc_("R","A",PETSC_NULL,&nconv,eps->Tl,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,&info);
166: #else
167: LAPACKtrevc_("R","A",PETSC_NULL,&nconv,eps->Tl,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,rwork,&info);
168: #endif
169: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);
171: /* AW = W * Z */
172: SlepcUpdateVectors(eps->nconv,eps->W,0,eps->nconv,Z,eps->nconv,PETSC_FALSE);
173: }
174:
175: PetscFree(Z);
176: PetscFree(work);
177: #if defined(PETSC_USE_COMPLEX)
178: PetscFree(rwork);
179: #endif
180: if (eps->ispositive) {
181: VecDestroy(w);
182: }
183: eps->evecsavailable = PETSC_TRUE;
184: return(0);
185: #endif
186: }
190: /*
191: EPSDefaultGetWork - Gets a number of work vectors.
193: Input Parameters:
194: + eps - eigensolver context
195: - nw - number of work vectors to allocate
197: Notes:
198: Call this only if no work vectors have been allocated.
200: */
201: PetscErrorCode EPSDefaultGetWork(EPS eps, PetscInt nw)
202: {
207: if (eps->nwork != nw) {
208: if (eps->nwork > 0) {
209: VecDestroyVecs(eps->work,eps->nwork);
210: }
211: eps->nwork = nw;
212: VecDuplicateVecs(eps->vec_initial,nw,&eps->work);
213: PetscLogObjectParents(eps,nw,eps->work);
214: }
215:
216: return(0);
217: }
221: /*
222: EPSDefaultFreeWork - Free work vectors.
224: Input Parameters:
225: . eps - eigensolver context
227: */
228: PetscErrorCode EPSDefaultFreeWork(EPS eps)
229: {
234: if (eps->work) {
235: VecDestroyVecs(eps->work,eps->nwork);
236: }
237: return(0);
238: }