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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 18*9566063dSJacob Faibussowitsch PetscCallMPI(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 */ 22*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&A)); 23*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,4,4,4,4)); 24*9566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSEQAIJ)); 25*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 26*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 27*9566063dSJacob Faibussowitsch row = 0; col=0; val=1.0; PetscCall(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES)); 28*9566063dSJacob Faibussowitsch row = 1; col=3; val=2.0; PetscCall(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES)); 29*9566063dSJacob Faibussowitsch row = 2; col=2; val=3.0; PetscCall(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES)); 30*9566063dSJacob Faibussowitsch row = 3; col=0; val=4.0; PetscCall(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES)); 31*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 32*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 33*9566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(A,"A_")); 34*9566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 35*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n")); 36c4762a1bSJed Brown 37c4762a1bSJed Brown /* Create seqaij R */ 38*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&R)); 39*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(R,2,4,2,4)); 40*9566063dSJacob Faibussowitsch PetscCall(MatSetType(R,MATSEQAIJ)); 41*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(R)); 42*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(R)); 43*9566063dSJacob Faibussowitsch row = 0; col=0; PetscCall(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES)); 44*9566063dSJacob Faibussowitsch row = 0; col=1; PetscCall(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES)); 45c4762a1bSJed Brown 46*9566063dSJacob Faibussowitsch row = 1; col=1; PetscCall(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES)); 47*9566063dSJacob Faibussowitsch row = 1; col=2; PetscCall(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES)); 48*9566063dSJacob Faibussowitsch row = 1; col=3; PetscCall(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES)); 49*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(R,MAT_FINAL_ASSEMBLY)); 50*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(R,MAT_FINAL_ASSEMBLY)); 51*9566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(R,"R_")); 52*9566063dSJacob Faibussowitsch PetscCall(MatView(R,PETSC_VIEWER_STDOUT_WORLD)); 53*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n")); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* C = A*R^T */ 56*9566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(A,R,MAT_INITIAL_MATRIX,2.0,&C)); 57*9566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(C,"ARt_")); 58*9566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 59*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n")); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* Create MatTransposeColoring from symbolic C=A*R^T */ 62*9566063dSJacob Faibussowitsch PetscCall(MatColoringCreate(C,&mc)); 63*9566063dSJacob Faibussowitsch PetscCall(MatColoringSetDistance(mc,2)); 64*9566063dSJacob Faibussowitsch /* PetscCall(MatColoringSetType(mc,MATCOLORINGSL)); */ 65*9566063dSJacob Faibussowitsch PetscCall(MatColoringSetFromOptions(mc)); 66*9566063dSJacob Faibussowitsch PetscCall(MatColoringApply(mc,&iscoloring)); 67*9566063dSJacob Faibussowitsch PetscCall(MatColoringDestroy(&mc)); 68*9566063dSJacob Faibussowitsch PetscCall(MatTransposeColoringCreate(C,iscoloring,&matcoloring)); 69*9566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* Create Rt_dense */ 72*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&Rt_dense)); 73*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Rt_dense,4,matcoloring->ncolors,PETSC_DECIDE,PETSC_DECIDE)); 74*9566063dSJacob Faibussowitsch PetscCall(MatSetType(Rt_dense,MATDENSE)); 75*9566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(Rt_dense,NULL)); 76*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Rt_dense,MAT_FINAL_ASSEMBLY)); 77*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Rt_dense,MAT_FINAL_ASSEMBLY)); 78*9566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Rt_dense,&m,&n)); 79*9566063dSJacob Faibussowitsch PetscCall(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 */ 82*9566063dSJacob Faibussowitsch PetscCall(MatTransColoringApplySpToDen(matcoloring,R,Rt_dense)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* C_dense = A*Rt_dense */ 85*9566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,Rt_dense,MAT_INITIAL_MATRIX,2.0,&C_dense)); 86*9566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(C_dense,"ARt_dense_")); 87*9566063dSJacob Faibussowitsch /*PetscCall(MatView(C_dense,PETSC_VIEWER_STDOUT_WORLD)); */ 88*9566063dSJacob Faibussowitsch /*PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n")); */ 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* Recover C from C_dense */ 91*9566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&C_sparse)); 92*9566063dSJacob Faibussowitsch PetscCall(MatTransColoringApplyDenToSp(matcoloring,C_dense,C_sparse)); 93*9566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(C_sparse,"ARt_color_")); 94*9566063dSJacob Faibussowitsch PetscCall(MatView(C_sparse,PETSC_VIEWER_STDOUT_WORLD)); 95*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n")); 96c4762a1bSJed Brown 97*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C_dense)); 98*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C_sparse)); 99*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Rt_dense)); 100*9566063dSJacob Faibussowitsch PetscCall(MatTransposeColoringDestroy(&matcoloring)); 101*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Test PtAP = P^T*A*P, P = R^T */ 104*9566063dSJacob Faibussowitsch PetscCall(MatTranspose(R,MAT_INITIAL_MATRIX,&P)); 105*9566063dSJacob Faibussowitsch PetscCall(MatPtAP(A,P,MAT_INITIAL_MATRIX,2.0,&PtAP)); 106*9566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(PtAP,"PtAP_")); 107*9566063dSJacob Faibussowitsch PetscCall(MatView(PtAP,PETSC_VIEWER_STDOUT_WORLD)); 108*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* Test C = RARt */ 111*9566063dSJacob Faibussowitsch PetscCall(MatRARt(A,R,MAT_INITIAL_MATRIX,2.0,&C)); 112*9566063dSJacob Faibussowitsch PetscCall(MatRARt(A,R,MAT_REUSE_MATRIX,2.0,&C)); 113*9566063dSJacob Faibussowitsch PetscCall(MatEqual(C,PtAP,&equal)); 11428b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: PtAP != RARt"); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* Free spaces */ 117*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 118*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 119*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&R)); 120*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&PtAP)); 121*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 122b122ec5aSJacob 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