| 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> ./ex22f [-n <n>] [-tau <tau>] [SLEPc opts] | ||
| 11 | ! | ||
| 12 | ! Description: Delay differential equation. Fortran90 equivalent of ex22.c | ||
| 13 | ! | ||
| 14 | ! The command line options are: | ||
| 15 | ! -n <n>, where <n> = number of grid subdivisions | ||
| 16 | ! -tau <tau>, where <tau> = delay parameter | ||
| 17 | ! | ||
| 18 | ! ---------------------------------------------------------------------- | ||
| 19 | ! Solve parabolic partial differential equation with time delay tau | ||
| 20 | ! | ||
| 21 | ! u_t = u_xx + aa*u(t) + bb*u(t-tau) | ||
| 22 | ! u(0,t) = u(pi,t) = 0 | ||
| 23 | ! | ||
| 24 | ! with aa = 20 and bb(x) = -4.1+x*(1-exp(x-pi)). | ||
| 25 | ! | ||
| 26 | ! Discretization leads to a DDE of dimension n | ||
| 27 | ! | ||
| 28 | ! -u' = A*u(t) + B*u(t-tau) | ||
| 29 | ! | ||
| 30 | ! which results in the nonlinear eigenproblem | ||
| 31 | ! | ||
| 32 | ! (-lambda*I + A + exp(-tau*lambda)*B)*u = 0 | ||
| 33 | ! ---------------------------------------------------------------------- | ||
| 34 | ! | ||
| 35 | #include <slepc/finclude/slepcnep.h> | ||
| 36 | 6 | program ex22f | |
| 37 | 6 | use slepcnep | |
| 38 | implicit none | ||
| 39 | |||
| 40 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 41 | ! Declarations | ||
| 42 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 43 | |||
| 44 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
24 | Mat :: Id, A, B, mats(3) ! problem matrices |
| 45 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
24 | FN :: f1, f2, f3, funs(3) ! functions to define the nonlinear operator |
| 46 | NEP :: nep ! nonlinear eigensolver context | ||
| 47 | NEPType :: tname | ||
| 48 | PetscScalar :: one, bb, coeffs(2), scal | ||
| 49 | PetscReal :: tau, h, aa, xi, tol | ||
| 50 | PetscInt :: n, i, k, nev, Istart, Iend | ||
| 51 | PetscMPIInt :: rank | ||
| 52 | PetscErrorCode :: ierr | ||
| 53 | PetscBool :: flg, terse | ||
| 54 | |||
| 55 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 56 | ! Beginning of program | ||
| 57 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 58 | |||
| 59 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(SlepcInitialize(PETSC_NULL_CHARACTER, ierr)) |
| 60 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr)) |
| 61 | 6 | n = 128 | |
| 62 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-n', n, flg, ierr)) |
| 63 | 6 | tau = 0.001 | |
| 64 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-tau', tau, flg, ierr)) |
| 65 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 66 | 6 | write (*, '(/a,i4,a,f6.3)') 'Delay Eigenproblem, n =', n, ', tau =', tau | |
| 67 | end if | ||
| 68 | |||
| 69 | 6 | one = 1.0 | |
| 70 | 6 | aa = 20.0 | |
| 71 | 6 | h = PETSC_PI/real(n + 1) | |
| 72 | |||
| 73 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 74 | ! Create problem matrices | ||
| 75 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 76 | |||
| 77 | ! ** Id is the identity matrix | ||
| 78 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatCreateConstantDiagonal(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, n, n, one, Id, ierr)) |
| 79 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetOption(Id, MAT_HERMITIAN, PETSC_TRUE, ierr)) |
| 80 | |||
| 81 | ! ** A = 1/h^2*tridiag(1,-2,1) + aa*I | ||
| 82 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatCreate(PETSC_COMM_WORLD, A, ierr)) |
| 83 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n, ierr)) |
| 84 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetFromOptions(A, ierr)) |
| 85 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatGetOwnershipRange(A, Istart, Iend, ierr)) |
| 86 | 6 | coeffs(1) = 1.0/(h*h) | |
| 87 | 6 | coeffs(2) = -2.0/(h*h) + aa | |
| 88 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
774 | do i = Istart, Iend - 1 |
| 89 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
768 | if (i > 0) then |
| 90 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
762 | PetscCallA(MatSetValue(A, i, i - 1, coeffs(1), INSERT_VALUES, ierr)) |
| 91 | end if | ||
| 92 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
768 | if (i < n - 1) then |
| 93 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
762 | PetscCallA(MatSetValue(A, i, i + 1, coeffs(1), INSERT_VALUES, ierr)) |
| 94 | end if | ||
| 95 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
774 | PetscCallA(MatSetValue(A, i, i, coeffs(2), INSERT_VALUES, ierr)) |
| 96 | end do | ||
| 97 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr)) |
| 98 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr)) |
| 99 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetOption(A, MAT_HERMITIAN, PETSC_TRUE, ierr)) |
| 100 | |||
| 101 | ! ** B = diag(bb(xi)) | ||
| 102 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatCreate(PETSC_COMM_WORLD, B, ierr)) |
| 103 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, n, n, ierr)) |
| 104 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetFromOptions(B, ierr)) |
| 105 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatGetOwnershipRange(B, Istart, Iend, ierr)) |
| 106 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
774 | do i = Istart, Iend - 1 |
| 107 | 768 | xi = (i + 1)*h | |
| 108 | 768 | bb = -4.1 + xi*(1.0 - exp(xi - PETSC_PI)) | |
| 109 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
774 | PetscCallA(MatSetValue(B, i, i, bb, INSERT_VALUES, ierr)) |
| 110 | end do | ||
| 111 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY, ierr)) |
| 112 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY, ierr)) |
| 113 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetOption(B, MAT_HERMITIAN, PETSC_TRUE, ierr)) |
| 114 | |||
| 115 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 116 | ! Create problem functions, f1=-lambda, f2=1.0, f3=exp(-tau*lambda) | ||
| 117 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 118 | |||
| 119 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNCreate(PETSC_COMM_WORLD, f1, ierr)) |
| 120 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNSetType(f1, FNRATIONAL, ierr)) |
| 121 | 6 | k = 2 | |
| 122 | 6 | coeffs(1) = -1.0 | |
| 123 | 6 | coeffs(2) = 0.0 | |
| 124 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNRationalSetNumerator(f1, k, coeffs, ierr)) |
| 125 | |||
| 126 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNCreate(PETSC_COMM_WORLD, f2, ierr)) |
| 127 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNSetType(f2, FNRATIONAL, ierr)) |
| 128 | 6 | k = 1 | |
| 129 | 6 | coeffs(1) = 1.0 | |
| 130 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNRationalSetNumerator(f2, k, coeffs, ierr)) |
| 131 | |||
| 132 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNCreate(PETSC_COMM_WORLD, f3, ierr)) |
| 133 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNSetType(f3, FNEXP, ierr)) |
| 134 | 6 | scal = -tau | |
| 135 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNSetScale(f3, scal, one, ierr)) |
| 136 | |||
| 137 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 138 | ! Create the eigensolver and set various options | ||
| 139 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 140 | |||
| 141 | ! ** Create eigensolver context | ||
| 142 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPCreate(PETSC_COMM_WORLD, nep, ierr)) |
| 143 | |||
| 144 | ! ** Set the split operator. Note that A is passed first so that | ||
| 145 | ! ** SUBSET_NONZERO_PATTERN can be used | ||
| 146 | 6 | k = 3 | |
| 147 | 6 | mats(1) = A | |
| 148 | 6 | mats(2) = Id | |
| 149 | 6 | mats(3) = B | |
| 150 | 6 | funs(1) = f2 | |
| 151 | 6 | funs(2) = f1 | |
| 152 | 6 | funs(3) = f3 | |
| 153 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetSplitOperator(nep, k, mats, funs, SUBSET_NONZERO_PATTERN, ierr)) |
| 154 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetProblemType(nep, NEP_GENERAL, ierr)) |
| 155 | |||
| 156 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 157 | ! Customize nonlinear solver; set runtime options | ||
| 158 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 159 | |||
| 160 | 6 | tol = 1e-9 | |
| 161 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetTolerances(nep, tol, PETSC_CURRENT_INTEGER, ierr)) |
| 162 | 6 | k = 1 | |
| 163 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetDimensions(nep, k, PETSC_DETERMINE_INTEGER, PETSC_DETERMINE_INTEGER, ierr)) |
| 164 | 6 | k = 0 | |
| 165 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPRIISetLagPreconditioner(nep, k, ierr)) |
| 166 | |||
| 167 | ! ** Set solver parameters at runtime | ||
| 168 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSetFromOptions(nep, ierr)) |
| 169 | |||
| 170 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 171 | ! Solve the eigensystem | ||
| 172 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 173 | |||
| 174 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPSolve(nep, ierr)) |
| 175 | |||
| 176 | ! ** Optional: Get some information from the solver and display it | ||
| 177 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetType(nep, tname, ierr)) |
| 178 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 179 | 6 | write (*, '(a,a)') ' Solution method: ', tname | |
| 180 | end if | ||
| 181 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPGetDimensions(nep, nev, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, ierr)) |
| 182 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (rank == 0) then |
| 183 | 6 | write (*, '(a,i4)') ' Number of requested eigenvalues:', nev | |
| 184 | end if | ||
| 185 | |||
| 186 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 187 | ! Display solution and clean up | ||
| 188 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 189 | |||
| 190 | ! ** show detailed info unless -terse option is given by user | ||
| 191 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-terse', terse, ierr)) |
| 192 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
6 | if (terse) then |
| 193 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPErrorView(nep, NEP_ERROR_RELATIVE, PETSC_NULL_VIEWER, ierr)) |
| 194 | else | ||
| 195 | ✗ | PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL, ierr)) | |
| 196 | ✗ | PetscCallA(NEPConvergedReasonView(nep, PETSC_VIEWER_STDOUT_WORLD, ierr)) | |
| 197 | ✗ | PetscCallA(NEPErrorView(nep, NEP_ERROR_RELATIVE, PETSC_VIEWER_STDOUT_WORLD, ierr)) | |
| 198 | ✗ | PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD, ierr)) | |
| 199 | end if | ||
| 200 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(NEPDestroy(nep, ierr)) |
| 201 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatDestroy(Id, ierr)) |
| 202 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatDestroy(A, ierr)) |
| 203 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatDestroy(B, ierr)) |
| 204 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNDestroy(f1, ierr)) |
| 205 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNDestroy(f2, ierr)) |
| 206 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(FNDestroy(f3, ierr)) |
| 207 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
|
6 | PetscCallA(SlepcFinalize(ierr)) |
| 208 | 1 | end program ex22f | |
| 209 | |||
| 210 | !/*TEST | ||
| 211 | ! | ||
| 212 | ! test: | ||
| 213 | ! suffix: 1 | ||
| 214 | ! args: -terse | ||
| 215 | ! requires: !single | ||
| 216 | ! filter: sed -e "s/[+-]0\.0*i//g" | ||
| 217 | ! | ||
| 218 | !TEST*/ | ||
| 219 |