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