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; 21c4762a1bSJed Brown PetscErrorCode ierr; 22c20d7725SJed Brown PetscBool equal; 23c4762a1bSJed Brown 24c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 25c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 26*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm,&size)); 27c20d7725SJed Brown 28c4762a1bSJed Brown /* 29c4762a1bSJed Brown Assemble the matrix for the five point stencil, YET AGAIN 30c4762a1bSJed Brown */ 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(comm,&A1)); 32c4762a1bSJed Brown m=2,n=2; 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A1,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A1)); 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A1)); 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A1,&istart,&iend)); 37c4762a1bSJed Brown for (ii=istart; ii<iend; ii++) { 38c4762a1bSJed Brown v = -1.0; i = ii/n; j = ii - i*n; 39*5f80ce2aSJacob Faibussowitsch if (i>0) {J = ii - n; CHKERRQ(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 40*5f80ce2aSJacob Faibussowitsch if (i<m-1) {J = ii + n; CHKERRQ(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 41*5f80ce2aSJacob Faibussowitsch if (j>0) {J = ii - 1; CHKERRQ(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 42*5f80ce2aSJacob Faibussowitsch if (j<n-1) {J = ii + 1; CHKERRQ(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 43*5f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(MatSetValues(A1,1,&ii,1,&ii,&v,INSERT_VALUES)); 44c4762a1bSJed Brown } 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A1,MAT_FINAL_ASSEMBLY)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A1,MAT_FINAL_ASSEMBLY)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A1,PETSC_VIEWER_STDOUT_WORLD)); 48c20d7725SJed Brown 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A1,MAT_COPY_VALUES,&A2)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A1,MAT_COPY_VALUES,&A3)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A1,MAT_COPY_VALUES,&A4)); 52c20d7725SJed Brown 53c4762a1bSJed Brown /*create a nest matrix */ 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(comm,&nest)); 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(nest,MATNEST)); 56c4762a1bSJed Brown mata[0]=A1,mata[1]=A2,mata[2]=A3,mata[3]=A4; 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNestSetSubMats(nest,2,NULL,2,NULL,mata)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(nest)); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(nest,MATAIJ,MAT_INITIAL_MATRIX,&aij)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(aij,PETSC_VIEWER_STDOUT_WORLD)); 61c20d7725SJed Brown 62c20d7725SJed Brown /* create a dense matrix */ 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(nest,&M,NULL)); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(nest,&m,NULL)); 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(comm,m,PETSC_DECIDE,M,K,NULL,&B)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(B,PETSC_NULL)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 69c20d7725SJed Brown 70c20d7725SJed Brown /* C = nest*B_dense */ 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(nest,B,C,10,&equal)); 742c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in C != nest*B_dense"); 753fbca975SStefano Zampini 763fbca975SStefano Zampini /* Test B = nest*C, reuse C and B with MatProductCreateWithMat() */ 773fbca975SStefano Zampini /* C has been obtained from nest*B. Clear internal data structures related to factors to prevent circular references */ 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductClear(C)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductCreateWithMat(nest,C,NULL,B)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetType(B,MATPRODUCT_AB)); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFromOptions(B)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSymbolic(B)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(B)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(nest,C,B,10,&equal)); 852c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in B != nest*C_dense"); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(nest,MATAIJ,MAT_INPLACE_MATRIX,&nest)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(nest,aij,&equal)); 882c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in aij != nest"); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&nest)); 90c20d7725SJed Brown 91c20d7725SJed Brown if (size > 1) { /* Do not know why this test fails for size = 1 */ 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateTranspose(A1,&A5)); /* A1 is symmetric */ 93c20d7725SJed Brown mata[0] = A5; 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(comm,&nest)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(nest,MATNEST)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNestSetSubMats(nest,2,NULL,2,NULL,mata)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(nest)); 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C1)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C1)); 100c20d7725SJed Brown 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(nest,B,C1,10,&equal)); 1022c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in C1 != C"); 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C1)); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A5)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&nest)); 106c20d7725SJed Brown } 107c20d7725SJed Brown 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&aij)); 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A1)); 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A2)); 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A3)); 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A4)); 115c4762a1bSJed Brown ierr = PetscFinalize(); 116c4762a1bSJed Brown return ierr; 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown /*TEST 120c4762a1bSJed Brown 121c4762a1bSJed Brown test: 122c4762a1bSJed Brown nsize: 2 123c4762a1bSJed Brown 1243fbca975SStefano Zampini test: 1253fbca975SStefano Zampini suffix: 2 1263fbca975SStefano Zampini 127c4762a1bSJed Brown TEST*/ 128