Actual source code: ex49.c
petsc-3.9.0 2018-04-07
2: static char help[] = "Tests SeqSBAIJ factorizations for different block sizes\n\n";
4: /*T
5: requires: !single
6: T*/
8: #include <petscksp.h>
10: int main(int argc,char **args)
11: {
12: Vec x,b,u;
13: Mat A; /* linear system matrix */
14: KSP ksp; /* linear solver context */
15: PetscRandom rctx; /* random number generator context */
16: PetscReal norm; /* norm of solution error */
17: PetscInt i,j,k,l,n = 27,its,bs = 2,Ii,J;
19: PetscScalar v;
20:
21: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
22: PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);
23: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
25: MatCreate(PETSC_COMM_WORLD,&A);
26: MatSetSizes(A,n*bs,n*bs,PETSC_DETERMINE,PETSC_DETERMINE);
27: MatSetBlockSize(A,bs);
28: MatSetType(A,MATSEQSBAIJ);
29: MatSetFromOptions(A);
30: MatSetUp(A);
32: /*
33: Don't bother to preallocate matrix
34: */
35: PetscRandomCreate(PETSC_COMM_SELF,&rctx);
36: for (i=0; i<n; i++) {
37: for (j=0; j<n; j++) {
38: PetscRandomGetValue(rctx,&v);
39: if (PetscRealPart(v) < .25 || i == j) {
40: for (k=0; k<bs; k++) {
41: for (l=0; l<bs; l++) {
42: PetscRandomGetValue(rctx,&v);
43: Ii = i*bs + k;
44: J = j*bs + l;
45: if (Ii == J) v += 10.;
46: MatSetValue(A,Ii,J,v,INSERT_VALUES);
47: }
48: }
49: }
50: }
51: }
53: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
54: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
56: VecCreate(PETSC_COMM_WORLD,&u);
57: VecSetSizes(u,PETSC_DECIDE,n*bs);
58: VecSetFromOptions(u);
59: VecDuplicate(u,&b);
60: VecDuplicate(b,&x);
62: VecSet(u,1.0);
63: MatMult(A,u,b);
65: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66: Create the linear solver and set various options
67: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69: /*
70: Create linear solver context
71: */
72: KSPCreate(PETSC_COMM_WORLD,&ksp);
74: /*
75: Set operators. Here the matrix that defines the linear system
76: also serves as the preconditioning matrix.
77: */
78: KSPSetOperators(ksp,A,A);
80: KSPSetFromOptions(ksp);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Solve the linear system
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: KSPSolve(ksp,b,x);
88: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89: Check solution and clean up
90: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92: /*
93: Check the error
94: */
95: VecAXPY(x,-1.0,u);
96: VecNorm(x,NORM_2,&norm);
97: KSPGetIterationNumber(ksp,&its);
99: /*
100: Print convergence information. PetscPrintf() produces a single
101: print statement from all processes that share a communicator.
102: An alternative is PetscFPrintf(), which prints to a file.
103: */
104: if (norm > 100*PETSC_SMALL) {
105: PetscPrintf(PETSC_COMM_WORLD,"Norm of residual %g iterations %D bs %D\n",(double)norm,its,bs);
106: }
108: /*
109: Free work space. All PETSc objects should be destroyed when they
110: are no longer needed.
111: */
112: KSPDestroy(&ksp);
113: VecDestroy(&u);
114: VecDestroy(&x);
115: VecDestroy(&b);
116: MatDestroy(&A);
117: PetscRandomDestroy(&rctx);
119: /*
120: Always call PetscFinalize() before exiting a program. This routine
121: - finalizes the PETSc libraries as well as MPI
122: - provides summary and diagnostic information if certain runtime
123: options are chosen (e.g., -log_view).
124: */
125: PetscFinalize();
126: return ierr;
127: }
130: /*TEST
132: test:
133: args: -bs {{1 2 3 4 5 6 7 8 9 10 11 12}} -pc_type cholesky
135: TEST*/