xref: /petsc/src/mat/tests/ex33.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
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 
11c4762a1bSJed Brown PetscErrorCode Print_memory(PetscLogDouble mem)
12c4762a1bSJed Brown {
13c4762a1bSJed Brown   PetscErrorCode ierr;
14c4762a1bSJed Brown   double         max_mem,min_mem;
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   PetscFunctionBeginUser;
17ffc4695bSBarry Smith   ierr = MPI_Reduce(&mem, &max_mem, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD);CHKERRMPI(ierr);
18ffc4695bSBarry Smith   ierr = MPI_Reduce(&mem, &min_mem, 1, MPI_DOUBLE, MPI_MIN, 0, MPI_COMM_WORLD);CHKERRMPI(ierr);
19c4762a1bSJed Brown   max_mem = max_mem / 1024.0 / 1024.0;
20c4762a1bSJed Brown   min_mem = min_mem / 1024.0 / 1024.0;
21c4762a1bSJed Brown   ierr = PetscPrintf(MPI_COMM_WORLD, " max and min memory across all processors %.4f Mb, %.4f Mb.\n", (double)max_mem,(double)min_mem);CHKERRQ(ierr);
22c4762a1bSJed Brown   PetscFunctionReturn(0);
23c4762a1bSJed Brown }
24c4762a1bSJed Brown 
25c4762a1bSJed Brown /*
26c4762a1bSJed Brown    Illustrate how to use MPI derived data types.
27c4762a1bSJed Brown    It would save memory significantly. See MatMPIDenseScatter()
28c4762a1bSJed Brown */
29c4762a1bSJed Brown PetscErrorCode TestMPIDerivedDataType()
30c4762a1bSJed Brown {
31c4762a1bSJed Brown   PetscErrorCode    ierr;
32c4762a1bSJed Brown   MPI_Datatype      type1, type2,rtype1,rtype2;
33c4762a1bSJed Brown   PetscInt          i,j;
34c4762a1bSJed Brown   PetscScalar       buffer[24]; /* An array of 4 rows, 6 cols */
35c4762a1bSJed Brown   MPI_Status        status;
36c4762a1bSJed Brown   PetscMPIInt       rank,size,disp[2];
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   PetscFunctionBeginUser;
39ffc4695bSBarry Smith   ierr = MPI_Comm_size(MPI_COMM_WORLD, &size);CHKERRMPI(ierr);
40c4762a1bSJed Brown   if (size < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must use at least 2 processors");
41ffc4695bSBarry Smith   ierr = MPI_Comm_rank(MPI_COMM_WORLD, &rank);CHKERRMPI(ierr);
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   if (rank == 0) {
44c4762a1bSJed Brown     /* proc[0] sends 2 rows to proc[1] */
45c4762a1bSJed Brown     for (i=0; i<24; i++) buffer[i] = (PetscScalar)i;
46c4762a1bSJed Brown 
47c4762a1bSJed Brown     disp[0] = 0;  disp[1] = 2;
48ffc4695bSBarry Smith     ierr = MPI_Type_create_indexed_block(2, 1, (const PetscMPIInt *)disp, MPIU_SCALAR, &type1);CHKERRMPI(ierr);
49c4762a1bSJed Brown     /* one column has 4 entries */
50ffc4695bSBarry Smith     ierr = MPI_Type_create_resized(type1,0,4*sizeof(PetscScalar),&type2);CHKERRMPI(ierr);
51ffc4695bSBarry Smith     ierr = MPI_Type_commit(&type2);CHKERRMPI(ierr);
52ffc4695bSBarry Smith     ierr = MPI_Send(buffer, 6, type2, 1, 123, MPI_COMM_WORLD);CHKERRMPI(ierr);
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   } else if (rank == 1) {
55c4762a1bSJed Brown     /* proc[1] receives 2 rows from proc[0], and put them into contiguous rows, starting at the row 1 (disp[0]) */
56c4762a1bSJed Brown     PetscInt blen = 2;
57c4762a1bSJed Brown     for (i=0; i<24; i++) buffer[i] = 0.0;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown     disp[0] = 1;
60ffc4695bSBarry Smith     ierr = MPI_Type_create_indexed_block(1, blen, (const PetscMPIInt *)disp, MPIU_SCALAR, &rtype1);CHKERRMPI(ierr);
61ffc4695bSBarry Smith     ierr = MPI_Type_create_resized(rtype1, 0, 4*sizeof(PetscScalar), &rtype2);CHKERRMPI(ierr);
62c4762a1bSJed Brown 
63ffc4695bSBarry Smith     ierr = MPI_Type_commit(&rtype2);CHKERRMPI(ierr);
64ffc4695bSBarry Smith     ierr = MPI_Recv(buffer, 6, rtype2, 0, 123, MPI_COMM_WORLD, &status);CHKERRMPI(ierr);
65c4762a1bSJed Brown     for (i=0; i<4; i++) {
66c4762a1bSJed Brown       for (j=0; j<6; j++) {
67*1e1ea65dSPierre Jolivet         ierr = PetscPrintf(MPI_COMM_SELF,"  %g", (double)PetscRealPart(buffer[i+j*4]));CHKERRQ(ierr);
68c4762a1bSJed Brown       }
69c4762a1bSJed Brown       ierr = PetscPrintf(MPI_COMM_SELF,"\n");CHKERRQ(ierr);
70c4762a1bSJed Brown     }
71c4762a1bSJed Brown   }
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   if (rank == 0) {
74ffc4695bSBarry Smith     ierr = MPI_Type_free(&type1);CHKERRMPI(ierr);
75ffc4695bSBarry Smith     ierr = MPI_Type_free(&type2);CHKERRMPI(ierr);
76c4762a1bSJed Brown   } else if (rank == 1) {
77ffc4695bSBarry Smith     ierr = MPI_Type_free(&rtype1);CHKERRMPI(ierr);
78ffc4695bSBarry Smith     ierr = MPI_Type_free(&rtype2);CHKERRMPI(ierr);
79c4762a1bSJed Brown   }
80ffc4695bSBarry Smith   ierr = MPI_Barrier(MPI_COMM_WORLD);CHKERRMPI(ierr);
81c4762a1bSJed Brown   PetscFunctionReturn(0);
82c4762a1bSJed Brown }
83c4762a1bSJed Brown 
84c4762a1bSJed Brown int main(int argc, char **args)
85c4762a1bSJed Brown {
86c4762a1bSJed Brown   PetscInt          mA = 2700,nX = 80,nz = 40;
87c4762a1bSJed Brown   /* PetscInt        mA=6,nX=5,nz=2; //small test */
88c4762a1bSJed Brown   PetscLogDouble    mem;
89c4762a1bSJed Brown   Mat               A,X,Y;
90c4762a1bSJed Brown   PetscErrorCode    ierr;
91c4762a1bSJed Brown   PetscBool         flg = PETSC_FALSE;
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
94c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-test_mpiderivedtype",&flg,NULL);CHKERRQ(ierr);
95c4762a1bSJed Brown   if (flg) {
96c4762a1bSJed Brown     ierr = TestMPIDerivedDataType();CHKERRQ(ierr);
97c4762a1bSJed Brown     ierr = PetscFinalize();
98c4762a1bSJed Brown     return ierr;
99c4762a1bSJed Brown   }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-mem_view",&flg,NULL);CHKERRQ(ierr);
102c4762a1bSJed Brown   ierr = PetscMemoryGetCurrentUsage(&mem);CHKERRQ(ierr);
103c4762a1bSJed Brown   if (flg) {
104c4762a1bSJed Brown     ierr = PetscPrintf(MPI_COMM_WORLD, "Before start,");CHKERRQ(ierr);
105c4762a1bSJed Brown     ierr = Print_memory(mem);CHKERRQ(ierr);
106c4762a1bSJed Brown   }
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   ierr = MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,mA,mA,nz,PETSC_NULL,nz,PETSC_NULL,&A);CHKERRQ(ierr);
109c4762a1bSJed Brown   ierr = MatSetRandom(A,PETSC_NULL);CHKERRQ(ierr);
110c4762a1bSJed Brown   ierr = PetscMemoryGetCurrentUsage(&mem);CHKERRQ(ierr);
111c4762a1bSJed Brown   if (flg) {
112c4762a1bSJed Brown     ierr = PetscPrintf(MPI_COMM_WORLD, "After creating A,");CHKERRQ(ierr);
113c4762a1bSJed Brown     ierr = Print_memory(mem);CHKERRQ(ierr);
114c4762a1bSJed Brown   }
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,mA,nX,PETSC_NULL,&X);CHKERRQ(ierr);
117c4762a1bSJed Brown   ierr = MatSetRandom(X,PETSC_NULL);CHKERRQ(ierr);
118c4762a1bSJed Brown   ierr = PetscMemoryGetCurrentUsage(&mem);CHKERRQ(ierr);
119c4762a1bSJed Brown   if (flg) {
120c4762a1bSJed Brown     ierr = PetscPrintf(MPI_COMM_WORLD, "After creating X,");CHKERRQ(ierr);
121c4762a1bSJed Brown     ierr = Print_memory(mem);CHKERRQ(ierr);
122c4762a1bSJed Brown   }
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   ierr = MatMatMult(A,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Y);CHKERRQ(ierr);
125c4762a1bSJed Brown   ierr = PetscMemoryGetCurrentUsage(&mem);CHKERRQ(ierr);
126c4762a1bSJed Brown   if (flg) {
127c4762a1bSJed Brown     ierr = PetscPrintf(MPI_COMM_WORLD, "After MatMatMult,");CHKERRQ(ierr);
128c4762a1bSJed Brown     ierr = Print_memory(mem);CHKERRQ(ierr);
129c4762a1bSJed Brown   }
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   /* Test reuse */
132c4762a1bSJed Brown   ierr = MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y);CHKERRQ(ierr);
133c4762a1bSJed Brown   ierr = PetscMemoryGetCurrentUsage(&mem);CHKERRQ(ierr);
134c4762a1bSJed Brown   if (flg) {
135c4762a1bSJed Brown     ierr = PetscPrintf(MPI_COMM_WORLD, "After reuse MatMatMult,");CHKERRQ(ierr);
136c4762a1bSJed Brown     ierr = Print_memory(mem);CHKERRQ(ierr);
137c4762a1bSJed Brown   }
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   /* Check accuracy */
140c4762a1bSJed Brown   ierr = MatMatMultEqual(A,X,Y,10,&flg);CHKERRQ(ierr);
141c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatMatMult()");
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
144c4762a1bSJed Brown   ierr = MatDestroy(&X);CHKERRQ(ierr);
145c4762a1bSJed Brown   ierr = MatDestroy(&Y);CHKERRQ(ierr);
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   ierr = PetscFinalize();
148c4762a1bSJed Brown   return ierr;
149c4762a1bSJed Brown }
150c4762a1bSJed Brown 
151c4762a1bSJed Brown /*TEST
152c4762a1bSJed Brown 
153c4762a1bSJed Brown    test:
154c4762a1bSJed Brown       suffix: 1
155c4762a1bSJed Brown       nsize: 4
156c4762a1bSJed Brown       output_file: output/ex33.out
157c4762a1bSJed Brown 
158c4762a1bSJed Brown    test:
159c4762a1bSJed Brown       suffix: 2
160c4762a1bSJed Brown       nsize: 8
161c4762a1bSJed Brown       output_file: output/ex33.out
162c4762a1bSJed Brown 
163c4762a1bSJed Brown    test:
164c4762a1bSJed Brown       suffix: 3
165c4762a1bSJed Brown       nsize: 2
166c4762a1bSJed Brown       args: -test_mpiderivedtype
167c4762a1bSJed Brown       output_file: output/ex33_3.out
168c4762a1bSJed Brown 
169c4762a1bSJed Brown TEST*/
170