xref: /petsc/src/dm/impls/plex/tutorials/ex1f90.F90 (revision ce78bad369055609e946c9d2c25ea67a45873e27)
1c4762a1bSJed Brown      program DMPlexTestField
2*ce78bad3SBarry Smith#include <petsc/finclude/petscdmplex.h>
3*ce78bad3SBarry Smith#include <petsc/finclude/petscdmlabel.h>
4*ce78bad3SBarry Smith      use petscdm
5c4762a1bSJed Brown      implicit none
6c4762a1bSJed Brown
7c4762a1bSJed Brown      DM :: dm
8c4762a1bSJed Brown      DMLabel :: label
9c4762a1bSJed Brown      Vec :: u
10c4762a1bSJed Brown      PetscViewer :: viewer
11c4762a1bSJed Brown      PetscSection :: section
12c4762a1bSJed Brown      PetscInt :: dim,numFields,numBC
13c4762a1bSJed Brown      PetscInt :: i,val
14c4762a1bSJed Brown      PetscInt, target, dimension(3) ::  numComp
15c4762a1bSJed Brown      PetscInt, pointer :: pNumComp(:)
16c4762a1bSJed Brown      PetscInt, target, dimension(12) ::  numDof
17c4762a1bSJed Brown      PetscInt, pointer :: pNumDof(:)
18c4762a1bSJed Brown      PetscInt, target, dimension(1) ::  bcField
19c4762a1bSJed Brown      PetscInt, pointer :: pBcField(:)
20c4762a1bSJed Brown      PetscInt, parameter :: zero = 0, one = 1, two = 2, eight = 8
212e3d3ef9SMatthew G. Knepley      PetscMPIInt :: size
22c4762a1bSJed Brown      IS, target, dimension(1) ::   bcCompIS
23c4762a1bSJed Brown      IS, target, dimension(1) ::   bcPointIS
24c4762a1bSJed Brown      IS, pointer :: pBcCompIS(:)
25c4762a1bSJed Brown      IS, pointer :: pBcPointIS(:)
26c4762a1bSJed Brown      PetscErrorCode :: ierr
27c4762a1bSJed Brown
28d8606c27SBarry Smith      PetscCallA(PetscInitialize(ierr))
29d8606c27SBarry Smith      PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, size, ierr))
30c4762a1bSJed Brown!     Create a mesh
31d8606c27SBarry Smith      PetscCallA(DMCreate(PETSC_COMM_WORLD, dm, ierr))
32d8606c27SBarry Smith      PetscCallA(DMSetType(dm, DMPLEX, ierr))
33d8606c27SBarry Smith      PetscCallA(DMSetFromOptions(dm, ierr))
34*ce78bad3SBarry Smith      PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OBJECT, '-dm_view', ierr))
35d8606c27SBarry Smith      PetscCallA(DMGetDimension(dm, dim, ierr))
36c4762a1bSJed Brown!     Create a scalar field u, a vector field v, and a surface vector field w
37c4762a1bSJed Brown      numFields  = 3
38c4762a1bSJed Brown      numComp(1) = 1
39c4762a1bSJed Brown      numComp(2) = dim
40c4762a1bSJed Brown      numComp(3) = dim-1
41c4762a1bSJed Brown      pNumComp => numComp
42c4762a1bSJed Brown      do i = 1, numFields*(dim+1)
43c4762a1bSJed Brown         numDof(i) = 0
44c4762a1bSJed Brown      end do
45c4762a1bSJed Brown!     Let u be defined on vertices
46c4762a1bSJed Brown      numDof(0*(dim+1)+1)     = 1
47c4762a1bSJed Brown!     Let v be defined on cells
48c4762a1bSJed Brown      numDof(1*(dim+1)+dim+1) = dim
49c4762a1bSJed Brown!     Let v be defined on faces
50c4762a1bSJed Brown      numDof(2*(dim+1)+dim)   = dim-1
51c4762a1bSJed Brown      pNumDof => numDof
52c4762a1bSJed Brown!     Setup boundary conditions
53c4762a1bSJed Brown      numBC = 1
54c4762a1bSJed Brown!     Test label retrieval
55d8606c27SBarry Smith      PetscCallA(DMGetLabel(dm, 'marker', label, ierr))
56d8606c27SBarry Smith      PetscCallA(DMLabelGetValue(label, zero, val, ierr))
57dcb3e689SBarry Smith      PetscCheckA(size .ne. 1 .or. val .eq. -1,PETSC_COMM_SELF,PETSC_ERR_PLIB,'Error in library')
58d8606c27SBarry Smith      PetscCallA(DMLabelGetValue(label, eight, val, ierr))
59dcb3e689SBarry Smith      PetscCheckA(size .ne. 1 .or. val .eq. 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,'Error in library')
60c4762a1bSJed Brown!     Prescribe a Dirichlet condition on u on the boundary
61c4762a1bSJed Brown!       Label "marker" is made by the mesh creation routine
62c4762a1bSJed Brown      bcField(1) = 0
63c4762a1bSJed Brown      pBcField => bcField
64d8606c27SBarry Smith      PetscCallA(ISCreateStride(PETSC_COMM_WORLD, one, zero, one, bcCompIS(1), ierr))
65c4762a1bSJed Brown      pBcCompIS => bcCompIS
66d8606c27SBarry Smith      PetscCallA(DMGetStratumIS(dm, 'marker', one, bcPointIS(1),ierr))
67c4762a1bSJed Brown      pBcPointIS => bcPointIS
68c4762a1bSJed Brown!     Create a PetscSection with this data layout
69d8606c27SBarry Smith      PetscCallA(DMSetNumFields(dm, numFields,ierr))
70*ce78bad3SBarry Smith      PetscCallA(DMPlexCreateSection(dm,PETSC_NULL_DMLABEL_ARRAY,pNumComp,pNumDof,numBC,pBcField,pBcCompIS,pBcPointIS,PETSC_NULL_IS,section,ierr))
71d8606c27SBarry Smith      PetscCallA(ISDestroy(bcCompIS(1), ierr))
72d8606c27SBarry Smith      PetscCallA(ISDestroy(bcPointIS(1), ierr))
73c4762a1bSJed Brown!     Name the Field variables
74d8606c27SBarry Smith      PetscCallA(PetscSectionSetFieldName(section, zero, 'u', ierr))
75d8606c27SBarry Smith      PetscCallA(PetscSectionSetFieldName(section, one,  'v', ierr))
76d8606c27SBarry Smith      PetscCallA(PetscSectionSetFieldName(section, two,  'w', ierr))
772e3d3ef9SMatthew G. Knepley      if (size .eq. 1) then
78d8606c27SBarry Smith        PetscCallA(PetscSectionView(section, PETSC_VIEWER_STDOUT_WORLD, ierr))
792e3d3ef9SMatthew G. Knepley      endif
80c4762a1bSJed Brown!     Tell the DM to use this data layout
81d8606c27SBarry Smith      PetscCallA(DMSetLocalSection(dm, section, ierr))
82c4762a1bSJed Brown!     Create a Vec with this layout and view it
83d8606c27SBarry Smith      PetscCallA(DMGetGlobalVector(dm, u, ierr))
84d8606c27SBarry Smith      PetscCallA(PetscViewerCreate(PETSC_COMM_WORLD, viewer, ierr))
85d8606c27SBarry Smith      PetscCallA(PetscViewerSetType(viewer, PETSCVIEWERVTK, ierr))
86d8606c27SBarry Smith      PetscCallA(PetscViewerFileSetName(viewer, 'sol.vtu', ierr))
87d8606c27SBarry Smith      PetscCallA(VecView(u, viewer, ierr))
88d8606c27SBarry Smith      PetscCallA(PetscViewerDestroy(viewer, ierr))
89d8606c27SBarry Smith      PetscCallA(DMRestoreGlobalVector(dm, u, ierr))
90c4762a1bSJed Brown!     Cleanup
91d8606c27SBarry Smith      PetscCallA(PetscSectionDestroy(section, ierr))
92d8606c27SBarry Smith      PetscCallA(DMDestroy(dm, ierr))
93c4762a1bSJed Brown
94d8606c27SBarry Smith      PetscCallA(PetscFinalize(ierr))
95c4762a1bSJed Brown      end program DMPlexTestField
96c4762a1bSJed Brown
97c4762a1bSJed Brown!/*TEST
98c4762a1bSJed Brown!  build:
99dfd57a17SPierre Jolivet!    requires: defined(PETSC_USING_F90FREEFORM)
100c4762a1bSJed Brown!
101c4762a1bSJed Brown!  test:
102c4762a1bSJed Brown!    suffix: 0
103c4762a1bSJed Brown!    requires: triangle
10430602db0SMatthew G. Knepley!    args: -info :~sys,mat:
105c4762a1bSJed Brown!
106c4762a1bSJed Brown!  test:
1072e3d3ef9SMatthew G. Knepley!    suffix: 0_2
1082e3d3ef9SMatthew G. Knepley!    requires: triangle
1092e3d3ef9SMatthew G. Knepley!    nsize: 2
1102e3d3ef9SMatthew G. Knepley!    args: -info :~sys,mat,partitioner:
1112e3d3ef9SMatthew G. Knepley!
1122e3d3ef9SMatthew G. Knepley!  test:
113c4762a1bSJed Brown!    suffix: 1
114c4762a1bSJed Brown!    requires: ctetgen
11530602db0SMatthew G. Knepley!    args: -dm_plex_dim 3 -info :~sys,mat:
116c4762a1bSJed Brown!
117c4762a1bSJed Brown!TEST*/
118