Actual source code: ex7.c

  1: static char help[] = "Create a Plex sphere from quads and create a P1 section\n\n";

  3: #include <petscdmplex.h>

  5: static PetscErrorCode SetupSection(DM dm)
  6: {
  7:   PetscSection   s;
  8:   PetscInt       vStart, vEnd, v;

 11:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
 12:   PetscSectionCreate(PetscObjectComm((PetscObject) dm), &s);
 13:   PetscSectionSetNumFields(s, 1);
 14:   PetscSectionSetFieldComponents(s, 0, 1);
 15:   PetscSectionSetChart(s, vStart, vEnd);
 16:   for (v = vStart; v < vEnd; ++v) {
 17:     PetscSectionSetDof(s, v, 1);
 18:     PetscSectionSetFieldDof(s, v, 0, 1);
 19:   }
 20:   PetscSectionSetUp(s);
 21:   DMSetLocalSection(dm, s);
 22:   PetscSectionDestroy(&s);
 23:   return 0;
 24: }

 26: int main(int argc, char **argv)
 27: {
 28:   DM             dm;
 29:   Vec            u;

 31:   PetscInitialize(&argc, &argv, NULL,help);
 32:   DMCreate(PETSC_COMM_WORLD, &dm);
 33:   DMSetType(dm, DMPLEX);
 34:   DMSetFromOptions(dm);
 35:   PetscObjectSetName((PetscObject) dm, "Sphere");
 36:   DMViewFromOptions(dm, NULL, "-dm_view");

 38:   SetupSection(dm);
 39:   DMGetGlobalVector(dm, &u);
 40:   VecSet(u, 2);
 41:   VecViewFromOptions(u, NULL, "-vec_view");
 42:   DMRestoreGlobalVector(dm, &u);
 43:   DMDestroy(&dm);
 44:   PetscFinalize();
 45:   return 0;
 46: }

 48: /*TEST

 50:   testset:
 51:     requires: !__float128
 52:     args: -dm_plex_shape sphere -dm_view

 54:     test:
 55:       suffix: 2d_quad
 56:       args: -dm_plex_simplex 0

 58:     test:
 59:       suffix: 2d_tri
 60:       args:

 62:     test:
 63:       suffix: 3d_tri
 64:       args: -dm_plex_dim 3

 66:   testset:
 67:     requires: !__float128
 68:     args: -dm_plex_shape sphere -petscpartitioner_type simple -dm_view

 70:     test:
 71:       suffix: 2d_quad_parallel
 72:       nsize: 2
 73:       args: -dm_plex_simplex 0

 75:     test:
 76:       suffix: 2d_tri_parallel
 77:       nsize: 2

 79:     test:
 80:       suffix: 3d_tri_parallel
 81:       nsize: 2
 82:       args: -dm_plex_dim 3

 84: TEST*/