Actual source code: test29.c

slepc-3.17.0 2022-03-31
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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: */

 11: static char help[] = "Illustrates the computation of left eigenvectors for generalized eigenproblems.\n\n"
 12:   "The command line options are:\n"
 13:   "  -f1 <filename> -f2 <filename>, PETSc binary files containing A and B\n\n";

 15: #include <slepceps.h>

 17: /*
 18:    User-defined routines
 19: */
 20: PetscErrorCode ComputeResidualNorm(Mat,Mat,PetscBool,PetscScalar,PetscScalar,Vec,Vec,Vec*,PetscReal*);

 22: int main(int argc,char **argv)
 23: {
 24:   Mat            A,B;
 25:   EPS            eps;
 26:   EPSType        type;
 27:   PetscInt       i,nconv;
 28:   PetscBool      twosided,flg;
 29:   PetscReal      nrmr,nrml=0.0,re,im,lev;
 30:   PetscScalar    *kr,*ki;
 31:   Vec            t,*xr,*xi,*yr,*yi,*z;
 32:   char           filename[PETSC_MAX_PATH_LEN];
 33:   PetscViewer    viewer;

 35:   SlepcInitialize(&argc,&argv,(char*)0,help);

 37:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 38:         Load the matrices that define the eigensystem, Ax=kBx
 39:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 41:   PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized eigenproblem stored in file.\n\n");
 42:   PetscOptionsGetString(NULL,NULL,"-f1",filename,sizeof(filename),&flg);

 45: #if defined(PETSC_USE_COMPLEX)
 46:   PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrices from binary files...\n");
 47: #else
 48:   PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrices from binary files...\n");
 49: #endif
 50:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
 51:   MatCreate(PETSC_COMM_WORLD,&A);
 52:   MatSetFromOptions(A);
 53:   MatLoad(A,viewer);
 54:   PetscViewerDestroy(&viewer);

 56:   PetscOptionsGetString(NULL,NULL,"-f2",filename,sizeof(filename),&flg);
 57:   if (flg) {
 58:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
 59:     MatCreate(PETSC_COMM_WORLD,&B);
 60:     MatSetFromOptions(B);
 61:     MatLoad(B,viewer);
 62:     PetscViewerDestroy(&viewer);
 63:   } else {
 64:     PetscPrintf(PETSC_COMM_WORLD," Matrix B was not provided, setting B=I\n\n");
 65:     B = NULL;
 66:   }
 67:   MatCreateVecs(A,NULL,&t);

 69:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 70:                 Create the eigensolver and set various options
 71:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 73:   EPSCreate(PETSC_COMM_WORLD,&eps);
 74:   EPSSetOperators(eps,A,B);

 76:   /* use a two-sided algorithm to compute left eigenvectors as well */
 77:   EPSSetTwoSided(eps,PETSC_TRUE);

 79:   /* allow user to change settings at run time */
 80:   EPSSetFromOptions(eps);
 81:   EPSGetTwoSided(eps,&twosided);

 83:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 84:                       Solve the eigensystem
 85:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 87:   EPSSolve(eps);

 89:   /*
 90:      Optional: Get some information from the solver and display it
 91:   */
 92:   EPSGetType(eps,&type);
 93:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);

 95:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 96:                     Display solution and clean up
 97:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 99:   /*
100:      Get number of converged approximate eigenpairs
101:   */
102:   EPSGetConverged(eps,&nconv);
103:   PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %" PetscInt_FMT "\n\n",nconv);
104:   PetscMalloc2(nconv,&kr,nconv,&ki);
105:   VecDuplicateVecs(t,3,&z);
106:   VecDuplicateVecs(t,nconv,&xr);
107:   VecDuplicateVecs(t,nconv,&xi);
108:   if (twosided) {
109:     VecDuplicateVecs(t,nconv,&yr);
110:     VecDuplicateVecs(t,nconv,&yi);
111:   }

113:   if (nconv>0) {
114:     /*
115:        Display eigenvalues and relative errors
116:     */
117:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,
118:          "           k            ||Ax-kBx||         ||y'A-y'Bk||\n"
119:          "   ---------------- ------------------ ------------------\n"));

121:     for (i=0;i<nconv;i++) {
122:       /*
123:         Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
124:         ki (imaginary part)
125:       */
126:       EPSGetEigenpair(eps,i,&kr[i],&ki[i],xr[i],xi[i]);
127:       if (twosided) EPSGetLeftEigenvector(eps,i,yr[i],yi[i]);
128:       /*
129:          Compute the residual norms associated to each eigenpair
130:       */
131:       ComputeResidualNorm(A,B,PETSC_FALSE,kr[i],ki[i],xr[i],xi[i],z,&nrmr);
132:       if (twosided) ComputeResidualNorm(A,B,PETSC_TRUE,kr[i],ki[i],yr[i],yi[i],z,&nrml);

134: #if defined(PETSC_USE_COMPLEX)
135:       re = PetscRealPart(kr[i]);
136:       im = PetscImaginaryPart(kr[i]);
137: #else
138:       re = kr[i];
139:       im = ki[i];
140: #endif
141:       if (im!=0.0) PetscPrintf(PETSC_COMM_WORLD," %8f%+8fi %12g       %12g\n",(double)re,(double)im,(double)nrmr,(double)nrml);
142:       else PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g       %12g\n",(double)re,(double)nrmr,(double)nrml);
143:     }
144:     PetscPrintf(PETSC_COMM_WORLD,"\n");
145:     /*
146:        Check bi-orthogonality of eigenvectors
147:     */
148:     if (twosided) {
149:       VecCheckOrthogonality(xr,nconv,yr,nconv,B,NULL,&lev);
150:       if (lev<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"  Level of bi-orthogonality of eigenvectors < 100*eps\n\n");
151:       else PetscPrintf(PETSC_COMM_WORLD,"  Level of bi-orthogonality of eigenvectors: %g\n\n",(double)lev);
152:     }
153:   }

155:   EPSDestroy(&eps);
156:   MatDestroy(&A);
157:   MatDestroy(&B);
158:   VecDestroy(&t);
159:   PetscFree2(kr,ki);
160:   VecDestroyVecs(3,&z);
161:   VecDestroyVecs(nconv,&xr);
162:   VecDestroyVecs(nconv,&xi);
163:   if (twosided) {
164:     VecDestroyVecs(nconv,&yr);
165:     VecDestroyVecs(nconv,&yi);
166:   }
167:   SlepcFinalize();
168:   return 0;
169: }

171: /*
172:    ComputeResidualNorm - Computes the norm of the residual vector
173:    associated with an eigenpair.

175:    Input Parameters:
176:      trans - whether A' must be used instead of A
177:      kr,ki - eigenvalue
178:      xr,xi - eigenvector
179:      z     - three work vectors (the second one not referenced in complex scalars)
180: */
181: PetscErrorCode ComputeResidualNorm(Mat A,Mat B,PetscBool trans,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,Vec *z,PetscReal *norm)
182: {
183:   Vec            u,w=NULL;
184:   PetscScalar    alpha;
185: #if !defined(PETSC_USE_COMPLEX)
186:   Vec            v;
187:   PetscReal      ni,nr;
188: #endif
189:   PetscErrorCode (*matmult)(Mat,Vec,Vec) = trans? MatMultHermitianTranspose: MatMult;

191:   u = z[0];
192:   if (B) w = z[2];

194: #if !defined(PETSC_USE_COMPLEX)
195:   v = z[1];
196:   if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
197: #endif
198:     (*matmult)(A,xr,u);                          /* u=A*x */
199:     if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
200:       if (B) (*matmult)(B,xr,w);             /* w=B*x */
201:       else w = xr;
202:       alpha = trans? -PetscConj(kr): -kr;
203:       VecAXPY(u,alpha,w);                        /* u=A*x-k*B*x */
204:     }
205:     VecNorm(u,NORM_2,norm);
206: #if !defined(PETSC_USE_COMPLEX)
207:   } else {
208:     (*matmult)(A,xr,u);                          /* u=A*xr */
209:     if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
210:       if (B) (*matmult)(B,xr,v);             /* v=B*xr */
211:       else VecCopy(xr,v);
212:       VecAXPY(u,-kr,v);                          /* u=A*xr-kr*B*xr */
213:       if (B) (*matmult)(B,xi,w);             /* w=B*xi */
214:       else w = xi;
215:       VecAXPY(u,trans?-ki:ki,w);                 /* u=A*xr-kr*B*xr+ki*B*xi */
216:     }
217:     VecNorm(u,NORM_2,&nr);
218:     (*matmult)(A,xi,u);                          /* u=A*xi */
219:     if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
220:       VecAXPY(u,-kr,w);                          /* u=A*xi-kr*B*xi */
221:       VecAXPY(u,trans?ki:-ki,v);                 /* u=A*xi-kr*B*xi-ki*B*xr */
222:     }
223:     VecNorm(u,NORM_2,&ni);
224:     *norm = SlepcAbsEigenvalue(nr,ni);
225:   }
226: #endif
227:   PetscFunctionReturn(0);
228: }

230: /*TEST

232:    testset:
233:       args: -f1 ${SLEPC_DIR}/share/slepc/datafiles/matrices/bfw62a.petsc -f2 ${SLEPC_DIR}/share/slepc/datafiles/matrices/bfw62b.petsc -eps_nev 4 -st_type sinvert -eps_target -190000
234:       filter: grep -v "method" | sed -e "s/[+-]0\.0*i//g" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
235:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
236:       test:
237:          suffix: 1
238:       test:
239:          suffix: 1_rqi
240:          args: -eps_type power -eps_power_shift_type rayleigh -eps_nev 2 -eps_target -2000
241:       test:
242:          suffix: 1_rqi_singular
243:          args: -eps_type power -eps_power_shift_type rayleigh -eps_nev 1 -eps_target -195500

245:    test:
246:       suffix: 2
247:       args: -f1 ${DATAFILESPATH}/matrices/complex/mhd1280a.petsc -f2 ${DATAFILESPATH}/matrices/complex/mhd1280b.petsc -eps_nev 6 -eps_tol 1e-11
248:       filter: sed -e "s/-892/+892/" | sed -e "s/-759/+759/" | sed -e "s/-674/+674/" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
249:       requires: double complex datafilespath !defined(PETSC_USE_64BIT_INDICES)
250:       timeoutfactor: 2

252: TEST*/