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