1 static char help[] = "Tests adaptive refinement using DMForest, and uses HDF5.\n\n"; 2 3 #include <petscdmforest.h> 4 #include <petscdmplex.h> 5 #include <petscviewerhdf5.h> 6 7 int main (int argc, char **argv) 8 { 9 DM base, forest, plex; 10 PetscSection s; 11 PetscViewer viewer; 12 Vec g = NULL, g2 = NULL; 13 PetscReal nrm; 14 PetscBool adapt = PETSC_FALSE, userSection = PETSC_FALSE; 15 PetscInt vStart, vEnd, v, i; 16 17 CHKERRQ(PetscInitialize(&argc, &argv, NULL, help)); 18 CHKERRQ(PetscOptionsGetBool(NULL, NULL, "-adapt", &adapt, NULL)); 19 CHKERRQ(PetscOptionsGetBool(NULL, NULL, "-user_section", &userSection, NULL)); 20 21 /* Create a base DMPlex mesh */ 22 CHKERRQ(DMCreate(PETSC_COMM_WORLD, &base)); 23 CHKERRQ(DMSetType(base, DMPLEX)); 24 CHKERRQ(DMSetFromOptions(base)); 25 CHKERRQ(DMViewFromOptions(base, NULL, "-dm_view")); 26 27 /* Covert Plex mesh to Forest and destroy base */ 28 CHKERRQ(DMCreate(PETSC_COMM_WORLD, &forest)); 29 CHKERRQ(DMSetType(forest, DMP4EST)); 30 CHKERRQ(DMForestSetBaseDM(forest, base)); 31 CHKERRQ(DMSetUp(forest)); 32 CHKERRQ(DMDestroy(&base)); 33 CHKERRQ(DMViewFromOptions(forest, NULL, "-dm_view")); 34 35 if (adapt) { 36 /* Adaptively refine the cell 0 of the mesh */ 37 for (i = 0; i < 3; ++i) { 38 DM postforest; 39 DMLabel adaptLabel = NULL; 40 41 CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 42 CHKERRQ(DMLabelSetValue(adaptLabel, 0, DM_ADAPT_REFINE)); 43 CHKERRQ(DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest)); 44 CHKERRQ(DMForestSetAdaptivityLabel(postforest, adaptLabel)); 45 CHKERRQ(DMLabelDestroy(&adaptLabel)); 46 CHKERRQ(DMSetUp(postforest)); 47 CHKERRQ(DMDestroy(&forest)); 48 forest = postforest; 49 } 50 } else { 51 /* Adaptively refine all cells of the mesh */ 52 PetscInt cStart, cEnd, c; 53 54 for (i = 0; i < 3; ++i) { 55 DM postforest; 56 DMLabel adaptLabel = NULL; 57 58 CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 59 60 CHKERRQ(DMForestGetCellChart(forest, &cStart, &cEnd)); 61 for (c = cStart; c < cEnd; ++c) { 62 CHKERRQ(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE)); 63 } 64 65 CHKERRQ(DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest)); 66 CHKERRQ(DMForestSetAdaptivityLabel(postforest, adaptLabel)); 67 CHKERRQ(DMLabelDestroy(&adaptLabel)); 68 CHKERRQ(DMSetUp(postforest)); 69 CHKERRQ(DMDestroy(&forest)); 70 forest = postforest; 71 } 72 } 73 CHKERRQ(DMViewFromOptions(forest, NULL, "-dm_view")); 74 75 /* Setup the section*/ 76 if (userSection) { 77 CHKERRQ(DMConvert(forest, DMPLEX, &plex)); 78 CHKERRQ(DMViewFromOptions(plex, NULL, "-dm_view")); 79 CHKERRQ(DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd)); 80 CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "Vertices [%D, %D)\n", vStart, vEnd)); 81 CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 82 CHKERRQ(PetscSectionCreate(PetscObjectComm((PetscObject) forest), &s)); 83 CHKERRQ(PetscSectionSetNumFields(s, 1)); 84 CHKERRQ(PetscSectionSetChart(s, vStart, vEnd)); 85 for (v = vStart; v < vEnd; ++v) { 86 CHKERRQ(PetscSectionSetDof(s, v, 1)); 87 CHKERRQ(PetscSectionSetFieldDof(s, v, 0, 1)); 88 } 89 CHKERRQ(PetscSectionSetUp(s)); 90 CHKERRQ(DMSetLocalSection(forest, s)); 91 CHKERRQ(PetscObjectViewFromOptions((PetscObject) s, NULL, "-my_section_view")); 92 CHKERRQ(PetscSectionDestroy(&s)); 93 CHKERRQ(DMDestroy(&plex)); 94 } else { 95 PetscFE fe; 96 PetscInt dim; 97 98 CHKERRQ(DMGetDimension(forest, &dim)); 99 CHKERRQ(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, 1, PETSC_DETERMINE, &fe)); 100 CHKERRQ(DMAddField(forest, NULL, (PetscObject) fe)); 101 CHKERRQ(PetscFEDestroy(&fe)); 102 CHKERRQ(DMCreateDS(forest)); 103 } 104 105 /* Create the global vector*/ 106 CHKERRQ(DMCreateGlobalVector(forest, &g)); 107 CHKERRQ(PetscObjectSetName((PetscObject) g, "g")); 108 CHKERRQ(VecSet(g, 1.0)); 109 110 /* Test global to local*/ 111 Vec l; 112 CHKERRQ(DMCreateLocalVector(forest, &l)); 113 CHKERRQ(VecZeroEntries(l)); 114 CHKERRQ(DMGlobalToLocal(forest, g, INSERT_VALUES, l)); 115 CHKERRQ(VecZeroEntries(g)); 116 CHKERRQ(DMLocalToGlobal(forest, l, INSERT_VALUES, g)); 117 CHKERRQ(VecDestroy(&l)); 118 119 /* Save a vector*/ 120 CHKERRQ(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_WRITE, &viewer)); 121 CHKERRQ(VecView(g, viewer)); 122 CHKERRQ(PetscViewerDestroy(&viewer)); 123 124 /* Load another vector to load into*/ 125 CHKERRQ(DMCreateGlobalVector(forest, &g2)); 126 CHKERRQ(PetscObjectSetName((PetscObject) g2, "g")); 127 CHKERRQ(VecZeroEntries(g2)); 128 129 /* Load a vector*/ 130 CHKERRQ(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_READ, &viewer)); 131 CHKERRQ(VecLoad(g2, viewer)); 132 CHKERRQ(PetscViewerDestroy(&viewer)); 133 134 /* Check if the data is the same*/ 135 CHKERRQ(VecAXPY(g2, -1.0, g)); 136 CHKERRQ(VecNorm(g2, NORM_INFINITY, &nrm)); 137 PetscCheckFalse(PetscAbsReal(nrm) > PETSC_SMALL,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Invalid difference norm %g", (double) nrm); 138 139 CHKERRQ(VecDestroy(&g)); 140 CHKERRQ(VecDestroy(&g2)); 141 CHKERRQ(DMDestroy(&forest)); 142 CHKERRQ(PetscFinalize()); 143 return 0; 144 } 145 146 /*TEST 147 148 build: 149 requires: hdf5 p4est 150 151 test: 152 suffix: 0 153 nsize: {{1 2 5}} 154 args: -adapt -dm_plex_simplex 0 -dm_plex_box_faces 2,2 155 156 TEST*/ 157