xref: /petsc/src/dm/impls/plex/tests/ex35.c (revision d40579cb7c56c9875be29dee4bb7a2221ffffb4a)
1 static char help[] = "Exhaustive memory tracking for DMPlex.\n\n\n";
2 
3 #include <petscdmplex.h>
4 
5 static PetscErrorCode EstimateMemory(DM dm, PetscLogDouble *est)
6 {
7   DMLabel        marker;
8   PetscInt       cdim, depth, d, pStart, pEnd, p, Nd[4] = {0, 0, 0, 0}, lsize = 0, rmem = 0, imem = 0;
9   PetscInt       coneSecMem = 0, coneMem = 0, supportSecMem = 0, supportMem = 0, labelMem = 0;
10 
11   PetscFunctionBeginUser;
12   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Memory Estimates\n"));
13   PetscCall(DMGetCoordinateDim(dm, &cdim));
14   PetscCall(DMPlexGetDepth(dm, &depth));
15   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
16   for (d = 0; d <= depth; ++d) {
17     PetscInt start, end;
18 
19     PetscCall(DMPlexGetDepthStratum(dm, d, &start, &end));
20     Nd[d] = end - start;
21   }
22   /* Coordinates: 3 Nv reals + 2*Nv + 2*Nv ints */
23   rmem += cdim*Nd[0];
24   imem += 2*Nd[0] + 2*Nd[0];
25   PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Coordinate mem:  %" PetscInt_FMT " %" PetscInt_FMT "\n", (PetscInt)(cdim*Nd[0]*sizeof(PetscReal)), (PetscInt)(4*Nd[0]*sizeof(PetscInt))));
26   /* Depth:       Nc+Nf+Ne+Nv ints */
27   for (d = 0; d <= depth; ++d) labelMem += Nd[d];
28   /* Cell Type:   Nc+Nf+Ne+Nv ints */
29   for (d = 0; d <= depth; ++d) labelMem += Nd[d];
30   /* Marker */
31   PetscCall(DMGetLabel(dm, "marker", &marker));
32   if (marker) PetscCall(DMLabelGetStratumSize(marker, 1, &lsize));
33   labelMem += lsize;
34   PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Label mem:       %" PetscInt_FMT "\n", (PetscInt)(labelMem*sizeof(PetscInt))));
35   //imem += labelMem;
36   /* Cones and Orientations:       4 Nc + 3 Nf + 2 Ne ints + (Nc+Nf+Ne) ints no separate orientation section */
37   for (d = 0; d <= depth; ++d) coneSecMem += 2*Nd[d];
38   for (p = pStart; p < pEnd; ++p) {
39     PetscInt csize;
40 
41     PetscCall(DMPlexGetConeSize(dm, p, &csize));
42     coneMem += csize;
43   }
44   PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Cone mem:        %" PetscInt_FMT " %" PetscInt_FMT " (%" PetscInt_FMT ")\n", (PetscInt)(coneMem*sizeof(PetscInt)), (PetscInt)(coneSecMem*sizeof(PetscInt)), (PetscInt)(coneMem*sizeof(PetscInt))));
45   imem += 2*coneMem + coneSecMem;
46   /* Supports:       4 Nc + 3 Nf + 2 Ne ints + Nc+Nf+Ne ints */
47   for (d = 0; d <= depth; ++d) supportSecMem += 2*Nd[d];
48   for (p = pStart; p < pEnd; ++p) {
49     PetscInt ssize;
50 
51     PetscCall(DMPlexGetSupportSize(dm, p, &ssize));
52     supportMem += ssize;
53   }
54   PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Support mem:     %" PetscInt_FMT " %" PetscInt_FMT "\n", (PetscInt)(supportMem*sizeof(PetscInt)), (PetscInt)(supportSecMem*sizeof(PetscInt))));
55   imem += supportMem + supportSecMem;
56   *est = ((PetscLogDouble) imem)*sizeof(PetscInt) + ((PetscLogDouble) rmem)*sizeof(PetscReal);
57   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  Estimated memory %" PetscInt_FMT "\n", (PetscInt) *est));
58   PetscFunctionReturn(0);
59 }
60 
61 int main(int argc, char **argv)
62 {
63   DM             dm;
64   PetscBool      trace = PETSC_FALSE, checkMemory = PETSC_TRUE, auxMemory = PETSC_FALSE;
65   PetscLogDouble before, after, est = 0, clean, max;
66 
67   PetscFunctionBeginUser;
68   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
69   PetscCall(PetscOptionsGetBool(NULL, NULL, "-trace", &trace, NULL));
70   PetscCall(PetscOptionsGetBool(NULL, NULL, "-check_memory", &checkMemory, NULL));
71   PetscCall(PetscOptionsGetBool(NULL, NULL, "-aux_memory", &auxMemory, NULL));
72   PetscCall(PetscMemorySetGetMaximumUsage());
73   PetscCall(PetscMallocGetCurrentUsage(&before));
74   if (trace) PetscCall(PetscMallocTraceSet(NULL, PETSC_TRUE, 5000.));
75   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
76   PetscCall(DMSetType(dm, DMPLEX));
77   PetscCall(DMSetFromOptions(dm));
78   if (trace) PetscCall(PetscMallocTraceSet(NULL, PETSC_FALSE, 5000));
79   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
80   PetscCall(PetscMallocGetCurrentUsage(&after));
81   PetscCall(PetscMemoryGetMaximumUsage(&max));
82   PetscCall(EstimateMemory(dm, &est));
83   PetscCall(DMDestroy(&dm));
84   PetscCall(PetscMallocGetCurrentUsage(&clean));
85   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Measured Memory\n"));
86   if (auxMemory) {
87     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  Initial memory         %" PetscInt_FMT "\n  Extra memory for build %" PetscInt_FMT "\n  Memory after destroy   %" PetscInt_FMT "\n",
88                           (PetscInt) before, (PetscInt) (max-after), (PetscInt) clean));
89   }
90   if (checkMemory) {
91     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  Memory for mesh  %" PetscInt_FMT "\n", (PetscInt) (after-before)));
92     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Discrepancy %" PetscInt_FMT "\n", (PetscInt) PetscAbsReal(after-before-est)));
93   }
94   PetscCall(PetscFinalize());
95   return 0;
96 }
97 
98 /*TEST
99   build:
100     requires: !defined(PETSC_USE_64BIT_INDICES) double !complex !defined(PETSCTEST_VALGRIND)
101 
102   # Memory checks cannot be included in tests because the allocated memory differs among environments
103   testset:
104     args: -malloc_requested_size -dm_plex_box_faces 5,5 -check_memory 0
105     test:
106       suffix: tri
107       requires: triangle
108       args: -dm_plex_simplex 1 -dm_plex_interpolate 0
109 
110     test:
111       suffix: tri_interp
112       requires: triangle
113       args: -dm_plex_simplex 1 -dm_plex_interpolate 1
114 
115     test:
116       suffix: quad
117       args: -dm_plex_simplex 0 -dm_plex_interpolate 0
118 
119     test:
120       suffix: quad_interp
121       args: -dm_plex_simplex 0 -dm_plex_interpolate 1
122 
123   # Memory checks cannot be included in tests because the allocated memory differs among environments
124   testset:
125     args: -malloc_requested_size -dm_plex_dim 3 -dm_plex_box_faces 5,5,5 -check_memory 0
126 
127     # Filter out label memory because tet mesher produce different surface meshes for different compilers
128     test:
129       suffix: tet
130       requires: ctetgen
131       filter: grep -v "Label mem:"
132       args: -dm_plex_simplex 1 -dm_plex_interpolate 0
133 
134     # Filter out label memory because tet mesher produce different surface meshes for different compilers
135     test:
136       suffix: tet_interp
137       requires: ctetgen
138       filter: grep -v "Label mem:"
139       args: -dm_plex_simplex 1 -dm_plex_interpolate 1
140 
141     test:
142       suffix: hex
143       args: -dm_plex_simplex 0 -dm_plex_interpolate 0
144 
145     test:
146       suffix: hex_interp
147       args: -dm_plex_simplex 0 -dm_plex_interpolate 1
148 TEST*/
149