1f52d5e2bSBlaise Bourdinprogram ex47f90 2f52d5e2bSBlaise Bourdin#include "petsc/finclude/petsc.h" 3f52d5e2bSBlaise Bourdin#include "petsc/finclude/petscvec.h" 4f52d5e2bSBlaise Bourdin use petsc 5f52d5e2bSBlaise Bourdin use petscvec 6f52d5e2bSBlaise Bourdin implicit none 7f52d5e2bSBlaise Bourdin 8f52d5e2bSBlaise Bourdin Type(tDM) :: dm 9f52d5e2bSBlaise Bourdin Type(tPetscSection) :: section 10f52d5e2bSBlaise Bourdin Character(len=PETSC_MAX_PATH_LEN) :: IOBuffer 11f52d5e2bSBlaise Bourdin PetscInt :: dof,p,pStart,pEnd,d 12f52d5e2bSBlaise Bourdin Type(tVec) :: v 13f52d5e2bSBlaise Bourdin PetscInt :: zero = 0 14f52d5e2bSBlaise Bourdin PetscInt :: one = 1 15f52d5e2bSBlaise Bourdin PetscInt :: two = 2 16f52d5e2bSBlaise Bourdin PetscScalar,Dimension(:),Pointer :: val 17f52d5e2bSBlaise Bourdin PetscScalar, pointer :: x(:) 18f52d5e2bSBlaise Bourdin PetscErrorCode :: ierr 19f52d5e2bSBlaise Bourdin 20f52d5e2bSBlaise Bourdin PetscCallA(PetscInitialize(ierr)) 21f52d5e2bSBlaise Bourdin 22f52d5e2bSBlaise Bourdin PetscCallA(DMCreate(PETSC_COMM_WORLD, dm, ierr)) 23f52d5e2bSBlaise Bourdin PetscCallA(DMSetType(dm, DMPLEX, ierr)) 24f52d5e2bSBlaise Bourdin PetscCallA(DMSetFromOptions(dm, ierr)) 25ce78bad3SBarry Smith PetscCallA(DMViewFromOptions(dm,PETSC_NULL_OBJECT,'-d_view',ierr)) 26f52d5e2bSBlaise Bourdin 27f52d5e2bSBlaise Bourdin PetscCallA(PetscSectionCreate(PETSC_COMM_WORLD,section,ierr)) 28f52d5e2bSBlaise Bourdin PetscCallA(DMPlexGetChart(dm,pStart,pEnd,ierr)) 29f52d5e2bSBlaise Bourdin PetscCallA(PetscSectionSetChart(section, pStart, pEnd,ierr)) 30f52d5e2bSBlaise Bourdin PetscCallA(DMPlexGetHeightStratum(dm,zero,pStart,pEnd,ierr)) 31f52d5e2bSBlaise Bourdin Do p = pStart,pEnd-1 32f52d5e2bSBlaise Bourdin PetscCallA(PetscSectionSetDof(section,p,one,ierr)) 33f52d5e2bSBlaise Bourdin End Do 34f52d5e2bSBlaise Bourdin PetscCallA(DMPlexGetDepthStratum(dm,zero,pStart,pEnd,ierr)) 35f52d5e2bSBlaise Bourdin Do p = pStart,pEnd-1 36f52d5e2bSBlaise Bourdin PetscCallA(PetscSectionSetDof(section,p,two,ierr)) 37f52d5e2bSBlaise Bourdin End Do 38f52d5e2bSBlaise Bourdin PetscCallA(PetscSectionSetUp(section,ierr)) 39f52d5e2bSBlaise Bourdin PetscCallA(DMSetLocalSection(dm, section,ierr)) 40ce78bad3SBarry Smith PetscCallA(PetscSectionViewFromOptions(section,PETSC_NULL_OBJECT,'-s_view',ierr)) 41f52d5e2bSBlaise Bourdin 42f52d5e2bSBlaise Bourdin PetscCallA(DMCreateGlobalVector(dm,v,ierr)) 43f52d5e2bSBlaise Bourdin 44f52d5e2bSBlaise Bourdin PetscCallA(DMPlexGetChart(dm,pStart,pEnd,ierr)) 45f52d5e2bSBlaise Bourdin Do p = pStart,pEnd-1 46f52d5e2bSBlaise Bourdin PetscCallA(PetscSectionGetDof(section,p,dof,ierr)) 47f52d5e2bSBlaise Bourdin Allocate(val(dof)) 48f52d5e2bSBlaise Bourdin Do d = 1,dof 49*ccfd86f1SBarry Smith val(d) = 100*p + d-1 50f52d5e2bSBlaise Bourdin End Do 51ce78bad3SBarry Smith PetscCallA(VecSetValuesSection(v,section,p,val,INSERT_VALUES,ierr)) 52f52d5e2bSBlaise Bourdin DeAllocate(val) 53f52d5e2bSBlaise Bourdin End Do 54f52d5e2bSBlaise Bourdin PetscCallA(VecView(v,PETSC_VIEWER_STDOUT_WORLD,ierr)) 55f52d5e2bSBlaise Bourdin 56f52d5e2bSBlaise Bourdin Do p = pStart,pEnd-1 57f52d5e2bSBlaise Bourdin PetscCallA(PetscSectionGetDof(section,p,dof,ierr)) 58ce78bad3SBarry Smith PetscCallA(VecGetValuesSection(v,section,p,x,ierr)) 59dcb3e689SBarry Smith write(IOBuffer,*) 'Point ',p,' dof ',dof,'\n' 60f52d5e2bSBlaise Bourdin PetscCallA(PetscPrintf(PETSC_COMM_SELF,IOBuffer,ierr)) 61ce78bad3SBarry Smith PetscCallA(VecRestoreValuesSection(v,section,p,x,ierr)) 62f52d5e2bSBlaise Bourdin End Do 63f52d5e2bSBlaise Bourdin 64f52d5e2bSBlaise Bourdin PetscCallA(PetscSectionDestroy(section,ierr)) 65f52d5e2bSBlaise Bourdin PetscCallA(VecDestroy(v,ierr)) 66f52d5e2bSBlaise Bourdin PetscCallA(DMDestroy(dm,ierr)) 67f52d5e2bSBlaise Bourdin PetscCallA(PetscFinalize(ierr)) 68f52d5e2bSBlaise Bourdinend program ex47f90 69f52d5e2bSBlaise Bourdin 70e2447916SStefano Zampini!/*TEST 71e2447916SStefano Zampini! 72e2447916SStefano Zampini! test: 73e2447916SStefano Zampini! suffix: 0 74e2447916SStefano Zampini! args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/quads-q2.msh 75e2447916SStefano Zampini! 76e2447916SStefano Zampini!TEST*/ 77