Actual source code: ex1f90.F90
1: #include <petsc/finclude/petscdmplex.h>
2: #include <petsc/finclude/petscdmlabel.h>
3: program DMPlexTestField
4: use petscdm
5: use petscdmplex
6: implicit none
8: DM :: dm
9: DMLabel :: label
10: Vec :: u
11: PetscViewer :: viewer
12: PetscSection :: section
13: PetscInt, parameter :: numFields = 3, numBC = 1
14: PetscInt :: dim, 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: PetscMPIInt :: size
22: IS, target, dimension(1) :: bcCompIS
23: IS, target, dimension(1) :: bcPointIS
24: IS, pointer :: pBcCompIS(:)
25: IS, pointer :: pBcPointIS(:)
26: PetscErrorCode :: ierr
28: PetscCallA(PetscInitialize(ierr))
29: PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, size, ierr))
30: ! Create a mesh
31: PetscCallA(DMCreate(PETSC_COMM_WORLD, dm, ierr))
32: PetscCallA(DMSetType(dm, DMPLEX, ierr))
33: PetscCallA(DMSetFromOptions(dm, ierr))
34: PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OBJECT, '-dm_view', ierr))
35: PetscCallA(DMGetDimension(dm, dim, ierr))
36: ! Create a scalar field u, a vector field v, and a surface vector field w
37: numComp = [1_PETSC_INT_KIND, dim, dim - 1_PETSC_INT_KIND]
38: pNumComp => numComp
39: numDof(1:numFields*(dim + 1)) = 0
40: ! Let u be defined on vertices
41: numDof(0*(dim + 1) + 1) = 1
42: ! Let v be defined on cells
43: numDof(1*(dim + 1) + dim + 1) = dim
44: ! Let v be defined on faces
45: numDof(2*(dim + 1) + dim) = dim - 1
46: pNumDof => numDof
47: ! Test label retrieval
48: PetscCallA(DMGetLabel(dm, 'marker', label, ierr))
49: PetscCallA(DMLabelGetValue(label, 0_PETSC_INT_KIND, val, ierr))
50: PetscCheckA(size /= 1 .or. val == -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Error in library')
51: PetscCallA(DMLabelGetValue(label, 8_PETSC_INT_KIND, val, ierr))
52: PetscCheckA(size /= 1 .or. val == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Error in library')
53: ! Prescribe a Dirichlet condition on u on the boundary
54: ! Label "marker" is made by the mesh creation routine
55: bcField(1) = 0
56: pBcField => bcField
57: PetscCallA(ISCreateStride(PETSC_COMM_WORLD, 1_PETSC_INT_KIND, 0_PETSC_INT_KIND, 1_PETSC_INT_KIND, bcCompIS(1), ierr))
58: pBcCompIS => bcCompIS
59: PetscCallA(DMGetStratumIS(dm, 'marker', 1_PETSC_INT_KIND, bcPointIS(1), ierr))
60: pBcPointIS => bcPointIS
61: ! Create a PetscSection with this data layout
62: PetscCallA(DMSetNumFields(dm, numFields, ierr))
63: PetscCallA(DMPlexCreateSection(dm, PETSC_NULL_DMLABEL_ARRAY, pNumComp, pNumDof, numBC, pBcField, pBcCompIS, pBcPointIS, PETSC_NULL_IS, section, ierr))
64: PetscCallA(ISDestroy(bcCompIS(1), ierr))
65: PetscCallA(ISDestroy(bcPointIS(1), ierr))
66: ! Name the Field variables
67: PetscCallA(PetscSectionSetFieldName(section, 0_PETSC_INT_KIND, 'u', ierr))
68: PetscCallA(PetscSectionSetFieldName(section, 1_PETSC_INT_KIND, 'v', ierr))
69: PetscCallA(PetscSectionSetFieldName(section, 2_PETSC_INT_KIND, 'w', ierr))
70: if (size == 1) then
71: PetscCallA(PetscSectionView(section, PETSC_VIEWER_STDOUT_WORLD, ierr))
72: end if
73: ! Tell the DM to use this data layout
74: PetscCallA(DMSetLocalSection(dm, section, ierr))
75: ! Create a Vec with this layout and view it
76: PetscCallA(DMGetGlobalVector(dm, u, ierr))
77: PetscCallA(PetscViewerCreate(PETSC_COMM_WORLD, viewer, ierr))
78: PetscCallA(PetscViewerSetType(viewer, PETSCVIEWERVTK, ierr))
79: PetscCallA(PetscViewerFileSetName(viewer, 'sol.vtu', ierr))
80: PetscCallA(VecView(u, viewer, ierr))
81: PetscCallA(PetscViewerDestroy(viewer, ierr))
82: PetscCallA(DMRestoreGlobalVector(dm, u, ierr))
83: ! Cleanup
84: PetscCallA(PetscSectionDestroy(section, ierr))
85: PetscCallA(DMDestroy(dm, ierr))
87: PetscCallA(PetscFinalize(ierr))
88: end program DMPlexTestField
90: !/*TEST
91: !
92: ! test:
93: ! suffix: 0
94: ! requires: triangle
95: ! args: -info :~sys,mat:
96: !
97: ! test:
98: ! suffix: 0_2
99: ! requires: triangle
100: ! nsize: 2
101: ! args: -info :~sys,mat,partitioner:
102: !
103: ! test:
104: ! suffix: 1
105: ! requires: ctetgen
106: ! args: -dm_plex_dim 3 -info :~sys,mat:
107: !
108: !TEST*/