| 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: Test Fortran interface of spectrum-slicing Krylov-Schur. | ||
| 11 | ! | ||
| 12 | ! ---------------------------------------------------------------------- | ||
| 13 | ! | ||
| 14 | #include <slepc/finclude/slepceps.h> | ||
| 15 | 20 | program test17f | |
| 16 | 12 | use slepceps | |
| 17 | implicit none | ||
| 18 | |||
| 19 | #define MAXSUB 16 | ||
| 20 | #define MAXSHI 16 | ||
| 21 | |||
| 22 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 23 | ! Declarations | ||
| 24 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 25 | Mat :: A, B, As, Bs, Au | ||
| 26 | EPS :: eps | ||
| 27 | ST :: st | ||
| 28 | KSP :: ksp | ||
| 29 | PC :: pc | ||
| 30 | Vec :: v | ||
| 31 | PetscScalar :: value | ||
| 32 | PetscInt :: n, m, i, j, k, Istart, Iend | ||
| 33 | PetscInt :: nev, ncv, mpd, nval | ||
| 34 | PetscInt :: row, col, nloc, nlocs, mlocs | ||
| 35 | PetscInt :: II, npart, inertias(MAXSHI) | ||
| 36 | PetscBool :: flg, lock | ||
| 37 | PetscMPIInt :: nprc, rank | ||
| 38 | PetscReal :: int0, int1, keep, subint(MAXSUB) | ||
| 39 | PetscReal :: shifts(MAXSHI) | ||
| 40 | PetscScalar :: eval, one, mone, zero | ||
| 41 | PetscErrorCode :: ierr | ||
| 42 | MPI_Comm :: comm | ||
| 43 | |||
| 44 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 45 | ! Beginning of program | ||
| 46 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 47 | |||
| 48 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(SlepcInitialize(PETSC_NULL_CHARACTER, ierr)) |
| 49 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, nprc, ierr)) |
| 50 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr)) |
| 51 | 12 | n = 35 | |
| 52 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-n', n, flg, ierr)) |
| 53 | 12 | m = n*n | |
| 54 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | if (rank == 0) then |
| 55 | 6 | write (*, '(/a,i3,a/)') 'Spectrum-slicing test, n =', n, ' (Fortran)' | |
| 56 | end if | ||
| 57 | |||
| 58 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatCreate(PETSC_COMM_WORLD, A, ierr)) |
| 59 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m, m, ierr)) |
| 60 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatSetFromOptions(A, ierr)) |
| 61 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatCreate(PETSC_COMM_WORLD, B, ierr)) |
| 62 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, m, ierr)) |
| 63 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatSetFromOptions(B, ierr)) |
| 64 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatGetOwnershipRange(A, Istart, Iend, ierr)) |
| 65 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
7362 | do II = Istart, Iend - 1 |
| 66 | 7350 | i = II/n | |
| 67 | 7350 | j = II - i*n | |
| 68 | 7350 | value = -1.0 | |
| 69 | 7350 | row = II | |
| 70 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
7350 | if (i > 0) then |
| 71 | 7140 | col = II - n | |
| 72 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
7140 | PetscCallA(MatSetValue(A, row, col, value, INSERT_VALUES, ierr)) |
| 73 | end if | ||
| 74 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
7350 | if (i < n - 1) then |
| 75 | 7140 | col = II + n | |
| 76 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
7140 | PetscCallA(MatSetValue(A, row, col, value, INSERT_VALUES, ierr)) |
| 77 | end if | ||
| 78 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
7350 | if (j > 0) then |
| 79 | 7140 | col = II - 1 | |
| 80 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
7140 | PetscCallA(MatSetValue(A, row, col, value, INSERT_VALUES, ierr)) |
| 81 | end if | ||
| 82 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
7350 | if (j < n - 1) then |
| 83 | 7140 | col = II + 1 | |
| 84 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
7140 | PetscCallA(MatSetValue(A, row, col, value, INSERT_VALUES, ierr)) |
| 85 | end if | ||
| 86 | 7350 | col = II | |
| 87 | 7350 | value = 4.0 | |
| 88 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
7350 | PetscCallA(MatSetValue(A, row, col, value, INSERT_VALUES, ierr)) |
| 89 | 7350 | value = 2.0 | |
| 90 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
7362 | PetscCallA(MatSetValue(B, row, col, value, INSERT_VALUES, ierr)) |
| 91 | end do | ||
| 92 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | if (Istart == 0) then |
| 93 | 6 | row = 0 | |
| 94 | 6 | col = 0 | |
| 95 | 6 | value = 6.0 | |
| 96 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetValue(B, row, col, value, INSERT_VALUES, ierr)) |
| 97 | 6 | row = 0 | |
| 98 | 6 | col = 1 | |
| 99 | 6 | value = -1.0 | |
| 100 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetValue(B, row, col, value, INSERT_VALUES, ierr)) |
| 101 | 6 | row = 1 | |
| 102 | 6 | col = 0 | |
| 103 | 6 | value = -1.0 | |
| 104 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetValue(B, row, col, value, INSERT_VALUES, ierr)) |
| 105 | 6 | row = 1 | |
| 106 | 6 | col = 1 | |
| 107 | 6 | value = 1.0 | |
| 108 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
6 | PetscCallA(MatSetValue(B, row, col, value, INSERT_VALUES, ierr)) |
| 109 | end if | ||
| 110 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr)) |
| 111 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr)) |
| 112 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY, ierr)) |
| 113 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY, ierr)) |
| 114 | |||
| 115 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 116 | ! Create eigensolver and set various options | ||
| 117 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 118 | |||
| 119 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSCreate(PETSC_COMM_WORLD, eps, ierr)) |
| 120 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSSetOperators(eps, A, B, ierr)) |
| 121 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSSetProblemType(eps, EPS_GHEP, ierr)) |
| 122 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSSetType(eps, EPSKRYLOVSCHUR, ierr)) |
| 123 | |||
| 124 | ! ** Set interval and other settings for spectrum slicing | ||
| 125 | |||
| 126 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSSetWhichEigenpairs(eps, EPS_ALL, ierr)) |
| 127 | 12 | int0 = 1.1 | |
| 128 | 12 | int1 = 1.3 | |
| 129 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSSetInterval(eps, int0, int1, ierr)) |
| 130 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSGetST(eps, st, ierr)) |
| 131 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(STSetType(st, STSINVERT, ierr)) |
| 132 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
12 | if (nprc > 0) then |
| 133 | 12 | npart = nprc | |
| 134 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurSetPartitions(eps, npart, ierr)) |
| 135 | end if | ||
| 136 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurGetKSP(eps, ksp, ierr)) |
| 137 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(KSPGetPC(ksp, pc, ierr)) |
| 138 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(KSPSetType(ksp, KSPPREONLY, ierr)) |
| 139 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(PCSetType(pc, PCCHOLESKY, ierr)) |
| 140 | |||
| 141 | ! ** Test interface functions of Krylov-Schur solver | ||
| 142 | |||
| 143 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurGetRestart(eps, keep, ierr)) |
| 144 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | if (rank == 0) then |
| 145 | 6 | write (*, '(a,f7.4)') ' Restart parameter before changing = ', keep | |
| 146 | end if | ||
| 147 | 12 | keep = 0.4 | |
| 148 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurSetRestart(eps, keep, ierr)) |
| 149 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurGetRestart(eps, keep, ierr)) |
| 150 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | if (rank == 0) then |
| 151 | 6 | write (*, '(a,f7.4)') ' ... changed to ', keep | |
| 152 | end if | ||
| 153 | |||
| 154 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurGetLocking(eps, lock, ierr)) |
| 155 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | if (rank == 0) then |
| 156 | 6 | write (*, '(a,l4)') ' Locking flag before changing = ', lock | |
| 157 | end if | ||
| 158 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurSetLocking(eps, PETSC_FALSE, ierr)) |
| 159 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurGetLocking(eps, lock, ierr)) |
| 160 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | if (rank == 0) then |
| 161 | 6 | write (*, '(a,l4)') ' ... changed to ', lock | |
| 162 | end if | ||
| 163 | |||
| 164 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurGetDimensions(eps, nev, ncv, mpd, ierr)) |
| 165 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | if (rank == 0) then |
| 166 | 6 | write (*, '(a,i2,a,i2,a,i2)') ' Sub-solve dimensions before changing: nev=', nev, ', ncv=', ncv, ', mpd=', mpd | |
| 167 | end if | ||
| 168 | 12 | nev = 30 | |
| 169 | 12 | ncv = 60 | |
| 170 | 12 | mpd = 60 | |
| 171 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurSetDimensions(eps, nev, ncv, mpd, ierr)) |
| 172 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurGetDimensions(eps, nev, ncv, mpd, ierr)) |
| 173 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | if (rank == 0) then |
| 174 | 6 | write (*, '(a,i2,a,i2,a,i2)') ' ... changed to: nev=', nev, ', ncv=', ncv, ', mpd=', mpd | |
| 175 | end if | ||
| 176 | |||
| 177 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
12 | if (nprc > 0) then |
| 178 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurGetPartitions(eps, npart, ierr)) |
| 179 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | if (rank == 0) then |
| 180 | 6 | write (*, '(a,i2,a)') ' Using ', npart, ' partitions' | |
| 181 | end if | ||
| 182 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCheckA(npart <= MAXSUB, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, 'Too many subintervals') |
| 183 | 12 | subint(1) = int0 | |
| 184 | 12 | subint(npart + 1) = int1 | |
| 185 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
24 | do i = 2, npart |
| 186 | 24 | subint(i) = int0 + (i - 1)*(int1 - int0)/npart | |
| 187 | end do | ||
| 188 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurSetSubintervals(eps, subint, ierr)) |
| 189 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurGetSubintervals(eps, subint, ierr)) |
| 190 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | if (rank == 0) then |
| 191 | 6 | write (*, *) 'Using sub-interval separations =' | |
| 192 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | do i = 2, npart |
| 193 | 12 | write (*, '(f7.4)') subint(i) | |
| 194 | end do | ||
| 195 | end if | ||
| 196 | end if | ||
| 197 | |||
| 198 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSSetFromOptions(eps, ierr)) |
| 199 | |||
| 200 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 201 | ! Compute all eigenvalues in interval and display info | ||
| 202 | ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 203 | |||
| 204 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSSetUp(eps, ierr)) |
| 205 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurGetInertias(eps, k, PETSC_NULL_REAL_ARRAY, PETSC_NULL_INTEGER_ARRAY, ierr)) |
| 206 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCheckA(k <= MAXSHI, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, 'Too many shifts') |
| 207 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurGetInertias(eps, k, shifts, inertias, ierr)) |
| 208 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | if (rank == 0) then |
| 209 | 6 | write (*, *) 'Inertias after EPSSetUp:' | |
| 210 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
24 | do i = 1, k |
| 211 | 24 | write (*, '(a,f4.1,a,i3,a)') ' .. ', shifts(i), ' (', inertias(i), ')' | |
| 212 | end do | ||
| 213 | end if | ||
| 214 | |||
| 215 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSSolve(eps, ierr)) |
| 216 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSGetDimensions(eps, nev, ncv, mpd, ierr)) |
| 217 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSGetInterval(eps, int0, int1, ierr)) |
| 218 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | if (rank == 0) then |
| 219 | 6 | write (*, '(a,i2,a,f7.4,a,f7.4,a)') ' Found ', nev, ' eigenvalues in interval [', int0, ',', int1, ']' | |
| 220 | end if | ||
| 221 | |||
| 222 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
12 | if (nprc > 0) then |
| 223 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurGetSubcommInfo(eps, k, nval, v, ierr)) |
| 224 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | if (rank == 0) then |
| 225 | 6 | write (*, '(a,i2,a,i2,a,i2,a)') ' Process ', rank, ' has worked in sub-interval ', k, ', containing ', nval, ' eigenvalues' | |
| 226 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
174 | do i = 0, nval - 1 |
| 227 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
168 | PetscCallA(EPSKrylovSchurGetSubcommPairs(eps, i, eval, v, ierr)) |
| 228 | 174 | write (*, '(f7.4)') PetscRealPart(eval) | |
| 229 | end do | ||
| 230 | end if | ||
| 231 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(VecDestroy(v, ierr)) |
| 232 | |||
| 233 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurGetSubcommMats(eps, As, Bs, ierr)) |
| 234 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatGetLocalSize(A, nloc, PETSC_NULL_INTEGER, ierr)) |
| 235 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatGetLocalSize(As, nlocs, mlocs, ierr)) |
| 236 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
12 | if (rank == 0) then |
| 237 | 6 | write (*, '(a,i2,a,i5,a,i5,a)') ' Process ', rank, ' owns ', nloc, ', rows of the global matrices, and ', nlocs, ' rows in the subcommunicator' | |
| 238 | end if | ||
| 239 | |||
| 240 | ! ** Modify A on subcommunicators | ||
| 241 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(PetscObjectGetComm(As, comm, ierr)) |
| 242 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatCreate(comm, Au, ierr)) |
| 243 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatSetSizes(Au, nlocs, mlocs, m, m, ierr)) |
| 244 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatSetFromOptions(Au, ierr)) |
| 245 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatGetOwnershipRange(Au, Istart, Iend, ierr)) |
| 246 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
|
14712 | do II = Istart, Iend - 1 |
| 247 | 14700 | value = 0.5 | |
| 248 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
14712 | PetscCallA(MatSetValue(Au, II, II, value, INSERT_VALUES, ierr)) |
| 249 | end do | ||
| 250 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatAssemblyBegin(Au, MAT_FINAL_ASSEMBLY, ierr)) |
| 251 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatAssemblyEnd(Au, MAT_FINAL_ASSEMBLY, ierr)) |
| 252 | 12 | one = 1.0 | |
| 253 | 12 | mone = -1.0 | |
| 254 | 12 | zero = 0.0 | |
| 255 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSKrylovSchurUpdateSubcommMats(eps, one, mone, Au, zero, zero, PETSC_NULL_MAT, DIFFERENT_NONZERO_PATTERN, PETSC_TRUE, ierr)) |
| 256 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatDestroy(Au, ierr)) |
| 257 | end if | ||
| 258 | |||
| 259 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(EPSDestroy(eps, ierr)) |
| 260 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatDestroy(A, ierr)) |
| 261 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
|
12 | PetscCallA(MatDestroy(B, ierr)) |
| 262 | |||
| 263 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
|
12 | PetscCallA(SlepcFinalize(ierr)) |
| 264 | 2 | end program test17f | |
| 265 | |||
| 266 | !/*TEST | ||
| 267 | ! | ||
| 268 | ! test: | ||
| 269 | ! suffix: 1 | ||
| 270 | ! nsize: 2 | ||
| 271 | ! | ||
| 272 | !TEST*/ | ||
| 273 |