| 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 | ! Program usage: mpiexec -n <np> ./ex10f [-help] [-n <n>] [all SLEPc options] | ||
| 11 | ! | ||
| 12 | ! Description: Illustrates the use of shell spectral transformations. | ||
| 13 | ! The problem to be solved is the same as ex1.c and corresponds to the | ||
| 14 | ! Laplacian operator in 1 dimension | ||
| 15 | ! | ||
| 16 | ! The command line options are: | ||
| 17 | ! nm <n>, where <n> is the number of grid subdivisions = matrix dimension | ||
| 18 | ! ---------------------------------------------------------------------- | ||
| 19 | |||
| 20 | ! Module contains data needed by shell ST | ||
| 21 | ! | ||
| 22 | #include <slepc/finclude/slepceps.h> | ||
| 23 | module ex10fmodule | ||
| 24 | use slepceps | ||
| 25 | implicit none | ||
| 26 | |||
| 27 | KSP myksp | ||
| 28 | |||
| 29 | contains | ||
| 30 | ! ------------------------------------------------------------------- | ||
| 31 | ! STApply_User - This routine demonstrates the use of a user-provided spectral | ||
| 32 | ! transformation. The transformation implemented in this code is just OP=A^-1. | ||
| 33 | ! | ||
| 34 | ! Input Parameters: | ||
| 35 | ! st - spectral transformation context | ||
| 36 | ! x - input vector | ||
| 37 | ! | ||
| 38 | ! Output Parameter: | ||
| 39 | ! y - output vector | ||
| 40 | ! | ||
| 41 | 240 | subroutine STApply_User(st, x, y, ierr) | |
| 42 | use slepceps | ||
| 43 | implicit none | ||
| 44 | |||
| 45 | ST :: st | ||
| 46 | Vec :: x, y | ||
| 47 | PetscErrorCode :: ierr | ||
| 48 | |||
| 49 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
240 | PetscCall(KSPSolve(myksp, x, y, ierr)) |
| 50 | end subroutine | ||
| 51 | |||
| 52 | ! ------------------------------------------------------------------- | ||
| 53 | ! STApplyTranspose_User - This is not required unless using a two-sided eigensolver | ||
| 54 | ! | ||
| 55 | ! Input Parameters: | ||
| 56 | ! st - spectral transformation context | ||
| 57 | ! x - input vector | ||
| 58 | ! | ||
| 59 | ! Output Parameter: | ||
| 60 | ! y - output vector | ||
| 61 | ! | ||
| 62 | 40 | subroutine STApplyTranspose_User(st, x, y, ierr) | |
| 63 | use slepceps | ||
| 64 | implicit none | ||
| 65 | |||
| 66 | ST :: st | ||
| 67 | Vec :: x, y | ||
| 68 | PetscErrorCode :: ierr | ||
| 69 | |||
| 70 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
40 | PetscCall(KSPSolveTranspose(myksp, x, y, ierr)) |
| 71 | end subroutine | ||
| 72 | |||
| 73 | #if defined(PETSC_USE_COMPLEX) | ||
| 74 | ! ------------------------------------------------------------------- | ||
| 75 | ! STApplyHermitianTranspose_User - This is not required unless using a two-sided eigensolver | ||
| 76 | ! in complex scalars | ||
| 77 | ! | ||
| 78 | ! Input Parameters: | ||
| 79 | ! st - spectral transformation context | ||
| 80 | ! x - input vector | ||
| 81 | ! | ||
| 82 | ! Output Parameter: | ||
| 83 | ! y - output vector | ||
| 84 | ! | ||
| 85 | 80 | subroutine STApplyHermitianTranspose_User(st, x, y, ierr) | |
| 86 | use slepceps | ||
| 87 | implicit none | ||
| 88 | |||
| 89 | ST :: st | ||
| 90 | Vec :: x, y, w | ||
| 91 | PetscErrorCode :: ierr | ||
| 92 | |||
| 93 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
|
80 | PetscCall(VecDuplicate(x, w, ierr)) |
| 94 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
|
80 | PetscCall(VecCopy(x, w, ierr)) |
| 95 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
|
80 | PetscCall(VecConjugate(w, ierr)) |
| 96 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
|
80 | PetscCall(KSPSolveTranspose(myksp, w, y, ierr)) |
| 97 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
|
80 | PetscCall(VecConjugate(y, ierr)) |
| 98 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
|
80 | PetscCall(VecDestroy(w, ierr)) |
| 99 | 80 | end subroutine | |
| 100 | #endif | ||
| 101 | |||
| 102 | ! ------------------------------------------------------------------- | ||
| 103 | ! STBackTransform_User - This routine demonstrates the use of a user-provided spectral | ||
| 104 | ! transformation | ||
| 105 | ! | ||
| 106 | ! Input Parameters: | ||
| 107 | ! st - spectral transformation context | ||
| 108 | ! n - number of eigenvalues to transform | ||
| 109 | ! | ||
| 110 | ! Output Parameters: | ||
| 111 | ! eigr - real part of eigenvalues | ||
| 112 | ! eigi - imaginary part of eigenvalues | ||
| 113 | ! | ||
| 114 | 3356 | subroutine STBackTransform_User(st, n, eigr, eigi, ierr) | |
| 115 | use slepceps | ||
| 116 | implicit none | ||
| 117 | |||
| 118 | ST :: st | ||
| 119 | PetscInt :: n, j | ||
| 120 | PetscScalar :: eigr(*), eigi(*) | ||
| 121 | PetscErrorCode :: ierr | ||
| 122 | |||
| 123 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
10128 | do j = 1, n |
| 124 | 10128 | eigr(j) = 1.0/eigr(j) | |
| 125 | end do | ||
| 126 | 3356 | ierr = 0 | |
| 127 | 3356 | end subroutine | |
| 128 | |||
| 129 | end module ex10fmodule | ||
| 130 | |||
| 131 | ! ---------------------------------------------------------------------- | ||
| 132 | |||
| 133 | 24 | program ex10f | |
| 134 | 18 | use slepceps | |
| 135 | use ex10fmodule | ||
| 136 | implicit none | ||
| 137 | |||
| 138 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 139 | ! Declarations | ||
| 140 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 141 | ! | ||
| 142 | Mat :: A ! operator matrix | ||
| 143 | EPS :: eps ! eigenproblem solver context | ||
| 144 | ST :: st | ||
| 145 | EPSType :: tname | ||
| 146 | PetscInt :: n, i, Istart, Iend, one, two, three | ||
| 147 | PetscInt :: nev, row(1), col(3) | ||
| 148 | PetscScalar :: val(3) | ||
| 149 | PetscBool :: flg, isShell, terse | ||
| 150 | PetscMPIInt :: rank | ||
| 151 | PetscErrorCode :: ierr | ||
| 152 | |||
| 153 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 154 | ! Beginning of program | ||
| 155 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 156 | |||
| 157 | 18 | one = 1 | |
| 158 | 18 | two = 2 | |
| 159 | 18 | three = 3 | |
| 160 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(SlepcInitialize(PETSC_NULL_CHARACTER, ierr)) |
| 161 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr)) |
| 162 | 18 | n = 30 | |
| 163 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-n', n, flg, ierr)) |
| 164 | |||
| 165 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
18 | if (rank == 0) then |
| 166 | 18 | write (*, '(/a,i6/)') '1-D Laplacian Eigenproblem (shell-enabled), n=', n | |
| 167 | end if | ||
| 168 | |||
| 169 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 170 | ! Compute the operator matrix that defines the eigensystem, Ax=kx | ||
| 171 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 172 | |||
| 173 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(MatCreate(PETSC_COMM_WORLD, A, ierr)) |
| 174 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n, ierr)) |
| 175 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(MatSetFromOptions(A, ierr)) |
| 176 | |||
| 177 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(MatGetOwnershipRange(A, Istart, Iend, ierr)) |
| 178 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
18 | if (Istart == 0) then |
| 179 | 18 | row(1) = 0 | |
| 180 | 18 | col(1) = 0 | |
| 181 | 18 | col(2) = 1 | |
| 182 | 18 | val(1) = 2.0 | |
| 183 | 18 | val(2) = -1.0 | |
| 184 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(MatSetValues(A, one, row, two, col, val, INSERT_VALUES, ierr)) |
| 185 | 18 | Istart = Istart + 1 | |
| 186 | end if | ||
| 187 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
18 | if (Iend == n) then |
| 188 | 18 | row(1) = n - 1 | |
| 189 | 18 | col(1) = n - 2 | |
| 190 | 18 | col(2) = n - 1 | |
| 191 | 18 | val(1) = -1.0 | |
| 192 | 18 | val(2) = 2.0 | |
| 193 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(MatSetValues(A, one, row, two, col, val, INSERT_VALUES, ierr)) |
| 194 | 18 | Iend = Iend - 1 | |
| 195 | end if | ||
| 196 | 18 | val(1) = -1.0 | |
| 197 | 18 | val(2) = 2.0 | |
| 198 | 18 | val(3) = -1.0 | |
| 199 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
522 | do i = Istart, Iend - 1 |
| 200 | 504 | row(1) = i | |
| 201 | 504 | col(1) = i - 1 | |
| 202 | 504 | col(2) = i | |
| 203 | 504 | col(3) = i + 1 | |
| 204 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
522 | PetscCallA(MatSetValues(A, one, row, three, col, val, INSERT_VALUES, ierr)) |
| 205 | end do | ||
| 206 | |||
| 207 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr)) |
| 208 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr)) |
| 209 | |||
| 210 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 211 | ! Create the eigensolver and set various options | ||
| 212 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 213 | |||
| 214 | ! ** Create eigensolver context | ||
| 215 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(EPSCreate(PETSC_COMM_WORLD, eps, ierr)) |
| 216 | |||
| 217 | ! ** Set operators. In this case, it is a standard eigenvalue problem | ||
| 218 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(EPSSetOperators(eps, A, PETSC_NULL_MAT, ierr)) |
| 219 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(EPSSetProblemType(eps, EPS_NHEP, ierr)) |
| 220 | |||
| 221 | ! ** Set solver parameters at runtime | ||
| 222 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(EPSSetFromOptions(eps, ierr)) |
| 223 | |||
| 224 | ! ** Initialize shell spectral transformation if selected by user | ||
| 225 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(EPSGetST(eps, st, ierr)) |
| 226 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(PetscObjectTypeCompare(st, STSHELL, isShell, ierr)) |
| 227 | |||
| 228 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
18 | if (isShell) then |
| 229 | ! ** Change sorting criterion since this ST example computes values | ||
| 230 | ! ** closest to 0 | ||
| 231 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSSetWhichEigenpairs(eps, EPS_SMALLEST_REAL, ierr)) |
| 232 | |||
| 233 | ! ** In Fortran, instead of a context for the user-defined spectral transform | ||
| 234 | ! ** we use a module containing any application-specific data, initialized here | ||
| 235 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(KSPCreate(PETSC_COMM_WORLD, myksp, ierr)) |
| 236 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(KSPAppendOptionsPrefix(myksp, "st_", ierr)) |
| 237 | |||
| 238 | ! ** (Required) Set the user-defined routine for applying the operator | ||
| 239 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(STShellSetApply(st, STApply_User, ierr)) |
| 240 | |||
| 241 | ! ** (Optional) Set the user-defined routine for applying the transposed operator | ||
| 242 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(STShellSetApplyTranspose(st, STApplyTranspose_User, ierr)) |
| 243 | |||
| 244 | #if defined(PETSC_USE_COMPLEX) | ||
| 245 | ! ** (Optional) Set the user-defined routine for applying the conjugate-transposed operator | ||
| 246 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
|
8 | PetscCallA(STShellSetApplyHermitianTranspose(st, STApplyHermitianTranspose_User, ierr)) |
| 247 | #endif | ||
| 248 | |||
| 249 | ! ** (Optional) Set the user-defined routine for back-transformation | ||
| 250 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(STShellSetBackTransform(st, STBackTransform_User, ierr)) |
| 251 | |||
| 252 | ! ** (Optional) Set a name for the transformation, used for STView() | ||
| 253 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(PetscObjectSetName(st, 'MyTransformation', ierr)) |
| 254 | |||
| 255 | ! ** (Optional) Do any setup required for the new transformation | ||
| 256 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(KSPSetOperators(myksp, A, A, ierr)) |
| 257 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(KSPSetFromOptions(myksp, ierr)) |
| 258 | end if | ||
| 259 | |||
| 260 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 261 | ! Solve the eigensystem | ||
| 262 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 263 | |||
| 264 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(EPSSolve(eps, ierr)) |
| 265 | |||
| 266 | ! ** Optional: Get some information from the solver and display it | ||
| 267 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(EPSGetType(eps, tname, ierr)) |
| 268 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
18 | if (rank == 0) then |
| 269 | 18 | write (*, '(a,a/)') ' Solution method: ', tname | |
| 270 | end if | ||
| 271 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(EPSGetDimensions(eps, nev, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, ierr)) |
| 272 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
18 | if (rank == 0) then |
| 273 | 18 | write (*, '(a,i2)') ' Number of requested eigenvalues:', nev | |
| 274 | end if | ||
| 275 | |||
| 276 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 277 | ! Display solution and clean up | ||
| 278 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 279 | |||
| 280 | ! ** show detailed info unless -terse option is given by user | ||
| 281 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-terse', terse, ierr)) |
| 282 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
18 | if (terse) then |
| 283 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(EPSErrorView(eps, EPS_ERROR_RELATIVE, PETSC_NULL_VIEWER, ierr)) |
| 284 | else | ||
| 285 | ✗ | PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL, ierr)) | |
| 286 | ✗ | PetscCallA(EPSConvergedReasonView(eps, PETSC_VIEWER_STDOUT_WORLD, ierr)) | |
| 287 | ✗ | PetscCallA(EPSErrorView(eps, EPS_ERROR_RELATIVE, PETSC_VIEWER_STDOUT_WORLD, ierr)) | |
| 288 | ✗ | PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD, ierr)) | |
| 289 | end if | ||
| 290 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
18 | if (isShell) then |
| 291 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(KSPDestroy(myksp, ierr)) |
| 292 | end if | ||
| 293 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(EPSDestroy(eps, ierr)) |
| 294 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
18 | PetscCallA(MatDestroy(A, ierr)) |
| 295 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
|
18 | PetscCallA(SlepcFinalize(ierr)) |
| 296 | 3 | end program ex10f | |
| 297 | |||
| 298 | !/*TEST | ||
| 299 | ! | ||
| 300 | ! testset: | ||
| 301 | ! args: -eps_nev 5 -eps_non_hermitian -terse | ||
| 302 | ! output_file: output/ex10_1.out | ||
| 303 | ! requires: !single | ||
| 304 | ! test: | ||
| 305 | ! suffix: 1_sinvert | ||
| 306 | ! args: -st_type sinvert | ||
| 307 | ! test: | ||
| 308 | ! suffix: 1_shell | ||
| 309 | ! args: -st_type shell -eps_two_sided {{0 1}} | ||
| 310 | ! | ||
| 311 | !TEST*/ | ||
| 312 |