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