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 6*9371c9d4SSatish Balay int main(int argc, char **argv) { 7c4762a1bSJed Brown Mat A, R, C, C_dense, C_sparse, Rt_dense, P, PtAP; 8c4762a1bSJed Brown PetscInt row, col, m, n; 9c4762a1bSJed Brown MatScalar one = 1.0, val; 10c4762a1bSJed Brown MatColoring mc; 11c4762a1bSJed Brown MatTransposeColoring matcoloring = 0; 12c4762a1bSJed Brown ISColoring iscoloring; 13c4762a1bSJed Brown PetscBool equal; 14c4762a1bSJed Brown PetscMPIInt size; 15c4762a1bSJed Brown 16327415f7SBarry Smith PetscFunctionBeginUser; 179566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 19be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* Create seqaij A */ 229566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &A)); 239566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, 4, 4, 4, 4)); 249566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSEQAIJ)); 259566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 269566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 27*9371c9d4SSatish Balay row = 0; 28*9371c9d4SSatish Balay col = 0; 29*9371c9d4SSatish Balay val = 1.0; 30*9371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &row, 1, &col, &val, ADD_VALUES)); 31*9371c9d4SSatish Balay row = 1; 32*9371c9d4SSatish Balay col = 3; 33*9371c9d4SSatish Balay val = 2.0; 34*9371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &row, 1, &col, &val, ADD_VALUES)); 35*9371c9d4SSatish Balay row = 2; 36*9371c9d4SSatish Balay col = 2; 37*9371c9d4SSatish Balay val = 3.0; 38*9371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &row, 1, &col, &val, ADD_VALUES)); 39*9371c9d4SSatish Balay row = 3; 40*9371c9d4SSatish Balay col = 0; 41*9371c9d4SSatish Balay val = 4.0; 42*9371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &row, 1, &col, &val, ADD_VALUES)); 439566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 449566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 459566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(A, "A_")); 469566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 479566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* Create seqaij R */ 509566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &R)); 519566063dSJacob Faibussowitsch PetscCall(MatSetSizes(R, 2, 4, 2, 4)); 529566063dSJacob Faibussowitsch PetscCall(MatSetType(R, MATSEQAIJ)); 539566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(R)); 549566063dSJacob Faibussowitsch PetscCall(MatSetUp(R)); 55*9371c9d4SSatish Balay row = 0; 56*9371c9d4SSatish Balay col = 0; 57*9371c9d4SSatish Balay PetscCall(MatSetValues(R, 1, &row, 1, &col, &one, ADD_VALUES)); 58*9371c9d4SSatish Balay row = 0; 59*9371c9d4SSatish Balay col = 1; 60*9371c9d4SSatish Balay PetscCall(MatSetValues(R, 1, &row, 1, &col, &one, ADD_VALUES)); 61c4762a1bSJed Brown 62*9371c9d4SSatish Balay row = 1; 63*9371c9d4SSatish Balay col = 1; 64*9371c9d4SSatish Balay PetscCall(MatSetValues(R, 1, &row, 1, &col, &one, ADD_VALUES)); 65*9371c9d4SSatish Balay row = 1; 66*9371c9d4SSatish Balay col = 2; 67*9371c9d4SSatish Balay PetscCall(MatSetValues(R, 1, &row, 1, &col, &one, ADD_VALUES)); 68*9371c9d4SSatish Balay row = 1; 69*9371c9d4SSatish Balay col = 3; 70*9371c9d4SSatish Balay PetscCall(MatSetValues(R, 1, &row, 1, &col, &one, ADD_VALUES)); 719566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(R, MAT_FINAL_ASSEMBLY)); 729566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(R, MAT_FINAL_ASSEMBLY)); 739566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(R, "R_")); 749566063dSJacob Faibussowitsch PetscCall(MatView(R, PETSC_VIEWER_STDOUT_WORLD)); 759566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* C = A*R^T */ 789566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(A, R, MAT_INITIAL_MATRIX, 2.0, &C)); 799566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(C, "ARt_")); 809566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 819566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* Create MatTransposeColoring from symbolic C=A*R^T */ 849566063dSJacob Faibussowitsch PetscCall(MatColoringCreate(C, &mc)); 859566063dSJacob Faibussowitsch PetscCall(MatColoringSetDistance(mc, 2)); 869566063dSJacob Faibussowitsch /* PetscCall(MatColoringSetType(mc,MATCOLORINGSL)); */ 879566063dSJacob Faibussowitsch PetscCall(MatColoringSetFromOptions(mc)); 889566063dSJacob Faibussowitsch PetscCall(MatColoringApply(mc, &iscoloring)); 899566063dSJacob Faibussowitsch PetscCall(MatColoringDestroy(&mc)); 909566063dSJacob Faibussowitsch PetscCall(MatTransposeColoringCreate(C, iscoloring, &matcoloring)); 919566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* Create Rt_dense */ 949566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Rt_dense)); 959566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Rt_dense, 4, matcoloring->ncolors, PETSC_DECIDE, PETSC_DECIDE)); 969566063dSJacob Faibussowitsch PetscCall(MatSetType(Rt_dense, MATDENSE)); 979566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(Rt_dense, NULL)); 989566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Rt_dense, MAT_FINAL_ASSEMBLY)); 999566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Rt_dense, MAT_FINAL_ASSEMBLY)); 1009566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Rt_dense, &m, &n)); 1019566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Rt_dense: %" PetscInt_FMT ",%" PetscInt_FMT "\n", m, n)); 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Get Rt_dense by Apply MatTransposeColoring to R */ 1049566063dSJacob Faibussowitsch PetscCall(MatTransColoringApplySpToDen(matcoloring, R, Rt_dense)); 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* C_dense = A*Rt_dense */ 1079566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, Rt_dense, MAT_INITIAL_MATRIX, 2.0, &C_dense)); 1089566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(C_dense, "ARt_dense_")); 1099566063dSJacob Faibussowitsch /*PetscCall(MatView(C_dense,PETSC_VIEWER_STDOUT_WORLD)); */ 1109566063dSJacob Faibussowitsch /*PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n")); */ 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* Recover C from C_dense */ 1139566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_DO_NOT_COPY_VALUES, &C_sparse)); 1149566063dSJacob Faibussowitsch PetscCall(MatTransColoringApplyDenToSp(matcoloring, C_dense, C_sparse)); 1159566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(C_sparse, "ARt_color_")); 1169566063dSJacob Faibussowitsch PetscCall(MatView(C_sparse, PETSC_VIEWER_STDOUT_WORLD)); 1179566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 118c4762a1bSJed Brown 1199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C_dense)); 1209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C_sparse)); 1219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Rt_dense)); 1229566063dSJacob Faibussowitsch PetscCall(MatTransposeColoringDestroy(&matcoloring)); 1239566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* Test PtAP = P^T*A*P, P = R^T */ 1269566063dSJacob Faibussowitsch PetscCall(MatTranspose(R, MAT_INITIAL_MATRIX, &P)); 1279566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_INITIAL_MATRIX, 2.0, &PtAP)); 1289566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(PtAP, "PtAP_")); 1299566063dSJacob Faibussowitsch PetscCall(MatView(PtAP, PETSC_VIEWER_STDOUT_WORLD)); 1309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* Test C = RARt */ 1339566063dSJacob Faibussowitsch PetscCall(MatRARt(A, R, MAT_INITIAL_MATRIX, 2.0, &C)); 1349566063dSJacob Faibussowitsch PetscCall(MatRARt(A, R, MAT_REUSE_MATRIX, 2.0, &C)); 1359566063dSJacob Faibussowitsch PetscCall(MatEqual(C, PtAP, &equal)); 13628b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: PtAP != RARt"); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* Free spaces */ 1399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&R)); 1429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&PtAP)); 1439566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 144b122ec5aSJacob Faibussowitsch return 0; 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /*TEST 148c4762a1bSJed Brown 149c4762a1bSJed Brown test: 150c4762a1bSJed Brown output_file: output/ex161.out 151c4762a1bSJed Brown test: 152c4762a1bSJed Brown suffix: 2 153c20d7725SJed Brown args: -matmattransmult_via color 154c4762a1bSJed Brown output_file: output/ex161.out 155c4762a1bSJed Brown 156c4762a1bSJed Brown test: 157c4762a1bSJed Brown suffix: 3 158c20d7725SJed Brown args: -matmattransmult_via color -matden2sp_brows 3 159c4762a1bSJed Brown output_file: output/ex161.out 160c4762a1bSJed Brown 161c4762a1bSJed Brown test: 162c4762a1bSJed Brown suffix: 4 163c20d7725SJed Brown args: -matmattransmult_via color -matrart_via r*art 164c4762a1bSJed Brown output_file: output/ex161.out 165c4762a1bSJed Brown 166c4762a1bSJed Brown test: 167c4762a1bSJed Brown suffix: 5 168c20d7725SJed Brown args: -matmattransmult_via color -matrart_via coloring_rart 169c4762a1bSJed Brown output_file: output/ex161.out 170c4762a1bSJed Brown 171c4762a1bSJed Brown TEST*/ 172