xref: /petsc/src/binding/petsc4py/demo/poisson2d/poisson2d.py (revision 552edb6364df478b294b3111f33a8f37ca096b20)
155a74a43SLisandro Dalcin# Poisson in 2D
255a74a43SLisandro Dalcin# =============
35808f684SSatish Balay#
455a74a43SLisandro Dalcin# Solve a constant coefficient Poisson problem on a regular grid. The source
555a74a43SLisandro Dalcin# code for this demo can be `downloaded here <../../_static/poisson2d.py>`__
65808f684SSatish Balay#
755a74a43SLisandro Dalcin# .. math::
85808f684SSatish Balay#
955a74a43SLisandro Dalcin#     - u_{xx} - u_{yy} = 1 \quad\textsf{in}\quad [0,1]^2\\
1055a74a43SLisandro Dalcin#     u = 0 \quad\textsf{on the boundary.}
115808f684SSatish Balay
1255a74a43SLisandro Dalcin# This is a naïve, parallel implementation, using :math:`n` interior grid
1355a74a43SLisandro Dalcin# points per dimension and a lexicographic ordering of the nodes.
145808f684SSatish Balay
1555a74a43SLisandro Dalcin# This code is kept as simple as possible. However, simplicity comes at a
1655a74a43SLisandro Dalcin# price. Here we use a naive decomposition that does not lead to an optimal
1755a74a43SLisandro Dalcin# communication complexity for the matrix-vector product. An optimal complexity
1855a74a43SLisandro Dalcin# decomposition of a structured grid could be achieved using `PETSc.DMDA`.
1955a74a43SLisandro Dalcin#
2055a74a43SLisandro Dalcin# This demo is structured as a script to be executed using:
2155a74a43SLisandro Dalcin#
2255a74a43SLisandro Dalcin# .. code-block:: console
2355a74a43SLisandro Dalcin#
2455a74a43SLisandro Dalcin#   $ python poisson2d.py
2555a74a43SLisandro Dalcin#
2655a74a43SLisandro Dalcin# potentially with additional options passed at the end of the command.
2755a74a43SLisandro Dalcin#
2855a74a43SLisandro Dalcin# At the start of your script, call `petsc4py.init` passing `sys.argv` so that
2955a74a43SLisandro Dalcin# command-line arguments to the script are passed through to PETSc.
3055a74a43SLisandro Dalcin
3155a74a43SLisandro Dalcinimport sys
3255a74a43SLisandro Dalcinimport petsc4py
3355a74a43SLisandro Dalcin
345808f684SSatish Balaypetsc4py.init(sys.argv)
355808f684SSatish Balay
3655a74a43SLisandro Dalcin# The full PETSc4py API is to be found in the `petsc4py.PETSc` module.
3755a74a43SLisandro Dalcin
385808f684SSatish Balayfrom petsc4py import PETSc
395808f684SSatish Balay
4055a74a43SLisandro Dalcin# PETSc is extensively programmable using the `PETSc.Options` database. For
4155a74a43SLisandro Dalcin# more information see `working with PETSc Options <petsc_options>`.
425808f684SSatish Balay
435808f684SSatish BalayOptDB = PETSc.Options()
445808f684SSatish Balay
4555a74a43SLisandro Dalcin# Grid size and spacing using a default value of ``5``. The user can specify a
4655a74a43SLisandro Dalcin# different number of points in each direction by passing the ``-n`` option to
4755a74a43SLisandro Dalcin# the script.
485808f684SSatish Balay
4955a74a43SLisandro Dalcinn = OptDB.getInt('n', 5)
5055a74a43SLisandro Dalcinh = 1.0 / (n + 1)
515808f684SSatish Balay
52*6f336411SStefano Zampini# Matrices are instances of the `PETSc.Mat` class.
53*6f336411SStefano Zampini
54*6f336411SStefano ZampiniA = PETSc.Mat()
55*6f336411SStefano Zampini
56*6f336411SStefano Zampini# Create the underlying PETSc C Mat object.
5755a74a43SLisandro Dalcin# You can omit the ``comm`` argument if your objects live on
5855a74a43SLisandro Dalcin# `PETSc.COMM_WORLD` but it is a dangerous choice to rely on default values
5955a74a43SLisandro Dalcin# for such important arguments.
605808f684SSatish Balay
6155a74a43SLisandro DalcinA.create(comm=PETSc.COMM_WORLD)
6255a74a43SLisandro Dalcin
6355a74a43SLisandro Dalcin# Specify global matrix shape with a tuple.
6455a74a43SLisandro Dalcin
6555a74a43SLisandro DalcinA.setSizes((n * n, n * n))
6655a74a43SLisandro Dalcin
6755a74a43SLisandro Dalcin# The call above implicitly assumes that we leave the parallel decomposition of
6855a74a43SLisandro Dalcin# the matrix rows to PETSc by using `PETSc.DECIDE` for local sizes.
6955a74a43SLisandro Dalcin# It is equivalent to:
7055a74a43SLisandro Dalcin#
7155a74a43SLisandro Dalcin# .. code-block:: python
7255a74a43SLisandro Dalcin#
7355a74a43SLisandro Dalcin#     A.setSizes(((PETSc.DECIDE, n * n), (PETSc.DECIDE, n * n)))
7455a74a43SLisandro Dalcin
75*6f336411SStefano Zampini# Here we use a sparse matrix of AIJ type
76*6f336411SStefano Zampini# Various `matrix formats <petsc4py.PETSc.Mat.Type>` can be selected:
7755a74a43SLisandro Dalcin
7855a74a43SLisandro DalcinA.setType(PETSc.Mat.Type.AIJ)
7955a74a43SLisandro Dalcin
8055a74a43SLisandro Dalcin# Finally we allow the user to set any options they want to on the matrix from
8155a74a43SLisandro Dalcin# the command line:
8255a74a43SLisandro Dalcin
8355a74a43SLisandro DalcinA.setFromOptions()
8455a74a43SLisandro Dalcin
8555a74a43SLisandro Dalcin# Insertion into some matrix types is vastly more efficient if we preallocate
8655a74a43SLisandro Dalcin# space rather than allow this to happen dynamically. Here we hint the number
8755a74a43SLisandro Dalcin# of nonzeros to be expected on each row.
8855a74a43SLisandro Dalcin
8955a74a43SLisandro DalcinA.setPreallocationNNZ(5)
9055a74a43SLisandro Dalcin
9155a74a43SLisandro Dalcin# We can now write out our finite difference matrix assembly using conventional
9255a74a43SLisandro Dalcin# Python syntax. `Mat.getOwnershipRange` is used to retrieve the range of rows
9355a74a43SLisandro Dalcin# local to this processor.
9455a74a43SLisandro Dalcin
9555a74a43SLisandro Dalcin
9655a74a43SLisandro Dalcindef index_to_grid(r):
9755a74a43SLisandro Dalcin    """Convert a row number into a grid point."""
9855a74a43SLisandro Dalcin    return (r // n, r % n)
9955a74a43SLisandro Dalcin
10055a74a43SLisandro Dalcin
10155a74a43SLisandro Dalcinrstart, rend = A.getOwnershipRange()
10255a74a43SLisandro Dalcinfor row in range(rstart, rend):
10355a74a43SLisandro Dalcin    i, j = index_to_grid(row)
10455a74a43SLisandro Dalcin    A[row, row] = 4.0 / h**2
10555a74a43SLisandro Dalcin    if i > 0:
10655a74a43SLisandro Dalcin        column = row - n
10755a74a43SLisandro Dalcin        A[row, column] = -1.0 / h**2
10855a74a43SLisandro Dalcin    if i < n - 1:
10955a74a43SLisandro Dalcin        column = row + n
11055a74a43SLisandro Dalcin        A[row, column] = -1.0 / h**2
11155a74a43SLisandro Dalcin    if j > 0:
11255a74a43SLisandro Dalcin        column = row - 1
11355a74a43SLisandro Dalcin        A[row, column] = -1.0 / h**2
11455a74a43SLisandro Dalcin    if j < n - 1:
11555a74a43SLisandro Dalcin        column = row + 1
11655a74a43SLisandro Dalcin        A[row, column] = -1.0 / h**2
11755a74a43SLisandro Dalcin
118*6f336411SStefano Zampini# At this stage, any exchange of information required in the matrix assembly
11955a74a43SLisandro Dalcin# process has not occurred. We achieve this by calling `Mat.assemblyBegin` and
12055a74a43SLisandro Dalcin# then `Mat.assemblyEnd`.
12155a74a43SLisandro Dalcin
12255a74a43SLisandro DalcinA.assemblyBegin()
12355a74a43SLisandro DalcinA.assemblyEnd()
12455a74a43SLisandro Dalcin
12555a74a43SLisandro Dalcin# We set up an additional option so that the user can print the matrix by
12655a74a43SLisandro Dalcin# passing ``-view_mat`` to the script.
12755a74a43SLisandro Dalcin
12855a74a43SLisandro DalcinA.viewFromOptions('-view_mat')
12955a74a43SLisandro Dalcin
13055a74a43SLisandro Dalcin# PETSc represents all linear solvers as preconditioned Krylov subspace methods
13155a74a43SLisandro Dalcin# of type `PETSc.KSP`. Here we create a KSP object for a conjugate gradient
13255a74a43SLisandro Dalcin# solver preconditioned with an algebraic multigrid method.
13355a74a43SLisandro Dalcin
13455a74a43SLisandro Dalcinksp = PETSc.KSP()
13555a74a43SLisandro Dalcinksp.create(comm=A.getComm())
13655a74a43SLisandro Dalcinksp.setType(PETSc.KSP.Type.CG)
13755a74a43SLisandro Dalcinksp.getPC().setType(PETSc.PC.Type.GAMG)
13855a74a43SLisandro Dalcin
13955a74a43SLisandro Dalcin# We set the matrix in our linear solver and allow the user to program the
14055a74a43SLisandro Dalcin# solver with options.
14155a74a43SLisandro Dalcin
1425808f684SSatish Balayksp.setOperators(A)
1435808f684SSatish Balayksp.setFromOptions()
1445808f684SSatish Balay
14555a74a43SLisandro Dalcin# Since the matrix knows its size and parallel distribution, we can retrieve
14655a74a43SLisandro Dalcin# appropriately-scaled vectors using `Mat.createVecs`. PETSc vectors are
147dd8e379bSPierre Jolivet# objects of type `PETSc.Vec`. Here we set the right-hand side of our system to
14855a74a43SLisandro Dalcin# a vector of ones, and then solve.
14955a74a43SLisandro Dalcin
15055a74a43SLisandro Dalcinx, b = A.createVecs()
15155a74a43SLisandro Dalcinb.set(1.0)
1525808f684SSatish Balayksp.solve(b, x)
1535808f684SSatish Balay
15455a74a43SLisandro Dalcin# Finally, allow the user to print the solution by passing ``-view_sol`` to the
15555a74a43SLisandro Dalcin# script.
1565808f684SSatish Balay
15755a74a43SLisandro Dalcinx.viewFromOptions('-view_sol')
15855a74a43SLisandro Dalcin
15955a74a43SLisandro Dalcin# Things to try
16055a74a43SLisandro Dalcin# -------------
16155a74a43SLisandro Dalcin#
16255a74a43SLisandro Dalcin# - Show the solution with ``-view_sol``.
16355a74a43SLisandro Dalcin# - Show the matrix with ``-view_mat``.
16455a74a43SLisandro Dalcin# - Change the resolution with ``-n``.
16555a74a43SLisandro Dalcin# - Use a direct solver by passing ``-ksp_type preonly -pc_type lu``.
16655a74a43SLisandro Dalcin# - Run in parallel on two processors using:
16755a74a43SLisandro Dalcin#
16855a74a43SLisandro Dalcin#   .. code-block:: console
16955a74a43SLisandro Dalcin#
17055a74a43SLisandro Dalcin#       mpiexec -n 2 python poisson2d.py
171