1c4762a1bSJed Brown static char help[] = "Tests MatPtAP() for MPIMAIJ and MPIAIJ \n ";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscdmda.h>
4c4762a1bSJed Brown
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
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
19327415f7SBarry 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));
249371c9d4SSatish Balay M = 10;
259371c9d4SSatish Balay N = 10;
269371c9d4SSatish Balay Z = 10;
27c4762a1bSJed Brown dof = 10;
28c4762a1bSJed Brown
299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_3D", &Test_3D, NULL));
309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL));
319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL));
329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-Z", &Z, NULL));
33c4762a1bSJed Brown /* Set up distributed array for fine grid */
34c4762a1bSJed Brown if (!Test_3D) {
359566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, M, N, PETSC_DECIDE, PETSC_DECIDE, dof, 1, NULL, NULL, &coarsedm));
36c4762a1bSJed Brown } else {
379566063dSJacob 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));
38c4762a1bSJed Brown }
399566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(coarsedm));
409566063dSJacob Faibussowitsch PetscCall(DMSetUp(coarsedm));
41c4762a1bSJed Brown
42c4762a1bSJed Brown /* This makes sure the coarse DMDA has the same partition as the fine DMDA */
439566063dSJacob Faibussowitsch PetscCall(DMRefine(coarsedm, PetscObjectComm((PetscObject)coarsedm), &finedm));
44c4762a1bSJed Brown
45c4762a1bSJed Brown /*------------------------------------------------------------*/
469566063dSJacob Faibussowitsch PetscCall(DMSetMatType(finedm, MATAIJ));
479566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(finedm, &A));
48c4762a1bSJed Brown
49c4762a1bSJed Brown /* set val=one to A */
50c4762a1bSJed Brown if (size == 1) {
519566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg));
52c4762a1bSJed Brown if (flg) {
539566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArray(A, &array));
54c4762a1bSJed Brown for (i = 0; i < ia[nrows]; i++) array[i] = one;
559566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArray(A, &array));
56c4762a1bSJed Brown }
579566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg));
58c4762a1bSJed Brown } else {
59c4762a1bSJed Brown Mat AA, AB;
609566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(A, &AA, &AB, NULL));
619566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(AA, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg));
62c4762a1bSJed Brown if (flg) {
639566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArray(AA, &array));
64c4762a1bSJed Brown for (i = 0; i < ia[nrows]; i++) array[i] = one;
659566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArray(AA, &array));
66c4762a1bSJed Brown }
679566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(AA, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg));
689566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(AB, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg));
69c4762a1bSJed Brown if (flg) {
709566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArray(AB, &array));
71c4762a1bSJed Brown for (i = 0; i < ia[nrows]; i++) array[i] = one;
729566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArray(AB, &array));
73c4762a1bSJed Brown }
749566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(AB, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg));
75c4762a1bSJed Brown }
76c4762a1bSJed Brown /* Create interpolation between the fine and coarse grids */
779566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(coarsedm, finedm, &P, NULL));
78c4762a1bSJed Brown
79c4762a1bSJed Brown /* Test P^T * A * P - MatPtAP() */
80c4762a1bSJed Brown /*------------------------------*/
81c20d7725SJed Brown /* (1) Developer API */
829566063dSJacob Faibussowitsch PetscCall(MatProductCreate(A, P, NULL, &C));
839566063dSJacob Faibussowitsch PetscCall(MatProductSetType(C, MATPRODUCT_PtAP));
849566063dSJacob Faibussowitsch PetscCall(MatProductSetAlgorithm(C, "allatonce"));
859566063dSJacob Faibussowitsch PetscCall(MatProductSetFill(C, PETSC_DEFAULT));
869566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(C));
879566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(C));
889566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(C));
899566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(C)); /* Test reuse of symbolic C */
90c20d7725SJed Brown
91c2f6fc52SHong Zhang { /* Test MatProductView() */
92c2f6fc52SHong Zhang PetscViewer viewer;
939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(comm, NULL, &viewer));
949566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
959566063dSJacob Faibussowitsch PetscCall(MatProductView(C, viewer));
969566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer));
979566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer));
98c2f6fc52SHong Zhang }
99c2f6fc52SHong Zhang
1009566063dSJacob Faibussowitsch PetscCall(MatPtAPMultEqual(A, P, C, 10, &flg));
10128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in MatProduct_PtAP");
1029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
103c20d7725SJed Brown
104c20d7725SJed Brown /* (2) User API */
1059566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_INITIAL_MATRIX, fill, &C));
106c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */
107c4762a1bSJed Brown alpha = 1.0;
108c4762a1bSJed Brown for (i = 0; i < 1; i++) {
109c4762a1bSJed Brown alpha -= 0.1;
1109566063dSJacob Faibussowitsch PetscCall(MatScale(A, alpha));
1119566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_REUSE_MATRIX, fill, &C));
112c4762a1bSJed Brown }
113c4762a1bSJed Brown
114c4762a1bSJed Brown /* Free intermediate data structures created for reuse of C=Pt*A*P */
1159566063dSJacob Faibussowitsch PetscCall(MatProductClear(C));
116c4762a1bSJed Brown
1179566063dSJacob Faibussowitsch PetscCall(MatPtAPMultEqual(A, P, C, 10, &flg));
11828b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in MatPtAP");
119c4762a1bSJed Brown
1209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
1219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
1229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P));
1239566063dSJacob Faibussowitsch PetscCall(DMDestroy(&finedm));
1249566063dSJacob Faibussowitsch PetscCall(DMDestroy(&coarsedm));
1259566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
126b122ec5aSJacob Faibussowitsch return 0;
127c4762a1bSJed Brown }
128c4762a1bSJed Brown
129c4762a1bSJed Brown /*TEST
130c4762a1bSJed Brown
131c4762a1bSJed Brown test:
132c4762a1bSJed Brown args: -M 10 -N 10 -Z 10
133*3886731fSPierre Jolivet output_file: output/empty.out
134c4762a1bSJed Brown
135c4762a1bSJed Brown test:
136c4762a1bSJed Brown suffix: allatonce
137c4762a1bSJed Brown nsize: 4
138c2f6fc52SHong Zhang args: -M 10 -N 10 -Z 10
139c2f6fc52SHong Zhang output_file: output/ex89_2.out
140c4762a1bSJed Brown
141c4762a1bSJed Brown test:
142c4762a1bSJed Brown suffix: allatonce_merged
143c4762a1bSJed Brown nsize: 4
1443e662e0bSHong Zhang args: -M 10 -M 5 -M 10 -mat_product_algorithm allatonce_merged
145c2f6fc52SHong Zhang output_file: output/ex89_3.out
146c4762a1bSJed Brown
147c4762a1bSJed Brown test:
148c2f6fc52SHong Zhang suffix: nonscalable_3D
149c4762a1bSJed Brown nsize: 4
1503e662e0bSHong Zhang args: -M 10 -M 5 -M 10 -test_3D 1 -mat_product_algorithm nonscalable
151c2f6fc52SHong Zhang output_file: output/ex89_4.out
152c4762a1bSJed Brown
153c4762a1bSJed Brown test:
154c4762a1bSJed Brown suffix: allatonce_merged_3D
155c4762a1bSJed Brown nsize: 4
1563e662e0bSHong Zhang args: -M 10 -M 5 -M 10 -test_3D 1 -mat_product_algorithm allatonce_merged
157c2f6fc52SHong Zhang output_file: output/ex89_3.out
158c2f6fc52SHong Zhang
159c2f6fc52SHong Zhang test:
160c2f6fc52SHong Zhang suffix: nonscalable
161c2f6fc52SHong Zhang nsize: 4
1623e662e0bSHong Zhang args: -M 10 -N 10 -Z 10 -mat_product_algorithm nonscalable
163c2f6fc52SHong Zhang output_file: output/ex89_5.out
164c4762a1bSJed Brown
165c4762a1bSJed Brown TEST*/
166