xref: /petsc/src/binding/petsc4py/demo/hpddm/saddle_point.py (revision 6c5693054f5123506dab0f5da2d352ed973d0e50)
1dc97da86SPierre Jolivet# This script demonstrates solving a symmetric saddle-point linear system using PETSc and HPDDM preconditioner
2dc97da86SPierre Jolivet# It must be run with exactly 4 MPI processes
3dc97da86SPierre Jolivet# Example run:
4*a8cf87e0SJunchao Zhang#   mpiexec -n 4 python3 saddle_point.py -ksp_monitor -ksp_type fgmres -ksp_max_it 10 -ksp_rtol 1e-4 -fieldsplit_ksp_max_it 100 -fieldsplit_pc_hpddm_levels_1_eps_nev 10 -fieldsplit_pc_hpddm_levels_1_st_share_sub_ksp -fieldsplit_pc_hpddm_has_neumann -fieldsplit_pc_hpddm_define_subdomains -fieldsplit_1_pc_hpddm_schur_precondition geneo -fieldsplit_pc_hpddm_coarse_pc_type cholesky -fieldsplit_pc_hpddm_levels_1_sub_pc_type lu -fieldsplit_ksp_type fgmres -fieldsplit_1_pc_hpddm_coarse_correction balanced -fieldsplit_1_pc_hpddm_levels_1_eps_gen_non_hermitian
5dc97da86SPierre Jolivet# For more options, see ${PETSC_DIR}/src/ksp/ksp/tutorials/ex87.c
6dc97da86SPierre Jolivet
7dc97da86SPierre Jolivetimport sys
8dc97da86SPierre Jolivetimport petsc4py
9dc97da86SPierre Jolivet# Initialize PETSc with command-line arguments
10dc97da86SPierre Jolivetpetsc4py.init(sys.argv)
11dc97da86SPierre Jolivetfrom petsc4py import PETSc
12dc97da86SPierre Jolivet
13dc97da86SPierre Jolivet# Function to load matrices and index sets from binary files
14dc97da86SPierre Jolivetdef mat_and_is_load(prefix, identifier, A, aux_IS, aux_Mat, rank, size):
15dc97da86SPierre Jolivet# Load an index set (IS) from binary file
16dc97da86SPierre Jolivet  sizes = PETSc.IS().load(PETSc.Viewer().createBinary(f"{prefix}{identifier}_sizes_{rank}_{size}.dat", "r", comm = PETSc.COMM_SELF))
17dc97da86SPierre Jolivet# Get indices from the loaded IS
18dc97da86SPierre Jolivet  idx   = sizes.getIndices()
19dc97da86SPierre Jolivet# Set the local and global sizes of the matrix
20dc97da86SPierre Jolivet  A.setSizes([[idx[0], idx[2]], [idx[1], idx[3]]])
21dc97da86SPierre Jolivet# Configure matrix using runtime options
22dc97da86SPierre Jolivet  A.setFromOptions()
23dc97da86SPierre Jolivet# Load matrix A from binary file
24dc97da86SPierre Jolivet  A = A.load(PETSc.Viewer().createBinary(f"{prefix}{identifier}.dat", "r", comm = PETSc.COMM_WORLD))
25dc97da86SPierre Jolivet
26dc97da86SPierre Jolivet# Load an index set (IS) from binary file
27dc97da86SPierre Jolivet  aux_IS.load(PETSc.Viewer().createBinary(f"{prefix}{identifier}_is_{rank}_{size}.dat", "r", comm = PETSc.COMM_SELF))
28dc97da86SPierre Jolivet# Load the Neumann matrix of the current process
29dc97da86SPierre Jolivet  aux_Mat.load(PETSc.Viewer().createBinary(f"{prefix}{identifier}_aux_{rank}_{size}.dat", "r", comm = PETSc.COMM_SELF))
30dc97da86SPierre Jolivet
31dc97da86SPierre Jolivet# Get the size of the communicator
32dc97da86SPierre Jolivetsize = PETSc.COMM_WORLD.getSize()
33dc97da86SPierre Jolivet# Get the rank of the current process
34dc97da86SPierre Jolivetrank = PETSc.COMM_WORLD.getRank()
35dc97da86SPierre Jolivetif size != 4:
36dc97da86SPierre Jolivet  if rank == 0:
37dc97da86SPierre Jolivet    print("This example requires 4 processes")
38dc97da86SPierre Jolivet  quit()
39dc97da86SPierre Jolivet
40dc97da86SPierre Jolivet# Problem type (either 'elasticity' or 'stokes')
41dc97da86SPierre Jolivetsystem_str = PETSc.Options().getString("system", "elasticity")
42dc97da86SPierre Jolivetid_sys = 0 if system_str == "elasticity" else 1
43dc97da86SPierre Jolivetempty_A11 = False
44dc97da86SPierre Jolivet# Lower-left (1,1) block is never zero when problem type is 'elasticity'
45dc97da86SPierre Jolivetif id_sys == 1:
46dc97da86SPierre Jolivet  empty_A11 = PETSc.Options().getBool("empty_A11", False)
47dc97da86SPierre Jolivet
48dc97da86SPierre Jolivet# 2-by-2 block structure
49dc97da86SPierre JolivetA = [None, None, None, None]
50dc97da86SPierre Jolivet# Auxiliary data only for the diagonal blocks
51dc97da86SPierre Jolivetaux_Mat = [None, None]
52dc97da86SPierre Jolivetaux_IS = [None, None]
53dc97da86SPierre Jolivet
54dc97da86SPierre Jolivet# Create placeholder objects for the diagonal blocks
55dc97da86SPierre JolivetA[0] = PETSc.Mat().create(comm = PETSc.COMM_WORLD)
56dc97da86SPierre JolivetA[0].setFromOptions()
57dc97da86SPierre Jolivetaux_IS[0] = PETSc.IS().create(comm = PETSc.COMM_SELF)
58dc97da86SPierre Jolivetaux_Mat[0] = PETSc.Mat().create(comm = PETSc.COMM_SELF)
59dc97da86SPierre JolivetA[3] = PETSc.Mat().create(comm = PETSc.COMM_WORLD)
60dc97da86SPierre JolivetA[3].setFromOptions()
61dc97da86SPierre Jolivetaux_IS[1] = PETSc.IS().create(comm = PETSc.COMM_SELF)
62dc97da86SPierre Jolivetaux_Mat[1] = PETSc.Mat().create(comm = PETSc.COMM_SELF)
63dc97da86SPierre Jolivet
64dc97da86SPierre Jolivet# Load directory for input data
65dc97da86SPierre Jolivetload_dir = PETSc.Options().getString("load_dir", "${DATAFILESPATH}/matrices/hpddm/GENEO")
66dc97da86SPierre Jolivet# Specific prefix for each problem
67dc97da86SPierre Jolivetprefix = f"{load_dir}/{ 'B' if id_sys == 1 else 'A' }"
68dc97da86SPierre Jolivet
69dc97da86SPierre Jolivet# Diagonal blocks and auxiliary data for PCHPDDM
70dc97da86SPierre Jolivetmat_and_is_load(prefix, "00", A[0], aux_IS[0], aux_Mat[0], rank, size)
71dc97da86SPierre Jolivetmat_and_is_load(prefix, "11", A[3], aux_IS[1], aux_Mat[1], rank, size)
72dc97da86SPierre Jolivet
73dc97da86SPierre Jolivet# Coherent off-diagonal (0,1) block
74dc97da86SPierre JolivetA[2] = PETSc.Mat().create(comm = PETSc.COMM_WORLD)
75dc97da86SPierre Jolivetn, _ = A[0].getLocalSize()
76dc97da86SPierre JolivetN, _ = A[0].getSize()
77dc97da86SPierre Jolivetm, _ = A[3].getLocalSize()
78dc97da86SPierre JolivetM, _ = A[3].getSize()
79dc97da86SPierre Jolivet# Set matrix sizes based on the sizes of (0,0) and (1,1) blocks
80dc97da86SPierre JolivetA[2].setSizes([[m, M], [n, N]])
81dc97da86SPierre JolivetA[2].setFromOptions()
82dc97da86SPierre JolivetA[2].load(PETSc.Viewer().createBinary(f"{load_dir}/{ 'B' if id_sys == 1 else 'A' }10.dat", "r", comm = PETSc.COMM_WORLD))
83dc97da86SPierre Jolivet# Create a matrix that behaves likes A[1]' without explicitly assembling it
84dc97da86SPierre JolivetA[1] = PETSc.Mat().createTranspose(A[2])
85dc97da86SPierre Jolivet
86dc97da86SPierre Jolivet# Global MatNest
87bd6c8f9dSPierre JolivetS = PETSc.Mat().createNest([[A[0], A[1]], [A[2], A[3] if not empty_A11 else None]])
88dc97da86SPierre Jolivet
89dc97da86SPierre Jolivetksp = PETSc.KSP().create(comm = PETSc.COMM_WORLD)
90dc97da86SPierre Jolivetksp.setOperators(S)
91dc97da86SPierre Jolivetpc = ksp.getPC()
92dc97da86SPierre Jolivet
93dc97da86SPierre Jolivet# Use FIELDSPLIT as the outer preconditioner
94dc97da86SPierre Jolivetpc.setType(PETSc.PC.Type.FIELDSPLIT)
95dc97da86SPierre Jolivet# Use SCHUR since there are only two fields
96dc97da86SPierre Jolivetpc.setFieldSplitType(PETSc.PC.CompositeType.SCHUR)
97dc97da86SPierre Jolivet# Use SELF because HPDDM deals with a Schur complement matrix
98dc97da86SPierre Jolivetpc.setFieldSplitSchurPreType(PETSc.PC.FieldSplitSchurPreType.SELF)
99dc97da86SPierre Jolivet# Apply any command-line options (which may override options from the source file)
100dc97da86SPierre Jolivetpc.setFromOptions()
101dc97da86SPierre Jolivet# Setup the outer preconditioner so that one can query the inner preconditioners
102dc97da86SPierre Jolivetpc.setUp()
103dc97da86SPierre Jolivet
104dc97da86SPierre Jolivet# Retrieve the inner solvers (associated to the diagonal blocks)
105dc97da86SPierre Jolivetksp0, ksp1 = pc.getFieldSplitSubKSP()
106dc97da86SPierre Jolivet
107dc97da86SPierre Jolivet# Upper-left (0,0) block
108dc97da86SPierre Jolivetpc0 = ksp0.getPC()
109dc97da86SPierre Jolivet# Use HPDDM as the preconditioner
110dc97da86SPierre Jolivetpc0.setType(PETSc.PC.Type.HPDDM)
111bd6c8f9dSPierre Jolivet# Set the index set (local-to-global numbering) and auxiliary matrix (first diagonal block)
112dc97da86SPierre Jolivetpc0.setHPDDMAuxiliaryMat(aux_IS[0], aux_Mat[0])
113dc97da86SPierre Jolivetpc0.setFromOptions()
114dc97da86SPierre Jolivet
115dc97da86SPierre Jolivet# Schur complement
116dc97da86SPierre Jolivetpc1 = ksp1.getPC()
117dc97da86SPierre Jolivet# Use HPDDM as the preconditioner
118dc97da86SPierre Jolivetpc1.setType(PETSc.PC.Type.HPDDM)
119bd6c8f9dSPierre Jolivetif not empty_A11:
120bd6c8f9dSPierre Jolivet# Set the index set (local-to-global numbering) and auxiliary matrix (second diagonal block)
121bd6c8f9dSPierre Jolivet# If there is no block (-empty_A11), then these are computed automatically by HPDDM
122dc97da86SPierre Jolivet  pc1.setHPDDMAuxiliaryMat(aux_IS[1], aux_Mat[1])
123dc97da86SPierre Jolivetpc1.setFromOptions()
124dc97da86SPierre Jolivet
125dc97da86SPierre Jolivet# Create RHS (b) and solution (x) vectors, load b from file, and solve the system
126dc97da86SPierre Jolivetb, x = S.createVecs()
127dc97da86SPierre Jolivetb.load(PETSc.Viewer().createBinary(f"{load_dir}/rhs_{ 'B' if id_sys == 1 else 'A' }.dat", "r", comm = PETSc.COMM_WORLD))
128dc97da86SPierre Jolivetksp.setFromOptions()
129dc97da86SPierre Jolivetksp.solve(b, x)
130