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