xref: /petsc/src/mat/tests/ex101.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1c4762a1bSJed Brown static char help[] = "Testing PtAP for SeqMAIJ matrix, P, with SeqAIJ matrix, A.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc,char **argv)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   Mat            pA,P,aijP;
8c4762a1bSJed Brown   PetscScalar    pa[]={1.,-1.,0.,0.,1.,-1.,0.,0.,1.};
9c4762a1bSJed Brown   PetscInt       i,pij[]={0,1,2};
10c4762a1bSJed Brown   PetscInt       aij[3][3]={{0,1,2},{3,4,5},{6,7,8}};
11c4762a1bSJed Brown   Mat            A,mC,C;
12c4762a1bSJed Brown   PetscScalar    one=1.;
13c4762a1bSJed Brown   PetscErrorCode ierr;
14c4762a1bSJed Brown   PetscMPIInt    size;
15c20d7725SJed Brown   PetscBool      flg;
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
18ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
19*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   /* Create MAIJ matrix, P */
22c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_SELF,&pA);CHKERRQ(ierr);
23c4762a1bSJed Brown   ierr = MatSetSizes(pA,3,3,3,3);CHKERRQ(ierr);
24c4762a1bSJed Brown   ierr = MatSetType(pA,MATSEQAIJ);CHKERRQ(ierr);
25c4762a1bSJed Brown   ierr = MatSetUp(pA);CHKERRQ(ierr);
26c4762a1bSJed Brown   ierr = MatSetOption(pA,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
27c4762a1bSJed Brown   ierr = MatSetValues(pA,3,pij,3,pij,pa,ADD_VALUES);CHKERRQ(ierr);
28c4762a1bSJed Brown   ierr = MatAssemblyBegin(pA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
29c4762a1bSJed Brown   ierr = MatAssemblyEnd(pA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
30c4762a1bSJed Brown   ierr = MatCreateMAIJ(pA,3,&P);CHKERRQ(ierr);
31c4762a1bSJed Brown   ierr = MatDestroy(&pA);CHKERRQ(ierr);
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   /* Create AIJ equivalent matrix, aijP, for comparison testing */
34c4762a1bSJed Brown   ierr = MatConvert(P,MATSEQAIJ,MAT_INITIAL_MATRIX,&aijP);CHKERRQ(ierr);
35c4762a1bSJed Brown 
36c20d7725SJed Brown   /* Create AIJ matrix A */
37c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_SELF,&A);CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = MatSetSizes(A,9,9,9,9);CHKERRQ(ierr);
39c4762a1bSJed Brown   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
40c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
41c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
42c4762a1bSJed Brown   ierr = MatSetValues(A,3,aij[0],3,aij[0],pa,ADD_VALUES);CHKERRQ(ierr);
43c4762a1bSJed Brown   ierr = MatSetValues(A,3,aij[1],3,aij[1],pa,ADD_VALUES);CHKERRQ(ierr);
44c4762a1bSJed Brown   ierr = MatSetValues(A,3,aij[2],3,aij[2],pa,ADD_VALUES);CHKERRQ(ierr);
45c4762a1bSJed Brown   for (i=0; i<9; i++) {
46c4762a1bSJed Brown     ierr = MatSetValue(A,i,i,one,ADD_VALUES);CHKERRQ(ierr);
47c4762a1bSJed Brown   }
48c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
49c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Perform PtAP_SeqAIJ_SeqAIJ for comparison testing */
52c4762a1bSJed Brown   ierr = MatPtAP(A,aijP,MAT_INITIAL_MATRIX,1.,&C);CHKERRQ(ierr);
53c20d7725SJed Brown 
54c20d7725SJed Brown   /* Perform PtAP_SeqAIJ_SeqMAIJ */
55c20d7725SJed Brown   /* Developer API */
56c20d7725SJed Brown   ierr = MatProductCreate(A,P,NULL,&mC);CHKERRQ(ierr);
57c20d7725SJed Brown   ierr = MatProductSetType(mC,MATPRODUCT_PtAP);CHKERRQ(ierr);
58c20d7725SJed Brown   ierr = MatProductSetAlgorithm(mC,"default");CHKERRQ(ierr);
59c20d7725SJed Brown   ierr = MatProductSetFill(mC,PETSC_DEFAULT);CHKERRQ(ierr);
60c20d7725SJed Brown   ierr = MatProductSetFromOptions(mC);CHKERRQ(ierr);
61c20d7725SJed Brown   ierr = MatProductSymbolic(mC);CHKERRQ(ierr);
62c20d7725SJed Brown   ierr = MatProductNumeric(mC);CHKERRQ(ierr);
63c20d7725SJed Brown   ierr = MatProductNumeric(mC);CHKERRQ(ierr);
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* Check mC = C */
66c20d7725SJed Brown   ierr = MatEqual(C,mC,&flg);CHKERRQ(ierr);
67*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"MatProduct C != mC");
68c20d7725SJed Brown   ierr = MatDestroy(&mC);CHKERRQ(ierr);
69c20d7725SJed Brown 
70c20d7725SJed Brown   /* User API */
71c20d7725SJed Brown   ierr = MatPtAP(A,P,MAT_INITIAL_MATRIX,1.,&mC);CHKERRQ(ierr);
72c20d7725SJed Brown   ierr = MatPtAP(A,P,MAT_REUSE_MATRIX,1.,&mC);CHKERRQ(ierr);
73c20d7725SJed Brown 
74c20d7725SJed Brown   /* Check mC = C */
75c20d7725SJed Brown   ierr = MatEqual(C,mC,&flg);CHKERRQ(ierr);
76*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"MatPtAP C != mC");
77c20d7725SJed Brown   ierr = MatDestroy(&mC);CHKERRQ(ierr);
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* Cleanup */
80c4762a1bSJed Brown   ierr = MatDestroy(&P);CHKERRQ(ierr);
81c4762a1bSJed Brown   ierr = MatDestroy(&aijP);CHKERRQ(ierr);
82c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
83c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
84c4762a1bSJed Brown   ierr = PetscFinalize();
85c4762a1bSJed Brown   return ierr;
86c4762a1bSJed Brown }
87c4762a1bSJed Brown 
88c4762a1bSJed Brown /*TEST
89c4762a1bSJed Brown 
90c4762a1bSJed Brown    test:
91c4762a1bSJed Brown       output_file: output/ex101.out
92c4762a1bSJed Brown 
93c4762a1bSJed Brown TEST*/
94