Actual source code: pepkrylov.c
slepc-3.7.0 2016-05-16
1: /*
2: Common subroutines for all Krylov-type PEP solvers.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2016, 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/pepimpl.h> /*I "slepcpep.h" I*/
25: #include <slepcblaslapack.h>
26: #include pepkrylov.h
30: PetscErrorCode PEPExtractVectors_TOAR(PEP pep)
31: {
33: PetscInt i,j,deg=pep->nmat-1,lds,idxcpy=0,ldds,k,ld;
34: PetscScalar *X,*er,*ei,*SS,*vals,*ivals,sone=1.0,szero=0.0,*yi,*yr,*tr,*ti,alpha,t,*S,*pS0;
35: PetscBLASInt k_,lds_,one=1,ldds_;
36: PetscBool flg;
37: PetscReal norm,max;
38: Vec xr,xi,w[4];
39: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
40: Mat S0;
43: S = ctx->S;
44: ld = ctx->ld;
45: k = pep->nconv;
46: if (k==0) return(0);
47: lds = deg*ld;
48: DSGetLeadingDimension(pep->ds,&ldds);
49: PetscMalloc5(k,&er,k,&ei,k*k,&SS,pep->nmat,&vals,pep->nmat,&ivals);
50: STGetTransform(pep->st,&flg);
51: for (i=0;i<k;i++) {
52: er[i] = pep->eigr[i];
53: ei[i] = pep->eigi[i];
54: }
55: if (flg) {
56: for (i=0;i<k;i++) {
57: er[i] = pep->sfactor*pep->eigr[i];
58: ei[i] = pep->sfactor*pep->eigi[i];
59: }
60: }
61: STBackTransform(pep->st,k,er,ei);
63: DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
64: DSGetArray(pep->ds,DS_MAT_X,&X);
66: PetscBLASIntCast(k,&k_);
67: PetscBLASIntCast(lds,&lds_);
68: PetscBLASIntCast(ldds,&ldds_);
70: if (pep->extract==PEP_EXTRACT_NONE || pep->refine==PEP_REFINE_MULTIPLE) {
71: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&sone,S,&lds_,X,&ldds_,&szero,SS,&k_));
72: } else {
73: switch (pep->extract) {
74: case PEP_EXTRACT_NONE:
75: break;
76: case PEP_EXTRACT_NORM:
77: for (i=0;i<k;i++) {
78: PEPEvaluateBasis(pep,er[i],ei[i],vals,ivals);
79: max = 1.0;
80: for (j=1;j<deg;j++) {
81: norm = SlepcAbsEigenvalue(vals[j],ivals[j]);
82: if (max < norm) { max = norm; idxcpy = j; }
83: }
84: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+idxcpy*ld,&lds_,X+i*ldds,&one,&szero,SS+i*k,&one));
85: #if !defined(PETSC_USE_COMPLEX)
86: if (PetscRealPart(ei[i])!=0.0) {
87: i++;
88: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+idxcpy*ld,&lds_,X+i*ldds,&one,&szero,SS+i*k,&one));
89: }
90: #endif
91: }
92: break;
93: case PEP_EXTRACT_RESIDUAL:
94: VecDuplicate(pep->work[0],&xr);
95: VecDuplicate(pep->work[0],&w[0]);
96: VecDuplicate(pep->work[0],&w[1]);
97: #if !defined(PETSC_USE_COMPLEX)
98: VecDuplicate(pep->work[0],&w[2]);
99: VecDuplicate(pep->work[0],&w[3]);
100: VecDuplicate(pep->work[0],&xi);
101: #else
102: xi = NULL;
103: #endif
104: for (i=0;i<k;i++) {
105: max = 0.0;
106: for (j=0;j<deg;j++) {
107: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+j*ld,&lds_,X+i*ldds,&one,&szero,SS+i*k,&one));
108: BVMultVec(pep->V,1.0,0.0,xr,SS+i*k);
109: #if !defined(PETSC_USE_COMPLEX)
110: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+j*ld,&lds_,X+(i+1)*ldds,&one,&szero,SS+i*k,&one));
111: BVMultVec(pep->V,1.0,0.0,xi,SS+i*k);
112: #endif
113: PEPComputeResidualNorm_Private(pep,er[i],ei[i],xr,xi,w,&norm);
114: if (norm>max) { max = norm; idxcpy=j; }
115: }
116: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+idxcpy*ld,&lds_,X+i*ldds,&one,&szero,SS+i*k,&one));
117: #if !defined(PETSC_USE_COMPLEX)
118: if (PetscRealPart(ei[i])!=0.0) {
119: i++;
120: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+idxcpy*ld,&lds_,X+i*ldds,&one,&szero,SS+i*k,&one));
121: }
122: #endif
123: }
124: VecDestroy(&xr);
125: VecDestroy(&w[0]);
126: VecDestroy(&w[1]);
127: #if !defined(PETSC_USE_COMPLEX)
128: VecDestroy(&w[2]);
129: VecDestroy(&w[3]);
130: VecDestroy(&xi);
131: #endif
132: break;
133: case PEP_EXTRACT_STRUCTURED:
134: PetscMalloc2(k,&tr,k,&ti);
135: for (i=0;i<k;i++) {
136: PetscMemzero(SS+i*k,k*sizeof(PetscScalar));
137: #if !defined(PETSC_USE_COMPLEX)
138: if (ei[i]!=0.0) {
139: PetscMemzero(SS+(i+1)*k,k*sizeof(PetscScalar));
140: }
141: #endif
142: t = 0.0;
143: PEPEvaluateBasis(pep,er[i],ei[i],vals,ivals);
144: yr = X+i*ldds; yi = NULL;
145: for (j=0;j<deg;j++) {
146: alpha = PetscConj(vals[j]);
147: #if !defined(PETSC_USE_COMPLEX)
148: if (ei[i]!=0.0) {
149: PetscMemzero(tr,k*sizeof(PetscScalar));
150: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k_,&vals[j],X+i*ldds,&one,tr,&one));
151: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k_,&ivals[j],X+(i+1)*ldds,&one,tr,&one));
152: yr = tr;
153: PetscMemzero(ti,k*sizeof(PetscScalar));
154: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k_,&vals[j],X+(i+1)*ldds,&one,ti,&one));
155: alpha = -ivals[j];
156: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k_,&alpha,X+i*ldds,&one,ti,&one));
157: yi = ti;
158: alpha = 1.0;
159: } else { yr = X+i*ldds; yi = NULL;}
160: #endif
161: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&alpha,S+j*ld,&lds_,yr,&one,&sone,SS+i*k,&one));
162: t += SlepcAbsEigenvalue(vals[j],ivals[j])*SlepcAbsEigenvalue(vals[j],ivals[j]);
163: if (yi) {
164: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&alpha,S+j*ld,&lds_,yi,&one,&sone,SS+(i+1)*k,&one));
165: }
166: }
167: t = 1.0/t;
168: PetscStackCallBLAS("BLASscal",BLASscal_(&k_,&t,SS+i*k,&one));
169: if (yi) {
170: PetscStackCallBLAS("BLASscal",BLASscal_(&k_,&t,SS+(i+1)*k,&one));
171: i++;
172: }
173: }
174: PetscFree2(tr,ti);
175: break;
176: default:
177: SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Extraction not implemented in this solver");
178: }
179: }
181: /* update vectors V = V*S */
182: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&S0);
183: MatDenseGetArray(S0,&pS0);
184: for (i=0;i<k;i++) {
185: PetscMemcpy(pS0+i*k,SS+i*k,k*sizeof(PetscScalar));
186: }
187: MatDenseRestoreArray(S0,&pS0);
188: BVSetActiveColumns(pep->V,0,k);
189: BVMultInPlace(pep->V,S0,0,k);
190: MatDestroy(&S0);
191: BVSetActiveColumns(pep->V,0,k);
192: PetscFree5(er,ei,SS,vals,ivals);
193: return(0);
194: }
198: PetscErrorCode PEPReset_TOAR(PEP pep)
199: {
200: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
204: PetscFree(ctx->S);
205: PetscFree(ctx->qB);
206: return(0);
207: }