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*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 24c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 255f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm,&size)); 26c20d7725SJed Brown 27c4762a1bSJed Brown /* 28c4762a1bSJed Brown Assemble the matrix for the five point stencil, YET AGAIN 29c4762a1bSJed Brown */ 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(comm,&A1)); 31c4762a1bSJed Brown m=2,n=2; 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A1,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A1)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A1)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A1,&istart,&iend)); 36c4762a1bSJed Brown for (ii=istart; ii<iend; ii++) { 37c4762a1bSJed Brown v = -1.0; i = ii/n; j = ii - i*n; 385f80ce2aSJacob Faibussowitsch if (i>0) {J = ii - n; CHKERRQ(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 395f80ce2aSJacob Faibussowitsch if (i<m-1) {J = ii + n; CHKERRQ(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 405f80ce2aSJacob Faibussowitsch if (j>0) {J = ii - 1; CHKERRQ(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 415f80ce2aSJacob Faibussowitsch if (j<n-1) {J = ii + 1; CHKERRQ(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 425f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(MatSetValues(A1,1,&ii,1,&ii,&v,INSERT_VALUES)); 43c4762a1bSJed Brown } 445f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A1,MAT_FINAL_ASSEMBLY)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A1,MAT_FINAL_ASSEMBLY)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A1,PETSC_VIEWER_STDOUT_WORLD)); 47c20d7725SJed Brown 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A1,MAT_COPY_VALUES,&A2)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A1,MAT_COPY_VALUES,&A3)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A1,MAT_COPY_VALUES,&A4)); 51c20d7725SJed Brown 52c4762a1bSJed Brown /*create a nest matrix */ 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(comm,&nest)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(nest,MATNEST)); 55c4762a1bSJed Brown mata[0]=A1,mata[1]=A2,mata[2]=A3,mata[3]=A4; 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatNestSetSubMats(nest,2,NULL,2,NULL,mata)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(nest)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(nest,MATAIJ,MAT_INITIAL_MATRIX,&aij)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(aij,PETSC_VIEWER_STDOUT_WORLD)); 60c20d7725SJed Brown 61c20d7725SJed Brown /* create a dense matrix */ 625f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(nest,&M,NULL)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(nest,&m,NULL)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(comm,m,PETSC_DECIDE,M,K,NULL,&B)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(B,PETSC_NULL)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 68c20d7725SJed Brown 69c20d7725SJed Brown /* C = nest*B_dense */ 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(nest,B,C,10,&equal)); 7328b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in C != nest*B_dense"); 743fbca975SStefano Zampini 753fbca975SStefano Zampini /* Test B = nest*C, reuse C and B with MatProductCreateWithMat() */ 763fbca975SStefano Zampini /* C has been obtained from nest*B. Clear internal data structures related to factors to prevent circular references */ 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductClear(C)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductCreateWithMat(nest,C,NULL,B)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetType(B,MATPRODUCT_AB)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFromOptions(B)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSymbolic(B)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(B)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(nest,C,B,10,&equal)); 8428b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in B != nest*C_dense"); 855f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(nest,MATAIJ,MAT_INPLACE_MATRIX,&nest)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(nest,aij,&equal)); 8728b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in aij != nest"); 885f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&nest)); 89c20d7725SJed Brown 90c20d7725SJed Brown if (size > 1) { /* Do not know why this test fails for size = 1 */ 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateTranspose(A1,&A5)); /* A1 is symmetric */ 92c20d7725SJed Brown mata[0] = A5; 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(comm,&nest)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(nest,MATNEST)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(MatNestSetSubMats(nest,2,NULL,2,NULL,mata)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(nest)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C1)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C1)); 99c20d7725SJed Brown 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(nest,B,C1,10,&equal)); 10128b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in C1 != C"); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C1)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A5)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&nest)); 105c20d7725SJed Brown } 106c20d7725SJed Brown 1075f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&aij)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A1)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A2)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A3)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A4)); 114*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 115*b122ec5aSJacob Faibussowitsch return 0; 116c4762a1bSJed Brown } 117c4762a1bSJed Brown 118c4762a1bSJed Brown /*TEST 119c4762a1bSJed Brown 120c4762a1bSJed Brown test: 121c4762a1bSJed Brown nsize: 2 122c4762a1bSJed Brown 1233fbca975SStefano Zampini test: 1243fbca975SStefano Zampini suffix: 2 1253fbca975SStefano Zampini 126c4762a1bSJed Brown TEST*/ 127