Actual source code: ex33.c
petsc-3.8.3 2017-12-09
1: static char help[] = "Test MatGetInertia().\n\n";
3: /*
4: Examples of command line options:
5: ./ex33 -sigma 2.0 -pc_factor_mat_solver_package mumps -mat_mumps_icntl_13 1
6: ./ex33 -sigma <shift> -fA <matrix_file>
7: */
9: #include <petscksp.h>
10: int main(int argc,char **args)
11: {
12: Mat A,B,F;
14: KSP ksp;
15: PC pc;
16: PetscInt N, n=10, m, Istart, Iend, II, J, i,j;
17: PetscInt nneg, nzero, npos;
18: PetscScalar v,sigma;
19: PetscBool flag,loadA=PETSC_FALSE,loadB=PETSC_FALSE;
20: char file[2][PETSC_MAX_PATH_LEN];
21: PetscViewer viewer;
22: PetscMPIInt rank;
24: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
25: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
26: Compute the matrices that define the eigensystem, Ax=kBx
27: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
28: PetscOptionsGetString(NULL,NULL,"-fA",file[0],PETSC_MAX_PATH_LEN,&loadA);
29: if (loadA) {
30: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[0],FILE_MODE_READ,&viewer);
31: MatCreate(PETSC_COMM_WORLD,&A);
32: MatSetType(A,MATSBAIJ);
33: MatLoad(A,viewer);
34: PetscViewerDestroy(&viewer);
36: PetscOptionsGetString(NULL,NULL,"-fB",file[1],PETSC_MAX_PATH_LEN,&loadB);
37: if (loadB) {
38: /* load B to get A = A + sigma*B */
39: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[1],FILE_MODE_READ,&viewer);
40: MatCreate(PETSC_COMM_WORLD,&B);
41: MatSetType(B,MATSBAIJ);
42: MatLoad(B,viewer);
43: PetscViewerDestroy(&viewer);
44: }
45: }
47: if (!loadA) { /* Matrix A is copied from slepc-3.0.0-p6/src/examples/ex13.c. */
48: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
49: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
50: if (!flag) m=n;
51: N = n*m;
52: MatCreate(PETSC_COMM_WORLD,&A);
53: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
54: MatSetType(A,MATSBAIJ);
55: MatSetFromOptions(A);
56: MatSetUp(A);
58: MatSetOption(A,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
59: MatGetOwnershipRange(A,&Istart,&Iend);
60: for (II=Istart; II<Iend; II++) {
61: v = -1.0; i = II/n; j = II-i*n;
62: if (i>0) { J=II-n; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES); }
63: if (i<m-1) { J=II+n; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES); }
64: if (j>0) { J=II-1; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES); }
65: if (j<n-1) { J=II+1; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES); }
66: v=4.0; MatSetValues(A,1,&II,1,&II,&v,INSERT_VALUES);
68: }
69: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
70: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
71: }
72: /* MatView(A,PETSC_VIEWER_STDOUT_WORLD); */
74: if (!loadB) {
75: MatGetLocalSize(A,&m,&n);
76: MatCreate(PETSC_COMM_WORLD,&B);
77: MatSetSizes(B,m,n,PETSC_DECIDE,PETSC_DECIDE);
78: MatSetType(B,MATSBAIJ);
79: MatSetFromOptions(B);
80: MatSetUp(B);
81: MatSetOption(B,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
82: MatGetOwnershipRange(A,&Istart,&Iend);
84: for (II=Istart; II<Iend; II++) {
85: /* v=4.0; MatSetValues(B,1,&II,1,&II,&v,INSERT_VALUES); */
86: v=1.0; MatSetValues(B,1,&II,1,&II,&v,INSERT_VALUES);
87: }
88: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
89: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
90: }
91: /* MatView(B,PETSC_VIEWER_STDOUT_WORLD); */
93: /* Set a shift: A = A - sigma*B */
94: PetscOptionsGetScalar(NULL,NULL,"-sigma",&sigma,&flag);
95: if (flag) {
96: sigma = -1.0 * sigma;
97: MatAXPY(A,sigma,B,DIFFERENT_NONZERO_PATTERN); /* A <- A - sigma*B */
98: /* MatView(A,PETSC_VIEWER_STDOUT_WORLD); */
99: }
101: /* Test MatGetInertia() */
102: KSPCreate(PETSC_COMM_WORLD,&ksp);
103: KSPSetType(ksp,KSPPREONLY);
104: KSPSetOperators(ksp,A,A);
106: KSPGetPC(ksp,&pc);
107: PCSetType(pc,PCCHOLESKY);
108: PCSetFromOptions(pc);
110: PCSetUp(pc);
111: PCFactorGetMatrix(pc,&F);
112: MatGetInertia(F,&nneg,&nzero,&npos);
113: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
114: if (!rank) {
115: PetscPrintf(PETSC_COMM_SELF," MatInertia: nneg: %D, nzero: %D, npos: %D\n",nneg,nzero,npos);
116: }
118: /* Destroy */
119: KSPDestroy(&ksp);
120: MatDestroy(&A);
121: MatDestroy(&B);
122: PetscFinalize();
123: return ierr;
124: }