xref: /petsc/src/dm/impls/plex/tests/ex35.c (revision d0609ced746bc51b019815ca91d747429db24893)
142565d0bSMatthew G. Knepley static char help[] = "Exhaustive memory tracking for DMPlex.\n\n\n";
242565d0bSMatthew G. Knepley 
342565d0bSMatthew G. Knepley #include <petscdmplex.h>
442565d0bSMatthew G. Knepley 
542565d0bSMatthew G. Knepley static PetscErrorCode EstimateMemory(DM dm, PetscLogDouble *est)
642565d0bSMatthew G. Knepley {
742565d0bSMatthew G. Knepley   DMLabel        marker;
842565d0bSMatthew G. Knepley   PetscInt       cdim, depth, d, pStart, pEnd, p, Nd[4] = {0, 0, 0, 0}, lsize = 0, rmem = 0, imem = 0;
942565d0bSMatthew G. Knepley   PetscInt       coneSecMem = 0, coneMem = 0, supportSecMem = 0, supportMem = 0, labelMem = 0;
1042565d0bSMatthew G. Knepley 
1142565d0bSMatthew G. Knepley   PetscFunctionBeginUser;
129566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Memory Estimates\n"));
139566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &cdim));
149566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
159566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1642565d0bSMatthew G. Knepley   for (d = 0; d <= depth; ++d) {
1742565d0bSMatthew G. Knepley     PetscInt start, end;
1842565d0bSMatthew G. Knepley 
199566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, d, &start, &end));
2042565d0bSMatthew G. Knepley     Nd[d] = end - start;
2142565d0bSMatthew G. Knepley   }
2242565d0bSMatthew G. Knepley   /* Coordinates: 3 Nv reals + 2*Nv + 2*Nv ints */
2342565d0bSMatthew G. Knepley   rmem += cdim*Nd[0];
2442565d0bSMatthew G. Knepley   imem += 2*Nd[0] + 2*Nd[0];
259566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Coordinate mem:  %D %D\n", cdim*Nd[0]*sizeof(PetscReal), 4*Nd[0]*sizeof(PetscInt)));
2642565d0bSMatthew G. Knepley   /* Depth:       Nc+Nf+Ne+Nv ints */
2742565d0bSMatthew G. Knepley   for (d = 0; d <= depth; ++d) labelMem += Nd[d];
2842565d0bSMatthew G. Knepley   /* Cell Type:   Nc+Nf+Ne+Nv ints */
2942565d0bSMatthew G. Knepley   for (d = 0; d <= depth; ++d) labelMem += Nd[d];
3042565d0bSMatthew G. Knepley   /* Marker */
319566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &marker));
329566063dSJacob Faibussowitsch   if (marker) PetscCall(DMLabelGetStratumSize(marker, 1, &lsize));
3342565d0bSMatthew G. Knepley   labelMem += lsize;
349566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Label mem:       %D\n", labelMem*sizeof(PetscInt)));
3542565d0bSMatthew G. Knepley   //imem += labelMem;
3642565d0bSMatthew G. Knepley   /* Cones and Orientations:       4 Nc + 3 Nf + 2 Ne ints + (Nc+Nf+Ne) ints no separate orientation section */
3742565d0bSMatthew G. Knepley   for (d = 0; d <= depth; ++d) coneSecMem += 2*Nd[d];
3842565d0bSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
3942565d0bSMatthew G. Knepley     PetscInt csize;
4042565d0bSMatthew G. Knepley 
419566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(dm, p, &csize));
4242565d0bSMatthew G. Knepley     coneMem += csize;
4342565d0bSMatthew G. Knepley   }
449566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Cone mem:        %D %D (%D)\n", coneMem*sizeof(PetscInt), coneSecMem*sizeof(PetscInt), coneMem*sizeof(PetscInt)));
4542565d0bSMatthew G. Knepley   imem += 2*coneMem + coneSecMem;
4642565d0bSMatthew G. Knepley   /* Supports:       4 Nc + 3 Nf + 2 Ne ints + Nc+Nf+Ne ints */
4742565d0bSMatthew G. Knepley   for (d = 0; d <= depth; ++d) supportSecMem += 2*Nd[d];
4842565d0bSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
4942565d0bSMatthew G. Knepley     PetscInt ssize;
5042565d0bSMatthew G. Knepley 
519566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupportSize(dm, p, &ssize));
5242565d0bSMatthew G. Knepley     supportMem += ssize;
5342565d0bSMatthew G. Knepley   }
549566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Support mem:     %D %D\n", supportMem*sizeof(PetscInt), supportSecMem*sizeof(PetscInt)));
5542565d0bSMatthew G. Knepley   imem += supportMem + supportSecMem;
5642565d0bSMatthew G. Knepley   *est = ((PetscLogDouble) imem)*sizeof(PetscInt) + ((PetscLogDouble) rmem)*sizeof(PetscReal);
579566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  Estimated memory %D\n", (PetscInt) *est));
5842565d0bSMatthew G. Knepley   PetscFunctionReturn(0);
5942565d0bSMatthew G. Knepley }
6042565d0bSMatthew G. Knepley 
6142565d0bSMatthew G. Knepley int main(int argc, char **argv)
6242565d0bSMatthew G. Knepley {
6342565d0bSMatthew G. Knepley   DM             dm;
6442565d0bSMatthew G. Knepley   PetscBool      trace = PETSC_FALSE, checkMemory = PETSC_TRUE, auxMemory = PETSC_FALSE;
6542565d0bSMatthew G. Knepley   PetscLogDouble before, after, est = 0, clean, max;
6642565d0bSMatthew G. Knepley 
679566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-trace", &trace, NULL));
699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-check_memory", &checkMemory, NULL));
709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-aux_memory", &auxMemory, NULL));
719566063dSJacob Faibussowitsch   PetscCall(PetscMemorySetGetMaximumUsage());
729566063dSJacob Faibussowitsch   PetscCall(PetscMallocGetCurrentUsage(&before));
739566063dSJacob Faibussowitsch   if (trace) PetscCall(PetscMallocTraceSet(NULL, PETSC_TRUE, 5000.));
749566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
759566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
769566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
779566063dSJacob Faibussowitsch   if (trace) PetscCall(PetscMallocTraceSet(NULL, PETSC_FALSE, 5000));
789566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
799566063dSJacob Faibussowitsch   PetscCall(PetscMallocGetCurrentUsage(&after));
809566063dSJacob Faibussowitsch   PetscCall(PetscMemoryGetMaximumUsage(&max));
819566063dSJacob Faibussowitsch   PetscCall(EstimateMemory(dm, &est));
829566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
839566063dSJacob Faibussowitsch   PetscCall(PetscMallocGetCurrentUsage(&clean));
849566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Measured Memory\n"));
8542565d0bSMatthew G. Knepley   if (auxMemory) {
86*d0609cedSBarry Smith     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  Initial memory         %D\n  Extra memory for build %D\n  Memory after destroy   %D\n",
87*d0609cedSBarry Smith                           (PetscInt) before, (PetscInt) (max-after), (PetscInt) clean));
8842565d0bSMatthew G. Knepley   }
8942565d0bSMatthew G. Knepley   if (checkMemory) {
909566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  Memory for mesh  %D\n", (PetscInt) (after-before)));
919566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Discrepancy %D\n", (PetscInt) PetscAbsReal(after-before-est)));
9242565d0bSMatthew G. Knepley   }
939566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
94b122ec5aSJacob Faibussowitsch   return 0;
9542565d0bSMatthew G. Knepley }
9642565d0bSMatthew G. Knepley 
9742565d0bSMatthew G. Knepley /*TEST
9842565d0bSMatthew G. Knepley   build:
99dfd57a17SPierre Jolivet     requires: !defined(PETSC_USE_64BIT_INDICES) double !complex !defined(PETSCTEST_VALGRIND)
10042565d0bSMatthew G. Knepley 
101a5b23f4aSJose E. Roman   # Memory checks cannot be included in tests because the allocated memory differs among environments
10242565d0bSMatthew G. Knepley   testset:
10342565d0bSMatthew G. Knepley     args: -malloc_requested_size -dm_plex_box_faces 5,5 -check_memory 0
10442565d0bSMatthew G. Knepley     test:
10542565d0bSMatthew G. Knepley       suffix: tri
10642565d0bSMatthew G. Knepley       requires: triangle
10730602db0SMatthew G. Knepley       args: -dm_plex_simplex 1 -dm_plex_interpolate 0
10842565d0bSMatthew G. Knepley 
10942565d0bSMatthew G. Knepley     test:
11042565d0bSMatthew G. Knepley       suffix: tri_interp
11142565d0bSMatthew G. Knepley       requires: triangle
11230602db0SMatthew G. Knepley       args: -dm_plex_simplex 1 -dm_plex_interpolate 1
11342565d0bSMatthew G. Knepley 
11442565d0bSMatthew G. Knepley     test:
11542565d0bSMatthew G. Knepley       suffix: quad
11630602db0SMatthew G. Knepley       args: -dm_plex_simplex 0 -dm_plex_interpolate 0
11742565d0bSMatthew G. Knepley 
11842565d0bSMatthew G. Knepley     test:
11942565d0bSMatthew G. Knepley       suffix: quad_interp
12030602db0SMatthew G. Knepley       args: -dm_plex_simplex 0 -dm_plex_interpolate 1
12142565d0bSMatthew G. Knepley 
12230602db0SMatthew G. Knepley   # Memory checks cannot be included in tests because the allocated memory differs among environments
12342565d0bSMatthew G. Knepley   testset:
12430602db0SMatthew G. Knepley     args: -malloc_requested_size -dm_plex_dim 3 -dm_plex_box_faces 5,5,5 -check_memory 0
12542565d0bSMatthew G. Knepley     test:
12642565d0bSMatthew G. Knepley       suffix: tet
12742565d0bSMatthew G. Knepley       requires: ctetgen
12830602db0SMatthew G. Knepley       args: -dm_plex_simplex 1 -dm_plex_interpolate 0
12942565d0bSMatthew G. Knepley 
13042565d0bSMatthew G. Knepley     test:
13142565d0bSMatthew G. Knepley       suffix: tet_interp
13242565d0bSMatthew G. Knepley       requires: ctetgen
13330602db0SMatthew G. Knepley       args: -dm_plex_simplex 1 -dm_plex_interpolate 1
13442565d0bSMatthew G. Knepley 
13542565d0bSMatthew G. Knepley     test:
13642565d0bSMatthew G. Knepley       suffix: hex
13730602db0SMatthew G. Knepley       args: -dm_plex_simplex 0 -dm_plex_interpolate 0
13842565d0bSMatthew G. Knepley 
13942565d0bSMatthew G. Knepley     test:
14042565d0bSMatthew G. Knepley       suffix: hex_interp
14130602db0SMatthew G. Knepley       args: -dm_plex_simplex 0 -dm_plex_interpolate 1
14242565d0bSMatthew G. Knepley TEST*/
143