| 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> ./ex16f [-help] [-n <n>] [-m <m>] [SLEPc opts] | ||
| 11 | ! | ||
| 12 | ! Description: Simple example that solves a quadratic eigensystem with PEP. | ||
| 13 | ! This is the Fortran90 equivalent to ex16.c | ||
| 14 | ! | ||
| 15 | ! The command line options are: | ||
| 16 | ! -n <n>, where <n> = number of grid subdivisions in x dimension | ||
| 17 | ! -m <m>, where <m> = number of grid subdivisions in y dimension | ||
| 18 | ! | ||
| 19 | ! ---------------------------------------------------------------------- | ||
| 20 | ! | ||
| 21 | #include <slepc/finclude/slepcpep.h> | ||
| 22 | 4 | program ex16f | |
| 23 | 4 | use slepcpep | |
| 24 | implicit none | ||
| 25 | |||
| 26 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 27 | ! Declarations | ||
| 28 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 29 | |||
| 30 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
16 | Mat :: M, C, K, A(3) ! problem matrices |
| 31 | PEP :: pep ! polynomial eigenproblem solver context | ||
| 32 | PEPType :: tname | ||
| 33 | PetscInt :: N, nx, ny, i, j, Istart, Iend, II | ||
| 34 | PetscInt :: nev | ||
| 35 | PetscMPIInt :: rank | ||
| 36 | PetscErrorCode :: ierr | ||
| 37 | PetscBool :: flg, terse | ||
| 38 | PetscScalar :: val | ||
| 39 | PetscScalar, parameter :: two = 2.0, four = 4.0, mone = -1.0 | ||
| 40 | |||
| 41 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 42 | ! Beginning of program | ||
| 43 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 44 | |||
| 45 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(SlepcInitialize(PETSC_NULL_CHARACTER, ierr)) |
| 46 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr)) |
| 47 | 4 | nx = 10 | |
| 48 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-n', nx, flg, ierr)) |
| 49 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-m', ny, flg, ierr)) |
| 50 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
4 | if (.not. flg) then |
| 51 | 4 | ny = nx | |
| 52 | end if | ||
| 53 | 4 | N = nx*ny | |
| 54 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
4 | if (rank == 0) then |
| 55 | 4 | write (*, '(/a,i6,a,i4,a,i4,a)') 'Quadratic Eigenproblem, N=', N, ' (', nx, 'x', ny, ' grid)' | |
| 56 | end if | ||
| 57 | |||
| 58 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 59 | ! Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0 | ||
| 60 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 61 | |||
| 62 | ! ** K is the 2-D Laplacian | ||
| 63 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatCreate(PETSC_COMM_WORLD, K, ierr)) |
| 64 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatSetSizes(K, PETSC_DECIDE, PETSC_DECIDE, N, N, ierr)) |
| 65 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatSetFromOptions(K, ierr)) |
| 66 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatGetOwnershipRange(K, Istart, Iend, ierr)) |
| 67 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
404 | do II = Istart, Iend - 1 |
| 68 | 400 | i = II/nx | |
| 69 | 400 | j = II - i*nx | |
| 70 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
400 | if (i > 0) then |
| 71 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
360 | PetscCallA(MatSetValue(K, II, II - nx, mone, INSERT_VALUES, ierr)) |
| 72 | end if | ||
| 73 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
400 | if (i < ny - 1) then |
| 74 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
360 | PetscCallA(MatSetValue(K, II, II + nx, mone, INSERT_VALUES, ierr)) |
| 75 | end if | ||
| 76 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
400 | if (j > 0) then |
| 77 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
360 | PetscCallA(MatSetValue(K, II, II - 1, mone, INSERT_VALUES, ierr)) |
| 78 | end if | ||
| 79 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
400 | if (j < nx - 1) then |
| 80 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
360 | PetscCallA(MatSetValue(K, II, II + 1, mone, INSERT_VALUES, ierr)) |
| 81 | end if | ||
| 82 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
404 | PetscCallA(MatSetValue(K, II, II, four, INSERT_VALUES, ierr)) |
| 83 | end do | ||
| 84 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatAssemblyBegin(K, MAT_FINAL_ASSEMBLY, ierr)) |
| 85 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatAssemblyEnd(K, MAT_FINAL_ASSEMBLY, ierr)) |
| 86 | |||
| 87 | ! ** C is the 1-D Laplacian on horizontal lines | ||
| 88 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatCreate(PETSC_COMM_WORLD, C, ierr)) |
| 89 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, N, N, ierr)) |
| 90 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatSetFromOptions(C, ierr)) |
| 91 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatGetOwnershipRange(C, Istart, Iend, ierr)) |
| 92 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
404 | do II = Istart, Iend - 1 |
| 93 | 400 | i = II/nx | |
| 94 | 400 | j = II - i*nx | |
| 95 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
400 | if (j > 0) then |
| 96 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
360 | PetscCallA(MatSetValue(C, II, II - 1, mone, INSERT_VALUES, ierr)) |
| 97 | end if | ||
| 98 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
400 | if (j < nx - 1) then |
| 99 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
360 | PetscCallA(MatSetValue(C, II, II + 1, mone, INSERT_VALUES, ierr)) |
| 100 | end if | ||
| 101 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
404 | PetscCallA(MatSetValue(C, II, II, two, INSERT_VALUES, ierr)) |
| 102 | end do | ||
| 103 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY, ierr)) |
| 104 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY, ierr)) |
| 105 | |||
| 106 | ! ** M is a diagonal matrix | ||
| 107 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatCreate(PETSC_COMM_WORLD, M, ierr)) |
| 108 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatSetSizes(M, PETSC_DECIDE, PETSC_DECIDE, N, N, ierr)) |
| 109 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatSetFromOptions(M, ierr)) |
| 110 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatGetOwnershipRange(M, Istart, Iend, ierr)) |
| 111 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
404 | do II = Istart, Iend - 1 |
| 112 | 400 | val = II + 1 | |
| 113 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
404 | PetscCallA(MatSetValue(M, II, II, val, INSERT_VALUES, ierr)) |
| 114 | end do | ||
| 115 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY, ierr)) |
| 116 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY, ierr)) |
| 117 | |||
| 118 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 119 | ! Create the eigensolver and set various options | ||
| 120 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 121 | |||
| 122 | ! ** Create eigensolver context | ||
| 123 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(PEPCreate(PETSC_COMM_WORLD, pep, ierr)) |
| 124 | |||
| 125 | ! ** Set matrices and problem type | ||
| 126 | 4 | A(1) = K | |
| 127 | 4 | A(2) = C | |
| 128 | 4 | A(3) = M | |
| 129 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(PEPSetOperators(pep, 3_PETSC_INT_KIND, A, ierr)) |
| 130 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(PEPSetProblemType(pep, PEP_GENERAL, ierr)) |
| 131 | |||
| 132 | ! ** Set solver parameters at runtime | ||
| 133 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(PEPSetFromOptions(pep, ierr)) |
| 134 | |||
| 135 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 136 | ! Solve the eigensystem | ||
| 137 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 138 | |||
| 139 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(PEPSolve(pep, ierr)) |
| 140 | |||
| 141 | ! ** Optional: Get some information from the solver and display it | ||
| 142 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(PEPGetType(pep, tname, ierr)) |
| 143 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
4 | if (rank == 0) then |
| 144 | 4 | write (*, '(a,a)') ' Solution method: ', tname | |
| 145 | end if | ||
| 146 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(PEPGetDimensions(pep, nev, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, ierr)) |
| 147 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
4 | if (rank == 0) then |
| 148 | 4 | write (*, '(a,i4)') ' Number of requested eigenvalues:', nev | |
| 149 | end if | ||
| 150 | |||
| 151 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 152 | ! Display solution and clean up | ||
| 153 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 154 | |||
| 155 | ! ** show detailed info unless -terse option is given by user | ||
| 156 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-terse', terse, ierr)) |
| 157 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
4 | if (terse) then |
| 158 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(PEPErrorView(pep, PEP_ERROR_BACKWARD, PETSC_NULL_VIEWER, ierr)) |
| 159 | else | ||
| 160 | ✗ | PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL, ierr)) | |
| 161 | ✗ | PetscCallA(PEPConvergedReasonView(pep, PETSC_VIEWER_STDOUT_WORLD, ierr)) | |
| 162 | ✗ | PetscCallA(PEPErrorView(pep, PEP_ERROR_BACKWARD, PETSC_VIEWER_STDOUT_WORLD, ierr)) | |
| 163 | ✗ | PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD, ierr)) | |
| 164 | end if | ||
| 165 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(PEPDestroy(pep, ierr)) |
| 166 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatDestroy(K, ierr)) |
| 167 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatDestroy(C, ierr)) |
| 168 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4 | PetscCallA(MatDestroy(M, ierr)) |
| 169 |
0/2✗ Branch 0 not taken.
✗ Branch 1 not taken.
|
4 | PetscCallA(SlepcFinalize(ierr)) |
| 170 | ✗ | end program ex16f | |
| 171 | |||
| 172 | !/*TEST | ||
| 173 | ! | ||
| 174 | ! test: | ||
| 175 | ! suffix: 1 | ||
| 176 | ! args: -pep_nev 4 -pep_ncv 19 -terse | ||
| 177 | ! requires: !complex | ||
| 178 | ! test: | ||
| 179 | ! suffix: 2 | ||
| 180 | ! args: -pep_nev 4 -pep_smallest_real -terse | ||
| 181 | ! output_file: output/ex16f_1.out | ||
| 182 | ! requires: !complex | ||
| 183 | ! | ||
| 184 | !TEST*/ | ||
| 185 |