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