Actual source code: ex1f90.F90
1: #include <petsc/finclude/petscdmlabel.h>
2: module ex1f90module
3: use petscdm
4: implicit none
6: contains
7: subroutine ViewLabels(dm, viewer, ierr)
8: type(tDM) :: dm
9: type(tPetscViewer) :: viewer
10: PetscErrorCode, intent(out) :: ierr
12: DMLabel :: label
13: type(tIS) :: labelIS
14: character(len=PETSC_MAX_PATH_LEN):: labelName, IObuffer
15: PetscInt :: numLabels, l
17: PetscCall(DMGetNumLabels(dm, numLabels, ierr))
18: write (IObuffer, *) 'Number of labels: ', numLabels, '\n'
19: PetscCall(PetscViewerASCIIPrintf(viewer, IObuffer, ierr))
20: do l = 0, numLabels - 1
21: PetscCall(DMGetLabelName(dm, l, labelName, ierr))
22: write (IObuffer, *) 'label ', l, ' name: ', trim(labelName), '\n'
23: PetscCall(PetscViewerASCIIPrintf(viewer, IObuffer, ierr))
25: PetscCall(PetscViewerASCIIPrintf(viewer, 'IS of values\n', ierr))
26: PetscCall(DMGetLabel(dm, labelName, label, ierr))
27: PetscCall(DMLabelGetValueIS(label, labelIS, ierr))
28: ! PetscCall(PetscViewerASCIIPushTab(viewer,ierr))
29: PetscCall(ISView(labelIS, viewer, ierr))
30: ! PetscCall(PetscViewerASCIIPopTab(viewer,ierr))
31: PetscCall(ISDestroy(labelIS, ierr))
32: PetscCall(PetscViewerASCIIPrintf(viewer, '\n', ierr))
33: end do
35: PetscCall(PetscViewerASCIIPrintf(viewer, '\n\nCell Set label IS\n', ierr))
36: PetscCall(DMGetLabel(dm, 'Cell Sets', label, ierr))
37: PetscCall(DMLabelGetValueIS(label, labelIS, ierr))
38: PetscCall(ISView(labelIS, viewer, ierr))
39: PetscCall(ISDestroy(labelIS, ierr))
40: end subroutine viewLabels
41: end module ex1f90module
43: program ex1f90
44: use petscdm
45: use ex1f90module
47: implicit none
48: type(tDM) :: dm, dmDist
49: character(len=PETSC_MAX_PATH_LEN) :: filename
50: PetscBool :: interpolate = PETSC_FALSE
51: PetscBool :: flg
52: PetscErrorCode :: ierr
54: PetscCallA(PetscInitialize(ierr))
55: PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-i', filename, flg, ierr))
56: PetscCallA(PetscOptionsGetBool(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-interpolate', interpolate, flg, ierr))
58: PetscCallA(DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, 'ex1f90_plex', interpolate, dm, ierr))
59: PetscCallA(DMPlexDistribute(dm, 0_PETSC_INT_KIND, PETSC_NULL_SF, dmDist, ierr))
60: if (.not. PetscObjectIsNull(dmDist)) then
61: PetscCallA(DMDestroy(dm, ierr))
62: dm = dmDist
63: end if
65: PetscCallA(ViewLabels(dm, PETSC_VIEWER_STDOUT_WORLD, ierr))
66: PetscCallA(DMDestroy(dm, ierr))
67: PetscCallA(PetscFinalize(ierr))
68: end program ex1F90
70: !/*TEST
71: !
72: ! test:
73: ! suffix: 0
74: ! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -interpolate
75: ! requires: exodusii
76: !
77: !TEST*/