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