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 12c4762a1bSJed Brown int main(int argc,char **args) 13c4762a1bSJed Brown { 14c20d7725SJed Brown Mat A1,A2,A3,A4,A5,B,C,C1,nest; 15c4762a1bSJed Brown Mat mata[4]; 16c4762a1bSJed Brown Mat aij; 17c4762a1bSJed Brown MPI_Comm comm; 18c4762a1bSJed Brown PetscInt m,M,n,istart,iend,ii,i,J,j,K=10; 19c4762a1bSJed Brown PetscScalar v; 20c4762a1bSJed Brown PetscMPIInt size; 21c20d7725SJed Brown PetscBool equal; 22c4762a1bSJed Brown 23*327415f7SBarry Smith PetscFunctionBeginUser; 249566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 25c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 27c20d7725SJed Brown 28c4762a1bSJed Brown /* 29c4762a1bSJed Brown Assemble the matrix for the five point stencil, YET AGAIN 30c4762a1bSJed Brown */ 319566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&A1)); 32c4762a1bSJed Brown m=2,n=2; 339566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A1,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 349566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A1)); 359566063dSJacob Faibussowitsch PetscCall(MatSetUp(A1)); 369566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A1,&istart,&iend)); 37c4762a1bSJed Brown for (ii=istart; ii<iend; ii++) { 38c4762a1bSJed Brown v = -1.0; i = ii/n; j = ii - i*n; 399566063dSJacob Faibussowitsch if (i>0) {J = ii - n; PetscCall(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 409566063dSJacob Faibussowitsch if (i<m-1) {J = ii + n; PetscCall(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 419566063dSJacob Faibussowitsch if (j>0) {J = ii - 1; PetscCall(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 429566063dSJacob Faibussowitsch if (j<n-1) {J = ii + 1; PetscCall(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 439566063dSJacob Faibussowitsch v = 4.0; PetscCall(MatSetValues(A1,1,&ii,1,&ii,&v,INSERT_VALUES)); 44c4762a1bSJed Brown } 459566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A1,MAT_FINAL_ASSEMBLY)); 469566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A1,MAT_FINAL_ASSEMBLY)); 479566063dSJacob Faibussowitsch PetscCall(MatView(A1,PETSC_VIEWER_STDOUT_WORLD)); 48c20d7725SJed Brown 499566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A1,MAT_COPY_VALUES,&A2)); 509566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A1,MAT_COPY_VALUES,&A3)); 519566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A1,MAT_COPY_VALUES,&A4)); 52c20d7725SJed Brown 53c4762a1bSJed Brown /*create a nest matrix */ 549566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&nest)); 559566063dSJacob Faibussowitsch PetscCall(MatSetType(nest,MATNEST)); 56c4762a1bSJed Brown mata[0]=A1,mata[1]=A2,mata[2]=A3,mata[3]=A4; 579566063dSJacob Faibussowitsch PetscCall(MatNestSetSubMats(nest,2,NULL,2,NULL,mata)); 589566063dSJacob Faibussowitsch PetscCall(MatSetUp(nest)); 599566063dSJacob Faibussowitsch PetscCall(MatConvert(nest,MATAIJ,MAT_INITIAL_MATRIX,&aij)); 609566063dSJacob Faibussowitsch PetscCall(MatView(aij,PETSC_VIEWER_STDOUT_WORLD)); 61c20d7725SJed Brown 62c20d7725SJed Brown /* create a dense matrix */ 639566063dSJacob Faibussowitsch PetscCall(MatGetSize(nest,&M,NULL)); 649566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(nest,&m,NULL)); 659566063dSJacob Faibussowitsch PetscCall(MatCreateDense(comm,m,PETSC_DECIDE,M,K,NULL,&B)); 669566063dSJacob Faibussowitsch PetscCall(MatSetRandom(B,PETSC_NULL)); 67c20d7725SJed Brown 68c20d7725SJed Brown /* C = nest*B_dense */ 699566063dSJacob Faibussowitsch PetscCall(MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 709566063dSJacob Faibussowitsch PetscCall(MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C)); 719566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(nest,B,C,10,&equal)); 7228b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in C != nest*B_dense"); 733fbca975SStefano Zampini 743fbca975SStefano Zampini /* Test B = nest*C, reuse C and B with MatProductCreateWithMat() */ 753fbca975SStefano Zampini /* C has been obtained from nest*B. Clear internal data structures related to factors to prevent circular references */ 769566063dSJacob Faibussowitsch PetscCall(MatProductClear(C)); 779566063dSJacob Faibussowitsch PetscCall(MatProductCreateWithMat(nest,C,NULL,B)); 789566063dSJacob Faibussowitsch PetscCall(MatProductSetType(B,MATPRODUCT_AB)); 799566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(B)); 809566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(B)); 819566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(B)); 829566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(nest,C,B,10,&equal)); 8328b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in B != nest*C_dense"); 849566063dSJacob Faibussowitsch PetscCall(MatConvert(nest,MATAIJ,MAT_INPLACE_MATRIX,&nest)); 859566063dSJacob Faibussowitsch PetscCall(MatEqual(nest,aij,&equal)); 8628b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in aij != nest"); 879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&nest)); 88c20d7725SJed Brown 89c20d7725SJed Brown if (size > 1) { /* Do not know why this test fails for size = 1 */ 909566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(A1,&A5)); /* A1 is symmetric */ 91c20d7725SJed Brown mata[0] = A5; 929566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&nest)); 939566063dSJacob Faibussowitsch PetscCall(MatSetType(nest,MATNEST)); 949566063dSJacob Faibussowitsch PetscCall(MatNestSetSubMats(nest,2,NULL,2,NULL,mata)); 959566063dSJacob Faibussowitsch PetscCall(MatSetUp(nest)); 969566063dSJacob Faibussowitsch PetscCall(MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C1)); 979566063dSJacob Faibussowitsch PetscCall(MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C1)); 98c20d7725SJed Brown 999566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(nest,B,C1,10,&equal)); 10028b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in C1 != C"); 1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1)); 1029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A5)); 1039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&nest)); 104c20d7725SJed Brown } 105c20d7725SJed Brown 1069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&aij)); 1099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A1)); 1109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 1119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A3)); 1129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A4)); 1139566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 114b122ec5aSJacob Faibussowitsch return 0; 115c4762a1bSJed Brown } 116c4762a1bSJed Brown 117c4762a1bSJed Brown /*TEST 118c4762a1bSJed Brown 119c4762a1bSJed Brown test: 120c4762a1bSJed Brown nsize: 2 121c4762a1bSJed Brown 1223fbca975SStefano Zampini test: 1233fbca975SStefano Zampini suffix: 2 1243fbca975SStefano Zampini 125c4762a1bSJed Brown TEST*/ 126