| 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 | ! Description: Simple example to test the PEP Fortran interface. | ||
| 11 | ! | ||
| 12 | ! ---------------------------------------------------------------------- | ||
| 13 | ! | ||
| 14 | #include <slepc/finclude/slepcpep.h> | ||
| 15 | 6 | program test3f | |
| 16 | 6 | use slepcpep | |
| 17 | implicit none | ||
| 18 | |||
| 19 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 20 | ! Declarations | ||
| 21 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 22 | |||
| 23 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
24 | Mat :: A(3), B |
| 24 | PEP :: pep | ||
| 25 | ST :: st | ||
| 26 | KSP :: ksp | ||
| 27 | DS :: ds | ||
| 28 | PetscReal :: tol, tolabs, alpha, lambda | ||
| 29 | PetscScalar :: tget, val | ||
| 30 | PetscInt :: n, i, its, Istart, Iend, nev, ncv, mpd, nmat, np | ||
| 31 | PEPWhich :: which | ||
| 32 | PEPConvergedReason :: reason | ||
| 33 | PEPType :: tname | ||
| 34 | PEPExtract :: extr | ||
| 35 | PEPBasis :: basis | ||
| 36 | PEPScale :: scal | ||
| 37 | PEPRefine :: refine | ||
| 38 | PEPRefineScheme :: rscheme | ||
| 39 | PEPConv :: conv | ||
| 40 | PEPStop :: stp | ||
| 41 | PEPProblemType :: ptype | ||
| 42 | PetscMPIInt :: rank | ||
| 43 | PetscErrorCode :: ierr | ||
| 44 | PetscViewerAndFormat :: vf | ||
| 45 | |||
| 46 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 47 | ! Beginning of program | ||
| 48 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 49 | |||
| 50 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(SlepcInitialize(PETSC_NULL_CHARACTER, ierr)) |
| 51 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr)) |
| 52 | 6 | n = 20 | |
| 53 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 54 | 6 | write (*, '(/a,i3,a)') 'Diagonal Quadratic Eigenproblem, n =', n, ' (Fortran)' | |
| 55 | end if | ||
| 56 | |||
| 57 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatCreate(PETSC_COMM_WORLD, A(1), ierr)) |
| 58 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetSizes(A(1), PETSC_DECIDE, PETSC_DECIDE, n, n, ierr)) |
| 59 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetFromOptions(A(1), ierr)) |
| 60 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatGetOwnershipRange(A(1), Istart, Iend, ierr)) |
| 61 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
126 | do i = Istart, Iend - 1 |
| 62 | 120 | val = i + 1 | |
| 63 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
126 | PetscCallA(MatSetValue(A(1), i, i, val, INSERT_VALUES, ierr)) |
| 64 | end do | ||
| 65 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyBegin(A(1), MAT_FINAL_ASSEMBLY, ierr)) |
| 66 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyEnd(A(1), MAT_FINAL_ASSEMBLY, ierr)) |
| 67 | |||
| 68 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatCreate(PETSC_COMM_WORLD, A(2), ierr)) |
| 69 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetSizes(A(2), PETSC_DECIDE, PETSC_DECIDE, n, n, ierr)) |
| 70 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetFromOptions(A(2), ierr)) |
| 71 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatGetOwnershipRange(A(2), Istart, Iend, ierr)) |
| 72 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
126 | do i = Istart, Iend - 1 |
| 73 | 120 | val = 1 | |
| 74 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
126 | PetscCallA(MatSetValue(A(2), i, i, val, INSERT_VALUES, ierr)) |
| 75 | end do | ||
| 76 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyBegin(A(2), MAT_FINAL_ASSEMBLY, ierr)) |
| 77 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyEnd(A(2), MAT_FINAL_ASSEMBLY, ierr)) |
| 78 | |||
| 79 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatCreate(PETSC_COMM_WORLD, A(3), ierr)) |
| 80 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetSizes(A(3), PETSC_DECIDE, PETSC_DECIDE, n, n, ierr)) |
| 81 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetFromOptions(A(3), ierr)) |
| 82 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatGetOwnershipRange(A(3), Istart, Iend, ierr)) |
| 83 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
126 | do i = Istart, Iend - 1 |
| 84 | 120 | val = real(n)/real(i + 1) | |
| 85 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
126 | PetscCallA(MatSetValue(A(3), i, i, val, INSERT_VALUES, ierr)) |
| 86 | end do | ||
| 87 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyBegin(A(3), MAT_FINAL_ASSEMBLY, ierr)) |
| 88 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyEnd(A(3), MAT_FINAL_ASSEMBLY, ierr)) |
| 89 | |||
| 90 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 91 | ! Create eigensolver and test interface functions | ||
| 92 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 93 | |||
| 94 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPCreate(PETSC_COMM_WORLD, pep, ierr)) |
| 95 | 6 | nmat = 3 | |
| 96 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPSetOperators(pep, nmat, A, ierr)) |
| 97 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetNumMatrices(pep, nmat, ierr)) |
| 98 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 99 | 6 | write (*, '(a,i2)') ' Polynomial of degree ', nmat - 1 | |
| 100 | end if | ||
| 101 | 6 | i = 0 | |
| 102 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetOperators(pep, i, B, ierr)) |
| 103 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatView(B, PETSC_NULL_VIEWER, ierr)) |
| 104 | |||
| 105 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPSetType(pep, PEPTOAR, ierr)) |
| 106 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetType(pep, tname, ierr)) |
| 107 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 108 | 6 | write (*, '(a,a)') ' Type set to ', tname | |
| 109 | end if | ||
| 110 | |||
| 111 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetProblemType(pep, ptype, ierr)) |
| 112 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 113 | 6 | write (*, '(a,i2)') ' Problem type before changing = ', ptype | |
| 114 | end if | ||
| 115 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPSetProblemType(pep, PEP_HERMITIAN, ierr)) |
| 116 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetProblemType(pep, ptype, ierr)) |
| 117 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 118 | 6 | write (*, '(a,i2)') ' ... changed to ', ptype | |
| 119 | end if | ||
| 120 | |||
| 121 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetExtract(pep, extr, ierr)) |
| 122 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 123 | 6 | write (*, '(a,i2)') ' Extraction before changing = ', extr | |
| 124 | end if | ||
| 125 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPSetExtract(pep, PEP_EXTRACT_STRUCTURED, ierr)) |
| 126 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetExtract(pep, extr, ierr)) |
| 127 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 128 | 6 | write (*, '(a,i2)') ' ... changed to ', extr | |
| 129 | end if | ||
| 130 | |||
| 131 | 6 | alpha = .1 | |
| 132 | 6 | its = 5 | |
| 133 | 6 | lambda = 1. | |
| 134 | 6 | scal = PEP_SCALE_SCALAR | |
| 135 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPSetScale(pep, scal, alpha, PETSC_NULL_VEC, PETSC_NULL_VEC, its, lambda, ierr)) |
| 136 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetScale(pep, scal, alpha, PETSC_NULL_VEC, PETSC_NULL_VEC, its, lambda, ierr)) |
| 137 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 138 | 6 | write (*, '(a,i2,a,f7.4,a,i2)') ' Scaling: ', scal, ', alpha=', alpha, ', its=', its | |
| 139 | end if | ||
| 140 | |||
| 141 | 6 | basis = PEP_BASIS_CHEBYSHEV1 | |
| 142 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPSetBasis(pep, basis, ierr)) |
| 143 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetBasis(pep, basis, ierr)) |
| 144 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 145 | 6 | write (*, '(a,i2)') ' Polynomial basis: ', basis | |
| 146 | end if | ||
| 147 | |||
| 148 | 6 | np = 1 | |
| 149 | 6 | tol = 1e-9 | |
| 150 | 6 | its = 2 | |
| 151 | 6 | refine = PEP_REFINE_SIMPLE | |
| 152 | 6 | rscheme = PEP_REFINE_SCHEME_SCHUR | |
| 153 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPSetRefine(pep, refine, np, tol, its, rscheme, ierr)) |
| 154 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetRefine(pep, refine, np, tol, its, rscheme, ierr)) |
| 155 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 156 | 6 | write (*, '(a,i2,a,f7.4,a,i2,a,i2)') ' Refinement: ', refine, ', tol=', tol, ', its=', its, ', schem=', rscheme | |
| 157 | end if | ||
| 158 | |||
| 159 | 6 | tget = 4.8 | |
| 160 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPSetTarget(pep, tget, ierr)) |
| 161 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetTarget(pep, tget, ierr)) |
| 162 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPSetWhichEigenpairs(pep, PEP_TARGET_MAGNITUDE, ierr)) |
| 163 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetWhichEigenpairs(pep, which, ierr)) |
| 164 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 165 | 6 | write (*, '(a,i2,a,f4.1)') ' Which = ', which, ', target = ', PetscRealPart(tget) | |
| 166 | end if | ||
| 167 | |||
| 168 | 6 | nev = 4 | |
| 169 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPSetDimensions(pep, nev, PETSC_DETERMINE_INTEGER, PETSC_DETERMINE_INTEGER, ierr)) |
| 170 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetDimensions(pep, nev, ncv, mpd, ierr)) |
| 171 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 172 | 6 | write (*, '(a,i2,a,i2,a,i2)') ' Dimensions: nev=', nev, ', ncv=', ncv, ', mpd=', mpd | |
| 173 | end if | ||
| 174 | |||
| 175 | 6 | tol = 2.2e-4 | |
| 176 | 6 | its = 200 | |
| 177 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPSetTolerances(pep, tol, its, ierr)) |
| 178 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetTolerances(pep, tol, its, ierr)) |
| 179 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 180 | 6 | write (*, '(a,f8.5,a,i4)') ' Tolerance =', tol, ', max_its =', its | |
| 181 | end if | ||
| 182 | |||
| 183 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPSetConvergenceTest(pep, PEP_CONV_ABS, ierr)) |
| 184 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetConvergenceTest(pep, conv, ierr)) |
| 185 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPSetStoppingTest(pep, PEP_STOP_BASIC, ierr)) |
| 186 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetStoppingTest(pep, stp, ierr)) |
| 187 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 188 | 6 | write (*, '(a,i2,a,i2)') ' Convergence test =', conv, ', stopping test =', stp | |
| 189 | end if | ||
| 190 | |||
| 191 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, vf, ierr)) |
| 192 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPMonitorSet(pep, PEPMONITORFIRST, vf, PetscViewerAndFormatDestroy, ierr)) |
| 193 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPMonitorConvergedCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, PETSC_NULL_VEC, vf, ierr)) |
| 194 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPMonitorSet(pep, PEPMONITORCONVERGED, vf, PEPMonitorConvergedDestroy, ierr)) |
| 195 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPMonitorCancel(pep, ierr)) |
| 196 | |||
| 197 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetST(pep, st, ierr)) |
| 198 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(STGetKSP(st, ksp, ierr)) |
| 199 | 6 | tol = 1.e-8 | |
| 200 | 6 | tolabs = 1.e-35 | |
| 201 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(KSPSetTolerances(ksp, tol, tolabs, PETSC_CURRENT_REAL, PETSC_CURRENT_INTEGER, ierr)) |
| 202 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(STView(st, PETSC_NULL_VIEWER, ierr)) |
| 203 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetDS(pep, ds, ierr)) |
| 204 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(DSView(ds, PETSC_NULL_VIEWER, ierr)) |
| 205 | |||
| 206 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPSetFromOptions(pep, ierr)) |
| 207 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPSolve(pep, ierr)) |
| 208 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPGetConvergedReason(pep, reason, ierr)) |
| 209 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 210 | 6 | write (*, '(a,i2)') ' Finished - converged reason =', reason | |
| 211 | end if | ||
| 212 | |||
| 213 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 214 | ! Display solution and clean up | ||
| 215 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 216 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPErrorView(pep, PEP_ERROR_RELATIVE, PETSC_NULL_VIEWER, ierr)) |
| 217 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PEPDestroy(pep, ierr)) |
| 218 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatDestroy(A(1), ierr)) |
| 219 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatDestroy(A(2), ierr)) |
| 220 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatDestroy(A(3), ierr)) |
| 221 | |||
| 222 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
|
6 | PetscCallA(SlepcFinalize(ierr)) |
| 223 | 1 | end program test3f | |
| 224 | |||
| 225 | !/*TEST | ||
| 226 | ! | ||
| 227 | ! test: | ||
| 228 | ! suffix: 1 | ||
| 229 | ! args: -pep_tol 1e-6 -pep_ncv 22 | ||
| 230 | ! filter: sed -e "s/[+-]0\.0*i//g" | ||
| 231 | ! | ||
| 232 | !TEST*/ | ||
| 233 |