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