xref: /petsc/src/mat/tests/ex33.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static char help[] = "Test memory scalability of MatMatMult() for AIJ and DENSE matrices. \n\
2c4762a1bSJed Brown Modified from the code contributed by Ian Lin <iancclin@umich.edu> \n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown Example:
6c4762a1bSJed Brown   mpiexec -n <np> ./ex33 -mem_view -matmatmult_Bbn <Bbn>
7c4762a1bSJed Brown */
8c4762a1bSJed Brown 
9c4762a1bSJed Brown #include <petsc.h>
10c4762a1bSJed Brown 
119371c9d4SSatish Balay PetscErrorCode Print_memory(PetscLogDouble mem) {
12c4762a1bSJed Brown   double max_mem, min_mem;
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   PetscFunctionBeginUser;
159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Reduce(&mem, &max_mem, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD));
169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Reduce(&mem, &min_mem, 1, MPI_DOUBLE, MPI_MIN, 0, MPI_COMM_WORLD));
17c4762a1bSJed Brown   max_mem = max_mem / 1024.0 / 1024.0;
18c4762a1bSJed Brown   min_mem = min_mem / 1024.0 / 1024.0;
199566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(MPI_COMM_WORLD, " max and min memory across all processors %.4f Mb, %.4f Mb.\n", (double)max_mem, (double)min_mem));
20c4762a1bSJed Brown   PetscFunctionReturn(0);
21c4762a1bSJed Brown }
22c4762a1bSJed Brown 
23c4762a1bSJed Brown /*
24c4762a1bSJed Brown    Illustrate how to use MPI derived data types.
25c4762a1bSJed Brown    It would save memory significantly. See MatMPIDenseScatter()
26c4762a1bSJed Brown */
279371c9d4SSatish Balay PetscErrorCode TestMPIDerivedDataType() {
28c4762a1bSJed Brown   MPI_Datatype type1, type2, rtype1, rtype2;
29c4762a1bSJed Brown   PetscInt     i, j;
30c4762a1bSJed Brown   PetscScalar  buffer[24]; /* An array of 4 rows, 6 cols */
31c4762a1bSJed Brown   MPI_Status   status;
32c4762a1bSJed Brown   PetscMPIInt  rank, size, disp[2];
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   PetscFunctionBeginUser;
359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &size));
36be096a46SBarry Smith   PetscCheck(size >= 2, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, "Must use at least 2 processors");
379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank));
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   if (rank == 0) {
40c4762a1bSJed Brown     /* proc[0] sends 2 rows to proc[1] */
41c4762a1bSJed Brown     for (i = 0; i < 24; i++) buffer[i] = (PetscScalar)i;
42c4762a1bSJed Brown 
439371c9d4SSatish Balay     disp[0] = 0;
449371c9d4SSatish Balay     disp[1] = 2;
459566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_indexed_block(2, 1, (const PetscMPIInt *)disp, MPIU_SCALAR, &type1));
46c4762a1bSJed Brown     /* one column has 4 entries */
479566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_resized(type1, 0, 4 * sizeof(PetscScalar), &type2));
489566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&type2));
499566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(buffer, 6, type2, 1, 123, MPI_COMM_WORLD));
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   } else if (rank == 1) {
52c4762a1bSJed Brown     /* proc[1] receives 2 rows from proc[0], and put them into contiguous rows, starting at the row 1 (disp[0]) */
53c4762a1bSJed Brown     PetscInt blen = 2;
54c4762a1bSJed Brown     for (i = 0; i < 24; i++) buffer[i] = 0.0;
55c4762a1bSJed Brown 
56c4762a1bSJed Brown     disp[0] = 1;
579566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_indexed_block(1, blen, (const PetscMPIInt *)disp, MPIU_SCALAR, &rtype1));
589566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_create_resized(rtype1, 0, 4 * sizeof(PetscScalar), &rtype2));
59c4762a1bSJed Brown 
609566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&rtype2));
619566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Recv(buffer, 6, rtype2, 0, 123, MPI_COMM_WORLD, &status));
62c4762a1bSJed Brown     for (i = 0; i < 4; i++) {
63*48a46eb9SPierre Jolivet       for (j = 0; j < 6; j++) PetscCall(PetscPrintf(MPI_COMM_SELF, "  %g", (double)PetscRealPart(buffer[i + j * 4])));
649566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(MPI_COMM_SELF, "\n"));
65c4762a1bSJed Brown     }
66c4762a1bSJed Brown   }
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   if (rank == 0) {
699566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&type1));
709566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&type2));
71c4762a1bSJed Brown   } else if (rank == 1) {
729566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&rtype1));
739566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&rtype2));
74c4762a1bSJed Brown   }
759566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD));
76c4762a1bSJed Brown   PetscFunctionReturn(0);
77c4762a1bSJed Brown }
78c4762a1bSJed Brown 
799371c9d4SSatish Balay int main(int argc, char **args) {
80c4762a1bSJed Brown   PetscInt       mA = 2700, nX = 80, nz = 40;
81c4762a1bSJed Brown   /* PetscInt        mA=6,nX=5,nz=2; //small test */
82c4762a1bSJed Brown   PetscLogDouble mem;
83c4762a1bSJed Brown   Mat            A, X, Y;
84c4762a1bSJed Brown   PetscBool      flg = PETSC_FALSE;
85c4762a1bSJed Brown 
86327415f7SBarry Smith   PetscFunctionBeginUser;
879566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_mpiderivedtype", &flg, NULL));
89c4762a1bSJed Brown   if (flg) {
909566063dSJacob Faibussowitsch     PetscCall(TestMPIDerivedDataType());
919566063dSJacob Faibussowitsch     PetscCall(PetscFinalize());
92b122ec5aSJacob Faibussowitsch     return 0;
93c4762a1bSJed Brown   }
94c4762a1bSJed Brown 
959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-mem_view", &flg, NULL));
969566063dSJacob Faibussowitsch   PetscCall(PetscMemoryGetCurrentUsage(&mem));
97c4762a1bSJed Brown   if (flg) {
989566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(MPI_COMM_WORLD, "Before start,"));
999566063dSJacob Faibussowitsch     PetscCall(Print_memory(mem));
100c4762a1bSJed Brown   }
101c4762a1bSJed Brown 
1029566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, mA, mA, nz, PETSC_NULL, nz, PETSC_NULL, &A));
1039566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(A, PETSC_NULL));
1049566063dSJacob Faibussowitsch   PetscCall(PetscMemoryGetCurrentUsage(&mem));
105c4762a1bSJed Brown   if (flg) {
1069566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(MPI_COMM_WORLD, "After creating A,"));
1079566063dSJacob Faibussowitsch     PetscCall(Print_memory(mem));
108c4762a1bSJed Brown   }
109c4762a1bSJed Brown 
1109566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, mA, nX, PETSC_NULL, &X));
1119566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(X, PETSC_NULL));
1129566063dSJacob Faibussowitsch   PetscCall(PetscMemoryGetCurrentUsage(&mem));
113c4762a1bSJed Brown   if (flg) {
1149566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(MPI_COMM_WORLD, "After creating X,"));
1159566063dSJacob Faibussowitsch     PetscCall(Print_memory(mem));
116c4762a1bSJed Brown   }
117c4762a1bSJed Brown 
1189566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Y));
1199566063dSJacob Faibussowitsch   PetscCall(PetscMemoryGetCurrentUsage(&mem));
120c4762a1bSJed Brown   if (flg) {
1219566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(MPI_COMM_WORLD, "After MatMatMult,"));
1229566063dSJacob Faibussowitsch     PetscCall(Print_memory(mem));
123c4762a1bSJed Brown   }
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* Test reuse */
1269566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A, X, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y));
1279566063dSJacob Faibussowitsch   PetscCall(PetscMemoryGetCurrentUsage(&mem));
128c4762a1bSJed Brown   if (flg) {
1299566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(MPI_COMM_WORLD, "After reuse MatMatMult,"));
1309566063dSJacob Faibussowitsch     PetscCall(Print_memory(mem));
131c4762a1bSJed Brown   }
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   /* Check accuracy */
1349566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(A, X, Y, 10, &flg));
13528b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Error in MatMatMult()");
136c4762a1bSJed Brown 
1379566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1389566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&X));
1399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Y));
140c4762a1bSJed Brown 
1419566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
142b122ec5aSJacob Faibussowitsch   return 0;
143c4762a1bSJed Brown }
144c4762a1bSJed Brown 
145c4762a1bSJed Brown /*TEST
146c4762a1bSJed Brown 
147c4762a1bSJed Brown    test:
148c4762a1bSJed Brown       suffix: 1
149c4762a1bSJed Brown       nsize: 4
150c4762a1bSJed Brown       output_file: output/ex33.out
151c4762a1bSJed Brown 
152c4762a1bSJed Brown    test:
153c4762a1bSJed Brown       suffix: 2
154c4762a1bSJed Brown       nsize: 8
155c4762a1bSJed Brown       output_file: output/ex33.out
156c4762a1bSJed Brown 
157c4762a1bSJed Brown    test:
158c4762a1bSJed Brown       suffix: 3
159c4762a1bSJed Brown       nsize: 2
160c4762a1bSJed Brown       args: -test_mpiderivedtype
161c4762a1bSJed Brown       output_file: output/ex33_3.out
162c4762a1bSJed Brown 
163c4762a1bSJed Brown TEST*/
164