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-2007, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc. See the README file for conditions of use
9: and additional information.
10: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
11: */
13: #include src/eps/epsimpl.h
14: #include slepcblaslapack.h
18: PetscErrorCode EPSDestroy_Default(EPS eps)
19: {
24: PetscFree(eps->data);
26: /* free work vectors */
27: EPSDefaultFreeWork(eps);
28: EPSFreeSolution(eps);
29: return(0);
30: }
34: PetscErrorCode EPSBackTransform_Default(EPS eps)
35: {
37: int i;
41: for (i=0;i<eps->nconv;i++) {
42: STBackTransform(eps->OP,&eps->eigr[i],&eps->eigi[i]);
43: }
44: return(0);
45: }
49: /*
50: EPSComputeVectors_Default - Compute eigenvectors from the vectors
51: provided by the eigensolver. This version just copies the vectors
52: and is intended for solvers such as power that provide the eigenvector.
53: */
54: PetscErrorCode EPSComputeVectors_Default(EPS eps)
55: {
57: int i;
60: for (i=0;i<eps->nconv;i++) {
61: VecCopy(eps->V[i],eps->AV[i]);
62: if (eps->solverclass == EPS_TWO_SIDE) {
63: VecCopy(eps->W[i],eps->AW[i]);
64: }
65: }
66: eps->evecsavailable = PETSC_TRUE;
67: return(0);
68: }
72: /*
73: EPSComputeVectors_Hermitian - Copies the Lanczos vectors as eigenvectors
74: using purification for generalized eigenproblems.
75: */
76: PetscErrorCode EPSComputeVectors_Hermitian(EPS eps)
77: {
79: int i;
80: PetscReal norm;
83: for (i=0;i<eps->nconv;i++) {
84: if (eps->isgeneralized) {
85: /* Purify eigenvectors */
86: STApply(eps->OP,eps->V[i],eps->AV[i]);
87: VecNormalize(eps->AV[i],&norm);
88: } else {
89: VecCopy(eps->V[i],eps->AV[i]);
90: }
91: if (eps->solverclass == EPS_TWO_SIDE) {
92: VecCopy(eps->W[i],eps->AW[i]);
93: }
94: }
95: eps->evecsavailable = PETSC_TRUE;
96: return(0);
97: }
101: /*
102: EPSComputeVectors_Schur - Compute eigenvectors from the vectors
103: provided by the eigensolver. This version is intended for solvers
104: that provide Schur vectors. Given the partial Schur decomposition
105: OP*V=V*T, the following steps are performed:
106: 1) compute eigenvectors of T: T*Z=Z*D
107: 2) compute eigenvectors of OP: X=V*Z
108: If left eigenvectors are required then also do Z'*Tl=D*Z', Y=W*Z
109: */
110: PetscErrorCode EPSComputeVectors_Schur(EPS eps)
111: {
112: #if defined(SLEPC_MISSING_LAPACK_TREVC)
113: SETERRQ(PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable.");
114: #else
116: int i,mout,info,nv=eps->nv;
117: PetscScalar *Z,*work;
118: #if defined(PETSC_USE_COMPLEX)
119: PetscReal *rwork;
120: #endif
121: PetscReal norm;
122: Vec w;
123:
125: if (eps->ishermitian) {
126: EPSComputeVectors_Hermitian(eps);
127: return(0);
128: }
129: if (eps->ispositive) {
130: VecDuplicate(eps->V[0],&w);
131: }
133: PetscMalloc(nv*nv*sizeof(PetscScalar),&Z);
134: PetscMalloc(3*nv*sizeof(PetscScalar),&work);
135: #if defined(PETSC_USE_COMPLEX)
136: PetscMalloc(nv*sizeof(PetscReal),&rwork);
137: #endif
139: /* right eigenvectors */
140: #if !defined(PETSC_USE_COMPLEX)
141: LAPACKtrevc_("R","A",PETSC_NULL,&nv,eps->T,&eps->ncv,PETSC_NULL,&nv,Z,&nv,&nv,&mout,work,&info,1,1);
142: #else
143: LAPACKtrevc_("R","A",PETSC_NULL,&nv,eps->T,&eps->ncv,PETSC_NULL,&nv,Z,&nv,&nv,&mout,work,rwork,&info,1,1);
144: #endif
145: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);
147: /* AV = V * Z */
148: for (i=0;i<eps->nconv;i++) {
149: if (eps->ispositive) {
150: /* Purify eigenvectors */
151: VecSet(w,0.0);
152: VecMAXPY(w,nv,Z+nv*i,eps->V);
153: STApply(eps->OP,w,eps->AV[i]);
154: VecNormalize(eps->AV[i],&norm);
155: } else {
156: VecSet(eps->AV[i],0.0);
157: VecMAXPY(eps->AV[i],nv,Z+nv*i,eps->V);
158: }
159: }
160:
161: /* left eigenvectors */
162: if (eps->solverclass == EPS_TWO_SIDE) {
163: #if !defined(PETSC_USE_COMPLEX)
164: LAPACKtrevc_("R","A",PETSC_NULL,&nv,eps->Tl,&eps->nv,PETSC_NULL,&nv,Z,&nv,&nv,&mout,work,&info,1,1);
165: #else
166: LAPACKtrevc_("R","A",PETSC_NULL,&nv,eps->Tl,&eps->nv,PETSC_NULL,&nv,Z,&nv,&nv,&mout,work,rwork,&info,1,1);
167: #endif
168: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);
170: /* AW = W * Z */
171: for (i=0;i<eps->nconv;i++) {
172: VecSet(eps->AW[i],0.0);
173: VecMAXPY(eps->AW[i],nv,Z+nv*i,eps->W);
174: }
175: }
176:
177: PetscFree(Z);
178: PetscFree(work);
179: #if defined(PETSC_USE_COMPLEX)
180: PetscFree(rwork);
181: #endif
182: if (eps->ispositive) {
183: VecDestroy(w);
184: }
185: eps->evecsavailable = PETSC_TRUE;
186: return(0);
187: #endif
188: }
192: /*
193: EPSDefaultGetWork - Gets a number of work vectors.
195: Input Parameters:
196: + eps - eigensolver context
197: - nw - number of work vectors to allocate
199: Notes:
200: Call this only if no work vectors have been allocated.
202: */
203: PetscErrorCode EPSDefaultGetWork(EPS eps, int nw)
204: {
209: if (eps->nwork != nw) {
210: if (eps->nwork > 0) {
211: VecDestroyVecs(eps->work,eps->nwork);
212: }
213: eps->nwork = nw;
214: VecDuplicateVecs(eps->vec_initial,nw,&eps->work);
215: PetscLogObjectParents(eps,nw,eps->work);
216: }
217:
218: return(0);
219: }
223: /*
224: EPSDefaultFreeWork - Free work vectors.
226: Input Parameters:
227: . eps - eigensolver context
229: */
230: PetscErrorCode EPSDefaultFreeWork(EPS eps)
231: {
236: if (eps->work) {
237: VecDestroyVecs(eps->work,eps->nwork);
238: }
239: return(0);
240: }