xref: /petsc/src/mat/tests/ex161.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1 static char help[] = "Test MatTransposeColoring for SeqAIJ matrices. Used for '-matmattransmult_color' on  MatMatTransposeMult \n\n";
2 
3 #include <petscmat.h>
4 #include <petsc/private/matimpl.h> /* Need struct _p_MatTransposeColoring for this test. */
5 
6 int main(int argc,char **argv)
7 {
8   Mat                  A,R,C,C_dense,C_sparse,Rt_dense,P,PtAP;
9   PetscInt             row,col,m,n;
10   PetscErrorCode       ierr;
11   MatScalar            one         =1.0,val;
12   MatColoring          mc;
13   MatTransposeColoring matcoloring = 0;
14   ISColoring           iscoloring;
15   PetscBool            equal;
16   PetscMPIInt          size;
17 
18   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
19   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
20   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
21 
22   /* Create seqaij A */
23   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A));
24   CHKERRQ(MatSetSizes(A,4,4,4,4));
25   CHKERRQ(MatSetType(A,MATSEQAIJ));
26   CHKERRQ(MatSetFromOptions(A));
27   CHKERRQ(MatSetUp(A));
28   row  = 0; col=0; val=1.0; CHKERRQ(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES));
29   row  = 1; col=3; val=2.0; CHKERRQ(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES));
30   row  = 2; col=2; val=3.0; CHKERRQ(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES));
31   row  = 3; col=0; val=4.0; CHKERRQ(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES));
32   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
33   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
34   CHKERRQ(MatSetOptionsPrefix(A,"A_"));
35   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
36   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n"));
37 
38   /* Create seqaij R */
39   CHKERRQ(MatCreate(PETSC_COMM_SELF,&R));
40   CHKERRQ(MatSetSizes(R,2,4,2,4));
41   CHKERRQ(MatSetType(R,MATSEQAIJ));
42   CHKERRQ(MatSetFromOptions(R));
43   CHKERRQ(MatSetUp(R));
44   row  = 0; col=0; CHKERRQ(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES));
45   row  = 0; col=1; CHKERRQ(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES));
46 
47   row  = 1; col=1; CHKERRQ(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES));
48   row  = 1; col=2; CHKERRQ(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES));
49   row  = 1; col=3; CHKERRQ(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES));
50   CHKERRQ(MatAssemblyBegin(R,MAT_FINAL_ASSEMBLY));
51   CHKERRQ(MatAssemblyEnd(R,MAT_FINAL_ASSEMBLY));
52   CHKERRQ(MatSetOptionsPrefix(R,"R_"));
53   CHKERRQ(MatView(R,PETSC_VIEWER_STDOUT_WORLD));
54   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n"));
55 
56   /* C = A*R^T */
57   CHKERRQ(MatMatTransposeMult(A,R,MAT_INITIAL_MATRIX,2.0,&C));
58   CHKERRQ(MatSetOptionsPrefix(C,"ARt_"));
59   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
60   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n"));
61 
62   /* Create MatTransposeColoring from symbolic C=A*R^T */
63   CHKERRQ(MatColoringCreate(C,&mc));
64   CHKERRQ(MatColoringSetDistance(mc,2));
65   /* CHKERRQ(MatColoringSetType(mc,MATCOLORINGSL)); */
66   CHKERRQ(MatColoringSetFromOptions(mc));
67   CHKERRQ(MatColoringApply(mc,&iscoloring));
68   CHKERRQ(MatColoringDestroy(&mc));
69   CHKERRQ(MatTransposeColoringCreate(C,iscoloring,&matcoloring));
70   CHKERRQ(ISColoringDestroy(&iscoloring));
71 
72   /* Create Rt_dense */
73   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&Rt_dense));
74   CHKERRQ(MatSetSizes(Rt_dense,4,matcoloring->ncolors,PETSC_DECIDE,PETSC_DECIDE));
75   CHKERRQ(MatSetType(Rt_dense,MATDENSE));
76   CHKERRQ(MatSeqDenseSetPreallocation(Rt_dense,NULL));
77   CHKERRQ(MatAssemblyBegin(Rt_dense,MAT_FINAL_ASSEMBLY));
78   CHKERRQ(MatAssemblyEnd(Rt_dense,MAT_FINAL_ASSEMBLY));
79   CHKERRQ(MatGetLocalSize(Rt_dense,&m,&n));
80   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Rt_dense: %" PetscInt_FMT ",%" PetscInt_FMT "\n",m,n));
81 
82   /* Get Rt_dense by Apply MatTransposeColoring to R */
83   CHKERRQ(MatTransColoringApplySpToDen(matcoloring,R,Rt_dense));
84 
85   /* C_dense = A*Rt_dense */
86   CHKERRQ(MatMatMult(A,Rt_dense,MAT_INITIAL_MATRIX,2.0,&C_dense));
87   CHKERRQ(MatSetOptionsPrefix(C_dense,"ARt_dense_"));
88   /*CHKERRQ(MatView(C_dense,PETSC_VIEWER_STDOUT_WORLD)); */
89   /*CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n")); */
90 
91   /* Recover C from C_dense */
92   CHKERRQ(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&C_sparse));
93   CHKERRQ(MatTransColoringApplyDenToSp(matcoloring,C_dense,C_sparse));
94   CHKERRQ(MatSetOptionsPrefix(C_sparse,"ARt_color_"));
95   CHKERRQ(MatView(C_sparse,PETSC_VIEWER_STDOUT_WORLD));
96   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n"));
97 
98   CHKERRQ(MatDestroy(&C_dense));
99   CHKERRQ(MatDestroy(&C_sparse));
100   CHKERRQ(MatDestroy(&Rt_dense));
101   CHKERRQ(MatTransposeColoringDestroy(&matcoloring));
102   CHKERRQ(MatDestroy(&C));
103 
104   /* Test PtAP = P^T*A*P, P = R^T */
105   CHKERRQ(MatTranspose(R,MAT_INITIAL_MATRIX,&P));
106   CHKERRQ(MatPtAP(A,P,MAT_INITIAL_MATRIX,2.0,&PtAP));
107   CHKERRQ(MatSetOptionsPrefix(PtAP,"PtAP_"));
108   CHKERRQ(MatView(PtAP,PETSC_VIEWER_STDOUT_WORLD));
109   CHKERRQ(MatDestroy(&P));
110 
111   /* Test C = RARt */
112   CHKERRQ(MatRARt(A,R,MAT_INITIAL_MATRIX,2.0,&C));
113   CHKERRQ(MatRARt(A,R,MAT_REUSE_MATRIX,2.0,&C));
114   CHKERRQ(MatEqual(C,PtAP,&equal));
115   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: PtAP != RARt");
116 
117   /* Free spaces */
118   CHKERRQ(MatDestroy(&C));
119   CHKERRQ(MatDestroy(&A));
120   CHKERRQ(MatDestroy(&R));
121   CHKERRQ(MatDestroy(&PtAP));
122   ierr = PetscFinalize();
123   return ierr;
124 }
125 
126 /*TEST
127 
128    test:
129       output_file: output/ex161.out
130    test:
131       suffix: 2
132       args: -matmattransmult_via color
133       output_file: output/ex161.out
134 
135    test:
136       suffix: 3
137       args: -matmattransmult_via color -matden2sp_brows 3
138       output_file: output/ex161.out
139 
140    test:
141       suffix: 4
142       args: -matmattransmult_via color -matrart_via r*art
143       output_file: output/ex161.out
144 
145    test:
146       suffix: 5
147       args: -matmattransmult_via color -matrart_via coloring_rart
148       output_file: output/ex161.out
149 
150 TEST*/
151