xref: /petsc/src/mat/tests/ex161.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Test MatTransposeColoring for SeqAIJ matrices. Used for '-matmattransmult_color' on  MatMatTransposeMult \n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown #include <petsc/private/matimpl.h> /* Need struct _p_MatTransposeColoring for this test. */
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **argv)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat                  A,R,C,C_dense,C_sparse,Rt_dense,P,PtAP;
9c4762a1bSJed Brown   PetscInt             row,col,m,n;
10c4762a1bSJed Brown   MatScalar            one         =1.0,val;
11c4762a1bSJed Brown   MatColoring          mc;
12c4762a1bSJed Brown   MatTransposeColoring matcoloring = 0;
13c4762a1bSJed Brown   ISColoring           iscoloring;
14c4762a1bSJed Brown   PetscBool            equal;
15c4762a1bSJed Brown   PetscMPIInt          size;
16c4762a1bSJed Brown 
17*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
185f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
192c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   /* Create seqaij A */
225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,4,4,4,4));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATSEQAIJ));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
275f80ce2aSJacob Faibussowitsch   row  = 0; col=0; val=1.0; CHKERRQ(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES));
285f80ce2aSJacob Faibussowitsch   row  = 1; col=3; val=2.0; CHKERRQ(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES));
295f80ce2aSJacob Faibussowitsch   row  = 2; col=2; val=3.0; CHKERRQ(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES));
305f80ce2aSJacob Faibussowitsch   row  = 3; col=0; val=4.0; CHKERRQ(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOptionsPrefix(A,"A_"));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n"));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* Create seqaij R */
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&R));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(R,2,4,2,4));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(R,MATSEQAIJ));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(R));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(R));
435f80ce2aSJacob Faibussowitsch   row  = 0; col=0; CHKERRQ(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES));
445f80ce2aSJacob Faibussowitsch   row  = 0; col=1; CHKERRQ(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES));
45c4762a1bSJed Brown 
465f80ce2aSJacob Faibussowitsch   row  = 1; col=1; CHKERRQ(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES));
475f80ce2aSJacob Faibussowitsch   row  = 1; col=2; CHKERRQ(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES));
485f80ce2aSJacob Faibussowitsch   row  = 1; col=3; CHKERRQ(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(R,MAT_FINAL_ASSEMBLY));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(R,MAT_FINAL_ASSEMBLY));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOptionsPrefix(R,"R_"));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(R,PETSC_VIEWER_STDOUT_WORLD));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n"));
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   /* C = A*R^T */
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatTransposeMult(A,R,MAT_INITIAL_MATRIX,2.0,&C));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOptionsPrefix(C,"ARt_"));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n"));
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* Create MatTransposeColoring from symbolic C=A*R^T */
625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatColoringCreate(C,&mc));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatColoringSetDistance(mc,2));
645f80ce2aSJacob Faibussowitsch   /* CHKERRQ(MatColoringSetType(mc,MATCOLORINGSL)); */
655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatColoringSetFromOptions(mc));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatColoringApply(mc,&iscoloring));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatColoringDestroy(&mc));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTransposeColoringCreate(C,iscoloring,&matcoloring));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(ISColoringDestroy(&iscoloring));
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /* Create Rt_dense */
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&Rt_dense));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(Rt_dense,4,matcoloring->ncolors,PETSC_DECIDE,PETSC_DECIDE));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(Rt_dense,MATDENSE));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqDenseSetPreallocation(Rt_dense,NULL));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(Rt_dense,MAT_FINAL_ASSEMBLY));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(Rt_dense,MAT_FINAL_ASSEMBLY));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(Rt_dense,&m,&n));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Rt_dense: %" PetscInt_FMT ",%" PetscInt_FMT "\n",m,n));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* Get Rt_dense by Apply MatTransposeColoring to R */
825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTransColoringApplySpToDen(matcoloring,R,Rt_dense));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /* C_dense = A*Rt_dense */
855f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,Rt_dense,MAT_INITIAL_MATRIX,2.0,&C_dense));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOptionsPrefix(C_dense,"ARt_dense_"));
875f80ce2aSJacob Faibussowitsch   /*CHKERRQ(MatView(C_dense,PETSC_VIEWER_STDOUT_WORLD)); */
885f80ce2aSJacob Faibussowitsch   /*CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n")); */
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   /* Recover C from C_dense */
915f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&C_sparse));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTransColoringApplyDenToSp(matcoloring,C_dense,C_sparse));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOptionsPrefix(C_sparse,"ARt_color_"));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C_sparse,PETSC_VIEWER_STDOUT_WORLD));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n"));
96c4762a1bSJed Brown 
975f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C_dense));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C_sparse));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Rt_dense));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTransposeColoringDestroy(&matcoloring));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   /* Test PtAP = P^T*A*P, P = R^T */
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(R,MAT_INITIAL_MATRIX,&P));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPtAP(A,P,MAT_INITIAL_MATRIX,2.0,&PtAP));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOptionsPrefix(PtAP,"PtAP_"));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(PtAP,PETSC_VIEWER_STDOUT_WORLD));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&P));
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /* Test C = RARt */
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(MatRARt(A,R,MAT_INITIAL_MATRIX,2.0,&C));
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(MatRARt(A,R,MAT_REUSE_MATRIX,2.0,&C));
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(MatEqual(C,PtAP,&equal));
11428b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: PtAP != RARt");
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* Free spaces */
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&R));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&PtAP));
121*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
122*b122ec5aSJacob Faibussowitsch   return 0;
123c4762a1bSJed Brown }
124c4762a1bSJed Brown 
125c4762a1bSJed Brown /*TEST
126c4762a1bSJed Brown 
127c4762a1bSJed Brown    test:
128c4762a1bSJed Brown       output_file: output/ex161.out
129c4762a1bSJed Brown    test:
130c4762a1bSJed Brown       suffix: 2
131c20d7725SJed Brown       args: -matmattransmult_via color
132c4762a1bSJed Brown       output_file: output/ex161.out
133c4762a1bSJed Brown 
134c4762a1bSJed Brown    test:
135c4762a1bSJed Brown       suffix: 3
136c20d7725SJed Brown       args: -matmattransmult_via color -matden2sp_brows 3
137c4762a1bSJed Brown       output_file: output/ex161.out
138c4762a1bSJed Brown 
139c4762a1bSJed Brown    test:
140c4762a1bSJed Brown       suffix: 4
141c20d7725SJed Brown       args: -matmattransmult_via color -matrart_via r*art
142c4762a1bSJed Brown       output_file: output/ex161.out
143c4762a1bSJed Brown 
144c4762a1bSJed Brown    test:
145c4762a1bSJed Brown       suffix: 5
146c20d7725SJed Brown       args: -matmattransmult_via color -matrart_via coloring_rart
147c4762a1bSJed Brown       output_file: output/ex161.out
148c4762a1bSJed Brown 
149c4762a1bSJed Brown TEST*/
150