Actual source code: pepkrylov.c

slepc-3.8.2 2017-12-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    Common subroutines for all Krylov-type PEP solvers
 12: */

 14: #include <slepc/private/pepimpl.h>    /*I "slepcpep.h" I*/
 15: #include <slepcblaslapack.h>
 16: #include "pepkrylov.h"

 18: PetscErrorCode PEPExtractVectors_TOAR(PEP pep)
 19: {
 21:   PetscInt       i,j,deg=pep->nmat-1,lds,idxcpy=0,ldds,k,ld;
 22:   PetscScalar    *X,*er,*ei,*SS,*vals,*ivals,sone=1.0,szero=0.0,*yi,*yr,*tr,*ti,alpha,t,*S,*pS0;
 23:   PetscBLASInt   k_,lds_,one=1,ldds_;
 24:   PetscBool      flg;
 25:   PetscReal      norm,max,factor=1.0;
 26:   Vec            xr,xi,w[4];
 27:   PEP_TOAR       *ctx = (PEP_TOAR*)pep->data;
 28:   Mat            S0;

 31:   S = ctx->S;
 32:   ld = ctx->ld;
 33:   k = pep->nconv;
 34:   if (k==0) return(0);
 35:   lds = deg*ld;
 36:   DSGetLeadingDimension(pep->ds,&ldds);
 37:   PetscCalloc5(k,&er,k,&ei,k*k,&SS,pep->nmat,&vals,pep->nmat,&ivals);
 38:   STGetTransform(pep->st,&flg);
 39:   if (flg) factor = pep->sfactor;
 40:   for (i=0;i<k;i++) {
 41:     er[i] = factor*pep->eigr[i];
 42:     ei[i] = factor*pep->eigi[i];
 43:   }
 44:   STBackTransform(pep->st,k,er,ei);

 46:   DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
 47:   DSGetArray(pep->ds,DS_MAT_X,&X);

 49:   PetscBLASIntCast(k,&k_);
 50:   PetscBLASIntCast(lds,&lds_);
 51:   PetscBLASIntCast(ldds,&ldds_);

 53:   if (pep->extract==PEP_EXTRACT_NONE || pep->refine==PEP_REFINE_MULTIPLE) {
 54:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&sone,S,&lds_,X,&ldds_,&szero,SS,&k_));
 55:   } else {
 56:     switch (pep->extract) {
 57:     case PEP_EXTRACT_NONE:
 58:       break;
 59:     case PEP_EXTRACT_NORM:
 60:       for (i=0;i<k;i++) {
 61:         PEPEvaluateBasis(pep,er[i],ei[i],vals,ivals);
 62:         max = 1.0;
 63:         for (j=1;j<deg;j++) {
 64:           norm = SlepcAbsEigenvalue(vals[j],ivals[j]);
 65:           if (max < norm) { max = norm; idxcpy = j; }
 66:         }
 67:         PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+idxcpy*ld,&lds_,X+i*ldds,&one,&szero,SS+i*k,&one));
 68: #if !defined(PETSC_USE_COMPLEX)
 69:         if (PetscRealPart(ei[i])!=0.0) {
 70:           i++;
 71:           PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+idxcpy*ld,&lds_,X+i*ldds,&one,&szero,SS+i*k,&one));
 72:         }
 73: #endif
 74:       }
 75:       break;
 76:     case PEP_EXTRACT_RESIDUAL:
 77:       VecDuplicate(pep->work[0],&xr);
 78:       VecDuplicate(pep->work[0],&w[0]);
 79:       VecDuplicate(pep->work[0],&w[1]);
 80: #if !defined(PETSC_USE_COMPLEX)
 81:       VecDuplicate(pep->work[0],&w[2]);
 82:       VecDuplicate(pep->work[0],&w[3]);
 83:       VecDuplicate(pep->work[0],&xi);
 84: #else
 85:       xi = NULL;
 86: #endif
 87:       for (i=0;i<k;i++) {
 88:         max = 0.0;
 89:         for (j=0;j<deg;j++) {
 90:           PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+j*ld,&lds_,X+i*ldds,&one,&szero,SS+i*k,&one));
 91:           BVMultVec(pep->V,1.0,0.0,xr,SS+i*k);
 92: #if !defined(PETSC_USE_COMPLEX)
 93:           PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+j*ld,&lds_,X+(i+1)*ldds,&one,&szero,SS+i*k,&one));
 94:           BVMultVec(pep->V,1.0,0.0,xi,SS+i*k);
 95: #endif
 96:           PEPComputeResidualNorm_Private(pep,er[i],ei[i],xr,xi,w,&norm);
 97:           if (norm>max) { max = norm; idxcpy=j; }
 98:         }
 99:         PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+idxcpy*ld,&lds_,X+i*ldds,&one,&szero,SS+i*k,&one));
100: #if !defined(PETSC_USE_COMPLEX)
101:         if (PetscRealPart(ei[i])!=0.0) {
102:           i++;
103:           PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S+idxcpy*ld,&lds_,X+i*ldds,&one,&szero,SS+i*k,&one));
104:         }
105: #endif
106:       }
107:       VecDestroy(&xr);
108:       VecDestroy(&w[0]);
109:       VecDestroy(&w[1]);
110: #if !defined(PETSC_USE_COMPLEX)
111:       VecDestroy(&w[2]);
112:       VecDestroy(&w[3]);
113:       VecDestroy(&xi);
114: #endif
115:       break;
116:     case PEP_EXTRACT_STRUCTURED:
117:       PetscMalloc2(k,&tr,k,&ti);
118:       for (i=0;i<k;i++) {
119:         t = 0.0;
120:         PEPEvaluateBasis(pep,er[i],ei[i],vals,ivals);
121:         yr = X+i*ldds; yi = NULL;
122:         for (j=0;j<deg;j++) {
123:           alpha = PetscConj(vals[j]);
124: #if !defined(PETSC_USE_COMPLEX)
125:           if (ei[i]!=0.0) {
126:             PetscMemzero(tr,k*sizeof(PetscScalar));
127:             PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k_,&vals[j],X+i*ldds,&one,tr,&one));
128:             PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k_,&ivals[j],X+(i+1)*ldds,&one,tr,&one));
129:             yr = tr;
130:             PetscMemzero(ti,k*sizeof(PetscScalar));
131:             PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k_,&vals[j],X+(i+1)*ldds,&one,ti,&one));
132:             alpha = -ivals[j];
133:             PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k_,&alpha,X+i*ldds,&one,ti,&one));
134:             yi = ti;
135:             alpha = 1.0;
136:           } else { yr = X+i*ldds; yi = NULL; }
137: #endif
138:           PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&alpha,S+j*ld,&lds_,yr,&one,&sone,SS+i*k,&one));
139:           t += SlepcAbsEigenvalue(vals[j],ivals[j])*SlepcAbsEigenvalue(vals[j],ivals[j]);
140:           if (yi) {
141:             PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&alpha,S+j*ld,&lds_,yi,&one,&sone,SS+(i+1)*k,&one));
142:           }
143:         }
144:         t = 1.0/t;
145:         PetscStackCallBLAS("BLASscal",BLASscal_(&k_,&t,SS+i*k,&one));
146:         if (yi) {
147:           PetscStackCallBLAS("BLASscal",BLASscal_(&k_,&t,SS+(i+1)*k,&one));
148:           i++;
149:         }
150:       }
151:       PetscFree2(tr,ti);
152:       break;
153:     }
154:   }

156:   /* update vectors V = V*S */
157:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&S0);
158:   MatDenseGetArray(S0,&pS0);
159:   for (i=0;i<k;i++) {
160:     PetscMemcpy(pS0+i*k,SS+i*k,k*sizeof(PetscScalar));
161:   }
162:   MatDenseRestoreArray(S0,&pS0);
163:   BVMultInPlace(pep->V,S0,0,k);
164:   MatDestroy(&S0);
165:   PetscFree5(er,ei,SS,vals,ivals);
166:   return(0);
167: }