Actual source code: ex98f90.F90

  1: #include "petsc/finclude/petsc.h"
  2: program ex98f90
  3:   use petsc
  4:   implicit none

  6:   type(tDM)                          :: dm, pdm
  7:   type(tPetscSection)                :: section
  8:   character(len=PETSC_MAX_PATH_LEN)  :: ifilename, iobuffer
  9:   PetscInt                           :: sdim, s, pStart, pEnd, p, numVS, numPoints
 10:   PetscInt, dimension(:), pointer    :: constraints
 11:   type(tIS)                          :: setIS, pointIS
 12:   PetscInt, dimension(:), pointer    :: setID, pointID
 13:   PetscErrorCode                     :: ierr
 14:   PetscBool                          :: flg
 15:   PetscMPIInt                        :: numProc
 16:   MPIU_Comm                          :: comm

 18:   PetscCallA(PetscInitialize(ierr))
 19:   PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, numProc, ierr))

 21:   PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-i', ifilename, flg, ierr))
 22:   PetscCheckA(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 'missing input file name -i ')

 24:   PetscCallA(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, PETSC_NULL_CHARACTER, PETSC_TRUE, dm, ierr))
 25:   PetscCallA(DMPlexDistributeSetDefault(dm, PETSC_FALSE, ierr))
 26:   PetscCallA(DMSetFromOptions(dm, ierr))

 28:   if (numproc > 1) then
 29:     PetscCallA(DMPlexDistribute(dm, 0_PETSC_INT_KIND, PETSC_NULL_SF, pdm, ierr))
 30:     PetscCallA(DMDestroy(dm, ierr))
 31:     dm = pdm
 32:   end if
 33:   PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OBJECT, '-dm_view', ierr))

 35:   PetscCallA(DMGetDimension(dm, sdim, ierr))
 36:   PetscCallA(PetscObjectGetComm(dm, comm, ierr))
 37:   PetscCallA(PetscSectionCreate(comm, section, ierr))
 38:   PetscCallA(PetscSectionSetNumFields(section, 1_PETSC_INT_KIND, ierr))
 39:   PetscCallA(PetscSectionSetFieldName(section, 0_PETSC_INT_KIND, 'U', ierr))
 40:   PetscCallA(PetscSectionSetFieldComponents(section, 0_PETSC_INT_KIND, sdim, ierr))
 41:   PetscCallA(DMPlexGetChart(dm, pStart, pEnd, ierr))
 42:   PetscCallA(PetscSectionSetChart(section, pStart, pEnd, ierr))

 44:   ! initialize the section storage for a P1 field
 45:   PetscCallA(DMPlexGetDepthStratum(dm, 0_PETSC_INT_KIND, pStart, pEnd, ierr))
 46:   do p = pStart, pEnd - 1
 47:     PetscCallA(PetscSectionSetDof(section, p, sdim, ierr))
 48:     PetscCallA(PetscSectionSetFieldDof(section, p, 0_PETSC_INT_KIND, sdim, ierr))
 49:   end do

 51:   ! add constraints at all vertices belonging to a vertex set:
 52:   ! first pass is to reserve space
 53:   PetscCallA(DMGetLabelSize(dm, 'Vertex Sets', numVS, ierr))
 54:   write (iobuffer, '("# Vertex set: ",i3,"\n")') numVS
 55:   PetscCallA(PetscPrintf(PETSC_COMM_WORLD, iobuffer, ierr))
 56:   PetscCallA(DMGetLabelIdIS(dm, 'Vertex Sets', setIS, ierr))
 57:   PetscCallA(ISGetIndices(setIS, setID, ierr))
 58:   do s = 1, numVS
 59:     PetscCallA(DMGetStratumIS(dm, 'Vertex Sets', setID(s), pointIS, ierr))
 60:     PetscCallA(DMGetStratumSize(dm, 'Vertex Sets', setID(s), numPoints, ierr))
 61:     write (iobuffer, '("set ",i3," size ",i3,"\n")') s, numPoints
 62:     PetscCallA(PetscPrintf(PETSC_COMM_WORLD, iobuffer, ierr))
 63:     PetscCallA(ISGetIndices(pointIS, pointID, ierr))
 64:     do p = 1, numPoints
 65:       write (iobuffer, '("   point ",i3,"\n")') pointID(p)
 66:       PetscCallA(PetscPrintf(PETSC_COMM_WORLD, iobuffer, ierr))
 67:       PetscCallA(PetscSectionSetConstraintDof(section, pointID(p), 1_PETSC_INT_KIND, ierr))
 68:       PetscCallA(PetscSectionSetFieldConstraintDof(section, pointID(p), 0_PETSC_INT_KIND, 1_PETSC_INT_KIND, ierr))
 69:     end do
 70:     PetscCallA(ISRestoreIndices(pointIS, pointID, ierr))
 71:     PetscCallA(ISDestroy(pointIS, ierr))
 72:   end do

 74:   PetscCallA(PetscSectionSetUp(section, ierr))

 76:   ! add constraints at all vertices belonging to a vertex set:
 77:   ! second pass is to assign constraints to a specific component / dof
 78:   allocate (constraints(1))
 79:   do s = 1, numVS
 80:     PetscCallA(DMGetStratumIS(dm, 'Vertex Sets', setID(s), pointIS, ierr))
 81:     PetscCallA(DMGetStratumSize(dm, 'Vertex Sets', setID(s), numPoints, ierr))
 82:     PetscCallA(ISGetIndices(pointIS, pointID, ierr))
 83:     do p = 1, numPoints
 84:       constraints(1) = mod(setID(s), sdim)
 85:       PetscCallA(PetscSectionSetConstraintIndices(section, pointID(p), constraints, ierr))
 86:       PetscCallA(PetscSectionSetFieldConstraintIndices(section, pointID(p), 0_PETSC_INT_KIND, constraints, ierr))
 87:     end do
 88:     PetscCallA(ISRestoreIndices(pointIS, pointID, ierr))
 89:     PetscCallA(ISDestroy(pointIS, ierr))
 90:   end do
 91:   deallocate (constraints)
 92:   PetscCallA(ISRestoreIndices(setIS, setID, ierr))
 93:   PetscCallA(ISDestroy(setIS, ierr))
 94:   PetscCallA(PetscObjectViewFromOptions(section, PETSC_NULL_OBJECT, '-dm_section_view', ierr))

 96:   PetscCallA(PetscSectionDestroy(section, ierr))
 97:   PetscCallA(DMDestroy(dm, ierr))

 99:   PetscCallA(PetscFinalize(ierr))
100: end program ex98f90

102: ! /*TEST
103: !   build:
104: !     requires: exodusii pnetcdf !complex
105: !   testset:
106: !     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/SquareFaceSet.exo -dm_view -dm_section_view
107: !     nsize: 1
108: !     test:
109: !       suffix: 0
110: !       args:
111: ! TEST*/