| 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[] = "Test some NLEIGS interface functions.\n\n" | ||
| 12 | "Based on ex27.c. The command line options are:\n" | ||
| 13 | " -n <n>, where <n> = matrix dimension.\n"; | ||
| 14 | |||
| 15 | /* | ||
| 16 | Solve T(lambda)x=0 using NLEIGS solver | ||
| 17 | with T(lambda) = -D+sqrt(lambda)*I | ||
| 18 | where D is the Laplacian operator in 1 dimension | ||
| 19 | and with the interpolation interval [.01,16] | ||
| 20 | */ | ||
| 21 | |||
| 22 | #include <slepcnep.h> | ||
| 23 | |||
| 24 | /* | ||
| 25 | User-defined routines | ||
| 26 | */ | ||
| 27 | PetscErrorCode ComputeSingularities(NEP,PetscInt*,PetscScalar*,void*); | ||
| 28 | |||
| 29 | 10 | int main(int argc,char **argv) | |
| 30 | { | ||
| 31 | 10 | NEP nep; /* nonlinear eigensolver context */ | |
| 32 | 10 | Mat A[2]; | |
| 33 | 10 | PetscInt n=100,Istart,Iend,i,ns,nsin; | |
| 34 | 10 | PetscBool terse,fb; | |
| 35 | 10 | RG rg; | |
| 36 | 10 | FN f[2]; | |
| 37 | 10 | PetscScalar coeffs,shifts[]={1.06,1.1,1.12,1.15},*rkshifts,val; | |
| 38 | 10 | PetscErrorCode (*fsing)(NEP,PetscInt*,PetscScalar*,void*); | |
| 39 | |||
| 40 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBeginUser; |
| 41 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
| 42 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
| 43 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSquare root eigenproblem, n=%" PetscInt_FMT "\n\n",n)); |
| 44 | |||
| 45 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 46 | Create nonlinear eigensolver and set some options | ||
| 47 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 48 | |||
| 49 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep)); |
| 50 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPSetType(nep,NEPNLEIGS)); |
| 51 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPNLEIGSSetSingularitiesFunction(nep,ComputeSingularities,NULL)); |
| 52 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPGetRG(nep,&rg)); |
| 53 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(RGSetType(rg,RGINTERVAL)); |
| 54 | #if defined(PETSC_USE_COMPLEX) | ||
| 55 |
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.
|
5 | PetscCall(RGIntervalSetEndpoints(rg,0.01,16.0,-0.001,0.001)); |
| 56 | #else | ||
| 57 |
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.
|
5 | PetscCall(RGIntervalSetEndpoints(rg,0.01,16.0,0,0)); |
| 58 | #endif | ||
| 59 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPSetTarget(nep,1.1)); |
| 60 | |||
| 61 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 62 | Define the nonlinear problem in split form | ||
| 63 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 64 | |||
| 65 | /* Create matrices */ | ||
| 66 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatCreate(PETSC_COMM_WORLD,&A[0])); |
| 67 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n)); |
| 68 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatSetFromOptions(A[0])); |
| 69 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend)); |
| 70 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1010 | for (i=Istart;i<Iend;i++) { |
| 71 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1000 | if (i>0) PetscCall(MatSetValue(A[0],i,i-1,1.0,INSERT_VALUES)); |
| 72 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1000 | if (i<n-1) PetscCall(MatSetValue(A[0],i,i+1,1.0,INSERT_VALUES)); |
| 73 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1000 | PetscCall(MatSetValue(A[0],i,i,-2.0,INSERT_VALUES)); |
| 74 | } | ||
| 75 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY)); |
| 76 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY)); |
| 77 | |||
| 78 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&A[1])); |
| 79 | |||
| 80 | /* Define functions */ | ||
| 81 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(FNCreate(PETSC_COMM_WORLD,&f[0])); |
| 82 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(FNSetType(f[0],FNRATIONAL)); |
| 83 | 10 | coeffs = 1.0; | |
| 84 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(FNRationalSetNumerator(f[0],1,&coeffs)); |
| 85 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(FNCreate(PETSC_COMM_WORLD,&f[1])); |
| 86 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(FNSetType(f[1],FNSQRT)); |
| 87 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPSetSplitOperator(nep,2,A,f,SUBSET_NONZERO_PATTERN)); |
| 88 | |||
| 89 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 90 | Set some options | ||
| 91 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 92 | |||
| 93 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPNLEIGSSetFullBasis(nep,PETSC_FALSE)); |
| 94 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPNLEIGSSetRKShifts(nep,4,shifts)); |
| 95 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPSetFromOptions(nep)); |
| 96 | |||
| 97 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPNLEIGSGetFullBasis(nep,&fb)); |
| 98 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
20 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Using full basis = %s\n",fb?"true":"false")); |
| 99 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPNLEIGSGetRKShifts(nep,&ns,&rkshifts)); |
| 100 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
10 | if (ns) { |
| 101 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Using %" PetscInt_FMT " RK shifts =",ns)); |
| 102 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
50 | for (i=0;i<ns;i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD," %g",(double)PetscRealPart(rkshifts[i]))); |
| 103 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n")); |
| 104 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
10 | PetscCall(PetscFree(rkshifts)); |
| 105 | } | ||
| 106 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPNLEIGSGetSingularitiesFunction(nep,&fsing,NULL)); |
| 107 | 10 | nsin = 1; | |
| 108 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall((*fsing)(nep,&nsin,&val,NULL)); |
| 109 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," First returned singularity = %g\n",(double)PetscRealPart(val))); |
| 110 | |||
| 111 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 112 | Solve the eigensystem | ||
| 113 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 114 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPSolve(nep)); |
| 115 | |||
| 116 | /* show detailed info unless -terse option is given by user */ | ||
| 117 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse)); |
| 118 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
10 | if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_BACKWARD,NULL)); |
| 119 | else { | ||
| 120 | ✗ | PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL)); | |
| 121 | ✗ | PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD)); | |
| 122 | ✗ | PetscCall(NEPErrorView(nep,NEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD)); | |
| 123 | ✗ | PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); | |
| 124 | } | ||
| 125 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPDestroy(&nep)); |
| 126 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatDestroy(&A[0])); |
| 127 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatDestroy(&A[1])); |
| 128 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(FNDestroy(&f[0])); |
| 129 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(FNDestroy(&f[1])); |
| 130 |
3/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
10 | PetscCall(SlepcFinalize()); |
| 131 | return 0; | ||
| 132 | } | ||
| 133 | |||
| 134 | /* ------------------------------------------------------------------- */ | ||
| 135 | /* | ||
| 136 | ComputeSingularities - Computes maxnp points (at most) in the complex plane where | ||
| 137 | the function T(.) is not analytic. | ||
| 138 | |||
| 139 | In this case, we discretize the singularity region (-inf,0)~(-1e+6,-1e-5) | ||
| 140 | */ | ||
| 141 | 20 | PetscErrorCode ComputeSingularities(NEP nep,PetscInt *maxnp,PetscScalar *xi,void *pt) | |
| 142 | { | ||
| 143 | 20 | PetscReal h; | |
| 144 | 20 | PetscInt i; | |
| 145 | |||
| 146 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | PetscFunctionBeginUser; |
| 147 | 20 | h = 11.0/(*maxnp-1); | |
| 148 | 20 | xi[0] = -1e-5; xi[*maxnp-1] = -1e+6; | |
| 149 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
100000 | for (i=1;i<*maxnp-1;i++) xi[i] = -PetscPowReal(10,-5+h*i); |
| 150 |
5/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 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
20 | PetscFunctionReturn(PETSC_SUCCESS); |
| 151 | } | ||
| 152 | |||
| 153 | /*TEST | ||
| 154 | |||
| 155 | test: | ||
| 156 | suffix: 1 | ||
| 157 | args: -nep_nev 3 -nep_nleigs_interpolation_degree 20 -terse -nep_view | ||
| 158 | requires: double | ||
| 159 | filter: grep -v tolerance | sed -e "s/[+-]0\.0*i//g" | ||
| 160 | |||
| 161 | TEST*/ | ||
| 162 |