Actual source code: ex48f90.F90
1: #include "petsc/finclude/petsc.h"
2: program ex47f90
3: use petsc
4: implicit none
6: type(tDM) :: dm
7: type(tPetscSection) :: section
8: character(len=PETSC_MAX_PATH_LEN) :: IOBuffer
9: PetscInt :: dof, p, pStart, pEnd, d
10: type(tVec) :: v
11: PetscScalar, dimension(:), pointer :: val
12: PetscScalar, pointer :: x(:)
13: PetscErrorCode :: ierr
15: PetscCallA(PetscInitialize(ierr))
17: PetscCallA(DMCreate(PETSC_COMM_WORLD, dm, ierr))
18: PetscCallA(DMSetType(dm, DMPLEX, ierr))
19: PetscCallA(DMSetFromOptions(dm, ierr))
20: PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OBJECT, '-d_view', ierr))
22: PetscCallA(PetscSectionCreate(PETSC_COMM_WORLD, section, ierr))
23: PetscCallA(DMPlexGetChart(dm, pStart, pEnd, ierr))
24: PetscCallA(PetscSectionSetChart(section, pStart, pEnd, ierr))
25: PetscCallA(DMPlexGetHeightStratum(dm, 0_PETSC_INT_KIND, pStart, pEnd, ierr))
26: do p = pStart, pEnd - 1
27: PetscCallA(PetscSectionSetDof(section, p, 1_PETSC_INT_KIND, ierr))
28: end do
29: PetscCallA(DMPlexGetDepthStratum(dm, 0_PETSC_INT_KIND, pStart, pEnd, ierr))
30: do p = pStart, pEnd - 1
31: PetscCallA(PetscSectionSetDof(section, p, 2_PETSC_INT_KIND, ierr))
32: end do
33: PetscCallA(PetscSectionSetUp(section, ierr))
34: PetscCallA(DMSetLocalSection(dm, section, ierr))
35: PetscCallA(PetscSectionViewFromOptions(section, PETSC_NULL_OBJECT, '-s_view', ierr))
37: PetscCallA(DMCreateGlobalVector(dm, v, ierr))
39: PetscCallA(DMPlexGetChart(dm, pStart, pEnd, ierr))
40: do p = pStart, pEnd - 1
41: PetscCallA(PetscSectionGetDof(section, p, dof, ierr))
42: allocate (val(dof))
43: do d = 1, dof
44: val(d) = 100*p + d - 1
45: end do
46: PetscCallA(VecSetValuesSection(v, section, p, val, INSERT_VALUES, ierr))
47: deallocate (val)
48: end do
49: PetscCallA(VecView(v, PETSC_VIEWER_STDOUT_WORLD, ierr))
51: do p = pStart, pEnd - 1
52: PetscCallA(PetscSectionGetDof(section, p, dof, ierr))
53: PetscCallA(VecGetValuesSection(v, section, p, x, ierr))
54: write (IOBuffer, *) 'Point ', p, ' dof ', dof, '\n'
55: PetscCallA(PetscPrintf(PETSC_COMM_SELF, IOBuffer, ierr))
56: PetscCallA(VecRestoreValuesSection(v, section, p, x, ierr))
57: end do
59: PetscCallA(PetscSectionDestroy(section, ierr))
60: PetscCallA(VecDestroy(v, ierr))
61: PetscCallA(DMDestroy(dm, ierr))
62: PetscCallA(PetscFinalize(ierr))
63: end program ex47f90
65: !/*TEST
66: !
67: ! test:
68: ! suffix: 0
69: ! args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/quads-q2.msh
70: !
71: !TEST*/