xref: /petsc/src/dm/impls/plex/tests/ex1f90.F90 (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1      program main
2#include <petsc/finclude/petscdmplex.h>
3      use petscdmplex
4      use petscsys
5      implicit none
6!
7!
8      DM dm
9      PetscInt, target, dimension(4) :: EC
10      PetscInt, pointer :: pEC(:)
11      PetscInt, pointer :: pES(:)
12      PetscInt c, firstCell, numCells
13      PetscInt v, numVertices, numPoints
14      PetscInt i0,i4
15      PetscErrorCode ierr
16
17      i0 = 0
18      i4 = 4
19
20      call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
21     if (ierr .ne. 0) then
22        print*,'Unable to initialize PETSc'
23        stop
24      endif
25      call DMPlexCreate(PETSC_COMM_WORLD, dm, ierr);CHKERRA(ierr)
26      firstCell = 0
27      numCells = 2
28      numVertices = 6
29      numPoints = numCells+numVertices
30      call DMPlexSetChart(dm, i0, numPoints, ierr);CHKERRA(ierr)
31      do c=firstCell,numCells-1
32         call DMPlexSetConeSize(dm, c, i4, ierr);CHKERRA(ierr)
33      end do
34      call DMSetUp(dm, ierr);CHKERRA(ierr)
35
36      EC(1) = 2
37      EC(2) = 3
38      EC(3) = 4
39      EC(4) = 5
40      pEC => EC
41      c = 0
42      write(*,1000) 'cell',c,pEC
43 1000 format (a,i4,50i4)
44      call DMPlexSetCone(dm, c , pEC, ierr);CHKERRA(ierr)
45      call DMPlexGetCone(dm, c , pEC, ierr);CHKERRA(ierr)
46      write(*,1000) 'cell',c,pEC
47      EC(1) = 4
48      EC(2) = 5
49      EC(3) = 6
50      EC(4) = 7
51      pEC => EC
52      c = 1
53      write(*,1000) 'cell',c,pEC
54      call DMPlexSetCone(dm, c , pEC, ierr);CHKERRA(ierr)
55      call DMPlexGetCone(dm, c , pEC, ierr);CHKERRA(ierr)
56      write(*,1000) 'cell',c,pEC
57      call DMPlexRestoreCone(dm, c , pEC, ierr);CHKERRA(ierr)
58
59      call DMPlexSymmetrize(dm, ierr);CHKERRA(ierr)
60      call DMPlexStratify(dm, ierr);CHKERRA(ierr)
61
62      v = 4
63      call DMPlexGetSupport(dm, v , pES, ierr);CHKERRA(ierr)
64      write(*,1000) 'vertex',v,pES
65      call DMPlexRestoreSupport(dm, v , pES, ierr);CHKERRA(ierr)
66
67      call DMDestroy(dm,ierr);CHKERRA(ierr)
68      call PetscFinalize(ierr)
69      end
70
71! /*TEST
72!
73! test:
74!   suffix: 0
75!
76! TEST*/
77