xref: /petsc/src/mat/tests/ex89.c (revision c20d77252dee0f9c80fc6f8b1a6f948e11175edb)
1c4762a1bSJed Brown static char help[] ="Tests MatPtAP() for MPIMAIJ and MPIAIJ \n ";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmda.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc,char **argv)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   PetscErrorCode ierr;
8c4762a1bSJed Brown   DM             coarsedm,finedm;
9c4762a1bSJed Brown   PetscMPIInt    size,rank;
10c4762a1bSJed Brown   PetscInt       M,N,Z,i,nrows;
11c4762a1bSJed Brown   PetscScalar    one = 1.0;
12c4762a1bSJed Brown   PetscReal      fill=2.0;
13c4762a1bSJed Brown   Mat            A,P,C;
14*c20d7725SJed Brown   PetscScalar    *array,alpha;
15c4762a1bSJed Brown   PetscBool      Test_3D=PETSC_FALSE,flg;
16c4762a1bSJed Brown   const PetscInt *ia,*ja;
17c4762a1bSJed Brown   PetscInt       dof;
18c4762a1bSJed Brown   MPI_Comm       comm;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
21c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
22c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
23c4762a1bSJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
24c4762a1bSJed Brown   M = 10; N = 10; Z = 10;
25c4762a1bSJed Brown   dof  = 10;
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-test_3D",&Test_3D,NULL);CHKERRQ(ierr);
28c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr);
29c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr);
30c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-Z",&Z,NULL);CHKERRQ(ierr);
31c4762a1bSJed Brown   /* Set up distributed array for fine grid */
32c4762a1bSJed Brown   if (!Test_3D) {
33c4762a1bSJed Brown     ierr = DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,&coarsedm);CHKERRQ(ierr);
34c4762a1bSJed Brown   } else {
35c4762a1bSJed Brown     ierr = DMDACreate3d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,Z,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,NULL,&coarsedm);CHKERRQ(ierr);
36c4762a1bSJed Brown   }
37c4762a1bSJed Brown   ierr = DMSetFromOptions(coarsedm);CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = DMSetUp(coarsedm);CHKERRQ(ierr);
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   /* This makes sure the coarse DMDA has the same partition as the fine DMDA */
41c4762a1bSJed Brown   ierr = DMRefine(coarsedm,PetscObjectComm((PetscObject)coarsedm),&finedm);CHKERRQ(ierr);
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /*------------------------------------------------------------*/
44c4762a1bSJed Brown   ierr = DMSetMatType(finedm,MATAIJ);CHKERRQ(ierr);
45c4762a1bSJed Brown   ierr = DMCreateMatrix(finedm,&A);CHKERRQ(ierr);
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* set val=one to A */
48c4762a1bSJed Brown   if (size == 1) {
49c4762a1bSJed Brown     ierr = MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
50c4762a1bSJed Brown     if (flg) {
51c4762a1bSJed Brown       ierr = MatSeqAIJGetArray(A,&array);CHKERRQ(ierr);
52c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
53c4762a1bSJed Brown       ierr = MatSeqAIJRestoreArray(A,&array);CHKERRQ(ierr);
54c4762a1bSJed Brown     }
55c4762a1bSJed Brown     ierr = MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
56c4762a1bSJed Brown   } else {
57c4762a1bSJed Brown     Mat AA,AB;
58c4762a1bSJed Brown     ierr = MatMPIAIJGetSeqAIJ(A,&AA,&AB,NULL);CHKERRQ(ierr);
59c4762a1bSJed Brown     ierr = MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
60c4762a1bSJed Brown     if (flg) {
61c4762a1bSJed Brown       ierr = MatSeqAIJGetArray(AA,&array);CHKERRQ(ierr);
62c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
63c4762a1bSJed Brown       ierr = MatSeqAIJRestoreArray(AA,&array);CHKERRQ(ierr);
64c4762a1bSJed Brown     }
65c4762a1bSJed Brown     ierr = MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
66c4762a1bSJed Brown     ierr = MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
67c4762a1bSJed Brown     if (flg) {
68c4762a1bSJed Brown       ierr = MatSeqAIJGetArray(AB,&array);CHKERRQ(ierr);
69c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
70c4762a1bSJed Brown       ierr = MatSeqAIJRestoreArray(AB,&array);CHKERRQ(ierr);
71c4762a1bSJed Brown     }
72c4762a1bSJed Brown     ierr = MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
73c4762a1bSJed Brown   }
74c4762a1bSJed Brown   /* Create interpolation between the fine and coarse grids */
75c4762a1bSJed Brown   ierr = DMCreateInterpolation(coarsedm,finedm,&P,NULL);CHKERRQ(ierr);
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /* Test P^T * A * P - MatPtAP() */
78c4762a1bSJed Brown   /*------------------------------*/
79*c20d7725SJed Brown   /* (1) Developer API */
80*c20d7725SJed Brown   ierr = MatProductCreate(A,P,NULL,&C);CHKERRQ(ierr);
81*c20d7725SJed Brown   ierr = MatProductSetType(C,MATPRODUCT_PtAP);CHKERRQ(ierr);
82*c20d7725SJed Brown   ierr = MatProductSetAlgorithm(C,"default");CHKERRQ(ierr);
83*c20d7725SJed Brown   ierr = MatProductSetFill(C,PETSC_DEFAULT);CHKERRQ(ierr);
84*c20d7725SJed Brown   ierr = MatProductSetFromOptions(C);CHKERRQ(ierr);
85*c20d7725SJed Brown   ierr = MatProductSymbolic(C);CHKERRQ(ierr);
86*c20d7725SJed Brown   ierr = MatProductNumeric(C);CHKERRQ(ierr);
87*c20d7725SJed Brown   ierr = MatProductNumeric(C);CHKERRQ(ierr); /* Test reuse of symbolic C */
88*c20d7725SJed Brown 
89*c20d7725SJed Brown   ierr = MatPtAPMultEqual(A,P,C,10,&flg);CHKERRQ(ierr);
90*c20d7725SJed Brown   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatProduct_PtAP");
91*c20d7725SJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
92*c20d7725SJed Brown 
93*c20d7725SJed Brown   /* (2) User API */
94c4762a1bSJed Brown   ierr = MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr);
95c4762a1bSJed Brown   /* Test MAT_REUSE_MATRIX - reuse symbolic C */
96c4762a1bSJed Brown   alpha=1.0;
97c4762a1bSJed Brown   for (i=0; i<1; i++) {
98c4762a1bSJed Brown     alpha -= 0.1;
99c4762a1bSJed Brown     ierr   = MatScale(A,alpha);CHKERRQ(ierr);
100c4762a1bSJed Brown     ierr   = MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr);
101c4762a1bSJed Brown   }
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   /* Free intermediate data structures created for reuse of C=Pt*A*P */
104c4762a1bSJed Brown   ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr);
105c4762a1bSJed Brown 
106*c20d7725SJed Brown   ierr = MatPtAPMultEqual(A,P,C,10,&flg);CHKERRQ(ierr);
107*c20d7725SJed Brown   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatPtAP");
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
110c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
111c4762a1bSJed Brown   ierr = MatDestroy(&P);CHKERRQ(ierr);
112c4762a1bSJed Brown   ierr = DMDestroy(&finedm);CHKERRQ(ierr);
113c4762a1bSJed Brown   ierr = DMDestroy(&coarsedm);CHKERRQ(ierr);
114c4762a1bSJed Brown   ierr = PetscFinalize();
115c4762a1bSJed Brown   return ierr;
116c4762a1bSJed Brown }
117c4762a1bSJed Brown 
118c4762a1bSJed Brown /*TEST
119c4762a1bSJed Brown 
120c4762a1bSJed Brown    test:
121c4762a1bSJed Brown       args: -M 10 -N 10 -Z 10
122c4762a1bSJed Brown       output_file: output/ex89_1.out
123c4762a1bSJed Brown 
124c4762a1bSJed Brown    test:
125c4762a1bSJed Brown       suffix: allatonce
126c4762a1bSJed Brown       nsize: 4
127*c20d7725SJed Brown       args: -M 10 -N 10 -Z 10 -matptap_via allatonce
128c4762a1bSJed Brown       output_file: output/ex89_1.out
129c4762a1bSJed Brown 
130c4762a1bSJed Brown    test:
131c4762a1bSJed Brown       suffix: allatonce_merged
132c4762a1bSJed Brown       nsize: 4
133*c20d7725SJed Brown       args: -M 10 -M 5 -M 10 -matptap_via allatonce_merged
134c4762a1bSJed Brown       output_file: output/ex96_1.out
135c4762a1bSJed Brown 
136c4762a1bSJed Brown    test:
137c4762a1bSJed Brown       suffix: allatonce_3D
138c4762a1bSJed Brown       nsize: 4
139*c20d7725SJed Brown       args: -M 10 -M 5 -M 10 -test_3D 1 -matptap_via allatonce
140c4762a1bSJed Brown       output_file: output/ex96_1.out
141c4762a1bSJed Brown 
142c4762a1bSJed Brown    test:
143c4762a1bSJed Brown       suffix: allatonce_merged_3D
144c4762a1bSJed Brown       nsize: 4
145*c20d7725SJed Brown       args: -M 10 -M 5 -M 10 -test_3D 1 -matptap_via allatonce_merged
146c4762a1bSJed Brown       output_file: output/ex96_1.out
147c4762a1bSJed Brown 
148c4762a1bSJed Brown TEST*/
149