This file is indexed.

/usr/share/doc/python-slepc4py-docs/usrman/_sources/tutorial.txt is in python-slepc4py-docs 3.7.0-3build1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
Tutorial
========

This tutorial is intended for basic use of slepc4py. For more advanced
use, the reader is referred to SLEPc tutorials as well as to slepc4py
reference documentation.

Commented source of a simple example
------------------------------------

In this section, we include the source code of example ``demo/ex1.py``
available in the slepc4py distribution, with comments inserted inline.

The first thing to do is initialize the libraries. This is normally
not required, as it is done automatically at import time. However, if
you want to gain access to the facilities for accesing command-line
options, the following lines must be executed by the main script prior
to any petsc4py or slepc4py calls::

    import sys, slepc4py
    slepc4py.init(sys.argv)

Next, we have to import the relevant modules. Normally, both PETSc and
SLEPc modules have to be imported in all slepc4py programs. It may be
useful to import NumPy as well::

    from petsc4py import PETSc
    from slepc4py import SLEPc
    import numpy

At this point, we can use any petsc4py and slepc4py operations. For
instance, the following lines allow the user to specify an integer
command-line argument ``n`` with a default value of 30 (see the next
section for example usage of command-line options)::

    opts = PETSc.Options()
    n = opts.getInt('n', 30)

It is necessary to build a matrix to define an eigenproblem (or two in
the case of generalized eigenproblems). The following fragment of code
creates the matrix object and then fills the non-zero elements one by
one. The matrix of this particular example is tridiagonal, with value
2 in the diagonal, and -1 in off-diagonal positions. See petsc4py
documentation for details about matrix objects::

    A = PETSc.Mat().create()
    A.setSizes([n, n])
    A.setFromOptions()
    A.setUp()

    rstart, rend = A.getOwnershipRange()

    # first row
    if rstart == 0:
        A[0, :2] = [2, -1]
        rstart += 1
    # last row
    if rend == n:
        A[n-1, -2:] = [-1, 2]
        rend -= 1
    # other rows
    for i in range(rstart, rend):
        A[i, i-1:i+2] = [-1, 2, -1]

    A.assemble()

The solver object is created in a similar way as other objects in
petsc4py::

    E = SLEPc.EPS(); E.create()

Once the object is created, the eigenvalue problem must be
specified. At least one matrix must be provided. The problem type must
be indicated as well, in this case it is HEP (Hermitian eigenvalue
problem). Apart from these, other settings could be provided here (for
instance, the tolerance for the computation). After all options have
been set, the user should call the ``setFromOptions()`` operation, so
that any options specified at run time in the command line are passed
to the solver object::

    E.setOperators(A)
    E.setProblemType(SLEPc.EPS.ProblemType.HEP)
    E.setFromOptions()

After that, the ``solve()`` method will run the selected eigensolver,
keeping the solution stored internally::

    E.solve()

Once the computation has finished, we are ready to print the results.
First, some informative data can be retrieved from the solver object::

    Print = PETSc.Sys.Print

    Print()
    Print("******************************")
    Print("*** SLEPc Solution Results ***")
    Print("******************************")
    Print()

    its = E.getIterationNumber()
    Print("Number of iterations of the method: %d" % its)

    eps_type = E.getType()
    Print("Solution method: %s" % eps_type)

    nev, ncv, mpd = E.getDimensions()
    Print("Number of requested eigenvalues: %d" % nev)

    tol, maxit = E.getTolerances()
    Print("Stopping condition: tol=%.4g, maxit=%d" % (tol, maxit))

For retrieving the solution, it is necessary to find out how many
eigenpairs have converged to the requested precision::

    nconv = E.getConverged()
    Print("Number of converged eigenpairs %d" % nconv)

For each of the ``nconv`` eigenpairs, we can retrieve the eigenvalue
``k``, and the eigenvector, which is represented by means of two
petsc4py vectors ``vr`` and ``vi`` (the real and imaginary part of the
eigenvector, since for real matrices the eigenvalue and eigenvector
may be complex).  We also compute the corresponding relative errors in
order to make sure that the computed solution is indeed correct::

    if nconv > 0:
        # Create the results vectors
        vr, wr = A.getVecs()
        vi, wi = A.getVecs()
        #
        Print()
        Print("        k          ||Ax-kx||/||kx|| ")
        Print("----------------- ------------------")
        for i in range(nconv):
            k = E.getEigenpair(i, vr, vi)
            error = E.computeError(i)
            if k.imag != 0.0:
                Print(" %9f%+9f j %12g" % (k.real, k.imag, error))
            else:
                Print(" %12f      %12g" % (k.real, error))
        Print()

Example of command-line usage
-----------------------------

Now we illustrate how to specify command-line options in order to
extract the full potential of slepc4py.

A simple execution of the ``demo/ex1.py`` script will result in the
following output::

    $ python demo/ex1.py

    ******************************
    *** SLEPc Solution Results ***
    ******************************

    Number of iterations of the method: 4
    Solution method: krylovschur
    Number of requested eigenvalues: 1
    Stopping condition: tol=1e-07, maxit=100
    Number of converged eigenpairs 4

        k          ||Ax-kx||/||kx||
    ----------------- ------------------
         3.989739        5.76012e-09
         3.959060        1.41957e-08
         3.908279        6.74118e-08
         3.837916        8.34269e-08

For specifying different setting for the solver parameters, we can use
SLEPc command-line options with the ``-eps`` prefix. For instance, to
change the number of requested eigenvalues and the tolerance::

    $ python demo/ex1.py -eps_nev 10 -eps_tol 1e-11

The method used by the solver object can also be set at run time::

    $ python demo/ex1.py -eps_type lanczos

All the above settings can also be change within the source code by
making use of the appropriate slepc4py method. Since options can be
set from within the code and the command-line, it is often useful to
view the particular settings that are currently being used::

    $ python demo/ex1.py -eps_view

    EPS Object: 1 MPI processes
      type: krylovschur
        Krylov-Schur: 50% of basis vectors kept after restart
      problem type: symmetric eigenvalue problem
      selected portion of the spectrum: largest eigenvalues in magnitude
      number of eigenvalues (nev): 1
      number of column vectors (ncv): 16
      maximum dimension of projected problem (mpd): 16
      maximum number of iterations: 100
      tolerance: 1e-08
      convergence test: relative to the eigenvalue
    BV Object: 1 MPI processes
      type: svec
      17 columns of global length 30
      orthogonalization method: classical Gram-Schmidt
      orthogonalization refinement: if needed (eta: 0.7071)
    DS Object: 1 MPI processes
      type: hep
      solving the problem with: Implicit QR method (_steqr)
    ST Object: 1 MPI processes
      type: shift
      shift: 0
      number of matrices: 1

Note that for computing eigenvalues of smallest magnitude we can use
the option ``-eps_smallest_magnitude``, but for interior eigenvalues
things are not so straightforward. One possibility is to try with
harmonic extraction, for instance to get the eigenvalues closest to
0.6::

    $ python demo/ex1.py -eps_harmonic -eps_target 0.6

Depending on the problem, harmonic extraction may fail to converge. In
those cases, it is necessary to specify a spectral transformation
other than the default. In the command-line, this is indicated with
the ``-st_`` prefix. For example, shift-and-invert with a value of the
shift equal to 0.6 would be::

    $ python demo/ex1.py -st_type sinvert -eps_target 0.6