/* Usage: mpirun ex2 [-help] [all PETSc options] */ static char help[] = "Solves a linear system in parallel with SLES.\n\ Input parameters include:\n\ -random_exact_sol : use a random exact solution vector\n\ -view_exact_sol : write exact solution vector to stdout\n\ -m : number of mesh points in x-direction\n\ -n : number of mesh points in y-direction\n\n"; /* Include "sles.h" so that we can use SLES solvers. Note that this file automatically includes: petsc.h - base PETSc routines vec.h - vectors sys.h - system routines mat.h - matrices is.h - index sets ksp.h - Krylov subspace methods viewer.h - viewers pc.h - preconditioners */ #include "sles.h" int main(int argc,char **args) { Vec x, b, u; /* approx solution, RHS, exact solution */ Mat A; /* linear system matrix */ SLES sles; /* linear solver context */ PetscRandom rctx; /* random number generator context */ double norm; /* norm of solution error */ int i, j, I, J, Istart, Iend, ierr, m = 8, n = 7, its, flg; Scalar v, one = 1.0, neg_one = -1.0; /* These variables are currently unused */ /* PC pc; */ /* preconditioner context */ /* KSP ksp; */ /* Krylov subspace method context */ PetscInitialize(&argc,&args,(char *)0,help); ierr = OptionsGetInt(PETSC_NULL,"-m",&m,&flg); CHKERRA(ierr); ierr = OptionsGetInt(PETSC_NULL,"-n",&n,&flg); CHKERRA(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Compute the matrix and right-hand-side vector that define the linear system, Ax = b. - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ /* Create parallel matrix, specifying only its global dimensions. When using MatCreate(), the matrix format can be specified at runtime. Also, the parallel partitioning of the matrix is determined by PETSc at runtime. */ ierr = MatCreate(PETSC_COMM_WORLD,m*n,m*n,&A); CHKERRA(ierr); /* Currently, all PETSc parallel matrix formats are partitioned by contiguous chunks of rows across the processors. Determine which rows of the matrix are locally owned. */ ierr = MatGetOwnershipRange(A,&Istart,&Iend); CHKERRA(ierr); /* Set matrix elements for the 2-D, five-point stencil in parallel. - Each processor needs to insert only elements that it owns locally (but any non-local elements will be sent to the appropriate processor during matrix assembly). - Always specify global rows and columns of matrix entries. */ for ( I=Istart; I0 ) {J = I - n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES); CHKERRA(ierr);} if ( i0 ) {J = I - 1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES); CHKERRA(ierr);} if ( j -pc_type -ksp_monitor -ksp_rtol These options will override those specified above as long as SLESSetFromOptions() is called _after_ any other customization routines. */ ierr = SLESSetFromOptions(sles); CHKERRA(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Solve the linear system - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = SLESSolve(sles,b,x,&its); CHKERRA(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Check solution and clean up - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ /* Check the error */ ierr = VecAXPY(&neg_one,u,x); CHKERRA(ierr); ierr = VecNorm(x,NORM_2,&norm); CHKERRA(ierr); /* Print convergence information. PetscPrintf() produces a single print statement from all processes that share a communicator. */ if (norm > 1.e-12) PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g iterations %d\n",norm,its); else PetscPrintf(PETSC_COMM_WORLD,"Norm of error < 1.e-12 Iterations %d\n",its); /* Free work space. All PETSc objects should be destroyed when they are no longer needed. */ ierr = SLESDestroy(sles); CHKERRA(ierr); ierr = VecDestroy(u); CHKERRA(ierr); ierr = VecDestroy(x); CHKERRA(ierr); ierr = VecDestroy(b); CHKERRA(ierr); ierr = MatDestroy(A); CHKERRA(ierr); /* Always call PetscFinalize() before exiting a program. This routine - finalizes the PETSc libraries as well as MPI - provides summary and diagnostic information if certain runtime options are chosen (e.g., -log_summary). */ PetscFinalize(); return 0; }