| 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> ./test7f [-help] [-n <n>] [all SLEPc options] | ||
| 11 | ! | ||
| 12 | ! Description: Simple example that solves an eigensystem with the EPS object. | ||
| 13 | ! Same problem as ex1f but with simplified output. | ||
| 14 | ! | ||
| 15 | ! The command line options are: | ||
| 16 | ! -n <n>, where <n> = number of grid points = matrix size | ||
| 17 | ! | ||
| 18 | ! ---------------------------------------------------------------------- | ||
| 19 | ! | ||
| 20 | #include <slepc/finclude/slepceps.h> | ||
| 21 | 6 | program test7f | |
| 22 | 6 | use slepceps | |
| 23 | implicit none | ||
| 24 | |||
| 25 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 26 | ! Declarations | ||
| 27 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 28 | |||
| 29 | Mat :: A ! operator matrix | ||
| 30 | EPS :: eps ! eigenproblem solver context | ||
| 31 | EPSType :: tname | ||
| 32 | PetscInt :: n, i, Istart, Iend | ||
| 33 | PetscInt :: nev, nini | ||
| 34 | PetscInt :: col(3) | ||
| 35 | PetscInt :: i0, i1, i2, i3 | ||
| 36 | PetscMPIInt :: rank | ||
| 37 | PetscErrorCode :: ierr | ||
| 38 | PetscBool :: flg | ||
| 39 | PetscScalar :: val(3), one | ||
| 40 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | Vec :: v(1) |
| 41 | |||
| 42 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 43 | ! Beginning of program | ||
| 44 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 45 | |||
| 46 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(SlepcInitialize(PETSC_NULL_CHARACTER, ierr)) |
| 47 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr)) |
| 48 | 6 | n = 30 | |
| 49 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-n', n, flg, ierr)) |
| 50 | |||
| 51 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 52 | 6 | write (*, '(/a,i3,a)') '1-D Laplacian Eigenproblem, n =', n, ' (Fortran)' | |
| 53 | end if | ||
| 54 | |||
| 55 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 56 | ! Compute the operator matrix that defines the eigensystem, Ax=kx | ||
| 57 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 58 | |||
| 59 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatCreate(PETSC_COMM_WORLD, A, ierr)) |
| 60 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n, ierr)) |
| 61 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetFromOptions(A, ierr)) |
| 62 | |||
| 63 | 6 | i0 = 0 | |
| 64 | 6 | i1 = 1 | |
| 65 | 6 | i2 = 2 | |
| 66 | 6 | i3 = 3 | |
| 67 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatGetOwnershipRange(A, Istart, Iend, ierr)) |
| 68 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (Istart == 0) then |
| 69 | 6 | i = 0 | |
| 70 | 6 | col(1) = 0 | |
| 71 | 6 | col(2) = 1 | |
| 72 | 6 | val(1) = 2.0 | |
| 73 | 6 | val(2) = -1.0 | |
| 74 |
3/4✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 6 times.
|
12 | PetscCallA(MatSetValues(A, i1, [i], i2, col, val, INSERT_VALUES, ierr)) |
| 75 | 6 | Istart = Istart + 1 | |
| 76 | end if | ||
| 77 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (Iend == n) then |
| 78 | 6 | i = n - 1 | |
| 79 | 6 | col(1) = n - 2 | |
| 80 | 6 | col(2) = n - 1 | |
| 81 | 6 | val(1) = -1.0 | |
| 82 | 6 | val(2) = 2.0 | |
| 83 |
3/4✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 6 times.
|
12 | PetscCallA(MatSetValues(A, i1, [i], i2, col, val, INSERT_VALUES, ierr)) |
| 84 | 6 | Iend = Iend - 1 | |
| 85 | end if | ||
| 86 | 6 | val(1) = -1.0 | |
| 87 | 6 | val(2) = 2.0 | |
| 88 | 6 | val(3) = -1.0 | |
| 89 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
174 | do i = Istart, Iend - 1 |
| 90 | 168 | col(1) = i - 1 | |
| 91 | 168 | col(2) = i | |
| 92 | 168 | col(3) = i + 1 | |
| 93 |
3/4✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
|
342 | PetscCallA(MatSetValues(A, i1, [i], i3, col, val, INSERT_VALUES, ierr)) |
| 94 | end do | ||
| 95 | |||
| 96 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr)) |
| 97 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr)) |
| 98 | |||
| 99 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatCreateVecs(A, v(1), PETSC_NULL_VEC, ierr)) |
| 100 | 6 | one = 1.0 | |
| 101 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | if (Istart == 0) then |
| 102 | ✗ | PetscCallA(VecSetValue(v(1), i0, one, INSERT_VALUES, ierr)) | |
| 103 | end if | ||
| 104 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(VecAssemblyBegin(v(1), ierr)) |
| 105 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(VecAssemblyEnd(v(1), ierr)) |
| 106 | |||
| 107 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 108 | ! Create the eigensolver and display info | ||
| 109 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 110 | |||
| 111 | ! ** Create eigensolver context | ||
| 112 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(EPSCreate(PETSC_COMM_WORLD, eps, ierr)) |
| 113 | |||
| 114 | ! ** Set operators. In this case, it is a standard eigenvalue problem | ||
| 115 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(EPSSetOperators(eps, A, PETSC_NULL_MAT, ierr)) |
| 116 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(EPSSetProblemType(eps, EPS_HEP, ierr)) |
| 117 | |||
| 118 | ! ** Set solver parameters at runtime | ||
| 119 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(EPSSetFromOptions(eps, ierr)) |
| 120 | |||
| 121 | ! ** Set initial vectors | ||
| 122 | 6 | nini = 1 | |
| 123 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(EPSSetInitialSpace(eps, nini, v, ierr)) |
| 124 | |||
| 125 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 126 | ! Solve the eigensystem | ||
| 127 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 128 | |||
| 129 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(EPSSolve(eps, ierr)) |
| 130 | |||
| 131 | ! ** Optional: Get some information from the solver and display it | ||
| 132 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(EPSGetType(eps, tname, ierr)) |
| 133 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 134 | 6 | write (*, '(a,a)') ' Solution method: ', tname | |
| 135 | end if | ||
| 136 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(EPSGetDimensions(eps, nev, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, ierr)) |
| 137 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 138 | 6 | write (*, '(a,i2)') ' Number of requested eigenvalues:', nev | |
| 139 | end if | ||
| 140 | |||
| 141 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 142 | ! Display solution and clean up | ||
| 143 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 144 | |||
| 145 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(EPSErrorView(eps, EPS_ERROR_RELATIVE, PETSC_NULL_VIEWER, ierr)) |
| 146 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(EPSDestroy(eps, ierr)) |
| 147 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatDestroy(A, ierr)) |
| 148 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(VecDestroy(v(1), ierr)) |
| 149 | |||
| 150 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
|
6 | PetscCallA(SlepcFinalize(ierr)) |
| 151 | 1 | end program test7f | |
| 152 | |||
| 153 | !/*TEST | ||
| 154 | ! | ||
| 155 | ! test: | ||
| 156 | ! suffix: 1 | ||
| 157 | ! args: -eps_nev 4 -eps_ncv 19 | ||
| 158 | ! filter: sed -e "s/83791/83792/" | ||
| 159 | ! | ||
| 160 | !TEST*/ | ||
| 161 |