155a74a43SLisandro Dalcin 255a74a43SLisandro Dalcin''' 355a74a43SLisandro DalcinEx23 from PETSc example files implemented for PETSc4py. 455a74a43SLisandro Dalcinhttps://petsc.org/release/src/ksp/ksp/tutorials/ex23.c.html 555a74a43SLisandro DalcinBy: Miguel Arriaga 655a74a43SLisandro Dalcin 755a74a43SLisandro DalcinSolves a tridiagonal linear system. 855a74a43SLisandro Dalcin 955a74a43SLisandro DalcinVec x, b, u; approx solution, RHS, exact solution 1055a74a43SLisandro DalcinMat A; linear system matrix 1155a74a43SLisandro DalcinKSP ksp; linear solver context 1255a74a43SLisandro DalcinPC pc; preconditioner context 1355a74a43SLisandro DalcinPetscReal norm; norm of solution error 1455a74a43SLisandro Dalcin 1555a74a43SLisandro Dalcin''' 1655a74a43SLisandro Dalcinimport sys 1755a74a43SLisandro Dalcinimport petsc4py 1855a74a43SLisandro Dalcinpetsc4py.init(sys.argv) 1955a74a43SLisandro Dalcinfrom petsc4py import PETSc 2055a74a43SLisandro Dalcin 2155a74a43SLisandro Dalcinimport numpy as np 2255a74a43SLisandro Dalcin 2355a74a43SLisandro Dalcincomm = PETSc.COMM_WORLD 2455a74a43SLisandro Dalcinsize = comm.getSize() 2555a74a43SLisandro Dalcinrank = comm.getRank() 2655a74a43SLisandro Dalcinn = 12 # Size of problem 2755a74a43SLisandro Dalcintol = 1E-11 # Tolerance of Result. tol=1000.*PETSC_MACHINE_EPSILON; 2855a74a43SLisandro Dalcin 2955a74a43SLisandro Dalcin''' 3055a74a43SLisandro Dalcin- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 3155a74a43SLisandro Dalcin Compute the matrix and right-hand-side vector that define 3255a74a43SLisandro Dalcin the linear system, Ax = b. 3355a74a43SLisandro Dalcin- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 3455a74a43SLisandro Dalcin 3555a74a43SLisandro DalcinCreate vectors. Note that we form 1 vector from scratch and 3655a74a43SLisandro Dalcinthen duplicate as needed. For this simple case let PETSc decide how 3755a74a43SLisandro Dalcinmany elements of the vector are stored on each processor. The second 3855a74a43SLisandro Dalcinargument to VecSetSizes() below causes PETSc to decide. 3955a74a43SLisandro Dalcin''' 4055a74a43SLisandro Dalcin 4155a74a43SLisandro Dalcinx = PETSc.Vec().create(comm=comm) 4255a74a43SLisandro Dalcinx.setSizes(n) 4355a74a43SLisandro Dalcinx.setFromOptions() 4455a74a43SLisandro Dalcin 4555a74a43SLisandro Dalcinb = x.duplicate() 4655a74a43SLisandro Dalcinu = x.duplicate() 4755a74a43SLisandro Dalcin 4855a74a43SLisandro Dalcin''' 4955a74a43SLisandro DalcinIdentify the starting and ending mesh points on each 5055a74a43SLisandro Dalcinprocessor for the interior part of the mesh. We let PETSc decide 5155a74a43SLisandro Dalcinabove. 5255a74a43SLisandro Dalcin''' 5355a74a43SLisandro Dalcin 5455a74a43SLisandro Dalcinrstart,rend = x.getOwnershipRange() 5555a74a43SLisandro Dalcinnlocal = x.getLocalSize() 5655a74a43SLisandro Dalcin 5755a74a43SLisandro Dalcin''' 5855a74a43SLisandro DalcinCreate matrix. When using MatCreate(), the matrix format can 5955a74a43SLisandro Dalcinbe specified at runtime. 6055a74a43SLisandro Dalcin 6155a74a43SLisandro DalcinPerformance tuning note: For problems of substantial size, 6255a74a43SLisandro Dalcinpreallocation of matrix memory is crucial for attaining good 6355a74a43SLisandro Dalcinperformance. See the matrix chapter of the users manual for details. 6455a74a43SLisandro Dalcin 6555a74a43SLisandro DalcinWe pass in nlocal as the "local" size of the matrix to force it 6655a74a43SLisandro Dalcinto have the same parallel layout as the vector created above. 6755a74a43SLisandro Dalcin''' 6855a74a43SLisandro Dalcin 6955a74a43SLisandro DalcinA = PETSc.Mat().create(comm=comm) 7055a74a43SLisandro DalcinA.setSizes(n,nlocal) 7155a74a43SLisandro DalcinA.setFromOptions() 7255a74a43SLisandro DalcinA.setUp() 7355a74a43SLisandro Dalcin 7455a74a43SLisandro Dalcin''' 7555a74a43SLisandro DalcinAssemble matrix. 7655a74a43SLisandro Dalcin 7755a74a43SLisandro DalcinThe linear system is distributed across the processors by 7855a74a43SLisandro Dalcinchunks of contiguous rows, which correspond to contiguous 7955a74a43SLisandro Dalcinsections of the mesh on which the problem is discretized. 8055a74a43SLisandro DalcinFor matrix assembly, each processor contributes entries for 8155a74a43SLisandro Dalcinthe part that it owns locally. 8255a74a43SLisandro Dalcin''' 8355a74a43SLisandro Dalcin 8455a74a43SLisandro Dalcincol = np.zeros(3,dtype=PETSc.IntType) 8555a74a43SLisandro Dalcinvalue = np.zeros(3,dtype=PETSc.ScalarType) 8655a74a43SLisandro Dalcin 8755a74a43SLisandro Dalcinif not rstart: 8855a74a43SLisandro Dalcin rstart = 1 8955a74a43SLisandro Dalcin i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0 9055a74a43SLisandro Dalcin A.setValues(i,col[0:2],value[0:2]) 9155a74a43SLisandro Dalcin 9255a74a43SLisandro Dalcinif rend == n: 9355a74a43SLisandro Dalcin rend = n-1 9455a74a43SLisandro Dalcin i = n-1; col[0] = n-2; col[1] = n-1; value[0] = -1.0; value[1] = 2.0 9555a74a43SLisandro Dalcin A.setValues(i,col[0:2],value[0:2]) 9655a74a43SLisandro Dalcin 9755a74a43SLisandro Dalcin 9855a74a43SLisandro Dalcin''' Set entries corresponding to the mesh interior ''' 9955a74a43SLisandro Dalcinvalue[0] = -1.0; value[1] = 2.0; value[2] = -1.0 10055a74a43SLisandro Dalcinfor i in range(rstart,rend): 10155a74a43SLisandro Dalcin col[0] = i-1; col[1] = i; col[2] = i+1 10255a74a43SLisandro Dalcin A.setValues(i,col,value) 10355a74a43SLisandro Dalcin 10455a74a43SLisandro Dalcin 10555a74a43SLisandro Dalcin''' Assemble the matrix ''' 10655a74a43SLisandro DalcinA.assemblyBegin(A.AssemblyType.FINAL) 10755a74a43SLisandro DalcinA.assemblyEnd(A.AssemblyType.FINAL) 10855a74a43SLisandro Dalcin 10955a74a43SLisandro Dalcin''' Set exact solution; then compute right-hand-side vector. ''' 11055a74a43SLisandro Dalcinu.set(1.0) 11155a74a43SLisandro Dalcinb = A(u) 11255a74a43SLisandro Dalcin 11355a74a43SLisandro Dalcin''' 11455a74a43SLisandro Dalcin- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 11555a74a43SLisandro Dalcin Create the linear solver and set various options 11655a74a43SLisandro Dalcin- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ''' 11755a74a43SLisandro Dalcin''' 11855a74a43SLisandro DalcinCreate linear solver context 11955a74a43SLisandro Dalcin''' 12055a74a43SLisandro Dalcinksp = PETSc.KSP().create() 12155a74a43SLisandro Dalcin 12255a74a43SLisandro Dalcin''' 12355a74a43SLisandro DalcinSet operators. Here the matrix that defines the linear system 124*7addb90fSBarry Smithalso serves as the matrix from which the preconditioner is constructed. 12555a74a43SLisandro Dalcin''' 12655a74a43SLisandro Dalcinksp.setOperators(A,A) 12755a74a43SLisandro Dalcin 12855a74a43SLisandro Dalcin''' 12955a74a43SLisandro DalcinSet linear solver defaults for this problem (optional). 13055a74a43SLisandro Dalcin - By extracting the KSP and PC contexts from the KSP context, 13155a74a43SLisandro Dalcin we can then directly call any KSP and PC routines to set 13255a74a43SLisandro Dalcin various options. 13355a74a43SLisandro Dalcin - The following four statements are optional; all of these 13455a74a43SLisandro Dalcin parameters could alternatively be specified at runtime via 13555a74a43SLisandro Dalcin KSPSetFromOptions(); 13655a74a43SLisandro Dalcin''' 13755a74a43SLisandro Dalcinpc = ksp.getPC() 13855a74a43SLisandro Dalcinpc.setType('jacobi') 13955a74a43SLisandro Dalcinksp.setTolerances(rtol=1.e-7) 14055a74a43SLisandro Dalcin 14155a74a43SLisandro Dalcin''' 14255a74a43SLisandro DalcinSet runtime options, e.g., 14355a74a43SLisandro Dalcin-ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol> 14455a74a43SLisandro DalcinThese options will override those specified above as long as 14555a74a43SLisandro DalcinKSPSetFromOptions() is called _after_ any other customization 14655a74a43SLisandro Dalcinroutines. 14755a74a43SLisandro Dalcin''' 14855a74a43SLisandro Dalcinksp.setFromOptions() 14955a74a43SLisandro Dalcin 15055a74a43SLisandro Dalcin''' 15155a74a43SLisandro Dalcin- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 15255a74a43SLisandro Dalcin Solve the linear system 15355a74a43SLisandro Dalcin- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ''' 15455a74a43SLisandro Dalcin 15555a74a43SLisandro Dalcin''' Solve linear system ''' 15655a74a43SLisandro Dalcinksp.solve(b,x) 15755a74a43SLisandro Dalcin 15855a74a43SLisandro Dalcin''' 15955a74a43SLisandro DalcinView solver info; we could instead use the option -ksp_view to 16055a74a43SLisandro Dalcinprint this info to the screen at the conclusion of KSPSolve(). 16155a74a43SLisandro Dalcin''' 16255a74a43SLisandro Dalcin# ksp.view() 16355a74a43SLisandro Dalcin 16455a74a43SLisandro Dalcin''' 16555a74a43SLisandro Dalcin- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 16655a74a43SLisandro Dalcin Check solution and clean up 16755a74a43SLisandro Dalcin- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ''' 16855a74a43SLisandro Dalcin 16955a74a43SLisandro Dalcin''' Check the error ''' 17055a74a43SLisandro Dalcinx = x - u # x.axpy(-1.0,u) 17155a74a43SLisandro Dalcinnorm = x.norm(PETSc.NormType.NORM_2) 17255a74a43SLisandro Dalcinits = ksp.getIterationNumber() 17355a74a43SLisandro Dalcinif norm > tol: 17455a74a43SLisandro Dalcin PETSc.Sys.Print("Norm of error {}, Iterations {}\n".format(norm,its),comm=comm) 17555a74a43SLisandro Dalcinelse: 17655a74a43SLisandro Dalcin if size==1: 17755a74a43SLisandro Dalcin PETSc.Sys.Print("- Serial OK",comm=comm) 17855a74a43SLisandro Dalcin else: 17955a74a43SLisandro Dalcin PETSc.Sys.Print("- Parallel OK",comm=comm) 180