xref: /petsc/src/dm/impls/plex/tests/ex2f90.F90 (revision 9b88ac225e01f016352a5f4cd90e158abe5f5675)
1c4762a1bSJed Brown#include <petsc/finclude/petscdmplex.h>
2*c5e229c2SMartin Diehlprogram main
3ce78bad3SBarry Smith  use petscdm
4c4762a1bSJed Brown  use petscdmplex
5c4762a1bSJed Brown  implicit none
6c4762a1bSJed Brown
7c4762a1bSJed Brown  DM dm
8c4762a1bSJed Brown  PetscInt, target, dimension(3) :: EC
9c4762a1bSJed Brown  PetscInt, target, dimension(2) :: VE
10c4762a1bSJed Brown  PetscInt, pointer :: pEC(:), pVE(:)
11c4762a1bSJed Brown  PetscInt, pointer :: nClosure(:)
12c4762a1bSJed Brown  PetscInt, pointer :: nJoin(:)
13c4762a1bSJed Brown  PetscInt, pointer :: nMeet(:)
14ce78bad3SBarry Smith  PetscInt dim, cell, size, nC
15d33816bfSBarry Smith  PetscInt i0, i1, i2, i3, i6, i7
16c4762a1bSJed Brown  PetscInt i8, i9, i10, i11
17c4762a1bSJed Brown  PetscErrorCode ierr
18c4762a1bSJed Brown
19c4762a1bSJed Brown  i0 = 0
20c4762a1bSJed Brown  i1 = 1
21c4762a1bSJed Brown  i2 = 2
22c4762a1bSJed Brown  i3 = 3
23c4762a1bSJed Brown  i6 = 6
24c4762a1bSJed Brown  i7 = 7
25c4762a1bSJed Brown  i8 = 8
26c4762a1bSJed Brown  i9 = 9
27c4762a1bSJed Brown  i10 = 10
28c4762a1bSJed Brown  i11 = 11
29c4762a1bSJed Brown
30d8606c27SBarry Smith  PetscCallA(PetscInitialize(ierr))
31c4762a1bSJed Brown
32d8606c27SBarry Smith  PetscCallA(DMPlexCreate(PETSC_COMM_WORLD, dm, ierr))
33d8606c27SBarry Smith  PetscCallA(PetscObjectSetName(dm, 'Mesh', ierr))
34c4762a1bSJed Brown  dim = 2
35d8606c27SBarry Smith  PetscCallA(DMSetDimension(dm, dim, ierr))
36c4762a1bSJed Brown
37c4762a1bSJed Brown! Make Doublet Mesh from Fig 2 of Flexible Representation of Computational Meshes,
38c4762a1bSJed Brown! except indexing is from 0 instead of 1 and we obey the new restrictions on
39c4762a1bSJed Brown! numbering: cells, vertices, faces, edges
40d8606c27SBarry Smith  PetscCallA(DMPlexSetChart(dm, i0, i11, ierr))
41c4762a1bSJed Brown!     cells
42d8606c27SBarry Smith  PetscCallA(DMPlexSetConeSize(dm, i0, i3, ierr))
43d8606c27SBarry Smith  PetscCallA(DMPlexSetConeSize(dm, i1, i3, ierr))
44c4762a1bSJed Brown!     edges
45d8606c27SBarry Smith  PetscCallA(DMPlexSetConeSize(dm, i6, i2, ierr))
46d8606c27SBarry Smith  PetscCallA(DMPlexSetConeSize(dm, i7, i2, ierr))
47d8606c27SBarry Smith  PetscCallA(DMPlexSetConeSize(dm, i8, i2, ierr))
48d8606c27SBarry Smith  PetscCallA(DMPlexSetConeSize(dm, i9, i2, ierr))
49d8606c27SBarry Smith  PetscCallA(DMPlexSetConeSize(dm, i10, i2, ierr))
50c4762a1bSJed Brown
51d8606c27SBarry Smith  PetscCallA(DMSetUp(dm, ierr))
52c4762a1bSJed Brown
53c4762a1bSJed Brown  EC(1) = 6
54c4762a1bSJed Brown  EC(2) = 7
55c4762a1bSJed Brown  EC(3) = 8
56c4762a1bSJed Brown  pEC => EC
57d8606c27SBarry Smith  PetscCallA(DMPlexSetCone(dm, i0, pEC, ierr))
58c4762a1bSJed Brown  EC(1) = 7
59c4762a1bSJed Brown  EC(2) = 9
60c4762a1bSJed Brown  EC(3) = 10
61c4762a1bSJed Brown  pEC => EC
62d8606c27SBarry Smith  PetscCallA(DMPlexSetCone(dm, i1, pEC, ierr))
63c4762a1bSJed Brown
64c4762a1bSJed Brown  VE(1) = 2
65c4762a1bSJed Brown  VE(2) = 3
66c4762a1bSJed Brown  pVE => VE
67d8606c27SBarry Smith  PetscCallA(DMPlexSetCone(dm, i6, pVE, ierr))
68c4762a1bSJed Brown  VE(1) = 3
69c4762a1bSJed Brown  VE(2) = 4
70c4762a1bSJed Brown  pVE => VE
71d8606c27SBarry Smith  PetscCallA(DMPlexSetCone(dm, i7, pVE, ierr))
72c4762a1bSJed Brown  VE(1) = 4
73c4762a1bSJed Brown  VE(2) = 2
74c4762a1bSJed Brown  pVE => VE
75d8606c27SBarry Smith  PetscCallA(DMPlexSetCone(dm, i8, pVE, ierr))
76c4762a1bSJed Brown  VE(1) = 3
77c4762a1bSJed Brown  VE(2) = 5
78c4762a1bSJed Brown  pVE => VE
79d8606c27SBarry Smith  PetscCallA(DMPlexSetCone(dm, i9, pVE, ierr))
80c4762a1bSJed Brown  VE(1) = 5
81c4762a1bSJed Brown  VE(2) = 4
82c4762a1bSJed Brown  pVE => VE
83d8606c27SBarry Smith  PetscCallA(DMPlexSetCone(dm, i10, pVE, ierr))
84c4762a1bSJed Brown
85d8606c27SBarry Smith  PetscCallA(DMPlexSymmetrize(dm, ierr))
86d8606c27SBarry Smith  PetscCallA(DMPlexStratify(dm, ierr))
87d8606c27SBarry Smith  PetscCallA(DMView(dm, PETSC_VIEWER_STDOUT_WORLD, ierr))
88c4762a1bSJed Brown
89c4762a1bSJed Brown!     Test Closure
90c4762a1bSJed Brown  do cell = 0, 1
91ce78bad3SBarry Smith    PetscCallA(DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, nC, nClosure, ierr))
92c4762a1bSJed Brown!     Different Fortran compilers print a different number of columns
93c4762a1bSJed Brown!     per row producing different outputs in the test runs hence
94c4762a1bSJed Brown!     do not print the nClosure
95c4762a1bSJed Brown    write (*, 1000) 'nClosure ', nClosure
96c4762a1bSJed Brown1000 format(a, 30i4)
97ce78bad3SBarry Smith    PetscCallA(DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, nC, nClosure, ierr))
98c4762a1bSJed Brown  end do
99c4762a1bSJed Brown
100c4762a1bSJed Brown!     Test Join
101c4762a1bSJed Brown  size = 2
102c4762a1bSJed Brown  VE(1) = 6
103c4762a1bSJed Brown  VE(2) = 7
104c4762a1bSJed Brown  pVE => VE
105ce78bad3SBarry Smith  PetscCallA(DMPlexGetJoin(dm, size, pVE, PETSC_NULL_INTEGER, nJoin, ierr))
106c4762a1bSJed Brown  write (*, 1001) 'Join of', pVE
107c4762a1bSJed Brown  write (*, 1002) '  is', nJoin
108ce78bad3SBarry Smith  PetscCallA(DMPlexRestoreJoin(dm, size, pVE, PETSC_NULL_INTEGER, nJoin, ierr))
109c4762a1bSJed Brown  size = 2
110c4762a1bSJed Brown  VE(1) = 9
111c4762a1bSJed Brown  VE(2) = 7
112c4762a1bSJed Brown  pVE => VE
113ce78bad3SBarry Smith  PetscCallA(DMPlexGetJoin(dm, size, pVE, PETSC_NULL_INTEGER, nJoin, ierr))
114c4762a1bSJed Brown  write (*, 1001) 'Join of', pVE
115c4762a1bSJed Brown1001 format(a, 10i5)
116c4762a1bSJed Brown  write (*, 1002) '  is', nJoin
117c4762a1bSJed Brown1002 format(a, 10i5)
118ce78bad3SBarry Smith  PetscCallA(DMPlexRestoreJoin(dm, size, pVE, PETSC_NULL_INTEGER, nJoin, ierr))
119c4762a1bSJed Brown!     Test Full Join
120c4762a1bSJed Brown  size = 3
121c4762a1bSJed Brown  EC(1) = 3
122c4762a1bSJed Brown  EC(2) = 4
123c4762a1bSJed Brown  EC(3) = 5
124c4762a1bSJed Brown  pEC => EC
125ce78bad3SBarry Smith  PetscCallA(DMPlexGetFullJoin(dm, size, pEC, PETSC_NULL_INTEGER, nJoin, ierr))
126c4762a1bSJed Brown  write (*, 1001) 'Full Join of', pEC
127c4762a1bSJed Brown  write (*, 1002) '  is', nJoin
128ce78bad3SBarry Smith  PetscCallA(DMPlexRestoreJoin(dm, size, pEC, PETSC_NULL_INTEGER, nJoin, ierr))
129c4762a1bSJed Brown!     Test Meet
130c4762a1bSJed Brown  size = 2
131c4762a1bSJed Brown  VE(1) = 0
132c4762a1bSJed Brown  VE(2) = 1
133c4762a1bSJed Brown  pVE => VE
134ce78bad3SBarry Smith  PetscCallA(DMPlexGetMeet(dm, size, pVE, PETSC_NULL_INTEGER, nMeet, ierr))
135c4762a1bSJed Brown  write (*, 1001) 'Meet of', pVE
136c4762a1bSJed Brown  write (*, 1002) '  is', nMeet
137ce78bad3SBarry Smith  PetscCallA(DMPlexRestoreMeet(dm, size, pVE, PETSC_NULL_INTEGER, nMeet, ierr))
138c4762a1bSJed Brown  size = 2
139c4762a1bSJed Brown  VE(1) = 6
140c4762a1bSJed Brown  VE(2) = 7
141c4762a1bSJed Brown  pVE => VE
142ce78bad3SBarry Smith  PetscCallA(DMPlexGetMeet(dm, size, pVE, PETSC_NULL_INTEGER, nMeet, ierr))
143c4762a1bSJed Brown  write (*, 1001) 'Meet of', pVE
144c4762a1bSJed Brown  write (*, 1002) '  is', nMeet
145ce78bad3SBarry Smith  PetscCallA(DMPlexRestoreMeet(dm, size, pVE, PETSC_NULL_INTEGER, nMeet, ierr))
146c4762a1bSJed Brown
147d8606c27SBarry Smith  PetscCallA(DMDestroy(dm, ierr))
148d8606c27SBarry Smith  PetscCallA(PetscFinalize(ierr))
149c4762a1bSJed Brownend
150c4762a1bSJed Brown!
151c4762a1bSJed Brown!/*TEST
152c4762a1bSJed Brown!
153c4762a1bSJed Brown!   test:
154c4762a1bSJed Brown!     suffix: 0
155c4762a1bSJed Brown!
156c4762a1bSJed Brown!TEST*/
157