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