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 5*9371c9d4SSatish Balay static PetscErrorCode EstimateMemory(DM dm, PetscLogDouble *est) { 642565d0bSMatthew G. Knepley DMLabel marker; 742565d0bSMatthew G. Knepley PetscInt cdim, depth, d, pStart, pEnd, p, Nd[4] = {0, 0, 0, 0}, lsize = 0, rmem = 0, imem = 0; 842565d0bSMatthew G. Knepley PetscInt coneSecMem = 0, coneMem = 0, supportSecMem = 0, supportMem = 0, labelMem = 0; 942565d0bSMatthew G. Knepley 1042565d0bSMatthew G. Knepley PetscFunctionBeginUser; 119566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Memory Estimates\n")); 129566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 139566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 149566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 1542565d0bSMatthew G. Knepley for (d = 0; d <= depth; ++d) { 1642565d0bSMatthew G. Knepley PetscInt start, end; 1742565d0bSMatthew G. Knepley 189566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, d, &start, &end)); 1942565d0bSMatthew G. Knepley Nd[d] = end - start; 2042565d0bSMatthew G. Knepley } 2142565d0bSMatthew G. Knepley /* Coordinates: 3 Nv reals + 2*Nv + 2*Nv ints */ 2242565d0bSMatthew G. Knepley rmem += cdim * Nd[0]; 2342565d0bSMatthew G. Knepley imem += 2 * Nd[0] + 2 * Nd[0]; 2463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Coordinate mem: %" PetscInt_FMT " %" PetscInt_FMT "\n", (PetscInt)(cdim * Nd[0] * sizeof(PetscReal)), (PetscInt)(4 * Nd[0] * sizeof(PetscInt)))); 2542565d0bSMatthew G. Knepley /* Depth: Nc+Nf+Ne+Nv ints */ 2642565d0bSMatthew G. Knepley for (d = 0; d <= depth; ++d) labelMem += Nd[d]; 2742565d0bSMatthew G. Knepley /* Cell Type: Nc+Nf+Ne+Nv ints */ 2842565d0bSMatthew G. Knepley for (d = 0; d <= depth; ++d) labelMem += Nd[d]; 2942565d0bSMatthew G. Knepley /* Marker */ 309566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &marker)); 319566063dSJacob Faibussowitsch if (marker) PetscCall(DMLabelGetStratumSize(marker, 1, &lsize)); 3242565d0bSMatthew G. Knepley labelMem += lsize; 3363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Label mem: %" PetscInt_FMT "\n", (PetscInt)(labelMem * sizeof(PetscInt)))); 3442565d0bSMatthew G. Knepley //imem += labelMem; 3542565d0bSMatthew G. Knepley /* Cones and Orientations: 4 Nc + 3 Nf + 2 Ne ints + (Nc+Nf+Ne) ints no separate orientation section */ 3642565d0bSMatthew G. Knepley for (d = 0; d <= depth; ++d) coneSecMem += 2 * Nd[d]; 3742565d0bSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 3842565d0bSMatthew G. Knepley PetscInt csize; 3942565d0bSMatthew G. Knepley 409566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, p, &csize)); 4142565d0bSMatthew G. Knepley coneMem += csize; 4242565d0bSMatthew G. Knepley } 4363a3b9bcSJacob Faibussowitsch 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)))); 4442565d0bSMatthew G. Knepley imem += 2 * coneMem + coneSecMem; 4542565d0bSMatthew G. Knepley /* Supports: 4 Nc + 3 Nf + 2 Ne ints + Nc+Nf+Ne ints */ 4642565d0bSMatthew G. Knepley for (d = 0; d <= depth; ++d) supportSecMem += 2 * Nd[d]; 4742565d0bSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 4842565d0bSMatthew G. Knepley PetscInt ssize; 4942565d0bSMatthew G. Knepley 509566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, p, &ssize)); 5142565d0bSMatthew G. Knepley supportMem += ssize; 5242565d0bSMatthew G. Knepley } 5363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Support mem: %" PetscInt_FMT " %" PetscInt_FMT "\n", (PetscInt)(supportMem * sizeof(PetscInt)), (PetscInt)(supportSecMem * sizeof(PetscInt)))); 5442565d0bSMatthew G. Knepley imem += supportMem + supportSecMem; 5542565d0bSMatthew G. Knepley *est = ((PetscLogDouble)imem) * sizeof(PetscInt) + ((PetscLogDouble)rmem) * sizeof(PetscReal); 5663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Estimated memory %" PetscInt_FMT "\n", (PetscInt)*est)); 5742565d0bSMatthew G. Knepley PetscFunctionReturn(0); 5842565d0bSMatthew G. Knepley } 5942565d0bSMatthew G. Knepley 60*9371c9d4SSatish Balay int main(int argc, char **argv) { 6142565d0bSMatthew G. Knepley DM dm; 6242565d0bSMatthew G. Knepley PetscBool trace = PETSC_FALSE, checkMemory = PETSC_TRUE, auxMemory = PETSC_FALSE; 6342565d0bSMatthew G. Knepley PetscLogDouble before, after, est = 0, clean, max; 6442565d0bSMatthew G. Knepley 65327415f7SBarry Smith PetscFunctionBeginUser; 669566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-trace", &trace, NULL)); 689566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-check_memory", &checkMemory, NULL)); 699566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-aux_memory", &auxMemory, NULL)); 709566063dSJacob Faibussowitsch PetscCall(PetscMemorySetGetMaximumUsage()); 719566063dSJacob Faibussowitsch PetscCall(PetscMallocGetCurrentUsage(&before)); 729566063dSJacob Faibussowitsch if (trace) PetscCall(PetscMallocTraceSet(NULL, PETSC_TRUE, 5000.)); 739566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 749566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 759566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 769566063dSJacob Faibussowitsch if (trace) PetscCall(PetscMallocTraceSet(NULL, PETSC_FALSE, 5000)); 779566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 789566063dSJacob Faibussowitsch PetscCall(PetscMallocGetCurrentUsage(&after)); 799566063dSJacob Faibussowitsch PetscCall(PetscMemoryGetMaximumUsage(&max)); 809566063dSJacob Faibussowitsch PetscCall(EstimateMemory(dm, &est)); 819566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 829566063dSJacob Faibussowitsch PetscCall(PetscMallocGetCurrentUsage(&clean)); 839566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Measured Memory\n")); 8442565d0bSMatthew G. Knepley if (auxMemory) { 85*9371c9d4SSatish Balay 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)); 8642565d0bSMatthew G. Knepley } 8742565d0bSMatthew G. Knepley if (checkMemory) { 8863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Memory for mesh %" PetscInt_FMT "\n", (PetscInt)(after - before))); 8963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Discrepancy %" PetscInt_FMT "\n", (PetscInt)PetscAbsReal(after - before - est))); 9042565d0bSMatthew G. Knepley } 919566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 92b122ec5aSJacob Faibussowitsch return 0; 9342565d0bSMatthew G. Knepley } 9442565d0bSMatthew G. Knepley 9542565d0bSMatthew G. Knepley /*TEST 9642565d0bSMatthew G. Knepley build: 97dfd57a17SPierre Jolivet requires: !defined(PETSC_USE_64BIT_INDICES) double !complex !defined(PETSCTEST_VALGRIND) 9842565d0bSMatthew G. Knepley 99a5b23f4aSJose E. Roman # Memory checks cannot be included in tests because the allocated memory differs among environments 10042565d0bSMatthew G. Knepley testset: 10142565d0bSMatthew G. Knepley args: -malloc_requested_size -dm_plex_box_faces 5,5 -check_memory 0 10242565d0bSMatthew G. Knepley test: 10342565d0bSMatthew G. Knepley suffix: tri 10442565d0bSMatthew G. Knepley requires: triangle 10530602db0SMatthew G. Knepley args: -dm_plex_simplex 1 -dm_plex_interpolate 0 10642565d0bSMatthew G. Knepley 10742565d0bSMatthew G. Knepley test: 10842565d0bSMatthew G. Knepley suffix: tri_interp 10942565d0bSMatthew G. Knepley requires: triangle 11030602db0SMatthew G. Knepley args: -dm_plex_simplex 1 -dm_plex_interpolate 1 11142565d0bSMatthew G. Knepley 11242565d0bSMatthew G. Knepley test: 11342565d0bSMatthew G. Knepley suffix: quad 11430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_interpolate 0 11542565d0bSMatthew G. Knepley 11642565d0bSMatthew G. Knepley test: 11742565d0bSMatthew G. Knepley suffix: quad_interp 11830602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_interpolate 1 11942565d0bSMatthew G. Knepley 12030602db0SMatthew G. Knepley # Memory checks cannot be included in tests because the allocated memory differs among environments 12142565d0bSMatthew G. Knepley testset: 12230602db0SMatthew G. Knepley args: -malloc_requested_size -dm_plex_dim 3 -dm_plex_box_faces 5,5,5 -check_memory 0 123da87d8e5SMatthew G. Knepley 124da87d8e5SMatthew G. Knepley # Filter out label memory because tet mesher produce different surface meshes for different compilers 12542565d0bSMatthew G. Knepley test: 12642565d0bSMatthew G. Knepley suffix: tet 12742565d0bSMatthew G. Knepley requires: ctetgen 128da87d8e5SMatthew G. Knepley filter: grep -v "Label mem:" 12930602db0SMatthew G. Knepley args: -dm_plex_simplex 1 -dm_plex_interpolate 0 13042565d0bSMatthew G. Knepley 131da87d8e5SMatthew G. Knepley # Filter out label memory because tet mesher produce different surface meshes for different compilers 13242565d0bSMatthew G. Knepley test: 13342565d0bSMatthew G. Knepley suffix: tet_interp 13442565d0bSMatthew G. Knepley requires: ctetgen 135da87d8e5SMatthew G. Knepley filter: grep -v "Label mem:" 13630602db0SMatthew G. Knepley args: -dm_plex_simplex 1 -dm_plex_interpolate 1 13742565d0bSMatthew G. Knepley 13842565d0bSMatthew G. Knepley test: 13942565d0bSMatthew G. Knepley suffix: hex 14030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_interpolate 0 14142565d0bSMatthew G. Knepley 14242565d0bSMatthew G. Knepley test: 14342565d0bSMatthew G. Knepley suffix: hex_interp 14430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_interpolate 1 14542565d0bSMatthew G. Knepley TEST*/ 146