xref: /petsc/src/mat/tests/ex117.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests Cholesky factorization for a SBAIJ matrix, (bs=2).\n";
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown   This code is modified from the code contributed by JUNWANG@uwm.edu on Apr 13, 2007
5c4762a1bSJed Brown */
6c4762a1bSJed Brown 
7c4762a1bSJed Brown #include <petscmat.h>
8c4762a1bSJed Brown 
9*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
10*d71ae5a4SJacob Faibussowitsch {
11c4762a1bSJed Brown   Mat           mat, fact, B;
12c4762a1bSJed Brown   PetscInt      ind1[2], ind2[2];
13c4762a1bSJed Brown   PetscScalar   temp[4];
14c4762a1bSJed Brown   PetscInt      nnz[3];
15c4762a1bSJed Brown   IS            perm, colp;
16c4762a1bSJed Brown   MatFactorInfo info;
17c4762a1bSJed Brown   PetscMPIInt   size;
18c4762a1bSJed Brown 
19327415f7SBarry Smith   PetscFunctionBeginUser;
209566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, 0, help));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
22be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
23c4762a1bSJed Brown 
249371c9d4SSatish Balay   nnz[0] = 2;
259371c9d4SSatish Balay   nnz[1] = 1;
269371c9d4SSatish Balay   nnz[2] = 1;
279566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqSBAIJ(PETSC_COMM_SELF, 2, 6, 6, 0, nnz, &mat));
28c4762a1bSJed Brown 
299371c9d4SSatish Balay   ind1[0] = 0;
309371c9d4SSatish Balay   ind1[1] = 1;
319371c9d4SSatish Balay   temp[0] = 3;
329371c9d4SSatish Balay   temp[1] = 2;
339371c9d4SSatish Balay   temp[2] = 0;
349371c9d4SSatish Balay   temp[3] = 3;
359566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 2, ind1, 2, ind1, temp, INSERT_VALUES));
369371c9d4SSatish Balay   ind2[0] = 4;
379371c9d4SSatish Balay   ind2[1] = 5;
389371c9d4SSatish Balay   temp[0] = 1;
399371c9d4SSatish Balay   temp[1] = 1;
409371c9d4SSatish Balay   temp[2] = 2;
419371c9d4SSatish Balay   temp[3] = 1;
429566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 2, ind1, 2, ind2, temp, INSERT_VALUES));
439371c9d4SSatish Balay   ind1[0] = 2;
449371c9d4SSatish Balay   ind1[1] = 3;
459371c9d4SSatish Balay   temp[0] = 4;
469371c9d4SSatish Balay   temp[1] = 1;
479371c9d4SSatish Balay   temp[2] = 1;
489371c9d4SSatish Balay   temp[3] = 5;
499566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 2, ind1, 2, ind1, temp, INSERT_VALUES));
509371c9d4SSatish Balay   ind1[0] = 4;
519371c9d4SSatish Balay   ind1[1] = 5;
529371c9d4SSatish Balay   temp[0] = 5;
539371c9d4SSatish Balay   temp[1] = 1;
549371c9d4SSatish Balay   temp[2] = 1;
559371c9d4SSatish Balay   temp[3] = 6;
569566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 2, ind1, 2, ind1, temp, INSERT_VALUES));
57c4762a1bSJed Brown 
589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
60c4762a1bSJed Brown 
619566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(mat, MAT_SHARE_NONZERO_PATTERN, &B));
629371c9d4SSatish Balay   ind1[0] = 0;
639371c9d4SSatish Balay   ind1[1] = 1;
649371c9d4SSatish Balay   temp[0] = 3;
659371c9d4SSatish Balay   temp[1] = 2;
669371c9d4SSatish Balay   temp[2] = 0;
679371c9d4SSatish Balay   temp[3] = 3;
689566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 2, ind1, 2, ind1, temp, INSERT_VALUES));
699371c9d4SSatish Balay   ind2[0] = 4;
709371c9d4SSatish Balay   ind2[1] = 5;
719371c9d4SSatish Balay   temp[0] = 1;
729371c9d4SSatish Balay   temp[1] = 1;
739371c9d4SSatish Balay   temp[2] = 2;
749371c9d4SSatish Balay   temp[3] = 1;
759566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 2, ind1, 2, ind2, temp, INSERT_VALUES));
769371c9d4SSatish Balay   ind1[0] = 2;
779371c9d4SSatish Balay   ind1[1] = 3;
789371c9d4SSatish Balay   temp[0] = 4;
799371c9d4SSatish Balay   temp[1] = 1;
809371c9d4SSatish Balay   temp[2] = 1;
819371c9d4SSatish Balay   temp[3] = 5;
829566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 2, ind1, 2, ind1, temp, INSERT_VALUES));
839371c9d4SSatish Balay   ind1[0] = 4;
849371c9d4SSatish Balay   ind1[1] = 5;
859371c9d4SSatish Balay   temp[0] = 5;
869371c9d4SSatish Balay   temp[1] = 1;
879371c9d4SSatish Balay   temp[2] = 1;
889371c9d4SSatish Balay   temp[3] = 6;
899566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 2, ind1, 2, ind1, temp, INSERT_VALUES));
90c4762a1bSJed Brown 
919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mat: \n"));
959566063dSJacob Faibussowitsch   PetscCall(MatView(mat, PETSC_VIEWER_STDOUT_SELF));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* begin cholesky factorization */
989566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(mat, MATORDERINGNATURAL, &perm, &colp));
999566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&colp));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   info.fill = 1.0;
1029566063dSJacob Faibussowitsch   PetscCall(MatGetFactor(mat, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &fact));
1039566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactorSymbolic(fact, mat, perm, &info));
1049566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactorNumeric(fact, mat, &info));
105c4762a1bSJed Brown 
1069566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
1079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
1089566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&fact));
1099566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1109566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
111b122ec5aSJacob Faibussowitsch   return 0;
112c4762a1bSJed Brown }
113