Actual source code: ex117.c
petsc-3.6.1 2015-08-06
2: static char help[] = "Tests Cholesky factorization for a SBAIJ matrix, (bs=2).\n";
3: /*
4: This code is modified from the code contributed by JUNWANG@uwm.edu on Apr 13, 2007
5: */
7: #include <petscmat.h>
11: int main(int argc,char **args)
12: {
14: Mat mat,fact,B;
15: PetscInt ind1[2],ind2[2];
16: PetscScalar temp[4];
17: PetscInt nnz[3];
18: IS perm,colp;
19: MatFactorInfo info;
20: PetscMPIInt size;
22: PetscInitialize(&argc,&args,0,help);
23: MPI_Comm_size(PETSC_COMM_WORLD,&size);
24: if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
26: nnz[0]=2;nnz[1]=1;nnz[2]=1;
27: MatCreateSeqSBAIJ(PETSC_COMM_SELF,2,6,6,0,nnz,&mat);
28:
29: ind1[0]=0;ind1[1]=1;
30: temp[0]=3;temp[1]=2;temp[2]=0;temp[3]=3;
31: MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);
32: ind2[0]=4;ind2[1]=5;
33: temp[0]=1;temp[1]=1;temp[2]=2;temp[3]=1;
34: MatSetValues(mat,2,ind1,2,ind2,temp,INSERT_VALUES);
35: ind1[0]=2;ind1[1]=3;
36: temp[0]=4;temp[1]=1;temp[2]=1;temp[3]=5;
37: MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);
38: ind1[0]=4;ind1[1]=5;
39: temp[0]=5;temp[1]=1;temp[2]=1;temp[3]=6;
40: MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);
42: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
43: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
45: MatDuplicate(mat,MAT_SHARE_NONZERO_PATTERN,&B);
46: ind1[0]=0;ind1[1]=1;
47: temp[0]=3;temp[1]=2;temp[2]=0;temp[3]=3;
48: MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);
49: ind2[0]=4;ind2[1]=5;
50: temp[0]=1;temp[1]=1;temp[2]=2;temp[3]=1;
51: MatSetValues(mat,2,ind1,2,ind2,temp,INSERT_VALUES);
52: ind1[0]=2;ind1[1]=3;
53: temp[0]=4;temp[1]=1;temp[2]=1;temp[3]=5;
54: MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);
55: ind1[0]=4;ind1[1]=5;
56: temp[0]=5;temp[1]=1;temp[2]=1;temp[3]=6;
57: MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);
59: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
60: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
62: PetscPrintf(PETSC_COMM_WORLD,"mat: \n");
63: MatView(mat,PETSC_VIEWER_STDOUT_SELF);
65: /* begin cholesky factorization */
66: MatGetOrdering(mat,MATORDERINGNATURAL,&perm,&colp);
67: ISDestroy(&colp);
69: info.fill=1.0;
70: MatGetFactor(mat,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&fact);
71: MatCholeskyFactorSymbolic(fact,mat,perm,&info);
72: MatCholeskyFactorNumeric(fact,mat,&info);
74: ISDestroy(&perm);
75: MatDestroy(&mat);
76: MatDestroy(&fact);
77: MatDestroy(&B);
78: PetscFinalize();
79: return 0;
80: }