xref: /petsc/src/mat/tests/ex89.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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   DM             coarsedm,finedm;
8c4762a1bSJed Brown   PetscMPIInt    size,rank;
9c4762a1bSJed Brown   PetscInt       M,N,Z,i,nrows;
10c4762a1bSJed Brown   PetscScalar    one = 1.0;
11c4762a1bSJed Brown   PetscReal      fill=2.0;
12c4762a1bSJed Brown   Mat            A,P,C;
13c20d7725SJed Brown   PetscScalar    *array,alpha;
14c4762a1bSJed Brown   PetscBool      Test_3D=PETSC_FALSE,flg;
15c4762a1bSJed Brown   const PetscInt *ia,*ja;
16c4762a1bSJed Brown   PetscInt       dof;
17c4762a1bSJed Brown   MPI_Comm       comm;
18c4762a1bSJed Brown 
19*327415f7SBarry Smith   PetscFunctionBeginUser;
209566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
21c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
24c4762a1bSJed Brown   M = 10; N = 10; Z = 10;
25c4762a1bSJed Brown   dof  = 10;
26c4762a1bSJed Brown 
279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_3D",&Test_3D,NULL));
289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL));
309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-Z",&Z,NULL));
31c4762a1bSJed Brown   /* Set up distributed array for fine grid */
32c4762a1bSJed Brown   if (!Test_3D) {
339566063dSJacob Faibussowitsch     PetscCall(DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,&coarsedm));
34c4762a1bSJed Brown   } else {
359566063dSJacob Faibussowitsch     PetscCall(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));
36c4762a1bSJed Brown   }
379566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(coarsedm));
389566063dSJacob Faibussowitsch   PetscCall(DMSetUp(coarsedm));
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   /* This makes sure the coarse DMDA has the same partition as the fine DMDA */
419566063dSJacob Faibussowitsch   PetscCall(DMRefine(coarsedm,PetscObjectComm((PetscObject)coarsedm),&finedm));
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /*------------------------------------------------------------*/
449566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(finedm,MATAIJ));
459566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(finedm,&A));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* set val=one to A */
48c4762a1bSJed Brown   if (size == 1) {
499566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
50c4762a1bSJed Brown     if (flg) {
519566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(A,&array));
52c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
539566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(A,&array));
54c4762a1bSJed Brown     }
559566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
56c4762a1bSJed Brown   } else {
57c4762a1bSJed Brown     Mat AA,AB;
589566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJGetSeqAIJ(A,&AA,&AB,NULL));
599566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
60c4762a1bSJed Brown     if (flg) {
619566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(AA,&array));
62c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
639566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(AA,&array));
64c4762a1bSJed Brown     }
659566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
669566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
67c4762a1bSJed Brown     if (flg) {
689566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(AB,&array));
69c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
709566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(AB,&array));
71c4762a1bSJed Brown     }
729566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
73c4762a1bSJed Brown   }
74c4762a1bSJed Brown   /* Create interpolation between the fine and coarse grids */
759566063dSJacob Faibussowitsch   PetscCall(DMCreateInterpolation(coarsedm,finedm,&P,NULL));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /* Test P^T * A * P - MatPtAP() */
78c4762a1bSJed Brown   /*------------------------------*/
79c20d7725SJed Brown   /* (1) Developer API */
809566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(A,P,NULL,&C));
819566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(C,MATPRODUCT_PtAP));
829566063dSJacob Faibussowitsch   PetscCall(MatProductSetAlgorithm(C,"allatonce"));
839566063dSJacob Faibussowitsch   PetscCall(MatProductSetFill(C,PETSC_DEFAULT));
849566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(C));
859566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(C));
869566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(C));
879566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(C)); /* Test reuse of symbolic C */
88c20d7725SJed Brown 
89c2f6fc52SHong Zhang   { /* Test MatProductView() */
90c2f6fc52SHong Zhang     PetscViewer viewer;
919566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(comm,NULL, &viewer));
929566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO));
939566063dSJacob Faibussowitsch     PetscCall(MatProductView(C,viewer));
949566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
959566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
96c2f6fc52SHong Zhang   }
97c2f6fc52SHong Zhang 
989566063dSJacob Faibussowitsch   PetscCall(MatPtAPMultEqual(A,P,C,10,&flg));
9928b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatProduct_PtAP");
1009566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
101c20d7725SJed Brown 
102c20d7725SJed Brown   /* (2) User API */
1039566063dSJacob Faibussowitsch   PetscCall(MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&C));
104c4762a1bSJed Brown   /* Test MAT_REUSE_MATRIX - reuse symbolic C */
105c4762a1bSJed Brown   alpha=1.0;
106c4762a1bSJed Brown   for (i=0; i<1; i++) {
107c4762a1bSJed Brown     alpha -= 0.1;
1089566063dSJacob Faibussowitsch     PetscCall(MatScale(A,alpha));
1099566063dSJacob Faibussowitsch     PetscCall(MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&C));
110c4762a1bSJed Brown   }
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /* Free intermediate data structures created for reuse of C=Pt*A*P */
1139566063dSJacob Faibussowitsch   PetscCall(MatProductClear(C));
114c4762a1bSJed Brown 
1159566063dSJacob Faibussowitsch   PetscCall(MatPtAPMultEqual(A,P,C,10,&flg));
11628b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatPtAP");
117c4762a1bSJed Brown 
1189566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1199566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1209566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&P));
1219566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&finedm));
1229566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&coarsedm));
1239566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
124b122ec5aSJacob Faibussowitsch   return 0;
125c4762a1bSJed Brown }
126c4762a1bSJed Brown 
127c4762a1bSJed Brown /*TEST
128c4762a1bSJed Brown 
129c4762a1bSJed Brown    test:
130c4762a1bSJed Brown       args: -M 10 -N 10 -Z 10
131c4762a1bSJed Brown       output_file: output/ex89_1.out
132c4762a1bSJed Brown 
133c4762a1bSJed Brown    test:
134c4762a1bSJed Brown       suffix: allatonce
135c4762a1bSJed Brown       nsize: 4
136c2f6fc52SHong Zhang       args: -M 10 -N 10 -Z 10
137c2f6fc52SHong Zhang       output_file: output/ex89_2.out
138c4762a1bSJed Brown 
139c4762a1bSJed Brown    test:
140c4762a1bSJed Brown       suffix: allatonce_merged
141c4762a1bSJed Brown       nsize: 4
1423e662e0bSHong Zhang       args: -M 10 -M 5 -M 10 -mat_product_algorithm allatonce_merged
143c2f6fc52SHong Zhang       output_file: output/ex89_3.out
144c4762a1bSJed Brown 
145c4762a1bSJed Brown    test:
146c2f6fc52SHong Zhang       suffix: nonscalable_3D
147c4762a1bSJed Brown       nsize: 4
1483e662e0bSHong Zhang       args: -M 10 -M 5 -M 10 -test_3D 1 -mat_product_algorithm nonscalable
149c2f6fc52SHong Zhang       output_file: output/ex89_4.out
150c4762a1bSJed Brown 
151c4762a1bSJed Brown    test:
152c4762a1bSJed Brown       suffix: allatonce_merged_3D
153c4762a1bSJed Brown       nsize: 4
1543e662e0bSHong Zhang       args: -M 10 -M 5 -M 10 -test_3D 1 -mat_product_algorithm allatonce_merged
155c2f6fc52SHong Zhang       output_file: output/ex89_3.out
156c2f6fc52SHong Zhang 
157c2f6fc52SHong Zhang    test:
158c2f6fc52SHong Zhang       suffix: nonscalable
159c2f6fc52SHong Zhang       nsize: 4
1603e662e0bSHong Zhang       args: -M 10 -N 10 -Z 10 -mat_product_algorithm nonscalable
161c2f6fc52SHong Zhang       output_file: output/ex89_5.out
162c4762a1bSJed Brown 
163c4762a1bSJed Brown TEST*/
164