Actual source code: ex1f90.F90
1: #include <petsc/finclude/petscdmplex.h>
2: program main
3: use petscdmplex
4: implicit none
5: !
6: !
7: DM dm
8: PetscInt, dimension(4) :: EC
9: PetscInt, pointer :: pEC(:), pES(:)
10: PetscInt, parameter :: firstCell = 0, numCells = 2, numVertices = 6, numPoints = numCells + numVertices
11: PetscInt c, v
12: PetscErrorCode ierr
14: PetscCallA(PetscInitialize(ierr))
16: PetscCallA(DMPlexCreate(PETSC_COMM_WORLD, dm, ierr))
17: PetscCallA(DMPlexSetChart(dm, 0_PETSC_INT_KIND, numPoints, ierr))
18: do c = firstCell, numCells - 1
19: PetscCallA(DMPlexSetConeSize(dm, c, 4_PETSC_INT_KIND, ierr))
20: end do
21: PetscCallA(DMSetUp(dm, ierr))
23: EC = [2, 3, 4, 5]
24: c = 0
25: write (*, 1000) 'cell EC 0', c, EC
26: 1000 format(a, i4, 50i4)
27: PetscCallA(DMPlexSetCone(dm, c, EC, ierr))
28: PetscCallA(DMPlexGetCone(dm, c, pEC, ierr))
29: write (*, 1000) 'cell pEC 0', c, pEC
30: PetscCallA(DMPlexRestoreCone(dm, c, pEC, ierr))
31: EC = [4, 5, 6, 7]
32: c = 1
33: write (*, 1000) 'cell EC 1', c, EC
34: PetscCallA(DMPlexSetCone(dm, c, EC, ierr))
35: PetscCallA(DMPlexGetCone(dm, c, pEC, ierr))
36: write (*, 1000) 'cell pEC 1', c, pEC
37: PetscCallA(DMPlexRestoreCone(dm, c, pEC, ierr))
38: CHKMEMQ
40: PetscCallA(DMPlexSymmetrize(dm, ierr))
41: PetscCallA(DMPlexStratify(dm, ierr))
42: PetscCallA(DMPlexGetCone(dm, c, pEC, ierr))
43: write (*, 1000) 'cell pEC 3', c, pEC
44: PetscCallA(DMPlexRestoreCone(dm, c, pEC, ierr))
46: v = 4
47: PetscCallA(DMPlexGetSupport(dm, v, pES, ierr))
48: write (*, 1000) 'vertex', v, pES
49: PetscCallA(DMPlexRestoreSupport(dm, v, pES, ierr))
51: PetscCallA(DMDestroy(dm, ierr))
52: PetscCallA(PetscFinalize(ierr))
53: end
55: ! /*TEST
56: !
57: ! test:
58: ! suffix: 0
59: !
60: ! TEST*/