Actual source code: ex28.c
petsc-3.6.4 2016-04-12
2: static char help[] = "Test procedural KSPSetFromOptions() or at runtime; Test PCREDUNDANT.\n\n";
4: /*T
5: Concepts: KSP^basic parallel example;
6: Processors: n
7: T*/
8: #include <petscksp.h>
12: int main(int argc,char **args)
13: {
14: Vec x, b, u; /* approx solution, RHS, exact solution */
15: Mat A; /* linear system matrix */
16: KSP ksp; /* linear solver context */
17: PC pc; /* preconditioner context */
18: PetscReal norm; /* norm of solution error */
20: PetscInt i,n = 10,col[3],its,rstart,rend,nlocal;
21: PetscScalar neg_one = -1.0,one = 1.0,value[3];
22: PetscBool TEST_PROCEDURAL=PETSC_FALSE;
24: PetscInitialize(&argc,&args,(char*)0,help);
25: PetscOptionsGetInt(NULL,"-n",&n,NULL);
26: PetscOptionsGetBool(NULL,"-procedural",&TEST_PROCEDURAL,NULL);
28: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
29: Compute the matrix and right-hand-side vector that define
30: the linear system, Ax = b.
31: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
33: /*
34: Create vectors. Note that we form 1 vector from scratch and
35: then duplicate as needed. For this simple case let PETSc decide how
36: many elements of the vector are stored on each processor. The second
37: argument to VecSetSizes() below causes PETSc to decide.
38: */
39: VecCreate(PETSC_COMM_WORLD,&x);
40: VecSetSizes(x,PETSC_DECIDE,n);
41: VecSetFromOptions(x);
42: VecDuplicate(x,&b);
43: VecDuplicate(x,&u);
45: /* Identify the starting and ending mesh points on each
46: processor for the interior part of the mesh. We let PETSc decide
47: above. */
49: VecGetOwnershipRange(x,&rstart,&rend);
50: VecGetLocalSize(x,&nlocal);
52: /* Create a tridiagonal matrix. See ../tutorials/ex23.c */
53: MatCreate(PETSC_COMM_WORLD,&A);
54: MatSetSizes(A,nlocal,nlocal,n,n);
55: MatSetFromOptions(A);
56: MatSetUp(A);
57: /* Assemble matrix */
58: if (!rstart) {
59: rstart = 1;
60: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
61: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
62: }
63: if (rend == n) {
64: rend = n-1;
65: i = n-1; col[0] = n-2; col[1] = n-1; value[0] = -1.0; value[1] = 2.0;
66: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
67: }
69: /* Set entries corresponding to the mesh interior */
70: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
71: for (i=rstart; i<rend; i++) {
72: col[0] = i-1; col[1] = i; col[2] = i+1;
73: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
74: }
75: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
76: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
78: /* Set exact solution; then compute right-hand-side vector. */
79: VecSet(u,one);
80: MatMult(A,u,b);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Create the linear solver and set various options
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85: KSPCreate(PETSC_COMM_WORLD,&ksp);
86: KSPSetOperators(ksp,A,A);
88: /*
89: Set linear solver defaults for this problem (optional).
90: - By extracting the KSP and PC contexts from the KSP context,
91: we can then directly call any KSP and PC routines to set
92: various options.
93: - The following statements are optional; all of these
94: parameters could alternatively be specified at runtime via
95: KSPSetFromOptions();
96: */
97: if (TEST_PROCEDURAL) {
98: /* Example of runtime options: '-pc_redundant_number 3 -redundant_ksp_type gmres -redundant_pc_type bjacobi' */
99: PetscMPIInt size,rank,subsize;
100: Mat A_redundant;
101: KSP innerksp;
102: PC innerpc;
103: MPI_Comm subcomm;
105: KSPGetPC(ksp,&pc);
106: PCSetType(pc,PCREDUNDANT);
107: MPI_Comm_size(PETSC_COMM_WORLD,&size);
108: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
109: if (size < 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "Num of processes %d must greater than 2",size);
110: PCRedundantSetNumber(pc,size-2);
111: KSPSetFromOptions(ksp);
113: /* Get subcommunicator and redundant matrix */
114: KSPSetUp(ksp);
115: PCRedundantGetKSP(pc,&innerksp);
116: KSPGetPC(innerksp,&innerpc);
117: PCGetOperators(innerpc,NULL,&A_redundant);
118: PetscObjectGetComm((PetscObject)A_redundant,&subcomm);
119: MPI_Comm_size(subcomm,&subsize);
120: if (subsize==1 && !rank) {
121: printf("A_redundant:\n");
122: MatView(A_redundant,PETSC_VIEWER_STDOUT_SELF);
123: }
124: } else {
125: KSPSetFromOptions(ksp);
126: }
127:
128: /* Solve linear system */
129: KSPSolve(ksp,b,x);
131: /* Check the error */
132: VecAXPY(x,neg_one,u);
133: VecNorm(x,NORM_2,&norm);
134: KSPGetIterationNumber(ksp,&its);
135: if (norm > 1.e-14) {
136: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);
137: }
139: /* Free work space. */
140: VecDestroy(&x); VecDestroy(&u);
141: VecDestroy(&b); MatDestroy(&A);
142: KSPDestroy(&ksp);
143: PetscFinalize();
144: return 0;
145: }