1c4762a1bSJed Brown /* 2c4762a1bSJed Brown * ex195.c 3c4762a1bSJed Brown * 4c4762a1bSJed Brown * Created on: Aug 24, 2015 5c4762a1bSJed Brown * Author: Fande Kong <fdkong.jd@gmail.com> 6c4762a1bSJed Brown */ 7c4762a1bSJed Brown 8c4762a1bSJed Brown static char help[] = " Demonstrate the use of MatConvert_Nest_AIJ\n"; 9c4762a1bSJed Brown 10c4762a1bSJed Brown #include <petscmat.h> 11c4762a1bSJed Brown 12d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 13d71ae5a4SJacob Faibussowitsch { 14c20d7725SJed Brown Mat A1, A2, A3, A4, A5, B, C, C1, nest; 15c4762a1bSJed Brown Mat aij; 16c4762a1bSJed Brown MPI_Comm comm; 17c4762a1bSJed Brown PetscInt m, M, n, istart, iend, ii, i, J, j, K = 10; 18c4762a1bSJed Brown PetscScalar v; 19c4762a1bSJed Brown PetscMPIInt size; 203536838dSStefano Zampini PetscBool equal, test = PETSC_FALSE, test_null = PETSC_FALSE; 21c4762a1bSJed Brown 22327415f7SBarry Smith PetscFunctionBeginUser; 23*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 24c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 26c20d7725SJed Brown 27c4762a1bSJed Brown /* 28c4762a1bSJed Brown Assemble the matrix for the five point stencil, YET AGAIN 29c4762a1bSJed Brown */ 309566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &A1)); 31c4762a1bSJed Brown m = 2, n = 2; 329566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A1, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n)); 339566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A1)); 349566063dSJacob Faibussowitsch PetscCall(MatSetUp(A1)); 359566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A1, &istart, &iend)); 36c4762a1bSJed Brown for (ii = istart; ii < iend; ii++) { 379371c9d4SSatish Balay v = -1.0; 389371c9d4SSatish Balay i = ii / n; 399371c9d4SSatish Balay j = ii - i * n; 409371c9d4SSatish Balay if (i > 0) { 419371c9d4SSatish Balay J = ii - n; 429371c9d4SSatish Balay PetscCall(MatSetValues(A1, 1, &ii, 1, &J, &v, INSERT_VALUES)); 439371c9d4SSatish Balay } 449371c9d4SSatish Balay if (i < m - 1) { 459371c9d4SSatish Balay J = ii + n; 469371c9d4SSatish Balay PetscCall(MatSetValues(A1, 1, &ii, 1, &J, &v, INSERT_VALUES)); 479371c9d4SSatish Balay } 489371c9d4SSatish Balay if (j > 0) { 499371c9d4SSatish Balay J = ii - 1; 509371c9d4SSatish Balay PetscCall(MatSetValues(A1, 1, &ii, 1, &J, &v, INSERT_VALUES)); 519371c9d4SSatish Balay } 529371c9d4SSatish Balay if (j < n - 1) { 539371c9d4SSatish Balay J = ii + 1; 549371c9d4SSatish Balay PetscCall(MatSetValues(A1, 1, &ii, 1, &J, &v, INSERT_VALUES)); 559371c9d4SSatish Balay } 569371c9d4SSatish Balay v = 4.0; 579371c9d4SSatish Balay PetscCall(MatSetValues(A1, 1, &ii, 1, &ii, &v, INSERT_VALUES)); 58c4762a1bSJed Brown } 599566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A1, MAT_FINAL_ASSEMBLY)); 609566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A1, MAT_FINAL_ASSEMBLY)); 619566063dSJacob Faibussowitsch PetscCall(MatView(A1, PETSC_VIEWER_STDOUT_WORLD)); 62c20d7725SJed Brown 639566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A1, MAT_COPY_VALUES, &A2)); 649566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A1, MAT_COPY_VALUES, &A3)); 659566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A1, MAT_COPY_VALUES, &A4)); 66c20d7725SJed Brown 67c4762a1bSJed Brown /* create a nest matrix */ 689566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &nest)); 699566063dSJacob Faibussowitsch PetscCall(MatSetType(nest, MATNEST)); 703536838dSStefano Zampini PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_null", &test_null, NULL)); 713536838dSStefano Zampini if (test_null) { 723536838dSStefano Zampini IS is[2]; 733536838dSStefano Zampini PetscInt st; 743536838dSStefano Zampini 753536838dSStefano Zampini PetscCall(MatGetLocalSize(A1, &m, NULL)); 763536838dSStefano Zampini PetscCall(MatGetOwnershipRange(A1, &st, NULL)); 773536838dSStefano Zampini PetscCall(ISCreateStride(comm, m, st, 2, &is[0])); 783536838dSStefano Zampini PetscCall(ISCreateStride(comm, m, st + 1, 2, &is[1])); 793536838dSStefano Zampini PetscCall(MatNestSetSubMats(nest, 2, is, 2, is, NULL)); 803536838dSStefano Zampini PetscCall(ISDestroy(&is[0])); 813536838dSStefano Zampini PetscCall(ISDestroy(&is[1])); 823536838dSStefano Zampini } else { 833536838dSStefano Zampini Mat mata[] = {A1, A2, A3, A4}; 843536838dSStefano Zampini 859566063dSJacob Faibussowitsch PetscCall(MatNestSetSubMats(nest, 2, NULL, 2, NULL, mata)); 863536838dSStefano Zampini } 879566063dSJacob Faibussowitsch PetscCall(MatSetUp(nest)); 88a0228903SHong Zhang 89a0228903SHong Zhang /* test matrix product error messages */ 90a0228903SHong Zhang PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_nest*nest", &test, NULL)); 91fb842aefSJose E. Roman if (test) PetscCall(MatMatMult(nest, nest, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &C)); 92a0228903SHong Zhang 93a0228903SHong Zhang PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_matproductsetfromoptions", &test, NULL)); 94a0228903SHong Zhang if (test) { 95a0228903SHong Zhang PetscCall(MatProductCreate(nest, nest, NULL, &C)); 96a0228903SHong Zhang PetscCall(MatProductSetType(C, MATPRODUCT_AB)); 97a0228903SHong Zhang PetscCall(MatProductSymbolic(C)); 98a0228903SHong Zhang } 99a0228903SHong Zhang 100a0228903SHong Zhang PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_matproductsymbolic", &test, NULL)); 101a0228903SHong Zhang if (test) { 102a0228903SHong Zhang PetscCall(MatProductCreate(nest, nest, NULL, &C)); 103a0228903SHong Zhang PetscCall(MatProductSetType(C, MATPRODUCT_AB)); 104a0228903SHong Zhang PetscCall(MatProductSetFromOptions(C)); 105a0228903SHong Zhang PetscCall(MatProductSymbolic(C)); 106a0228903SHong Zhang } 107a0228903SHong Zhang 1089566063dSJacob Faibussowitsch PetscCall(MatConvert(nest, MATAIJ, MAT_INITIAL_MATRIX, &aij)); 1099566063dSJacob Faibussowitsch PetscCall(MatView(aij, PETSC_VIEWER_STDOUT_WORLD)); 110c20d7725SJed Brown 111c20d7725SJed Brown /* create a dense matrix */ 1129566063dSJacob Faibussowitsch PetscCall(MatGetSize(nest, &M, NULL)); 1139566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(nest, &m, NULL)); 1149566063dSJacob Faibussowitsch PetscCall(MatCreateDense(comm, m, PETSC_DECIDE, M, K, NULL, &B)); 115f3fa974cSJacob Faibussowitsch PetscCall(MatSetRandom(B, NULL)); 116c20d7725SJed Brown 117c20d7725SJed Brown /* C = nest*B_dense */ 118fb842aefSJose E. Roman PetscCall(MatMatMult(nest, B, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &C)); 119fb842aefSJose E. Roman PetscCall(MatMatMult(nest, B, MAT_REUSE_MATRIX, PETSC_DETERMINE, &C)); 1209566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(nest, B, C, 10, &equal)); 12128b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in C != nest*B_dense"); 1223fbca975SStefano Zampini 1233fbca975SStefano Zampini /* Test B = nest*C, reuse C and B with MatProductCreateWithMat() */ 1243fbca975SStefano Zampini /* C has been obtained from nest*B. Clear internal data structures related to factors to prevent circular references */ 1259566063dSJacob Faibussowitsch PetscCall(MatProductClear(C)); 1269566063dSJacob Faibussowitsch PetscCall(MatProductCreateWithMat(nest, C, NULL, B)); 1279566063dSJacob Faibussowitsch PetscCall(MatProductSetType(B, MATPRODUCT_AB)); 1289566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(B)); 1299566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(B)); 1309566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(B)); 1319566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(nest, C, B, 10, &equal)); 13228b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in B != nest*C_dense"); 1339566063dSJacob Faibussowitsch PetscCall(MatConvert(nest, MATAIJ, MAT_INPLACE_MATRIX, &nest)); 1349566063dSJacob Faibussowitsch PetscCall(MatEqual(nest, aij, &equal)); 13528b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in aij != nest"); 136c20d7725SJed Brown 1373536838dSStefano Zampini /* Test with virtual block */ 1389566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(A1, &A5)); /* A1 is symmetric */ 1393536838dSStefano Zampini PetscCall(MatNestSetSubMat(nest, 0, 0, A5)); 140fb842aefSJose E. Roman PetscCall(MatMatMult(nest, B, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &C1)); 141fb842aefSJose E. Roman PetscCall(MatMatMult(nest, B, MAT_REUSE_MATRIX, PETSC_DETERMINE, &C1)); 1429566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(nest, B, C1, 10, &equal)); 14328b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in C1 != C"); 1449566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1)); 1459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A5)); 146c20d7725SJed Brown 1473536838dSStefano Zampini PetscCall(MatDestroy(&nest)); 1489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&aij)); 1519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A1)); 1529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 1539566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A3)); 1549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A4)); 1559566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 156b122ec5aSJacob Faibussowitsch return 0; 157c4762a1bSJed Brown } 158c4762a1bSJed Brown 159c4762a1bSJed Brown /*TEST 160c4762a1bSJed Brown 161c4762a1bSJed Brown test: 1623536838dSStefano Zampini diff_args: -j 163c4762a1bSJed Brown nsize: 2 1643536838dSStefano Zampini suffix: 1 165c4762a1bSJed Brown 1663fbca975SStefano Zampini test: 1673536838dSStefano Zampini diff_args: -j 1683536838dSStefano Zampini nsize: 2 1693536838dSStefano Zampini suffix: 1_null 1703536838dSStefano Zampini args: -test_null 1713536838dSStefano Zampini 1723536838dSStefano Zampini test: 1733536838dSStefano Zampini diff_args: -j 1743fbca975SStefano Zampini suffix: 2 1753fbca975SStefano Zampini 1763536838dSStefano Zampini test: 1773536838dSStefano Zampini diff_args: -j 1783536838dSStefano Zampini suffix: 2_null 1793536838dSStefano Zampini args: -test_null 1803536838dSStefano Zampini 181c4762a1bSJed Brown TEST*/ 182