xref: /petsc/src/mat/tests/ex195.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
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,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 
24   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
25   comm = PETSC_COMM_WORLD;
26   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
27   /*
28      Assemble the matrix for the five point stencil, YET AGAIN
29   */
30   ierr = MatCreate(comm,&A1);CHKERRQ(ierr);
31   m=2,n=2;
32   ierr = MatSetSizes(A1,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr);
33   ierr = MatSetFromOptions(A1);CHKERRQ(ierr);
34   ierr = MatSetUp(A1);CHKERRQ(ierr);
35   ierr = MatGetOwnershipRange(A1,&istart,&iend);CHKERRQ(ierr);
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; ierr = MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
39     if (i<m-1) {J = ii + n; ierr = MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
40     if (j>0)   {J = ii - 1; ierr = MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
41     if (j<n-1) {J = ii + 1; ierr = MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
42     v = 4.0; ierr = MatSetValues(A1,1,&ii,1,&ii,&v,INSERT_VALUES);CHKERRQ(ierr);
43   }
44   ierr = MatAssemblyBegin(A1,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
45   ierr = MatAssemblyEnd(A1,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
46   ierr = MatView(A1,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
47   ierr = MatDuplicate(A1,MAT_COPY_VALUES,&A2);CHKERRQ(ierr);
48   ierr = MatDuplicate(A1,MAT_COPY_VALUES,&A3);CHKERRQ(ierr);
49   ierr = MatDuplicate(A1,MAT_COPY_VALUES,&A4);CHKERRQ(ierr);
50   /*create a nest matrix */
51   ierr = MatCreate(comm,&nest);CHKERRQ(ierr);
52   ierr = MatSetType(nest,MATNEST);CHKERRQ(ierr);
53   mata[0]=A1,mata[1]=A2,mata[2]=A3,mata[3]=A4;
54   ierr = MatNestSetSubMats(nest,2,NULL,2,NULL,mata);CHKERRQ(ierr);
55   ierr = MatSetUp(nest);CHKERRQ(ierr);
56   ierr = MatConvert(nest,MATAIJ,MAT_INITIAL_MATRIX,&aij);CHKERRQ(ierr);
57   ierr = MatView(aij,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
58   ierr = MatGetSize(nest,&M,NULL);CHKERRQ(ierr);
59   ierr = MatGetLocalSize(nest,&m,NULL);CHKERRQ(ierr);
60   ierr = MatCreateDense(comm,m,PETSC_DECIDE,M,K,NULL,&B);CHKERRQ(ierr);
61   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
62   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
63   ierr = MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
64   ierr = MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
65   ierr = MatDestroy(&C);CHKERRQ(ierr);
66   ierr = MatCreateDense(comm,m,PETSC_DECIDE,M,K,NULL,&C);CHKERRQ(ierr);
67   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
69   ierr = MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
70   ierr = MatDestroy(&C);CHKERRQ(ierr);
71   ierr = MatDestroy(&B);CHKERRQ(ierr);
72   ierr = MatDestroy(&nest);CHKERRQ(ierr);
73   ierr = MatCreateTranspose(A3,&A5);CHKERRQ(ierr);
74   mata[3]=A5;
75   ierr = MatCreate(comm,&nest);CHKERRQ(ierr);
76   ierr = MatSetType(nest,MATNEST);CHKERRQ(ierr);
77   ierr = MatNestSetSubMats(nest,2,NULL,2,NULL,mata);CHKERRQ(ierr);
78   ierr = MatSetUp(nest);CHKERRQ(ierr);
79   ierr = MatGetSize(nest,&M,NULL);CHKERRQ(ierr);
80   ierr = MatGetLocalSize(nest,&m,NULL);CHKERRQ(ierr);
81   ierr = MatCreateDense(comm,m,PETSC_DECIDE,M,K,NULL,&B);CHKERRQ(ierr);
82   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
83   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
84   ierr = MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
85   ierr = MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
86   ierr = MatDestroy(&C);CHKERRQ(ierr);
87   ierr = MatCreateDense(comm,m,PETSC_DECIDE,M,K,NULL,&C);CHKERRQ(ierr);
88   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90   ierr = MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
91   ierr = MatDestroy(&C);CHKERRQ(ierr);
92   ierr = MatDestroy(&B);CHKERRQ(ierr);
93   ierr = MatDestroy(&nest);CHKERRQ(ierr);
94   ierr = MatDestroy(&aij);CHKERRQ(ierr);
95   ierr = MatDestroy(&A1);CHKERRQ(ierr);
96   ierr = MatDestroy(&A2);CHKERRQ(ierr);
97   ierr = MatDestroy(&A3);CHKERRQ(ierr);
98   ierr = MatDestroy(&A4);CHKERRQ(ierr);
99   ierr = MatDestroy(&A5);CHKERRQ(ierr);
100   ierr = PetscFinalize();
101   return ierr;
102 }
103 
104 
105 /*TEST
106 
107    test:
108       nsize: 2
109 
110 TEST*/
111