| 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 NEP Fortran interface. | ||
| 11 | ! | ||
| 12 | ! ---------------------------------------------------------------------- | ||
| 13 | ! | ||
| 14 | #include <slepc/finclude/slepcnep.h> | ||
| 15 | 6 | program test2f | |
| 16 | 6 | use slepcnep | |
| 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 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
24 | FN :: f(3), g |
| 25 | NEP :: nep | ||
| 26 | DS :: ds | ||
| 27 | RG :: rg | ||
| 28 | PetscReal :: tol | ||
| 29 | PetscScalar :: coeffs(2), tget, val | ||
| 30 | PetscInt :: n, i, its, Istart, Iend | ||
| 31 | PetscInt :: nev, ncv, mpd, nterm | ||
| 32 | PetscInt :: nc, np | ||
| 33 | NEPWhich :: which | ||
| 34 | NEPConvergedReason :: reason | ||
| 35 | NEPType :: tname | ||
| 36 | NEPRefine :: refine | ||
| 37 | NEPRefineScheme :: rscheme | ||
| 38 | NEPConv :: conv | ||
| 39 | NEPStop :: stp | ||
| 40 | NEPProblemType :: ptype | ||
| 41 | MatStructure :: mstr | ||
| 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 Nonlinear Eigenproblem, n =', n, ' (Fortran)' | |
| 55 | end if | ||
| 56 | |||
| 57 | ! Matrices | ||
| 58 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatCreate(PETSC_COMM_WORLD, A(1), ierr)) |
| 59 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetSizes(A(1), PETSC_DECIDE, PETSC_DECIDE, n, n, ierr)) |
| 60 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetFromOptions(A(1), ierr)) |
| 61 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatGetOwnershipRange(A(1), Istart, Iend, ierr)) |
| 62 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
126 | do i = Istart, Iend - 1 |
| 63 | 120 | val = i + 1 | |
| 64 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
126 | PetscCallA(MatSetValue(A(1), i, i, val, INSERT_VALUES, ierr)) |
| 65 | end do | ||
| 66 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyBegin(A(1), MAT_FINAL_ASSEMBLY, ierr)) |
| 67 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyEnd(A(1), MAT_FINAL_ASSEMBLY, ierr)) |
| 68 | |||
| 69 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatCreate(PETSC_COMM_WORLD, A(2), ierr)) |
| 70 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetSizes(A(2), PETSC_DECIDE, PETSC_DECIDE, n, n, ierr)) |
| 71 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetFromOptions(A(2), ierr)) |
| 72 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatGetOwnershipRange(A(2), Istart, Iend, ierr)) |
| 73 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
126 | do i = Istart, Iend - 1 |
| 74 | 120 | val = 1 | |
| 75 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
126 | PetscCallA(MatSetValue(A(2), i, i, val, INSERT_VALUES, ierr)) |
| 76 | end do | ||
| 77 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyBegin(A(2), MAT_FINAL_ASSEMBLY, ierr)) |
| 78 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyEnd(A(2), MAT_FINAL_ASSEMBLY, ierr)) |
| 79 | |||
| 80 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatCreate(PETSC_COMM_WORLD, A(3), ierr)) |
| 81 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetSizes(A(3), PETSC_DECIDE, PETSC_DECIDE, n, n, ierr)) |
| 82 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetFromOptions(A(3), ierr)) |
| 83 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatGetOwnershipRange(A(3), Istart, Iend, ierr)) |
| 84 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
126 | do i = Istart, Iend - 1 |
| 85 | 120 | val = real(n)/real(i + 1) | |
| 86 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
126 | PetscCallA(MatSetValue(A(3), i, i, val, INSERT_VALUES, ierr)) |
| 87 | end do | ||
| 88 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyBegin(A(3), MAT_FINAL_ASSEMBLY, ierr)) |
| 89 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyEnd(A(3), MAT_FINAL_ASSEMBLY, ierr)) |
| 90 | |||
| 91 | ! Functions: f0=-lambda, f1=1.0, f2=sqrt(lambda) | ||
| 92 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNCreate(PETSC_COMM_WORLD, f(1), ierr)) |
| 93 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNSetType(f(1), FNRATIONAL, ierr)) |
| 94 | 6 | nc = 2 | |
| 95 | 6 | coeffs(1) = -1.0 | |
| 96 | 6 | coeffs(2) = 0.0 | |
| 97 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNRationalSetNumerator(f(1), nc, coeffs, ierr)) |
| 98 | |||
| 99 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNCreate(PETSC_COMM_WORLD, f(2), ierr)) |
| 100 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNSetType(f(2), FNRATIONAL, ierr)) |
| 101 | 6 | nc = 1 | |
| 102 | 6 | coeffs(1) = 1.0 | |
| 103 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNRationalSetNumerator(f(2), nc, coeffs, ierr)) |
| 104 | |||
| 105 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNCreate(PETSC_COMM_WORLD, f(3), ierr)) |
| 106 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNSetType(f(3), FNSQRT, ierr)) |
| 107 | |||
| 108 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 109 | ! Create eigensolver and test interface functions | ||
| 110 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 111 | |||
| 112 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPCreate(PETSC_COMM_WORLD, nep, ierr)) |
| 113 | 6 | nterm = 3 | |
| 114 | 6 | mstr = SAME_NONZERO_PATTERN | |
| 115 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetSplitOperator(nep, nterm, A, f, mstr, ierr)) |
| 116 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetSplitOperatorInfo(nep, nterm, mstr, ierr)) |
| 117 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 118 | 6 | write (*, '(a,i2,a)') ' Nonlinear function with ', nterm, ' terms' | |
| 119 | end if | ||
| 120 | 6 | i = 0 | |
| 121 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetSplitOperatorTerm(nep, i, B, g, ierr)) |
| 122 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatView(B, PETSC_NULL_VIEWER, ierr)) |
| 123 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNView(g, PETSC_NULL_VIEWER, ierr)) |
| 124 | |||
| 125 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetType(nep, NEPRII, ierr)) |
| 126 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetType(nep, tname, ierr)) |
| 127 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 128 | 6 | write (*, '(a,a)') ' Type set to ', tname | |
| 129 | end if | ||
| 130 | |||
| 131 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetProblemType(nep, ptype, ierr)) |
| 132 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 133 | 6 | write (*, '(a,i2)') ' Problem type before changing = ', ptype | |
| 134 | end if | ||
| 135 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetProblemType(nep, NEP_RATIONAL, ierr)) |
| 136 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetProblemType(nep, ptype, ierr)) |
| 137 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 138 | 6 | write (*, '(a,i2)') ' ... changed to ', ptype | |
| 139 | end if | ||
| 140 | |||
| 141 | 6 | np = 1 | |
| 142 | 6 | tol = 1e-9 | |
| 143 | 6 | its = 2 | |
| 144 | 6 | refine = NEP_REFINE_SIMPLE | |
| 145 | 6 | rscheme = NEP_REFINE_SCHEME_EXPLICIT | |
| 146 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetRefine(nep, refine, np, tol, its, rscheme, ierr)) |
| 147 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetRefine(nep, refine, np, tol, its, rscheme, ierr)) |
| 148 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 149 | 6 | write (*, '(a,i2,a,f12.9,a,i2,a,i2)') ' Refinement: ', refine, ', tol=', tol, ', its=', its, ', scheme=', rscheme | |
| 150 | end if | ||
| 151 | |||
| 152 | 6 | tget = 1.1 | |
| 153 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetTarget(nep, tget, ierr)) |
| 154 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetTarget(nep, tget, ierr)) |
| 155 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetWhichEigenpairs(nep, NEP_TARGET_MAGNITUDE, ierr)) |
| 156 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetWhichEigenpairs(nep, which, ierr)) |
| 157 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 158 | 6 | write (*, '(a,i2,a,f4.1)') ' Which = ', which, ', target = ', PetscRealPart(tget) | |
| 159 | end if | ||
| 160 | |||
| 161 | 6 | nev = 1 | |
| 162 | 6 | ncv = 12 | |
| 163 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetDimensions(nep, nev, ncv, PETSC_DETERMINE_INTEGER, ierr)) |
| 164 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetDimensions(nep, nev, ncv, mpd, ierr)) |
| 165 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 166 | 6 | write (*, '(a,i2,a,i2,a,i2)') ' Dimensions: nev=', nev, ', ncv=', ncv, ', mpd=', mpd | |
| 167 | end if | ||
| 168 | |||
| 169 | 6 | tol = 1.0e-6 | |
| 170 | 6 | its = 200 | |
| 171 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetTolerances(nep, tol, its, ierr)) |
| 172 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetTolerances(nep, tol, its, ierr)) |
| 173 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 174 | 6 | write (*, '(a,f9.6,a,i4)') ' Tolerance =', tol, ', max_its =', its | |
| 175 | end if | ||
| 176 | |||
| 177 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetConvergenceTest(nep, NEP_CONV_ABS, ierr)) |
| 178 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetConvergenceTest(nep, conv, ierr)) |
| 179 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetStoppingTest(nep, NEP_STOP_BASIC, ierr)) |
| 180 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetStoppingTest(nep, stp, ierr)) |
| 181 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 182 | 6 | write (*, '(a,i2,a,i2)') ' Convergence test =', conv, ', stopping test =', stp | |
| 183 | end if | ||
| 184 | |||
| 185 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, vf, ierr)) |
| 186 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPMonitorSet(nep, NEPMONITORFIRST, vf, PetscViewerAndFormatDestroy, ierr)) |
| 187 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPMonitorConvergedCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, PETSC_NULL_VEC, vf, ierr)) |
| 188 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPMonitorSet(nep, NEPMONITORCONVERGED, vf, NEPMonitorConvergedDestroy, ierr)) |
| 189 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPMonitorCancel(nep, ierr)) |
| 190 | |||
| 191 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetDS(nep, ds, ierr)) |
| 192 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(DSView(ds, PETSC_NULL_VIEWER, ierr)) |
| 193 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetFromOptions(nep, ierr)) |
| 194 | |||
| 195 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetRG(nep, rg, ierr)) |
| 196 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(RGView(rg, PETSC_NULL_VIEWER, ierr)) |
| 197 | |||
| 198 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSolve(nep, ierr)) |
| 199 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetConvergedReason(nep, reason, ierr)) |
| 200 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 201 | 6 | write (*, '(a,i2)') ' Finished - converged reason =', reason | |
| 202 | end if | ||
| 203 | |||
| 204 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 205 | ! Display solution and clean up | ||
| 206 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 207 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPErrorView(nep, NEP_ERROR_RELATIVE, PETSC_NULL_VIEWER, ierr)) |
| 208 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPDestroy(nep, ierr)) |
| 209 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatDestroy(A(1), ierr)) |
| 210 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatDestroy(A(2), ierr)) |
| 211 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatDestroy(A(3), ierr)) |
| 212 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNDestroy(f(1), ierr)) |
| 213 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNDestroy(f(2), ierr)) |
| 214 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNDestroy(f(3), ierr)) |
| 215 | |||
| 216 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
|
6 | PetscCallA(SlepcFinalize(ierr)) |
| 217 | 1 | end program test2f | |
| 218 | |||
| 219 | !/*TEST | ||
| 220 | ! | ||
| 221 | ! test: | ||
| 222 | ! suffix: 1 | ||
| 223 | ! requires: !single | ||
| 224 | ! | ||
| 225 | !TEST*/ | ||
| 226 |