xref: /petsc/src/mat/tests/ex195.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown /*
2*c4762a1bSJed Brown  * ex195.c
3*c4762a1bSJed Brown  *
4*c4762a1bSJed Brown  *  Created on: Aug 24, 2015
5*c4762a1bSJed Brown  *      Author: Fande Kong <fdkong.jd@gmail.com>
6*c4762a1bSJed Brown  */
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown static char help[] = " Demonstrate the use of MatConvert_Nest_AIJ\n";
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown #include <petscmat.h>
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown int main(int argc,char **args)
14*c4762a1bSJed Brown {
15*c4762a1bSJed Brown   Mat                 A1,A2,A3,A4,A5,B,C,nest;
16*c4762a1bSJed Brown   Mat                 mata[4];
17*c4762a1bSJed Brown   Mat                 aij;
18*c4762a1bSJed Brown   MPI_Comm            comm;
19*c4762a1bSJed Brown   PetscInt            m,M,n,istart,iend,ii,i,J,j,K=10;
20*c4762a1bSJed Brown   PetscScalar         v;
21*c4762a1bSJed Brown   PetscMPIInt         size;
22*c4762a1bSJed Brown   PetscErrorCode      ierr;
23*c4762a1bSJed Brown 
24*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
25*c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
26*c4762a1bSJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
27*c4762a1bSJed Brown   /*
28*c4762a1bSJed Brown      Assemble the matrix for the five point stencil, YET AGAIN
29*c4762a1bSJed Brown   */
30*c4762a1bSJed Brown   ierr = MatCreate(comm,&A1);CHKERRQ(ierr);
31*c4762a1bSJed Brown   m=2,n=2;
32*c4762a1bSJed Brown   ierr = MatSetSizes(A1,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = MatSetFromOptions(A1);CHKERRQ(ierr);
34*c4762a1bSJed Brown   ierr = MatSetUp(A1);CHKERRQ(ierr);
35*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A1,&istart,&iend);CHKERRQ(ierr);
36*c4762a1bSJed Brown   for (ii=istart; ii<iend; ii++) {
37*c4762a1bSJed Brown     v = -1.0; i = ii/n; j = ii - i*n;
38*c4762a1bSJed Brown     if (i>0)   {J = ii - n; ierr = MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
39*c4762a1bSJed Brown     if (i<m-1) {J = ii + n; ierr = MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
40*c4762a1bSJed Brown     if (j>0)   {J = ii - 1; ierr = MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
41*c4762a1bSJed Brown     if (j<n-1) {J = ii + 1; ierr = MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
42*c4762a1bSJed Brown     v = 4.0; ierr = MatSetValues(A1,1,&ii,1,&ii,&v,INSERT_VALUES);CHKERRQ(ierr);
43*c4762a1bSJed Brown   }
44*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A1,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
45*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A1,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = MatView(A1,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
47*c4762a1bSJed Brown   ierr = MatDuplicate(A1,MAT_COPY_VALUES,&A2);CHKERRQ(ierr);
48*c4762a1bSJed Brown   ierr = MatDuplicate(A1,MAT_COPY_VALUES,&A3);CHKERRQ(ierr);
49*c4762a1bSJed Brown   ierr = MatDuplicate(A1,MAT_COPY_VALUES,&A4);CHKERRQ(ierr);
50*c4762a1bSJed Brown   /*create a nest matrix */
51*c4762a1bSJed Brown   ierr = MatCreate(comm,&nest);CHKERRQ(ierr);
52*c4762a1bSJed Brown   ierr = MatSetType(nest,MATNEST);CHKERRQ(ierr);
53*c4762a1bSJed Brown   mata[0]=A1,mata[1]=A2,mata[2]=A3,mata[3]=A4;
54*c4762a1bSJed Brown   ierr = MatNestSetSubMats(nest,2,NULL,2,NULL,mata);CHKERRQ(ierr);
55*c4762a1bSJed Brown   ierr = MatSetUp(nest);CHKERRQ(ierr);
56*c4762a1bSJed Brown   ierr = MatConvert(nest,MATAIJ,MAT_INITIAL_MATRIX,&aij);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = MatView(aij,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
58*c4762a1bSJed Brown   ierr = MatGetSize(nest,&M,NULL);CHKERRQ(ierr);
59*c4762a1bSJed Brown   ierr = MatGetLocalSize(nest,&m,NULL);CHKERRQ(ierr);
60*c4762a1bSJed Brown   ierr = MatCreateDense(comm,m,PETSC_DECIDE,M,K,NULL,&B);CHKERRQ(ierr);
61*c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr = MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
64*c4762a1bSJed Brown   ierr = MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
65*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
66*c4762a1bSJed Brown   ierr = MatCreateDense(comm,m,PETSC_DECIDE,M,K,NULL,&C);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68*c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
69*c4762a1bSJed Brown   ierr = MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
70*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
71*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
72*c4762a1bSJed Brown   ierr = MatDestroy(&nest);CHKERRQ(ierr);
73*c4762a1bSJed Brown   ierr = MatCreateTranspose(A3,&A5);CHKERRQ(ierr);
74*c4762a1bSJed Brown   mata[3]=A5;
75*c4762a1bSJed Brown   ierr = MatCreate(comm,&nest);CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = MatSetType(nest,MATNEST);CHKERRQ(ierr);
77*c4762a1bSJed Brown   ierr = MatNestSetSubMats(nest,2,NULL,2,NULL,mata);CHKERRQ(ierr);
78*c4762a1bSJed Brown   ierr = MatSetUp(nest);CHKERRQ(ierr);
79*c4762a1bSJed Brown   ierr = MatGetSize(nest,&M,NULL);CHKERRQ(ierr);
80*c4762a1bSJed Brown   ierr = MatGetLocalSize(nest,&m,NULL);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ierr = MatCreateDense(comm,m,PETSC_DECIDE,M,K,NULL,&B);CHKERRQ(ierr);
82*c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr = MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
86*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = MatCreateDense(comm,m,PETSC_DECIDE,M,K,NULL,&C);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
92*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = MatDestroy(&nest);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = MatDestroy(&aij);CHKERRQ(ierr);
95*c4762a1bSJed Brown   ierr = MatDestroy(&A1);CHKERRQ(ierr);
96*c4762a1bSJed Brown   ierr = MatDestroy(&A2);CHKERRQ(ierr);
97*c4762a1bSJed Brown   ierr = MatDestroy(&A3);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = MatDestroy(&A4);CHKERRQ(ierr);
99*c4762a1bSJed Brown   ierr = MatDestroy(&A5);CHKERRQ(ierr);
100*c4762a1bSJed Brown   ierr = PetscFinalize();
101*c4762a1bSJed Brown   return ierr;
102*c4762a1bSJed Brown }
103*c4762a1bSJed Brown 
104*c4762a1bSJed Brown 
105*c4762a1bSJed Brown /*TEST
106*c4762a1bSJed Brown 
107*c4762a1bSJed Brown    test:
108*c4762a1bSJed Brown       nsize: 2
109*c4762a1bSJed Brown 
110*c4762a1bSJed Brown TEST*/
111