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