| Line | Branch | Exec | Source |
|---|---|---|---|
| 1 | /* | ||
| 2 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 3 | SLEPc - Scalable Library for Eigenvalue Problem Computations | ||
| 4 | Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain | ||
| 5 | |||
| 6 | This file is part of SLEPc. | ||
| 7 | SLEPc is distributed under a 2-clause BSD license (see LICENSE). | ||
| 8 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 9 | */ | ||
| 10 | |||
| 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"; | ||
| 14 | |||
| 15 | #include <slepceps.h> | ||
| 16 | |||
| 17 | /* | ||
| 18 | User-defined routines | ||
| 19 | */ | ||
| 20 | PetscErrorCode ComputeResidualNorm(Mat,Mat,PetscBool,PetscScalar,PetscScalar,Vec,Vec,Vec*,PetscReal*); | ||
| 21 | PetscErrorCode CheckBiorthogonality(Vec*,Vec*,PetscScalar*,PetscInt,Mat,PetscReal*); | ||
| 22 | |||
| 23 | 28 | int main(int argc,char **argv) | |
| 24 | { | ||
| 25 | 28 | Mat A,B; | |
| 26 | 28 | EPS eps; | |
| 27 | 28 | EPSType type; | |
| 28 | 28 | PetscInt i,nconv; | |
| 29 | 28 | PetscBool twosided,flg; | |
| 30 | 28 | PetscReal nrmr,nrml=0.0,re,im,lev; | |
| 31 | 28 | PetscScalar *kr,*ki; | |
| 32 | 28 | Vec t,*xr,*xi,*yr,*yi,*z; | |
| 33 | 28 | char filename[PETSC_MAX_PATH_LEN]; | |
| 34 | 28 | PetscViewer viewer; | |
| 35 | |||
| 36 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
28 | PetscFunctionBeginUser; |
| 37 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
| 38 | |||
| 39 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 40 | Load the matrices that define the eigensystem, Ax=kBx | ||
| 41 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 42 | |||
| 43 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized eigenproblem stored in file.\n\n")); |
| 44 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscOptionsGetString(NULL,NULL,"-f1",filename,sizeof(filename),&flg)); |
| 45 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
28 | PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Must indicate a file name for matrix A with the -f1 option"); |
| 46 | |||
| 47 | #if defined(PETSC_USE_COMPLEX) | ||
| 48 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
3 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrices from binary files...\n")); |
| 49 | #else | ||
| 50 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
25 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrices from binary files...\n")); |
| 51 | #endif | ||
| 52 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer)); |
| 53 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); |
| 54 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(MatSetFromOptions(A)); |
| 55 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(MatLoad(A,viewer)); |
| 56 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscViewerDestroy(&viewer)); |
| 57 | |||
| 58 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscOptionsGetString(NULL,NULL,"-f2",filename,sizeof(filename),&flg)); |
| 59 |
1/2✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
|
28 | if (flg) { |
| 60 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer)); |
| 61 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); |
| 62 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(MatSetFromOptions(B)); |
| 63 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(MatLoad(B,viewer)); |
| 64 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscViewerDestroy(&viewer)); |
| 65 | } else { | ||
| 66 | ✗ | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Matrix B was not provided, setting B=I\n\n")); | |
| 67 | ✗ | B = NULL; | |
| 68 | } | ||
| 69 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(MatCreateVecs(A,NULL,&t)); |
| 70 | |||
| 71 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 72 | Create the eigensolver and set various options | ||
| 73 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 74 | |||
| 75 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps)); |
| 76 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(EPSSetOperators(eps,A,B)); |
| 77 | |||
| 78 | /* use a two-sided algorithm to compute left eigenvectors as well */ | ||
| 79 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(EPSSetTwoSided(eps,PETSC_TRUE)); |
| 80 | |||
| 81 | /* allow user to change settings at run time */ | ||
| 82 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(EPSSetFromOptions(eps)); |
| 83 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(EPSGetTwoSided(eps,&twosided)); |
| 84 | |||
| 85 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 86 | Solve the eigensystem | ||
| 87 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 88 | |||
| 89 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(EPSSolve(eps)); |
| 90 | |||
| 91 | /* | ||
| 92 | Optional: Get some information from the solver and display it | ||
| 93 | */ | ||
| 94 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(EPSGetType(eps,&type)); |
| 95 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type)); |
| 96 | |||
| 97 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 98 | Display solution and clean up | ||
| 99 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 100 | |||
| 101 | /* | ||
| 102 | Get number of converged approximate eigenpairs | ||
| 103 | */ | ||
| 104 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(EPSGetConverged(eps,&nconv)); |
| 105 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %" PetscInt_FMT "\n\n",nconv)); |
| 106 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscMalloc2(nconv,&kr,nconv,&ki)); |
| 107 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(VecDuplicateVecs(t,3,&z)); |
| 108 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(VecDuplicateVecs(t,nconv,&xr)); |
| 109 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(VecDuplicateVecs(t,nconv,&xi)); |
| 110 |
1/2✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
|
28 | if (twosided) { |
| 111 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(VecDuplicateVecs(t,nconv,&yr)); |
| 112 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(VecDuplicateVecs(t,nconv,&yi)); |
| 113 | } | ||
| 114 | |||
| 115 |
1/2✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
|
28 | if (nconv>0) { |
| 116 | /* | ||
| 117 | Display eigenvalues and relative errors | ||
| 118 | */ | ||
| 119 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscPrintf(PETSC_COMM_WORLD, |
| 120 | " k ||Ax-kBx|| ||y'A-y'Bk||\n" | ||
| 121 | " ---------------- ------------------ ------------------\n")); | ||
| 122 | |||
| 123 |
2/2✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
|
141 | for (i=0;i<nconv;i++) { |
| 124 | /* | ||
| 125 | Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and | ||
| 126 | ki (imaginary part) | ||
| 127 | */ | ||
| 128 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
113 | PetscCall(EPSGetEigenpair(eps,i,&kr[i],&ki[i],xr[i],xi[i])); |
| 129 |
5/8✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
113 | if (twosided) PetscCall(EPSGetLeftEigenvector(eps,i,yr[i],yi[i])); |
| 130 | /* | ||
| 131 | Compute the residual norms associated to each eigenpair | ||
| 132 | */ | ||
| 133 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
113 | PetscCall(ComputeResidualNorm(A,B,PETSC_FALSE,kr[i],ki[i],xr[i],xi[i],z,&nrmr)); |
| 134 |
5/8✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
113 | if (twosided) PetscCall(ComputeResidualNorm(A,B,PETSC_TRUE,kr[i],ki[i],yr[i],yi[i],z,&nrml)); |
| 135 | |||
| 136 | #if defined(PETSC_USE_COMPLEX) | ||
| 137 | 18 | re = PetscRealPart(kr[i]); | |
| 138 | 18 | im = PetscImaginaryPart(kr[i]); | |
| 139 | #else | ||
| 140 | 95 | re = kr[i]; | |
| 141 | 95 | im = ki[i]; | |
| 142 | #endif | ||
| 143 |
6/8✓ Branch 0 taken 8 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
113 | if (im!=0.0) PetscCall(PetscPrintf(PETSC_COMM_WORLD," %8f%+8fi %12g %12g\n",(double)re,(double)im,(double)nrmr,(double)nrml)); |
| 144 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
113 | else PetscCall(PetscPrintf(PETSC_COMM_WORLD," %12f %12g %12g\n",(double)re,(double)nrmr,(double)nrml)); |
| 145 | } | ||
| 146 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n")); |
| 147 | /* | ||
| 148 | Check bi-orthogonality of eigenvectors | ||
| 149 | */ | ||
| 150 |
1/2✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
|
28 | if (twosided) { |
| 151 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(CheckBiorthogonality(xr,yr,ki,nconv,B,&lev)); |
| 152 |
5/8✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
28 | if (lev<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Level of bi-orthogonality of eigenvectors < 100*eps\n\n")); |
| 153 | ✗ | else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Level of bi-orthogonality of eigenvectors: %g\n\n",(double)lev)); | |
| 154 | } | ||
| 155 | } | ||
| 156 | |||
| 157 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(EPSDestroy(&eps)); |
| 158 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(MatDestroy(&A)); |
| 159 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(MatDestroy(&B)); |
| 160 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(VecDestroy(&t)); |
| 161 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscFree2(kr,ki)); |
| 162 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(VecDestroyVecs(3,&z)); |
| 163 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(VecDestroyVecs(nconv,&xr)); |
| 164 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(VecDestroyVecs(nconv,&xi)); |
| 165 |
1/2✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
|
28 | if (twosided) { |
| 166 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(VecDestroyVecs(nconv,&yr)); |
| 167 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(VecDestroyVecs(nconv,&yi)); |
| 168 | } | ||
| 169 |
3/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 1 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
28 | PetscCall(SlepcFinalize()); |
| 170 | return 0; | ||
| 171 | } | ||
| 172 | |||
| 173 | /* | ||
| 174 | ComputeResidualNorm - Computes the norm of the residual vector | ||
| 175 | associated with an eigenpair. | ||
| 176 | |||
| 177 | Input Parameters: | ||
| 178 | trans - whether A' must be used instead of A | ||
| 179 | kr,ki - eigenvalue | ||
| 180 | xr,xi - eigenvector | ||
| 181 | z - three work vectors (the second one not referenced in complex scalars) | ||
| 182 | */ | ||
| 183 | 226 | PetscErrorCode ComputeResidualNorm(Mat A,Mat B,PetscBool trans,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,Vec *z,PetscReal *norm) | |
| 184 | { | ||
| 185 | 226 | Vec u,w=NULL; | |
| 186 | 226 | PetscScalar alpha; | |
| 187 | #if !defined(PETSC_USE_COMPLEX) | ||
| 188 | 190 | Vec v; | |
| 189 | 190 | PetscReal ni,nr; | |
| 190 | #endif | ||
| 191 |
2/2✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
|
226 | PetscErrorCode (*matmult)(Mat,Vec,Vec) = trans? MatMultHermitianTranspose: MatMult; |
| 192 | |||
| 193 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
226 | PetscFunctionBegin; |
| 194 | 226 | u = z[0]; | |
| 195 |
1/2✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
|
226 | if (B) w = z[2]; |
| 196 | |||
| 197 | #if !defined(PETSC_USE_COMPLEX) | ||
| 198 | 190 | v = z[1]; | |
| 199 |
3/4✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 5 times.
|
190 | if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) { |
| 200 | #endif | ||
| 201 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
186 | PetscCall((*matmult)(A,xr,u)); /* u=A*x */ |
| 202 |
1/2✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
|
186 | if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) { |
| 203 |
5/8✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
186 | if (B) PetscCall((*matmult)(B,xr,w)); /* w=B*x */ |
| 204 | else w = xr; | ||
| 205 |
2/2✓ Branch 0 taken 3 times.
✓ Branch 1 taken 3 times.
|
186 | alpha = trans? -PetscConj(kr): -kr; |
| 206 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
186 | PetscCall(VecAXPY(u,alpha,w)); /* u=A*x-k*B*x */ |
| 207 | } | ||
| 208 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
186 | PetscCall(VecNorm(u,NORM_2,norm)); |
| 209 | #if !defined(PETSC_USE_COMPLEX) | ||
| 210 | } else { | ||
| 211 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
40 | PetscCall((*matmult)(A,xr,u)); /* u=A*xr */ |
| 212 |
1/2✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
|
40 | if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) { |
| 213 |
5/8✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
|
40 | if (B) PetscCall((*matmult)(B,xr,v)); /* v=B*xr */ |
| 214 | ✗ | else PetscCall(VecCopy(xr,v)); | |
| 215 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
40 | PetscCall(VecAXPY(u,-kr,v)); /* u=A*xr-kr*B*xr */ |
| 216 |
5/8✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
|
40 | if (B) PetscCall((*matmult)(B,xi,w)); /* w=B*xi */ |
| 217 | else w = xi; | ||
| 218 |
7/8✓ Branch 0 taken 5 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
|
40 | PetscCall(VecAXPY(u,trans?-ki:ki,w)); /* u=A*xr-kr*B*xr+ki*B*xi */ |
| 219 | } | ||
| 220 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
40 | PetscCall(VecNorm(u,NORM_2,&nr)); |
| 221 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
40 | PetscCall((*matmult)(A,xi,u)); /* u=A*xi */ |
| 222 |
1/2✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
|
40 | if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) { |
| 223 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
40 | PetscCall(VecAXPY(u,-kr,w)); /* u=A*xi-kr*B*xi */ |
| 224 |
7/8✓ Branch 0 taken 5 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
|
40 | PetscCall(VecAXPY(u,trans?ki:-ki,v)); /* u=A*xi-kr*B*xi-ki*B*xr */ |
| 225 | } | ||
| 226 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
40 | PetscCall(VecNorm(u,NORM_2,&ni)); |
| 227 | 40 | *norm = SlepcAbsEigenvalue(nr,ni); | |
| 228 | } | ||
| 229 | #endif | ||
| 230 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
50 | PetscFunctionReturn(PETSC_SUCCESS); |
| 231 | } | ||
| 232 | |||
| 233 | /* | ||
| 234 | CheckOrthogonality - Similar to VecCheckOrthogonality but taking into account complex eigenvalues. | ||
| 235 | |||
| 236 | Input Parameters: | ||
| 237 | V - right eigenvectors | ||
| 238 | W - left eigenvectors | ||
| 239 | eigi - imaginary part of eigenvalues (only to check if the corresponding vector is complex) | ||
| 240 | n - number of V,W vectors | ||
| 241 | B - matrix defining the inner product | ||
| 242 | Output Parameter: | ||
| 243 | lev - level of bi-orthogonality | ||
| 244 | */ | ||
| 245 | 28 | PetscErrorCode CheckBiorthogonality(Vec *V,Vec *W,PetscScalar *ki,PetscInt n,Mat B,PetscReal *lev) | |
| 246 | { | ||
| 247 | 28 | PetscInt i,j; | |
| 248 | 28 | PetscScalar val1; | |
| 249 | 28 | PetscBool wcmplx; | |
| 250 | 28 | Vec wr; | |
| 251 | #if !defined(PETSC_USE_COMPLEX) | ||
| 252 | 25 | Vec wi; | |
| 253 | 25 | PetscScalar val2,dotr,doti; | |
| 254 | #endif | ||
| 255 | |||
| 256 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
28 | PetscFunctionBegin; |
| 257 |
2/14✓ Branch 0 taken 6 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
28 | if (n<=0) PetscFunctionReturn(PETSC_SUCCESS); |
| 258 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(VecDuplicate(V[0],&wr)); |
| 259 | #if !defined(PETSC_USE_COMPLEX) | ||
| 260 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
25 | PetscCall(VecDuplicate(V[0],&wi)); |
| 261 | #endif | ||
| 262 | 28 | *lev = 0.0; | |
| 263 |
2/2✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
|
131 | for (i=0;i<n;i++) { |
| 264 | 103 | wcmplx = PETSC_FALSE; | |
| 265 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
103 | PetscCall(MatMult(B,W[i],wr)); |
| 266 | #if !defined(PETSC_USE_COMPLEX) | ||
| 267 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
85 | if (ki[i] != 0.0) { |
| 268 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
10 | PetscCall(MatMult(B,W[i+1],wi)); |
| 269 | wcmplx = PETSC_TRUE; | ||
| 270 | } | ||
| 271 | #endif | ||
| 272 |
2/2✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
|
576 | for (j=0;j<n;j++) { |
| 273 |
6/6✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
✓ Branch 4 taken 5 times.
✓ Branch 5 taken 5 times.
|
473 | if (j==i || (wcmplx && j==i+1)) continue; |
| 274 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
360 | PetscCall(VecDot(wr,V[j],&val1)); |
| 275 | #if !defined(PETSC_USE_COMPLEX) | ||
| 276 | 270 | dotr = val1; | |
| 277 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
270 | if (ki[j] != 0.0) { /* V complex */ |
| 278 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
40 | if (wcmplx) { |
| 279 | ✗ | PetscCall(VecDot(wi,V[j+1],&val2)); | |
| 280 | ✗ | dotr = val1+val2; | |
| 281 | ✗ | PetscCall(VecDot(wi,V[j],&val1)); | |
| 282 | ✗ | PetscCall(VecDot(wr,V[j+1],&val2)); | |
| 283 | ✗ | doti = val1-val2; | |
| 284 | } else { | ||
| 285 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
40 | PetscCall(VecDot(wr,V[j+1],&val1)); |
| 286 | 40 | doti = -val1; | |
| 287 | } | ||
| 288 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
230 | } else if (wcmplx) { |
| 289 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
40 | PetscCall(VecDot(wi,V[j],&doti)); |
| 290 | 190 | } else doti = 0.0; | |
| 291 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
270 | *lev = PetscMax(*lev,SlepcAbsEigenvalue(dotr,doti)); |
| 292 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
270 | if (ki[j] != 0.0) j++; |
| 293 | #else | ||
| 294 |
2/2✓ Branch 0 taken 3 times.
✓ Branch 1 taken 3 times.
|
171 | *lev = PetscMax(*lev,PetscAbsScalar(val1)); |
| 295 | #endif | ||
| 296 | } | ||
| 297 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
103 | if (wcmplx) i++; |
| 298 | } | ||
| 299 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(VecDestroy(&wr)); |
| 300 | #if !defined(PETSC_USE_COMPLEX) | ||
| 301 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
25 | PetscCall(VecDestroy(&wi)); |
| 302 | #endif | ||
| 303 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
6 | PetscFunctionReturn(PETSC_SUCCESS); |
| 304 | } | ||
| 305 | |||
| 306 | /*TEST | ||
| 307 | |||
| 308 | testset: | ||
| 309 | 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 | ||
| 310 | filter: grep -v "method" | sed -e "s/[+-]0\.0*i//g" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g" | ||
| 311 | requires: double !complex !defined(PETSC_USE_64BIT_INDICES) | ||
| 312 | test: | ||
| 313 | suffix: 1 | ||
| 314 | test: | ||
| 315 | suffix: 1_cmplxvals | ||
| 316 | args: -eps_target -250000 | ||
| 317 | test: | ||
| 318 | suffix: 1_rqi | ||
| 319 | args: -eps_type power -eps_power_shift_type rayleigh -eps_nev 2 -eps_target -2000 | ||
| 320 | test: | ||
| 321 | suffix: 1_rqi_singular | ||
| 322 | args: -eps_type power -eps_power_shift_type rayleigh -eps_nev 1 -eps_target -195500 | ||
| 323 | |||
| 324 | test: | ||
| 325 | suffix: 2 | ||
| 326 | args: -f1 ${DATAFILESPATH}/matrices/complex/mhd1280a.petsc -f2 ${DATAFILESPATH}/matrices/complex/mhd1280b.petsc -eps_nev 6 -eps_tol 1e-11 | ||
| 327 | 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" | ||
| 328 | requires: double complex datafilespath !defined(PETSC_USE_64BIT_INDICES) | ||
| 329 | timeoutfactor: 2 | ||
| 330 | |||
| 331 | test: | ||
| 332 | suffix: 3 | ||
| 333 | 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 -250000 | ||
| 334 | filter: grep -v "method" | sed -e "s/[+-]0\.0*i//g" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g" | grep -v bi-orthog | ||
| 335 | requires: double !complex !defined(PETSC_USE_64BIT_INDICES) | ||
| 336 | |||
| 337 | TEST*/ | ||
| 338 |