Actual source code: ex1f90.F90
petsc-3.10.2 2018-10-09
1: program DMPlexTestField
2: #include "petsc/finclude/petscdmplex.h"
3: #include "petsc/finclude/petscdmlabel.h"
4: use petscdmplex
5: use petscsys
6: implicit none
8: DM :: dm
9: DMLabel :: label
10: Vec :: u
11: PetscViewer :: viewer
12: PetscSection :: section
13: PetscInt :: dim,numFields,numBC
14: PetscInt :: i,val
15: PetscInt, target, dimension(3) :: numComp
16: PetscInt, pointer :: pNumComp(:)
17: PetscInt, target, dimension(12) :: numDof
18: PetscInt, pointer :: pNumDof(:)
19: PetscInt, target, dimension(1) :: bcField
20: PetscInt, pointer :: pBcField(:)
21: PetscInt :: zero,eight
22: IS, target, dimension(1) :: bcCompIS
23: IS, target, dimension(1) :: bcPointIS
24: IS, pointer :: pBcCompIS(:)
25: IS, pointer :: pBcPointIS(:)
26: PetscBool :: interpolate,flg
27: PetscErrorCode :: ierr
29: call PetscInitialize(PETSC_NULL_CHARACTER, ierr)
30: if (ierr .ne. 0) then
31: print*,'Unable to initialize PETSc'
32: stop
33: endif
34: dim = 2
35: call PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-dim', dim,flg, ierr);CHKERRA(ierr)
36: interpolate = PETSC_TRUE
37: ! Create a mesh
38: call DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_TRUE, PETSC_NULL_INTEGER, PETSC_NULL_REAL, PETSC_NULL_REAL, PETSC_NULL_INTEGER, interpolate, dm, ierr);CHKERRA(ierr)
39: ! Create a scalar field u, a vector field v, and a surface vector field w
40: numFields = 3
41: numComp(1) = 1
42: numComp(2) = dim
43: numComp(3) = dim-1
44: pNumComp => numComp
45: do i = 1, numFields*(dim+1)
46: numDof(i) = 0
47: end do
48: ! Let u be defined on vertices
49: numDof(0*(dim+1)+1) = 1
50: ! Let v be defined on cells
51: numDof(1*(dim+1)+dim+1) = dim
52: ! Let v be defined on faces
53: numDof(2*(dim+1)+dim) = dim-1
54: pNumDof => numDof
55: ! Setup boundary conditions
56: numBC = 1
57: ! Test label retrieval
58: call DMGetLabel(dm, 'marker', label, ierr);CHKERRA(ierr)
59: zero = 0
60: call DMLabelGetValue(label, zero, val, ierr);CHKERRA(ierr)
61: if (val .ne. -1) then
62: CHKERRA(1)
63: endif
64: eight = 8
65: call DMLabelGetValue(label, eight, val, ierr);CHKERRA(ierr)
66: if (val .ne. 1) then
67: CHKERRA(1)
68: endif
69: ! Prescribe a Dirichlet condition on u on the boundary
70: ! Label "marker" is made by the mesh creation routine
71: bcField(1) = 0
72: pBcField => bcField
73: call ISCreateStride(PETSC_COMM_WORLD, 1, 0, 1, bcCompIS(1), ierr);CHKERRA(ierr)
74: pBcCompIS => bcCompIS
75: call DMGetStratumIS(dm, 'marker', 1, bcPointIS(1),ierr);CHKERRA(ierr)
76: pBcPointIS => bcPointIS
77: ! Create a PetscSection with this data layout
78: call DMPlexCreateSection(dm,dim,numFields,pNumComp,pNumDof,numBC,pBcField,pBcCompIS,pBcPointIS,PETSC_NULL_IS,section,ierr)
79: CHKERRA(ierr)
80: call ISDestroy(bcCompIS(1), ierr);CHKERRA(ierr)
81: call ISDestroy(bcPointIS(1), ierr);CHKERRA(ierr)
82: ! Name the Field variables
83: call PetscSectionSetFieldName(section, 0, 'u', ierr);CHKERRA(ierr)
84: call PetscSectionSetFieldName(section, 1, 'v', ierr);CHKERRA(ierr)
85: call PetscSectionSetFieldName(section, 2, 'w', ierr);CHKERRA(ierr)
86: call PetscSectionView(section, PETSC_VIEWER_STDOUT_WORLD, ierr);CHKERRA(ierr)
87: ! Tell the DM to use this data layout
88: call DMSetSection(dm, section, ierr);CHKERRA(ierr)
89: ! Create a Vec with this layout and view it
90: call DMGetGlobalVector(dm, u, ierr);CHKERRA(ierr)
91: call PetscViewerCreate(PETSC_COMM_WORLD, viewer, ierr);CHKERRA(ierr)
92: call PetscViewerSetType(viewer, PETSCVIEWERVTK, ierr);CHKERRA(ierr)
93: call PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_VTK, ierr);CHKERRA(ierr)
94: call PetscViewerFileSetName(viewer, 'sol.vtk', ierr);CHKERRA(ierr)
95: call VecView(u, viewer, ierr);CHKERRA(ierr)
96: call PetscViewerDestroy(viewer, ierr);CHKERRA(ierr)
97: call DMRestoreGlobalVector(dm, u, ierr);CHKERRA(ierr)
98: ! Cleanup
99: call PetscSectionDestroy(section, ierr);CHKERRA(ierr)
100: call DMDestroy(dm, ierr);CHKERRA(ierr)
102: call PetscFinalize(ierr)
103: end program DMPlexTestField
105: !/*TEST
106: ! build:
107: ! requires: define(PETSC_USING_F90FREEFORM)
108: !
109: ! test:
110: ! suffix: 0
111: ! requires: triangle
112: !
113: ! test:
114: ! suffix: 1
115: ! requires: ctetgen
116: ! args: -dim 3
117: !
118: !TEST*/