1 /* 2 * ex195.c 3 * 4 * Created on: Aug 24, 2015 5 * Author: Fande Kong <fdkong.jd@gmail.com> 6 */ 7 8 static char help[] = " Demonstrate the use of MatConvert_Nest_AIJ\n"; 9 10 #include <petscmat.h> 11 12 int main(int argc,char **args) 13 { 14 Mat A1,A2,A3,A4,A5,B,C,C1,nest; 15 Mat mata[4]; 16 Mat aij; 17 MPI_Comm comm; 18 PetscInt m,M,n,istart,iend,ii,i,J,j,K=10; 19 PetscScalar v; 20 PetscMPIInt size; 21 PetscBool equal; 22 23 PetscFunctionBeginUser; 24 PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 25 comm = PETSC_COMM_WORLD; 26 PetscCallMPI(MPI_Comm_size(comm,&size)); 27 28 /* 29 Assemble the matrix for the five point stencil, YET AGAIN 30 */ 31 PetscCall(MatCreate(comm,&A1)); 32 m=2,n=2; 33 PetscCall(MatSetSizes(A1,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 34 PetscCall(MatSetFromOptions(A1)); 35 PetscCall(MatSetUp(A1)); 36 PetscCall(MatGetOwnershipRange(A1,&istart,&iend)); 37 for (ii=istart; ii<iend; ii++) { 38 v = -1.0; i = ii/n; j = ii - i*n; 39 if (i>0) {J = ii - n; PetscCall(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 40 if (i<m-1) {J = ii + n; PetscCall(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 41 if (j>0) {J = ii - 1; PetscCall(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 42 if (j<n-1) {J = ii + 1; PetscCall(MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES));} 43 v = 4.0; PetscCall(MatSetValues(A1,1,&ii,1,&ii,&v,INSERT_VALUES)); 44 } 45 PetscCall(MatAssemblyBegin(A1,MAT_FINAL_ASSEMBLY)); 46 PetscCall(MatAssemblyEnd(A1,MAT_FINAL_ASSEMBLY)); 47 PetscCall(MatView(A1,PETSC_VIEWER_STDOUT_WORLD)); 48 49 PetscCall(MatDuplicate(A1,MAT_COPY_VALUES,&A2)); 50 PetscCall(MatDuplicate(A1,MAT_COPY_VALUES,&A3)); 51 PetscCall(MatDuplicate(A1,MAT_COPY_VALUES,&A4)); 52 53 /*create a nest matrix */ 54 PetscCall(MatCreate(comm,&nest)); 55 PetscCall(MatSetType(nest,MATNEST)); 56 mata[0]=A1,mata[1]=A2,mata[2]=A3,mata[3]=A4; 57 PetscCall(MatNestSetSubMats(nest,2,NULL,2,NULL,mata)); 58 PetscCall(MatSetUp(nest)); 59 PetscCall(MatConvert(nest,MATAIJ,MAT_INITIAL_MATRIX,&aij)); 60 PetscCall(MatView(aij,PETSC_VIEWER_STDOUT_WORLD)); 61 62 /* create a dense matrix */ 63 PetscCall(MatGetSize(nest,&M,NULL)); 64 PetscCall(MatGetLocalSize(nest,&m,NULL)); 65 PetscCall(MatCreateDense(comm,m,PETSC_DECIDE,M,K,NULL,&B)); 66 PetscCall(MatSetRandom(B,PETSC_NULL)); 67 68 /* C = nest*B_dense */ 69 PetscCall(MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 70 PetscCall(MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C)); 71 PetscCall(MatMatMultEqual(nest,B,C,10,&equal)); 72 PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in C != nest*B_dense"); 73 74 /* Test B = nest*C, reuse C and B with MatProductCreateWithMat() */ 75 /* C has been obtained from nest*B. Clear internal data structures related to factors to prevent circular references */ 76 PetscCall(MatProductClear(C)); 77 PetscCall(MatProductCreateWithMat(nest,C,NULL,B)); 78 PetscCall(MatProductSetType(B,MATPRODUCT_AB)); 79 PetscCall(MatProductSetFromOptions(B)); 80 PetscCall(MatProductSymbolic(B)); 81 PetscCall(MatProductNumeric(B)); 82 PetscCall(MatMatMultEqual(nest,C,B,10,&equal)); 83 PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in B != nest*C_dense"); 84 PetscCall(MatConvert(nest,MATAIJ,MAT_INPLACE_MATRIX,&nest)); 85 PetscCall(MatEqual(nest,aij,&equal)); 86 PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in aij != nest"); 87 PetscCall(MatDestroy(&nest)); 88 89 if (size > 1) { /* Do not know why this test fails for size = 1 */ 90 PetscCall(MatCreateTranspose(A1,&A5)); /* A1 is symmetric */ 91 mata[0] = A5; 92 PetscCall(MatCreate(comm,&nest)); 93 PetscCall(MatSetType(nest,MATNEST)); 94 PetscCall(MatNestSetSubMats(nest,2,NULL,2,NULL,mata)); 95 PetscCall(MatSetUp(nest)); 96 PetscCall(MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C1)); 97 PetscCall(MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C1)); 98 99 PetscCall(MatMatMultEqual(nest,B,C1,10,&equal)); 100 PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in C1 != C"); 101 PetscCall(MatDestroy(&C1)); 102 PetscCall(MatDestroy(&A5)); 103 PetscCall(MatDestroy(&nest)); 104 } 105 106 PetscCall(MatDestroy(&C)); 107 PetscCall(MatDestroy(&B)); 108 PetscCall(MatDestroy(&aij)); 109 PetscCall(MatDestroy(&A1)); 110 PetscCall(MatDestroy(&A2)); 111 PetscCall(MatDestroy(&A3)); 112 PetscCall(MatDestroy(&A4)); 113 PetscCall(PetscFinalize()); 114 return 0; 115 } 116 117 /*TEST 118 119 test: 120 nsize: 2 121 122 test: 123 suffix: 2 124 125 TEST*/ 126