1fa58d8a1SBlaise Bourdin#include "petsc/finclude/petsc.h" 2c5e229c2SMartin Diehlprogram ex98f90 3fa58d8a1SBlaise Bourdin use petsc 4fa58d8a1SBlaise Bourdin implicit none 5fa58d8a1SBlaise Bourdin 6fa58d8a1SBlaise Bourdin ! Get the fortran kind associated with PetscInt and PetscReal so that we can use literal constants. 7fa58d8a1SBlaise Bourdin PetscInt :: dummyPetscInt 8fa58d8a1SBlaise Bourdin PetscReal :: dummyPetscreal 9fa58d8a1SBlaise Bourdin integer, parameter :: kPI = kind(dummyPetscInt) 10fa58d8a1SBlaise Bourdin integer, parameter :: kPR = kind(dummyPetscReal) 11fa58d8a1SBlaise Bourdin 12fa58d8a1SBlaise Bourdin type(tDM) :: dm, pdm 13fa58d8a1SBlaise Bourdin type(tPetscSection) :: section 14fa58d8a1SBlaise Bourdin character(len=PETSC_MAX_PATH_LEN) :: ifilename, iobuffer 15fa58d8a1SBlaise Bourdin PetscInt :: sdim, s, pStart, pEnd, p, numVS, numPoints 16fa58d8a1SBlaise Bourdin PetscInt, dimension(:), pointer :: constraints 17fa58d8a1SBlaise Bourdin type(tIS) :: setIS, pointIS 18fa58d8a1SBlaise Bourdin PetscInt, dimension(:), pointer :: setID, pointID 19fa58d8a1SBlaise Bourdin PetscErrorCode :: ierr 20fa58d8a1SBlaise Bourdin PetscBool :: flg 21fa58d8a1SBlaise Bourdin PetscMPIInt :: numProc 22*b06eb4cdSBarry Smith MPIU_Comm :: comm 23fa58d8a1SBlaise Bourdin 24fa58d8a1SBlaise Bourdin PetscCallA(PetscInitialize(ierr)) 25fa58d8a1SBlaise Bourdin PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, numProc, ierr)) 26fa58d8a1SBlaise Bourdin 27dcb3e689SBarry Smith PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-i', ifilename, flg, ierr)) 28dcb3e689SBarry Smith PetscCheckA(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 'missing input file name -i <input file name>') 29fa58d8a1SBlaise Bourdin 30fa58d8a1SBlaise Bourdin PetscCallA(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, PETSC_NULL_CHARACTER, PETSC_TRUE, dm, ierr)) 31fa58d8a1SBlaise Bourdin PetscCallA(DMPlexDistributeSetDefault(dm, PETSC_FALSE, ierr)) 32fa58d8a1SBlaise Bourdin PetscCallA(DMSetFromOptions(dm, ierr)) 33fa58d8a1SBlaise Bourdin 34fa58d8a1SBlaise Bourdin if (numproc > 1) then 35fa58d8a1SBlaise Bourdin PetscCallA(DMPlexDistribute(dm, 0_kPI, PETSC_NULL_SF, pdm, ierr)) 36fa58d8a1SBlaise Bourdin PetscCallA(DMDestroy(dm, ierr)) 37fa58d8a1SBlaise Bourdin dm = pdm 38fa58d8a1SBlaise Bourdin end if 39ce78bad3SBarry Smith PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OBJECT, '-dm_view', ierr)) 40fa58d8a1SBlaise Bourdin 41fa58d8a1SBlaise Bourdin PetscCallA(DMGetDimension(dm, sdim, ierr)) 42fa58d8a1SBlaise Bourdin PetscCallA(PetscObjectGetComm(dm, comm, ierr)) 43fa58d8a1SBlaise Bourdin PetscCallA(PetscSectionCreate(comm, section, ierr)) 44fa58d8a1SBlaise Bourdin PetscCallA(PetscSectionSetNumFields(section, 1_kPI, ierr)) 45dcb3e689SBarry Smith PetscCallA(PetscSectionSetFieldName(section, 0_kPI, 'U', ierr)) 46fa58d8a1SBlaise Bourdin PetscCallA(PetscSectionSetFieldComponents(section, 0_kPI, sdim, ierr)) 47fa58d8a1SBlaise Bourdin PetscCallA(DMPlexGetChart(dm, pStart, pEnd, ierr)) 48fa58d8a1SBlaise Bourdin PetscCallA(PetscSectionSetChart(section, pStart, pEnd, ierr)) 49fa58d8a1SBlaise Bourdin 50fa58d8a1SBlaise Bourdin ! initialize the section storage for a P1 field 51fa58d8a1SBlaise Bourdin PetscCallA(DMPlexGetDepthStratum(dm, 0_kPI, pStart, pEnd, ierr)) 52fa58d8a1SBlaise Bourdin do p = pStart, pEnd - 1 53fa58d8a1SBlaise Bourdin PetscCallA(PetscSectionSetDof(section, p, sdim, ierr)) 54fa58d8a1SBlaise Bourdin PetscCallA(PetscSectionSetFieldDof(section, p, 0_kPI, sdim, ierr)) 55fa58d8a1SBlaise Bourdin end do 56fa58d8a1SBlaise Bourdin 57fa58d8a1SBlaise Bourdin ! add constraints at all vertices belonging to a vertex set: 58fa58d8a1SBlaise Bourdin ! first pass is to reserve space 59dcb3e689SBarry Smith PetscCallA(DMGetLabelSize(dm, 'Vertex Sets', numVS, ierr)) 60fa58d8a1SBlaise Bourdin write (iobuffer, '("# Vertex set: ",i3,"\n")') numVS 61fa58d8a1SBlaise Bourdin PetscCallA(PetscPrintf(PETSC_COMM_WORLD, iobuffer, ierr)) 62dcb3e689SBarry Smith PetscCallA(DMGetLabelIdIS(dm, 'Vertex Sets', setIS, ierr)) 63ce78bad3SBarry Smith PetscCallA(ISGetIndices(setIS, setID, ierr)) 64fa58d8a1SBlaise Bourdin do s = 1, numVS 65dcb3e689SBarry Smith PetscCallA(DMGetStratumIS(dm, 'Vertex Sets', setID(s), pointIS, ierr)) 66dcb3e689SBarry Smith PetscCallA(DMGetStratumSize(dm, 'Vertex Sets', setID(s), numPoints, ierr)) 67fa58d8a1SBlaise Bourdin write (iobuffer, '("set ",i3," size ",i3,"\n")') s, numPoints 68fa58d8a1SBlaise Bourdin PetscCallA(PetscPrintf(PETSC_COMM_WORLD, iobuffer, ierr)) 69ce78bad3SBarry Smith PetscCallA(ISGetIndices(pointIS, pointID, ierr)) 70fa58d8a1SBlaise Bourdin do p = 1, numPoints 71fa58d8a1SBlaise Bourdin write (iobuffer, '(" point ",i3,"\n")') pointID(p) 72fa58d8a1SBlaise Bourdin PetscCallA(PetscPrintf(PETSC_COMM_WORLD, iobuffer, ierr)) 73fa58d8a1SBlaise Bourdin PetscCallA(PetscSectionSetConstraintDof(section, pointID(p), 1_kPI, ierr)) 74fa58d8a1SBlaise Bourdin PetscCallA(PetscSectionSetFieldConstraintDof(section, pointID(p), 0_kPI, 1_kPI, ierr)) 75fa58d8a1SBlaise Bourdin end do 76ce78bad3SBarry Smith PetscCallA(ISRestoreIndices(pointIS, pointID, ierr)) 77fa58d8a1SBlaise Bourdin PetscCallA(ISDestroy(pointIS, ierr)) 78fa58d8a1SBlaise Bourdin end do 79fa58d8a1SBlaise Bourdin 80fa58d8a1SBlaise Bourdin PetscCallA(PetscSectionSetUp(section, ierr)) 81fa58d8a1SBlaise Bourdin 82fa58d8a1SBlaise Bourdin ! add constraints at all vertices belonging to a vertex set: 83fa58d8a1SBlaise Bourdin ! second pass is to assign constraints to a specific component / dof 84fa58d8a1SBlaise Bourdin allocate (constraints(1)) 85fa58d8a1SBlaise Bourdin do s = 1, numVS 86dcb3e689SBarry Smith PetscCallA(DMGetStratumIS(dm, 'Vertex Sets', setID(s), pointIS, ierr)) 87dcb3e689SBarry Smith PetscCallA(DMGetStratumSize(dm, 'Vertex Sets', setID(s), numPoints, ierr)) 88ce78bad3SBarry Smith PetscCallA(ISGetIndices(pointIS, pointID, ierr)) 89fa58d8a1SBlaise Bourdin do p = 1, numPoints 90fa58d8a1SBlaise Bourdin constraints(1) = mod(setID(s), sdim) 91ce78bad3SBarry Smith PetscCallA(PetscSectionSetConstraintIndices(section, pointID(p), constraints, ierr)) 92ce78bad3SBarry Smith PetscCallA(PetscSectionSetFieldConstraintIndices(section, pointID(p), 0_kPI, constraints, ierr)) 93fa58d8a1SBlaise Bourdin end do 94ce78bad3SBarry Smith PetscCallA(ISRestoreIndices(pointIS, pointID, ierr)) 95fa58d8a1SBlaise Bourdin PetscCallA(ISDestroy(pointIS, ierr)) 96fa58d8a1SBlaise Bourdin end do 97fa58d8a1SBlaise Bourdin deallocate (constraints) 98ce78bad3SBarry Smith PetscCallA(ISRestoreIndices(setIS, setID, ierr)) 99fa58d8a1SBlaise Bourdin PetscCallA(ISDestroy(setIS, ierr)) 100ce78bad3SBarry Smith PetscCallA(PetscObjectViewFromOptions(section, PETSC_NULL_OBJECT, '-dm_section_view', ierr)) 101fa58d8a1SBlaise Bourdin 102fa58d8a1SBlaise Bourdin PetscCallA(PetscSectionDestroy(section, ierr)) 103fa58d8a1SBlaise Bourdin PetscCallA(DMDestroy(dm, ierr)) 104fa58d8a1SBlaise Bourdin 105fa58d8a1SBlaise Bourdin PetscCallA(PetscFinalize(ierr)) 106fa58d8a1SBlaise Bourdinend program ex98f90 107fa58d8a1SBlaise Bourdin 108fa58d8a1SBlaise Bourdin! /*TEST 109fa58d8a1SBlaise Bourdin! build: 110fa58d8a1SBlaise Bourdin! requires: exodusii pnetcdf !complex 111fa58d8a1SBlaise Bourdin! testset: 11246ac1a18SMatthew G. Knepley! args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/SquareFaceSet.exo -dm_view -dm_section_view 113fa58d8a1SBlaise Bourdin! nsize: 1 114fa58d8a1SBlaise Bourdin! test: 115fa58d8a1SBlaise Bourdin! suffix: 0 116fa58d8a1SBlaise Bourdin! args: 117fa58d8a1SBlaise Bourdin! TEST*/ 118