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 12*9371c9d4SSatish Balay int main(int argc, char **args) { 13c20d7725SJed Brown Mat A1, A2, A3, A4, A5, B, C, C1, nest; 14c4762a1bSJed Brown Mat mata[4]; 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; 20c20d7725SJed Brown PetscBool equal; 21c4762a1bSJed Brown 22327415f7SBarry Smith PetscFunctionBeginUser; 239566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, 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++) { 37*9371c9d4SSatish Balay v = -1.0; 38*9371c9d4SSatish Balay i = ii / n; 39*9371c9d4SSatish Balay j = ii - i * n; 40*9371c9d4SSatish Balay if (i > 0) { 41*9371c9d4SSatish Balay J = ii - n; 42*9371c9d4SSatish Balay PetscCall(MatSetValues(A1, 1, &ii, 1, &J, &v, INSERT_VALUES)); 43*9371c9d4SSatish Balay } 44*9371c9d4SSatish Balay if (i < m - 1) { 45*9371c9d4SSatish Balay J = ii + n; 46*9371c9d4SSatish Balay PetscCall(MatSetValues(A1, 1, &ii, 1, &J, &v, INSERT_VALUES)); 47*9371c9d4SSatish Balay } 48*9371c9d4SSatish Balay if (j > 0) { 49*9371c9d4SSatish Balay J = ii - 1; 50*9371c9d4SSatish Balay PetscCall(MatSetValues(A1, 1, &ii, 1, &J, &v, INSERT_VALUES)); 51*9371c9d4SSatish Balay } 52*9371c9d4SSatish Balay if (j < n - 1) { 53*9371c9d4SSatish Balay J = ii + 1; 54*9371c9d4SSatish Balay PetscCall(MatSetValues(A1, 1, &ii, 1, &J, &v, INSERT_VALUES)); 55*9371c9d4SSatish Balay } 56*9371c9d4SSatish Balay v = 4.0; 57*9371c9d4SSatish 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)); 70c4762a1bSJed Brown mata[0] = A1, mata[1] = A2, mata[2] = A3, mata[3] = A4; 719566063dSJacob Faibussowitsch PetscCall(MatNestSetSubMats(nest, 2, NULL, 2, NULL, mata)); 729566063dSJacob Faibussowitsch PetscCall(MatSetUp(nest)); 739566063dSJacob Faibussowitsch PetscCall(MatConvert(nest, MATAIJ, MAT_INITIAL_MATRIX, &aij)); 749566063dSJacob Faibussowitsch PetscCall(MatView(aij, PETSC_VIEWER_STDOUT_WORLD)); 75c20d7725SJed Brown 76c20d7725SJed Brown /* create a dense matrix */ 779566063dSJacob Faibussowitsch PetscCall(MatGetSize(nest, &M, NULL)); 789566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(nest, &m, NULL)); 799566063dSJacob Faibussowitsch PetscCall(MatCreateDense(comm, m, PETSC_DECIDE, M, K, NULL, &B)); 809566063dSJacob Faibussowitsch PetscCall(MatSetRandom(B, PETSC_NULL)); 81c20d7725SJed Brown 82c20d7725SJed Brown /* C = nest*B_dense */ 839566063dSJacob Faibussowitsch PetscCall(MatMatMult(nest, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C)); 849566063dSJacob Faibussowitsch PetscCall(MatMatMult(nest, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C)); 859566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(nest, B, C, 10, &equal)); 8628b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in C != nest*B_dense"); 873fbca975SStefano Zampini 883fbca975SStefano Zampini /* Test B = nest*C, reuse C and B with MatProductCreateWithMat() */ 893fbca975SStefano Zampini /* C has been obtained from nest*B. Clear internal data structures related to factors to prevent circular references */ 909566063dSJacob Faibussowitsch PetscCall(MatProductClear(C)); 919566063dSJacob Faibussowitsch PetscCall(MatProductCreateWithMat(nest, C, NULL, B)); 929566063dSJacob Faibussowitsch PetscCall(MatProductSetType(B, MATPRODUCT_AB)); 939566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(B)); 949566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(B)); 959566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(B)); 969566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(nest, C, B, 10, &equal)); 9728b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in B != nest*C_dense"); 989566063dSJacob Faibussowitsch PetscCall(MatConvert(nest, MATAIJ, MAT_INPLACE_MATRIX, &nest)); 999566063dSJacob Faibussowitsch PetscCall(MatEqual(nest, aij, &equal)); 10028b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in aij != nest"); 1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&nest)); 102c20d7725SJed Brown 103c20d7725SJed Brown if (size > 1) { /* Do not know why this test fails for size = 1 */ 1049566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(A1, &A5)); /* A1 is symmetric */ 105c20d7725SJed Brown mata[0] = A5; 1069566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &nest)); 1079566063dSJacob Faibussowitsch PetscCall(MatSetType(nest, MATNEST)); 1089566063dSJacob Faibussowitsch PetscCall(MatNestSetSubMats(nest, 2, NULL, 2, NULL, mata)); 1099566063dSJacob Faibussowitsch PetscCall(MatSetUp(nest)); 1109566063dSJacob Faibussowitsch PetscCall(MatMatMult(nest, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C1)); 1119566063dSJacob Faibussowitsch PetscCall(MatMatMult(nest, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C1)); 112c20d7725SJed Brown 1139566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(nest, B, C1, 10, &equal)); 11428b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in C1 != C"); 1159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1)); 1169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A5)); 1179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&nest)); 118c20d7725SJed Brown } 119c20d7725SJed Brown 1209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&aij)); 1239566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A1)); 1249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 1259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A3)); 1269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A4)); 1279566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 128b122ec5aSJacob Faibussowitsch return 0; 129c4762a1bSJed Brown } 130c4762a1bSJed Brown 131c4762a1bSJed Brown /*TEST 132c4762a1bSJed Brown 133c4762a1bSJed Brown test: 134c4762a1bSJed Brown nsize: 2 135c4762a1bSJed Brown 1363fbca975SStefano Zampini test: 1373fbca975SStefano Zampini suffix: 2 1383fbca975SStefano Zampini 139c4762a1bSJed Brown TEST*/ 140