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