Actual source code: pepkrylov.c
slepc-3.7.4 2017-05-17
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,factor=1.0;
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: PetscCalloc5(k,&er,k,&ei,k*k,&SS,pep->nmat,&vals,pep->nmat,&ivals);
50: STGetTransform(pep->st,&flg);
51: if (flg) factor = pep->sfactor;
52: for (i=0;i<k;i++) {
53: er[i] = factor*pep->eigr[i];
54: ei[i] = factor*pep->eigi[i];
55: }
56: STBackTransform(pep->st,k,er,ei);
58: DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
59: DSGetArray(pep->ds,DS_MAT_X,&X);
61: PetscBLASIntCast(k,&k_);
62: PetscBLASIntCast(lds,&lds_);
63: PetscBLASIntCast(ldds,&ldds_);
65: if (pep->extract==PEP_EXTRACT_NONE || pep->refine==PEP_REFINE_MULTIPLE) {
66: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&sone,S,&lds_,X,&ldds_,&szero,SS,&k_));
67: } else {
68: switch (pep->extract) {
69: case PEP_EXTRACT_NONE:
70: break;
71: case PEP_EXTRACT_NORM:
72: for (i=0;i<k;i++) {
73: PEPEvaluateBasis(pep,er[i],ei[i],vals,ivals);
74: max = 1.0;
75: for (j=1;j<deg;j++) {
76: norm = SlepcAbsEigenvalue(vals[j],ivals[j]);
77: if (max < norm) { max = norm; idxcpy = j; }
78: }
79: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+idxcpy*ld,&lds_,X+i*ldds,&one,&szero,SS+i*k,&one));
80: #if !defined(PETSC_USE_COMPLEX)
81: if (PetscRealPart(ei[i])!=0.0) {
82: i++;
83: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+idxcpy*ld,&lds_,X+i*ldds,&one,&szero,SS+i*k,&one));
84: }
85: #endif
86: }
87: break;
88: case PEP_EXTRACT_RESIDUAL:
89: VecDuplicate(pep->work[0],&xr);
90: VecDuplicate(pep->work[0],&w[0]);
91: VecDuplicate(pep->work[0],&w[1]);
92: #if !defined(PETSC_USE_COMPLEX)
93: VecDuplicate(pep->work[0],&w[2]);
94: VecDuplicate(pep->work[0],&w[3]);
95: VecDuplicate(pep->work[0],&xi);
96: #else
97: xi = NULL;
98: #endif
99: for (i=0;i<k;i++) {
100: max = 0.0;
101: for (j=0;j<deg;j++) {
102: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+j*ld,&lds_,X+i*ldds,&one,&szero,SS+i*k,&one));
103: BVMultVec(pep->V,1.0,0.0,xr,SS+i*k);
104: #if !defined(PETSC_USE_COMPLEX)
105: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+j*ld,&lds_,X+(i+1)*ldds,&one,&szero,SS+i*k,&one));
106: BVMultVec(pep->V,1.0,0.0,xi,SS+i*k);
107: #endif
108: PEPComputeResidualNorm_Private(pep,er[i],ei[i],xr,xi,w,&norm);
109: if (norm>max) { max = norm; idxcpy=j; }
110: }
111: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+idxcpy*ld,&lds_,X+i*ldds,&one,&szero,SS+i*k,&one));
112: #if !defined(PETSC_USE_COMPLEX)
113: if (PetscRealPart(ei[i])!=0.0) {
114: i++;
115: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+idxcpy*ld,&lds_,X+i*ldds,&one,&szero,SS+i*k,&one));
116: }
117: #endif
118: }
119: VecDestroy(&xr);
120: VecDestroy(&w[0]);
121: VecDestroy(&w[1]);
122: #if !defined(PETSC_USE_COMPLEX)
123: VecDestroy(&w[2]);
124: VecDestroy(&w[3]);
125: VecDestroy(&xi);
126: #endif
127: break;
128: case PEP_EXTRACT_STRUCTURED:
129: PetscMalloc2(k,&tr,k,&ti);
130: for (i=0;i<k;i++) {
131: t = 0.0;
132: PEPEvaluateBasis(pep,er[i],ei[i],vals,ivals);
133: yr = X+i*ldds; yi = NULL;
134: for (j=0;j<deg;j++) {
135: alpha = PetscConj(vals[j]);
136: #if !defined(PETSC_USE_COMPLEX)
137: if (ei[i]!=0.0) {
138: PetscMemzero(tr,k*sizeof(PetscScalar));
139: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k_,&vals[j],X+i*ldds,&one,tr,&one));
140: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k_,&ivals[j],X+(i+1)*ldds,&one,tr,&one));
141: yr = tr;
142: PetscMemzero(ti,k*sizeof(PetscScalar));
143: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k_,&vals[j],X+(i+1)*ldds,&one,ti,&one));
144: alpha = -ivals[j];
145: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k_,&alpha,X+i*ldds,&one,ti,&one));
146: yi = ti;
147: alpha = 1.0;
148: } else { yr = X+i*ldds; yi = NULL; }
149: #endif
150: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&alpha,S+j*ld,&lds_,yr,&one,&sone,SS+i*k,&one));
151: t += SlepcAbsEigenvalue(vals[j],ivals[j])*SlepcAbsEigenvalue(vals[j],ivals[j]);
152: if (yi) {
153: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&alpha,S+j*ld,&lds_,yi,&one,&sone,SS+(i+1)*k,&one));
154: }
155: }
156: t = 1.0/t;
157: PetscStackCallBLAS("BLASscal",BLASscal_(&k_,&t,SS+i*k,&one));
158: if (yi) {
159: PetscStackCallBLAS("BLASscal",BLASscal_(&k_,&t,SS+(i+1)*k,&one));
160: i++;
161: }
162: }
163: PetscFree2(tr,ti);
164: break;
165: default:
166: SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Extraction not implemented in this solver");
167: }
168: }
170: /* update vectors V = V*S */
171: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&S0);
172: MatDenseGetArray(S0,&pS0);
173: for (i=0;i<k;i++) {
174: PetscMemcpy(pS0+i*k,SS+i*k,k*sizeof(PetscScalar));
175: }
176: MatDenseRestoreArray(S0,&pS0);
177: BVMultInPlace(pep->V,S0,0,k);
178: MatDestroy(&S0);
179: PetscFree5(er,ei,SS,vals,ivals);
180: return(0);
181: }
185: PetscErrorCode PEPReset_TOAR(PEP pep)
186: {
187: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
191: PetscFree(ctx->S);
192: PetscFree(ctx->qB);
193: return(0);
194: }