xref: /petsc/src/dm/impls/forest/tests/ex3.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
19b5c7ca5SMatthew G. Knepley static char help[] = "Tests adaptive refinement using DMForest, and uses HDF5.\n\n";
29b5c7ca5SMatthew G. Knepley 
39b5c7ca5SMatthew G. Knepley #include <petscdmforest.h>
49b5c7ca5SMatthew G. Knepley #include <petscdmplex.h>
59b5c7ca5SMatthew G. Knepley #include <petscviewerhdf5.h>
69b5c7ca5SMatthew G. Knepley 
79b5c7ca5SMatthew G. Knepley int main (int argc, char **argv)
89b5c7ca5SMatthew G. Knepley {
99b5c7ca5SMatthew G. Knepley   DM             base, forest, plex;
109b5c7ca5SMatthew G. Knepley   PetscSection   s;
119b5c7ca5SMatthew G. Knepley   PetscViewer    viewer;
129b5c7ca5SMatthew G. Knepley   Vec            g = NULL, g2 = NULL;
139b5c7ca5SMatthew G. Knepley   PetscReal      nrm;
149b5c7ca5SMatthew G. Knepley   PetscBool      adapt = PETSC_FALSE, userSection = PETSC_FALSE;
159b5c7ca5SMatthew G. Knepley   PetscInt       vStart, vEnd, v, i;
169b5c7ca5SMatthew G. Knepley   PetscErrorCode ierr;
179b5c7ca5SMatthew G. Knepley 
189b5c7ca5SMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
199b5c7ca5SMatthew G. Knepley   ierr = PetscOptionsGetBool(NULL, NULL, "-adapt", &adapt, NULL);CHKERRQ(ierr);
209b5c7ca5SMatthew G. Knepley   ierr = PetscOptionsGetBool(NULL, NULL, "-user_section", &userSection, NULL);CHKERRQ(ierr);
219b5c7ca5SMatthew G. Knepley 
229b5c7ca5SMatthew G. Knepley   /* Create a base DMPlex mesh */
2330602db0SMatthew G. Knepley   ierr = DMCreate(PETSC_COMM_WORLD, &base);CHKERRQ(ierr);
2430602db0SMatthew G. Knepley   ierr = DMSetType(base, DMPLEX);CHKERRQ(ierr);
259b5c7ca5SMatthew G. Knepley   ierr = DMSetFromOptions(base);CHKERRQ(ierr);
269b5c7ca5SMatthew G. Knepley   ierr = DMViewFromOptions(base, NULL, "-dm_view");CHKERRQ(ierr);
279b5c7ca5SMatthew G. Knepley 
289b5c7ca5SMatthew G. Knepley   /* Covert Plex mesh to Forest and destroy base */
299b5c7ca5SMatthew G. Knepley   ierr = DMCreate(PETSC_COMM_WORLD, &forest);CHKERRQ(ierr);
309b5c7ca5SMatthew G. Knepley   ierr = DMSetType(forest, DMP4EST);CHKERRQ(ierr);
319b5c7ca5SMatthew G. Knepley   ierr = DMForestSetBaseDM(forest, base);CHKERRQ(ierr);
329b5c7ca5SMatthew G. Knepley   ierr = DMSetUp(forest);CHKERRQ(ierr);
339b5c7ca5SMatthew G. Knepley   ierr = DMDestroy(&base);CHKERRQ(ierr);
349b5c7ca5SMatthew G. Knepley   ierr = DMViewFromOptions(forest, NULL, "-dm_view");CHKERRQ(ierr);
359b5c7ca5SMatthew G. Knepley 
369b5c7ca5SMatthew G. Knepley   if (adapt) {
379b5c7ca5SMatthew G. Knepley     /* Adaptively refine the cell 0 of the mesh */
389b5c7ca5SMatthew G. Knepley     for (i = 0; i < 3; ++i) {
399b5c7ca5SMatthew G. Knepley       DM      postforest;
409b5c7ca5SMatthew G. Knepley       DMLabel adaptLabel = NULL;
419b5c7ca5SMatthew G. Knepley 
429b5c7ca5SMatthew G. Knepley       ierr = DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel);CHKERRQ(ierr);
439b5c7ca5SMatthew G. Knepley       ierr = DMLabelSetValue(adaptLabel, 0, DM_ADAPT_REFINE);CHKERRQ(ierr);
449b5c7ca5SMatthew G. Knepley       ierr = DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest);CHKERRQ(ierr);
459b5c7ca5SMatthew G. Knepley       ierr = DMForestSetAdaptivityLabel(postforest, adaptLabel);CHKERRQ(ierr);
469b5c7ca5SMatthew G. Knepley       ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr);
479b5c7ca5SMatthew G. Knepley       ierr = DMSetUp(postforest);CHKERRQ(ierr);
489b5c7ca5SMatthew G. Knepley       ierr = DMDestroy(&forest);CHKERRQ(ierr);
499b5c7ca5SMatthew G. Knepley       forest = postforest;
509b5c7ca5SMatthew G. Knepley     }
519b5c7ca5SMatthew G. Knepley   } else {
529b5c7ca5SMatthew G. Knepley     /* Adaptively refine all cells of the mesh */
539b5c7ca5SMatthew G. Knepley     PetscInt cStart, cEnd, c;
549b5c7ca5SMatthew G. Knepley 
559b5c7ca5SMatthew G. Knepley     for (i = 0; i < 3; ++i) {
569b5c7ca5SMatthew G. Knepley       DM      postforest;
579b5c7ca5SMatthew G. Knepley       DMLabel adaptLabel = NULL;
589b5c7ca5SMatthew G. Knepley 
599b5c7ca5SMatthew G. Knepley       ierr = DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel);CHKERRQ(ierr);
609b5c7ca5SMatthew G. Knepley 
619b5c7ca5SMatthew G. Knepley       ierr = DMForestGetCellChart(forest, &cStart, &cEnd);CHKERRQ(ierr);
629b5c7ca5SMatthew G. Knepley       for (c = cStart; c < cEnd; ++c) {
639b5c7ca5SMatthew G. Knepley         ierr = DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE);CHKERRQ(ierr);
649b5c7ca5SMatthew G. Knepley       }
659b5c7ca5SMatthew G. Knepley 
669b5c7ca5SMatthew G. Knepley       ierr = DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest);CHKERRQ(ierr);
679b5c7ca5SMatthew G. Knepley       ierr = DMForestSetAdaptivityLabel(postforest, adaptLabel);CHKERRQ(ierr);
689b5c7ca5SMatthew G. Knepley       ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr);
699b5c7ca5SMatthew G. Knepley       ierr = DMSetUp(postforest);CHKERRQ(ierr);
709b5c7ca5SMatthew G. Knepley       ierr = DMDestroy(&forest);CHKERRQ(ierr);
719b5c7ca5SMatthew G. Knepley       forest = postforest;
729b5c7ca5SMatthew G. Knepley     }
739b5c7ca5SMatthew G. Knepley   }
749b5c7ca5SMatthew G. Knepley   ierr = DMViewFromOptions(forest, NULL, "-dm_view");CHKERRQ(ierr);
759b5c7ca5SMatthew G. Knepley 
769b5c7ca5SMatthew G. Knepley   /*  Setup the section*/
779b5c7ca5SMatthew G. Knepley   if (userSection) {
789b5c7ca5SMatthew G. Knepley     ierr = DMConvert(forest, DMPLEX, &plex);CHKERRQ(ierr);
799b5c7ca5SMatthew G. Knepley     ierr = DMViewFromOptions(plex, NULL, "-dm_view");CHKERRQ(ierr);
809b5c7ca5SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd);CHKERRQ(ierr);
819b5c7ca5SMatthew G. Knepley     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "Vertices [%D, %D)\n", vStart, vEnd);CHKERRQ(ierr);
829b5c7ca5SMatthew G. Knepley     ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL);CHKERRQ(ierr);
839b5c7ca5SMatthew G. Knepley     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) forest), &s);CHKERRQ(ierr);
849b5c7ca5SMatthew G. Knepley     ierr = PetscSectionSetNumFields(s, 1);CHKERRQ(ierr);
859b5c7ca5SMatthew G. Knepley     ierr = PetscSectionSetChart(s, vStart, vEnd);CHKERRQ(ierr);
869b5c7ca5SMatthew G. Knepley     for (v = vStart; v < vEnd; ++v) {
879b5c7ca5SMatthew G. Knepley       ierr = PetscSectionSetDof(s, v, 1);CHKERRQ(ierr);
889b5c7ca5SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(s, v, 0, 1);CHKERRQ(ierr);
899b5c7ca5SMatthew G. Knepley     }
909b5c7ca5SMatthew G. Knepley     ierr = PetscSectionSetUp(s);CHKERRQ(ierr);
919b5c7ca5SMatthew G. Knepley     ierr = DMSetLocalSection(forest, s);CHKERRQ(ierr);
929b5c7ca5SMatthew G. Knepley     ierr = PetscObjectViewFromOptions((PetscObject) s, NULL, "-my_section_view");CHKERRQ(ierr);
939b5c7ca5SMatthew G. Knepley     ierr = PetscSectionDestroy(&s);CHKERRQ(ierr);
949b5c7ca5SMatthew G. Knepley     ierr = DMDestroy(&plex);CHKERRQ(ierr);
959b5c7ca5SMatthew G. Knepley   } else {
969b5c7ca5SMatthew G. Knepley     PetscFE  fe;
979b5c7ca5SMatthew G. Knepley     PetscInt dim;
989b5c7ca5SMatthew G. Knepley 
999b5c7ca5SMatthew G. Knepley     ierr = DMGetDimension(forest, &dim);CHKERRQ(ierr);
1009b5c7ca5SMatthew G. Knepley     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
1019b5c7ca5SMatthew G. Knepley     ierr = DMAddField(forest, NULL, (PetscObject) fe);CHKERRQ(ierr);
1029b5c7ca5SMatthew G. Knepley     ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
1039b5c7ca5SMatthew G. Knepley     ierr = DMCreateDS(forest);CHKERRQ(ierr);
1049b5c7ca5SMatthew G. Knepley   }
1059b5c7ca5SMatthew G. Knepley 
1069b5c7ca5SMatthew G. Knepley   /* Create the global vector*/
1079b5c7ca5SMatthew G. Knepley   ierr = DMCreateGlobalVector(forest, &g);CHKERRQ(ierr);
1089b5c7ca5SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) g, "g");CHKERRQ(ierr);
1099b5c7ca5SMatthew G. Knepley   ierr = VecSet(g, 1.0);CHKERRQ(ierr);
1109b5c7ca5SMatthew G. Knepley 
1119b5c7ca5SMatthew G. Knepley   /* Test global to local*/
1129b5c7ca5SMatthew G. Knepley   Vec l;
1139b5c7ca5SMatthew G. Knepley   ierr = DMCreateLocalVector(forest, &l);CHKERRQ(ierr);
1149b5c7ca5SMatthew G. Knepley   ierr = VecZeroEntries(l);CHKERRQ(ierr);
1159b5c7ca5SMatthew G. Knepley   ierr = DMGlobalToLocal(forest, g, INSERT_VALUES, l);CHKERRQ(ierr);
1169b5c7ca5SMatthew G. Knepley   ierr = VecZeroEntries(g);CHKERRQ(ierr);
1179b5c7ca5SMatthew G. Knepley   ierr = DMLocalToGlobal(forest, l, INSERT_VALUES, g);CHKERRQ(ierr);
1189b5c7ca5SMatthew G. Knepley   ierr = VecDestroy(&l);CHKERRQ(ierr);
1199b5c7ca5SMatthew G. Knepley 
1209b5c7ca5SMatthew G. Knepley   /*  Save a vector*/
1211e1ea65dSPierre Jolivet   ierr = PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_WRITE, &viewer);CHKERRQ(ierr);
1229b5c7ca5SMatthew G. Knepley   ierr = VecView(g, viewer);CHKERRQ(ierr);
1239b5c7ca5SMatthew G. Knepley   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1249b5c7ca5SMatthew G. Knepley 
1259b5c7ca5SMatthew G. Knepley   /* Load another vector to load into*/
1269b5c7ca5SMatthew G. Knepley   ierr = DMCreateGlobalVector(forest, &g2);CHKERRQ(ierr);
1279b5c7ca5SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) g2, "g");CHKERRQ(ierr);
1289b5c7ca5SMatthew G. Knepley   ierr = VecZeroEntries(g2);CHKERRQ(ierr);
1299b5c7ca5SMatthew G. Knepley 
1309b5c7ca5SMatthew G. Knepley   /*  Load a vector*/
1311e1ea65dSPierre Jolivet   ierr = PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_READ, &viewer);CHKERRQ(ierr);
1329b5c7ca5SMatthew G. Knepley   ierr = VecLoad(g2, viewer);CHKERRQ(ierr);
1339b5c7ca5SMatthew G. Knepley   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1349b5c7ca5SMatthew G. Knepley 
1359b5c7ca5SMatthew G. Knepley   /*  Check if the data is the same*/
1369b5c7ca5SMatthew G. Knepley   ierr = VecAXPY(g2, -1.0, g);CHKERRQ(ierr);
1379b5c7ca5SMatthew G. Knepley   ierr = VecNorm(g2, NORM_INFINITY, &nrm);CHKERRQ(ierr);
138*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(PetscAbsReal(nrm) > PETSC_SMALL,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Invalid difference norm %g", (double) nrm);
1399b5c7ca5SMatthew G. Knepley 
1409b5c7ca5SMatthew G. Knepley   ierr = VecDestroy(&g);CHKERRQ(ierr);
1419b5c7ca5SMatthew G. Knepley   ierr = VecDestroy(&g2);CHKERRQ(ierr);
1429b5c7ca5SMatthew G. Knepley   ierr = DMDestroy(&forest);CHKERRQ(ierr);
1439b5c7ca5SMatthew G. Knepley   ierr = PetscFinalize();
1449b5c7ca5SMatthew G. Knepley   return ierr;
1459b5c7ca5SMatthew G. Knepley }
1469b5c7ca5SMatthew G. Knepley 
1479b5c7ca5SMatthew G. Knepley /*TEST
1489b5c7ca5SMatthew G. Knepley 
1499b5c7ca5SMatthew G. Knepley   build:
1509b5c7ca5SMatthew G. Knepley     requires: hdf5 p4est
1519b5c7ca5SMatthew G. Knepley 
1529b5c7ca5SMatthew G. Knepley   test:
1539b5c7ca5SMatthew G. Knepley     suffix: 0
1549b5c7ca5SMatthew G. Knepley     nsize: {{1 2 5}}
15530602db0SMatthew G. Knepley     args: -adapt -dm_plex_simplex 0 -dm_plex_box_faces 2,2
1569b5c7ca5SMatthew G. Knepley 
1579b5c7ca5SMatthew G. Knepley TEST*/
158