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