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