Actual source code: ex6.c
petsc-3.13.6 2020-09-29
2: static char help[] = "Solves a tridiagonal linear system with KSP. \n\
3: It illustrates how to do one symbolic factorization and multiple numeric factorizations using same matrix structure. \n\n";
5: #include <petscksp.h>
6: int main(int argc,char **args)
7: {
8: Vec x, b, u; /* approx solution, RHS, exact solution */
9: Mat A; /* linear system matrix */
10: KSP ksp; /* linear solver context */
11: PC pc; /* preconditioner context */
12: PetscReal norm; /* norm of solution error */
14: PetscInt i,col[3],its,rstart,rend,N=10,num_numfac;
15: PetscScalar value[3];
17: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
18: PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);
20: /* Create and assemble matrix. */
21: MatCreate(PETSC_COMM_WORLD,&A);
22: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
23: MatSetFromOptions(A);
24: MatSetUp(A);
25: MatGetOwnershipRange(A,&rstart,&rend);
27: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
28: for (i=rstart; i<rend; i++) {
29: col[0] = i-1; col[1] = i; col[2] = i+1;
30: if (i == 0) {
31: MatSetValues(A,1,&i,2,col+1,value+1,INSERT_VALUES);
32: } else if (i == N-1) {
33: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
34: } else {
35: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
36: }
37: }
38: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
39: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
40: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
42: /* Create vectors */
43: MatCreateVecs(A,&x,&b);
44: VecDuplicate(x,&u);
46: /* Set exact solution; then compute right-hand-side vector. */
47: VecSet(u,1.0);
48: MatMult(A,u,b);
50: /* Create the linear solver and set various options. */
51: KSPCreate(PETSC_COMM_WORLD,&ksp);
52: KSPGetPC(ksp,&pc);
53: PCSetType(pc,PCJACOBI);
54: KSPSetTolerances(ksp,1.e-5,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
55: KSPSetOperators(ksp,A,A);
56: KSPSetFromOptions(ksp);
58: num_numfac = 1;
59: PetscOptionsGetInt(NULL,NULL,"-num_numfac",&num_numfac,NULL);
60: while (num_numfac--) {
61: /* An example on how to update matrix A for repeated numerical factorization and solve. */
62: PetscScalar one=1.0;
63: PetscInt i = 0;
64: MatSetValues(A,1,&i,1,&i,&one,ADD_VALUES);
65: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
66: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
67: /* Update b */
68: MatMult(A,u,b);
70: /* Solve the linear system */
71: KSPSolve(ksp,b,x);
73: /* Check the solution and clean up */
74: VecAXPY(x,-1.0,u);
75: VecNorm(x,NORM_2,&norm);
76: KSPGetIterationNumber(ksp,&its);
77: if (norm > 100*PETSC_MACHINE_EPSILON) {
78: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);
79: }
80: }
82: /* Free work space. */
83: VecDestroy(&x); VecDestroy(&u);
84: VecDestroy(&b); MatDestroy(&A);
85: KSPDestroy(&ksp);
87: PetscFinalize();
88: return ierr;
89: }
91: /*TEST
93: test:
94: args: -num_numfac 2 -pc_type lu
96: test:
97: suffix: 2
98: args: -num_numfac 2 -pc_type lu -pc_factor_mat_solver_type mumps
99: requires: mumps
101: test:
102: suffix: 3
103: nsize: 3
104: args: -num_numfac 2 -pc_type lu -pc_factor_mat_solver_type mumps
105: requires: mumps
107: TEST*/