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 6d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 7d71ae5a4SJacob Faibussowitsch { 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 17327415f7SBarry Smith PetscFunctionBeginUser; 18*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 20be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* Create seqaij A */ 239566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &A)); 249566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, 4, 4, 4, 4)); 259566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSEQAIJ)); 269566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 279566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 289371c9d4SSatish Balay row = 0; 299371c9d4SSatish Balay col = 0; 309371c9d4SSatish Balay val = 1.0; 319371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &row, 1, &col, &val, ADD_VALUES)); 329371c9d4SSatish Balay row = 1; 339371c9d4SSatish Balay col = 3; 349371c9d4SSatish Balay val = 2.0; 359371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &row, 1, &col, &val, ADD_VALUES)); 369371c9d4SSatish Balay row = 2; 379371c9d4SSatish Balay col = 2; 389371c9d4SSatish Balay val = 3.0; 399371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &row, 1, &col, &val, ADD_VALUES)); 409371c9d4SSatish Balay row = 3; 419371c9d4SSatish Balay col = 0; 429371c9d4SSatish Balay val = 4.0; 439371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &row, 1, &col, &val, ADD_VALUES)); 449566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 459566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 469566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(A, "A_")); 479566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 489566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* Create seqaij R */ 519566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &R)); 529566063dSJacob Faibussowitsch PetscCall(MatSetSizes(R, 2, 4, 2, 4)); 539566063dSJacob Faibussowitsch PetscCall(MatSetType(R, MATSEQAIJ)); 549566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(R)); 559566063dSJacob Faibussowitsch PetscCall(MatSetUp(R)); 569371c9d4SSatish Balay row = 0; 579371c9d4SSatish Balay col = 0; 589371c9d4SSatish Balay PetscCall(MatSetValues(R, 1, &row, 1, &col, &one, ADD_VALUES)); 599371c9d4SSatish Balay row = 0; 609371c9d4SSatish Balay col = 1; 619371c9d4SSatish Balay PetscCall(MatSetValues(R, 1, &row, 1, &col, &one, ADD_VALUES)); 62c4762a1bSJed Brown 639371c9d4SSatish Balay row = 1; 649371c9d4SSatish Balay col = 1; 659371c9d4SSatish Balay PetscCall(MatSetValues(R, 1, &row, 1, &col, &one, ADD_VALUES)); 669371c9d4SSatish Balay row = 1; 679371c9d4SSatish Balay col = 2; 689371c9d4SSatish Balay PetscCall(MatSetValues(R, 1, &row, 1, &col, &one, ADD_VALUES)); 699371c9d4SSatish Balay row = 1; 709371c9d4SSatish Balay col = 3; 719371c9d4SSatish Balay PetscCall(MatSetValues(R, 1, &row, 1, &col, &one, ADD_VALUES)); 729566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(R, MAT_FINAL_ASSEMBLY)); 739566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(R, MAT_FINAL_ASSEMBLY)); 749566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(R, "R_")); 759566063dSJacob Faibussowitsch PetscCall(MatView(R, PETSC_VIEWER_STDOUT_WORLD)); 769566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* C = A*R^T */ 799566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(A, R, MAT_INITIAL_MATRIX, 2.0, &C)); 809566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(C, "ARt_")); 819566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 829566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* Create MatTransposeColoring from symbolic C=A*R^T */ 859566063dSJacob Faibussowitsch PetscCall(MatColoringCreate(C, &mc)); 869566063dSJacob Faibussowitsch PetscCall(MatColoringSetDistance(mc, 2)); 879566063dSJacob Faibussowitsch /* PetscCall(MatColoringSetType(mc,MATCOLORINGSL)); */ 889566063dSJacob Faibussowitsch PetscCall(MatColoringSetFromOptions(mc)); 899566063dSJacob Faibussowitsch PetscCall(MatColoringApply(mc, &iscoloring)); 909566063dSJacob Faibussowitsch PetscCall(MatColoringDestroy(&mc)); 919566063dSJacob Faibussowitsch PetscCall(MatTransposeColoringCreate(C, iscoloring, &matcoloring)); 929566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* Create Rt_dense */ 959566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Rt_dense)); 969566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Rt_dense, 4, matcoloring->ncolors, PETSC_DECIDE, PETSC_DECIDE)); 979566063dSJacob Faibussowitsch PetscCall(MatSetType(Rt_dense, MATDENSE)); 989566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(Rt_dense, NULL)); 999566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Rt_dense, MAT_FINAL_ASSEMBLY)); 1009566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Rt_dense, MAT_FINAL_ASSEMBLY)); 1019566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Rt_dense, &m, &n)); 1029566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Rt_dense: %" PetscInt_FMT ",%" PetscInt_FMT "\n", m, n)); 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* Get Rt_dense by Apply MatTransposeColoring to R */ 1059566063dSJacob Faibussowitsch PetscCall(MatTransColoringApplySpToDen(matcoloring, R, Rt_dense)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* C_dense = A*Rt_dense */ 1089566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, Rt_dense, MAT_INITIAL_MATRIX, 2.0, &C_dense)); 1099566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(C_dense, "ARt_dense_")); 1109566063dSJacob Faibussowitsch /*PetscCall(MatView(C_dense,PETSC_VIEWER_STDOUT_WORLD)); */ 1119566063dSJacob Faibussowitsch /*PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n")); */ 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* Recover C from C_dense */ 1149566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_DO_NOT_COPY_VALUES, &C_sparse)); 1159566063dSJacob Faibussowitsch PetscCall(MatTransColoringApplyDenToSp(matcoloring, C_dense, C_sparse)); 1169566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(C_sparse, "ARt_color_")); 1179566063dSJacob Faibussowitsch PetscCall(MatView(C_sparse, PETSC_VIEWER_STDOUT_WORLD)); 1189566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 119c4762a1bSJed Brown 1209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C_dense)); 1219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C_sparse)); 1229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Rt_dense)); 1239566063dSJacob Faibussowitsch PetscCall(MatTransposeColoringDestroy(&matcoloring)); 1249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* Test PtAP = P^T*A*P, P = R^T */ 1279566063dSJacob Faibussowitsch PetscCall(MatTranspose(R, MAT_INITIAL_MATRIX, &P)); 1289566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_INITIAL_MATRIX, 2.0, &PtAP)); 1299566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(PtAP, "PtAP_")); 1309566063dSJacob Faibussowitsch PetscCall(MatView(PtAP, PETSC_VIEWER_STDOUT_WORLD)); 1319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* Test C = RARt */ 1349566063dSJacob Faibussowitsch PetscCall(MatRARt(A, R, MAT_INITIAL_MATRIX, 2.0, &C)); 1359566063dSJacob Faibussowitsch PetscCall(MatRARt(A, R, MAT_REUSE_MATRIX, 2.0, &C)); 1369566063dSJacob Faibussowitsch PetscCall(MatEqual(C, PtAP, &equal)); 13728b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: PtAP != RARt"); 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* Free spaces */ 1409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&R)); 1439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&PtAP)); 1449566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 145b122ec5aSJacob Faibussowitsch return 0; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown 148c4762a1bSJed Brown /*TEST 149c4762a1bSJed Brown 150c4762a1bSJed Brown test: 151c4762a1bSJed Brown output_file: output/ex161.out 152c4762a1bSJed Brown test: 153c4762a1bSJed Brown suffix: 2 154c20d7725SJed Brown args: -matmattransmult_via color 155c4762a1bSJed Brown output_file: output/ex161.out 156c4762a1bSJed Brown 157c4762a1bSJed Brown test: 158c4762a1bSJed Brown suffix: 3 159c20d7725SJed Brown args: -matmattransmult_via color -matden2sp_brows 3 160c4762a1bSJed Brown output_file: output/ex161.out 161c4762a1bSJed Brown 162c4762a1bSJed Brown test: 163c4762a1bSJed Brown suffix: 4 164c20d7725SJed Brown args: -matmattransmult_via color -matrart_via r*art 165c4762a1bSJed Brown output_file: output/ex161.out 166c4762a1bSJed Brown 167c4762a1bSJed Brown test: 168c4762a1bSJed Brown suffix: 5 169c20d7725SJed Brown args: -matmattransmult_via color -matrart_via coloring_rart 170c4762a1bSJed Brown output_file: output/ex161.out 171c4762a1bSJed Brown 172c4762a1bSJed Brown TEST*/ 173