xref: /petsc/src/binding/petsc4py/demo/legacy/petsc-examples/ksp/ex2.py (revision 226f8a8a5081bc6ad7227cd631662400f0d6e2a0)
155a74a43SLisandro Dalcin
255a74a43SLisandro Dalcin'''
355a74a43SLisandro DalcinEx2 from PETSc example files implemented for PETSc4py.
455a74a43SLisandro Dalcinhttps://petsc.org/release/src/ksp/ksp/tutorials/ex2.c.html
555a74a43SLisandro DalcinBy: Miguel Arriaga
655a74a43SLisandro Dalcin
755a74a43SLisandro Dalcin
855a74a43SLisandro DalcinSolves a linear system in parallel with KSP.
955a74a43SLisandro DalcinInput parameters include:
1055a74a43SLisandro Dalcin    -view_exact_sol   : write exact solution vector to stdout
1155a74a43SLisandro Dalcin    -m <mesh_x>       : number of mesh points in x-direction
1255a74a43SLisandro Dalcin    -n <mesh_n>       : number of mesh points in y-direction
1355a74a43SLisandro Dalcin
1455a74a43SLisandro Dalcin
1555a74a43SLisandro DalcinVec            x,b,u;    # approx solution, RHS, exact solution
1655a74a43SLisandro DalcinMat            A;        # linear system matrix
1755a74a43SLisandro DalcinKSP            ksp;      # linear solver context
1855a74a43SLisandro DalcinPetscReal      norm;     # norm of solution error
1955a74a43SLisandro Dalcin'''
2055a74a43SLisandro Dalcinimport sys
2155a74a43SLisandro Dalcinimport petsc4py
2255a74a43SLisandro Dalcinpetsc4py.init(sys.argv)
2355a74a43SLisandro Dalcinfrom petsc4py import PETSc
2455a74a43SLisandro Dalcin
2555a74a43SLisandro Dalcinimport numpy as np
2655a74a43SLisandro Dalcin
2755a74a43SLisandro Dalcincomm = PETSc.COMM_WORLD
2855a74a43SLisandro Dalcinsize = comm.getSize()
2955a74a43SLisandro Dalcinrank = comm.getRank()
3055a74a43SLisandro Dalcin
3155a74a43SLisandro DalcinOptDB = PETSc.Options()
3255a74a43SLisandro Dalcinm  = OptDB.getInt('m', 8)
3355a74a43SLisandro Dalcinn  = OptDB.getInt('n', 7)
3455a74a43SLisandro Dalcin
3555a74a43SLisandro Dalcin''' - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3655a74a43SLisandro Dalcin        Compute the matrix and right-hand-side vector that define
3755a74a43SLisandro Dalcin        the linear system, Ax = b.
3855a74a43SLisandro Dalcin    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - '''
3955a74a43SLisandro Dalcin'''
4055a74a43SLisandro Dalcin    Create parallel matrix, specifying only its global dimensions.
4155a74a43SLisandro Dalcin    When using MatCreate(), the matrix format can be specified at
4255a74a43SLisandro Dalcin    runtime. Also, the parallel partitioning of the matrix is
4355a74a43SLisandro Dalcin    determined by PETSc at runtime.
4455a74a43SLisandro Dalcin
4555a74a43SLisandro Dalcin    Performance tuning note:  For problems of substantial size,
4655a74a43SLisandro Dalcin    preallocation of matrix memory is crucial for attaining good
4755a74a43SLisandro Dalcin    performance. See the matrix chapter of the users manual for details.
4855a74a43SLisandro Dalcin'''
4955a74a43SLisandro Dalcin
5055a74a43SLisandro DalcinA = PETSc.Mat().create(comm=comm)
5155a74a43SLisandro DalcinA.setSizes((m*n,m*n))
5255a74a43SLisandro DalcinA.setFromOptions()
5355a74a43SLisandro DalcinA.setPreallocationNNZ((5,5))
5455a74a43SLisandro Dalcin
5555a74a43SLisandro Dalcin'''
5655a74a43SLisandro Dalcin    Currently, all PETSc parallel matrix formats are partitioned by
5755a74a43SLisandro Dalcin    contiguous chunks of rows across the processors.  Determine which
5855a74a43SLisandro Dalcin    rows of the matrix are locally owned.
5955a74a43SLisandro Dalcin'''
6055a74a43SLisandro DalcinIstart,Iend = A.getOwnershipRange()
6155a74a43SLisandro Dalcin
6255a74a43SLisandro Dalcin'''
6355a74a43SLisandro Dalcin    Set matrix elements for the 2-D, five-point stencil in parallel.
6455a74a43SLisandro Dalcin    - Each processor needs to insert only elements that it owns
6555a74a43SLisandro Dalcin    locally (but any non-local elements will be sent to the
6655a74a43SLisandro Dalcin    appropriate processor during matrix assembly).
6755a74a43SLisandro Dalcin    - Always specify global rows and columns of matrix entries.
6855a74a43SLisandro Dalcin
6955a74a43SLisandro Dalcin    Note: this uses the less common natural ordering that orders first
7055a74a43SLisandro Dalcin    all the unknowns for x = h then for x = 2h etc; Hence you see J = Ii +- n
7155a74a43SLisandro Dalcin    instead of J = I +- m as you might expect. The more standard ordering
7255a74a43SLisandro Dalcin    would first do all variables for y = h, then y = 2h etc.
7355a74a43SLisandro Dalcin'''
7455a74a43SLisandro Dalcin
7555a74a43SLisandro Dalcinfor Ii in range(Istart,Iend):
7655a74a43SLisandro Dalcin    v = -1.0
7755a74a43SLisandro Dalcin    i = int(Ii/n)
7855a74a43SLisandro Dalcin    j = int(Ii - i*n)
7955a74a43SLisandro Dalcin
8055a74a43SLisandro Dalcin    if (i>0):
8155a74a43SLisandro Dalcin        J = Ii - n
8255a74a43SLisandro Dalcin        A.setValues(Ii,J,v,addv=True)
8355a74a43SLisandro Dalcin    if (i<m-1):
8455a74a43SLisandro Dalcin        J = Ii + n
8555a74a43SLisandro Dalcin        A.setValues(Ii,J,v,addv=True)
8655a74a43SLisandro Dalcin    if (j>0):
8755a74a43SLisandro Dalcin        J = Ii - 1
8855a74a43SLisandro Dalcin        A.setValues(Ii,J,v,addv=True)
8955a74a43SLisandro Dalcin    if (j<n-1):
9055a74a43SLisandro Dalcin        J = Ii + 1
9155a74a43SLisandro Dalcin        A.setValues(Ii,J,v,addv=True)
9255a74a43SLisandro Dalcin
9355a74a43SLisandro Dalcin    v = 4.0
9455a74a43SLisandro Dalcin    A.setValues(Ii,Ii,v,addv=True)
9555a74a43SLisandro Dalcin
9655a74a43SLisandro Dalcin'''
9755a74a43SLisandro Dalcin    Assemble matrix, using the 2-step process:
9855a74a43SLisandro Dalcin    MatAssemblyBegin(), MatAssemblyEnd()
9955a74a43SLisandro Dalcin    Computations can be done while messages are in transition
10055a74a43SLisandro Dalcin    by placing code between these two statements.
10155a74a43SLisandro Dalcin'''
10255a74a43SLisandro Dalcin
10355a74a43SLisandro DalcinA.assemblyBegin(A.AssemblyType.FINAL)
10455a74a43SLisandro DalcinA.assemblyEnd(A.AssemblyType.FINAL)
10555a74a43SLisandro Dalcin''' A is symmetric. Set symmetric flag to enable ICC/Cholesky preconditioner '''
10655a74a43SLisandro Dalcin
10755a74a43SLisandro DalcinA.setOption(A.Option.SYMMETRIC,True)
10855a74a43SLisandro Dalcin
10955a74a43SLisandro Dalcin'''
11055a74a43SLisandro Dalcin    Create parallel vectors.
11155a74a43SLisandro Dalcin    - We form 1 vector from scratch and then duplicate as needed.
11255a74a43SLisandro Dalcin    - When using VecCreate(), VecSetSizes and VecSetFromOptions()
11355a74a43SLisandro Dalcin    in this example, we specify only the
11455a74a43SLisandro Dalcin    vector's global dimension; the parallel partitioning is determined
11555a74a43SLisandro Dalcin    at runtime.
11655a74a43SLisandro Dalcin    - When solving a linear system, the vectors and matrices MUST
11755a74a43SLisandro Dalcin    be partitioned accordingly.  PETSc automatically generates
11855a74a43SLisandro Dalcin    appropriately partitioned matrices and vectors when MatCreate()
11955a74a43SLisandro Dalcin    and VecCreate() are used with the same communicator.
12055a74a43SLisandro Dalcin    - The user can alternatively specify the local vector and matrix
12155a74a43SLisandro Dalcin    dimensions when more sophisticated partitioning is needed
12255a74a43SLisandro Dalcin    (replacing the PETSC_DECIDE argument in the VecSetSizes() statement
12355a74a43SLisandro Dalcin    below).
12455a74a43SLisandro Dalcin'''
12555a74a43SLisandro Dalcin
12655a74a43SLisandro Dalcinu = PETSc.Vec().create(comm=comm)
12755a74a43SLisandro Dalcinu.setSizes(m*n)
12855a74a43SLisandro Dalcinu.setFromOptions()
12955a74a43SLisandro Dalcin
13055a74a43SLisandro Dalcinb = u.duplicate()
13155a74a43SLisandro Dalcinx = b.duplicate()
13255a74a43SLisandro Dalcin
13355a74a43SLisandro Dalcin'''
13455a74a43SLisandro Dalcin    Set exact solution; then compute right-hand-side vector.
13555a74a43SLisandro Dalcin    By default we use an exact solution of a vector with all
13655a74a43SLisandro Dalcin    elements of 1.0;
13755a74a43SLisandro Dalcin'''
13855a74a43SLisandro Dalcinu.set(1.0)
13955a74a43SLisandro Dalcinb = A(u)
14055a74a43SLisandro Dalcin
14155a74a43SLisandro Dalcin'''
14255a74a43SLisandro Dalcin    View the exact solution vector if desired
14355a74a43SLisandro Dalcin'''
14455a74a43SLisandro Dalcinflg = OptDB.getBool('view_exact_sol', False)
14555a74a43SLisandro Dalcinif flg:
14655a74a43SLisandro Dalcin    u.view()
14755a74a43SLisandro Dalcin
14855a74a43SLisandro Dalcin''' - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
14955a74a43SLisandro Dalcin            Create the linear solver and set various options
15055a74a43SLisandro Dalcin    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - '''
15155a74a43SLisandro Dalcinksp = PETSc.KSP().create(comm=comm)
15255a74a43SLisandro Dalcin
15355a74a43SLisandro Dalcin'''
15455a74a43SLisandro Dalcin    Set operators. Here the matrix that defines the linear system
155*7addb90fSBarry Smith    also serves as the matrix from which the preconditioner is constructed.
15655a74a43SLisandro Dalcin'''
15755a74a43SLisandro Dalcinksp.setOperators(A,A)
15855a74a43SLisandro Dalcin
15955a74a43SLisandro Dalcin'''
16055a74a43SLisandro Dalcin    Set linear solver defaults for this problem (optional).
16155a74a43SLisandro Dalcin    - By extracting the KSP and PC contexts from the KSP context,
16255a74a43SLisandro Dalcin    we can then directly call any KSP and PC routines to set
16355a74a43SLisandro Dalcin    various options.
16455a74a43SLisandro Dalcin    - The following two statements are optional; all of these
16555a74a43SLisandro Dalcin    parameters could alternatively be specified at runtime via
16655a74a43SLisandro Dalcin    KSPSetFromOptions().  All of these defaults can be
16755a74a43SLisandro Dalcin    overridden at runtime, as indicated below.
16855a74a43SLisandro Dalcin'''
16955a74a43SLisandro Dalcinrtol=1.e-2/((m+1)*(n+1))
17055a74a43SLisandro Dalcinksp.setTolerances(rtol=rtol,atol=1.e-50)
17155a74a43SLisandro Dalcin
17255a74a43SLisandro Dalcin'''
17355a74a43SLisandro DalcinSet runtime options, e.g.,
17455a74a43SLisandro Dalcin    -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
17555a74a43SLisandro DalcinThese options will override those specified above as long as
17655a74a43SLisandro DalcinKSPSetFromOptions() is called _after_ any other customization
17755a74a43SLisandro Dalcinroutines.
17855a74a43SLisandro Dalcin'''
17955a74a43SLisandro Dalcinksp.setFromOptions()
18055a74a43SLisandro Dalcin
18155a74a43SLisandro Dalcin''' - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
18255a74a43SLisandro Dalcin                    Solve the linear system
18355a74a43SLisandro Dalcin    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - '''
18455a74a43SLisandro Dalcin
18555a74a43SLisandro Dalcinksp.solve(b,x)
18655a74a43SLisandro Dalcin
18755a74a43SLisandro Dalcin''' - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
18855a74a43SLisandro Dalcin                    Check the solution and clean up
18955a74a43SLisandro Dalcin    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - '''
19055a74a43SLisandro Dalcinx = x - u # x.axpy(-1.0,u)
19155a74a43SLisandro Dalcinnorm = x.norm(PETSc.NormType.NORM_2)
19255a74a43SLisandro Dalcinits = ksp.getIterationNumber()
19355a74a43SLisandro Dalcin
19455a74a43SLisandro Dalcin'''
19555a74a43SLisandro Dalcin    Print convergence information.  PetscPrintf() produces a single
19655a74a43SLisandro Dalcin    print statement from all processes that share a communicator.
19755a74a43SLisandro Dalcin    An alternative is PetscFPrintf(), which prints to a file.
19855a74a43SLisandro Dalcin'''
19955a74a43SLisandro Dalcinif norm > rtol*10:
20055a74a43SLisandro Dalcin    PETSc.Sys.Print("Norm of error {}, Iterations {}".format(norm,its),comm=comm)
20155a74a43SLisandro Dalcinelse:
20255a74a43SLisandro Dalcin    if size==1:
20355a74a43SLisandro Dalcin        PETSc.Sys.Print("- Serial OK",comm=comm)
20455a74a43SLisandro Dalcin    else:
20555a74a43SLisandro Dalcin        PETSc.Sys.Print("- Parallel OK",comm=comm)
206