xref: /petsc/src/mat/tests/ex195.c (revision 030f984af8d8bb4c203755d35bded3c05b3d83ce)
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   PetscErrorCode      ierr;
22   PetscBool           equal;
23 
24   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
25   comm = PETSC_COMM_WORLD;
26   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
27 
28   /*
29      Assemble the matrix for the five point stencil, YET AGAIN
30   */
31   ierr = MatCreate(comm,&A1);CHKERRQ(ierr);
32   m=2,n=2;
33   ierr = MatSetSizes(A1,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr);
34   ierr = MatSetFromOptions(A1);CHKERRQ(ierr);
35   ierr = MatSetUp(A1);CHKERRQ(ierr);
36   ierr = MatGetOwnershipRange(A1,&istart,&iend);CHKERRQ(ierr);
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; ierr = MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
40     if (i<m-1) {J = ii + n; ierr = MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
41     if (j>0)   {J = ii - 1; ierr = MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
42     if (j<n-1) {J = ii + 1; ierr = MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
43     v = 4.0; ierr = MatSetValues(A1,1,&ii,1,&ii,&v,INSERT_VALUES);CHKERRQ(ierr);
44   }
45   ierr = MatAssemblyBegin(A1,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
46   ierr = MatAssemblyEnd(A1,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
47   ierr = MatView(A1,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
48 
49   ierr = MatDuplicate(A1,MAT_COPY_VALUES,&A2);CHKERRQ(ierr);
50   ierr = MatDuplicate(A1,MAT_COPY_VALUES,&A3);CHKERRQ(ierr);
51   ierr = MatDuplicate(A1,MAT_COPY_VALUES,&A4);CHKERRQ(ierr);
52 
53   /*create a nest matrix */
54   ierr = MatCreate(comm,&nest);CHKERRQ(ierr);
55   ierr = MatSetType(nest,MATNEST);CHKERRQ(ierr);
56   mata[0]=A1,mata[1]=A2,mata[2]=A3,mata[3]=A4;
57   ierr = MatNestSetSubMats(nest,2,NULL,2,NULL,mata);CHKERRQ(ierr);
58   ierr = MatSetUp(nest);CHKERRQ(ierr);
59   ierr = MatConvert(nest,MATAIJ,MAT_INITIAL_MATRIX,&aij);CHKERRQ(ierr);
60   ierr = MatView(aij,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
61 
62   /* create a dense matrix */
63   ierr = MatGetSize(nest,&M,NULL);CHKERRQ(ierr);
64   ierr = MatGetLocalSize(nest,&m,NULL);CHKERRQ(ierr);
65   ierr = MatCreateDense(comm,m,PETSC_DECIDE,M,K,NULL,&B);CHKERRQ(ierr);
66   ierr = MatSetRandom(B,PETSC_NULL);CHKERRQ(ierr);
67   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
69 
70   /* C = nest*B_dense */
71   ierr = MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
72   ierr = MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
73   ierr = MatMatMultEqual(nest,B,C,10,&equal);CHKERRQ(ierr);
74   if (!equal) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in C != nest*B_dense");
75 
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   ierr = MatProductClear(C);CHKERRQ(ierr);
79   ierr = MatProductCreateWithMat(nest,C,NULL,B);CHKERRQ(ierr);
80   ierr = MatProductSetType(B,MATPRODUCT_AB);CHKERRQ(ierr);
81   ierr = MatProductSetFromOptions(B);CHKERRQ(ierr);
82   ierr = MatProductSymbolic(B);CHKERRQ(ierr);
83   ierr = MatProductNumeric(B);CHKERRQ(ierr);
84   ierr = MatMatMultEqual(nest,C,B,10,&equal);CHKERRQ(ierr);
85   if (!equal) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in B != nest*C_dense");
86   ierr = MatConvert(nest,MATAIJ,MAT_INPLACE_MATRIX,&nest);CHKERRQ(ierr);
87   ierr = MatEqual(nest,aij,&equal);CHKERRQ(ierr);
88   if (!equal) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in aij != nest");
89   ierr = MatDestroy(&nest);CHKERRQ(ierr);
90 
91   if (size > 1) { /* Do not know why this test fails for size = 1 */
92     ierr = MatCreateTranspose(A1,&A5);CHKERRQ(ierr); /* A1 is symmetric */
93     mata[0] = A5;
94     ierr = MatCreate(comm,&nest);CHKERRQ(ierr);
95     ierr = MatSetType(nest,MATNEST);CHKERRQ(ierr);
96     ierr = MatNestSetSubMats(nest,2,NULL,2,NULL,mata);CHKERRQ(ierr);
97     ierr = MatSetUp(nest);CHKERRQ(ierr);
98     ierr = MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C1);CHKERRQ(ierr);
99     ierr = MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C1);CHKERRQ(ierr);
100 
101     ierr = MatMatMultEqual(nest,B,C1,10,&equal);CHKERRQ(ierr);
102     if (!equal) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in C1 != C");
103     ierr = MatDestroy(&C1);CHKERRQ(ierr);
104     ierr = MatDestroy(&A5);CHKERRQ(ierr);
105     ierr = MatDestroy(&nest);CHKERRQ(ierr);
106   }
107 
108   ierr = MatDestroy(&C);CHKERRQ(ierr);
109   ierr = MatDestroy(&B);CHKERRQ(ierr);
110   ierr = MatDestroy(&aij);CHKERRQ(ierr);
111   ierr = MatDestroy(&A1);CHKERRQ(ierr);
112   ierr = MatDestroy(&A2);CHKERRQ(ierr);
113   ierr = MatDestroy(&A3);CHKERRQ(ierr);
114   ierr = MatDestroy(&A4);CHKERRQ(ierr);
115   ierr = PetscFinalize();
116   return ierr;
117 }
118 
119 /*TEST
120 
121    test:
122       nsize: 2
123 
124    test:
125       suffix: 2
126 
127 TEST*/
128