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