xref: /petsc/src/dm/impls/plex/tests/ex48f90.F90 (revision d52a580b706c59ca78066c1e38754e45b6b56e2b)
1#include "petsc/finclude/petsc.h"
2program ex47f90
3  use petsc
4  implicit none
5
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  PetscInt                          :: zero = 0
12  PetscInt                          :: one = 1
13  PetscInt                          :: two = 2
14  PetscScalar, dimension(:), pointer  :: val
15  PetscScalar, pointer              :: x(:)
16  PetscErrorCode                    :: ierr
17
18  PetscCallA(PetscInitialize(ierr))
19
20  PetscCallA(DMCreate(PETSC_COMM_WORLD, dm, ierr))
21  PetscCallA(DMSetType(dm, DMPLEX, ierr))
22  PetscCallA(DMSetFromOptions(dm, ierr))
23  PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OBJECT, '-d_view', ierr))
24
25  PetscCallA(PetscSectionCreate(PETSC_COMM_WORLD, section, ierr))
26  PetscCallA(DMPlexGetChart(dm, pStart, pEnd, ierr))
27  PetscCallA(PetscSectionSetChart(section, pStart, pEnd, ierr))
28  PetscCallA(DMPlexGetHeightStratum(dm, zero, pStart, pEnd, ierr))
29  do p = pStart, pEnd - 1
30    PetscCallA(PetscSectionSetDof(section, p, one, ierr))
31  end do
32  PetscCallA(DMPlexGetDepthStratum(dm, zero, pStart, pEnd, ierr))
33  do p = pStart, pEnd - 1
34    PetscCallA(PetscSectionSetDof(section, p, two, ierr))
35  end do
36  PetscCallA(PetscSectionSetUp(section, ierr))
37  PetscCallA(DMSetLocalSection(dm, section, ierr))
38  PetscCallA(PetscSectionViewFromOptions(section, PETSC_NULL_OBJECT, '-s_view', ierr))
39
40  PetscCallA(DMCreateGlobalVector(dm, v, ierr))
41
42  PetscCallA(DMPlexGetChart(dm, pStart, pEnd, ierr))
43  do p = pStart, pEnd - 1
44    PetscCallA(PetscSectionGetDof(section, p, dof, ierr))
45    allocate (val(dof))
46    do d = 1, dof
47      val(d) = 100*p + d - 1
48    end do
49    PetscCallA(VecSetValuesSection(v, section, p, val, INSERT_VALUES, ierr))
50    deallocate (val)
51  end do
52  PetscCallA(VecView(v, PETSC_VIEWER_STDOUT_WORLD, ierr))
53
54  do p = pStart, pEnd - 1
55    PetscCallA(PetscSectionGetDof(section, p, dof, ierr))
56    PetscCallA(VecGetValuesSection(v, section, p, x, ierr))
57    write (IOBuffer, *) 'Point ', p, ' dof ', dof, '\n'
58    PetscCallA(PetscPrintf(PETSC_COMM_SELF, IOBuffer, ierr))
59    PetscCallA(VecRestoreValuesSection(v, section, p, x, ierr))
60  end do
61
62  PetscCallA(PetscSectionDestroy(section, ierr))
63  PetscCallA(VecDestroy(v, ierr))
64  PetscCallA(DMDestroy(dm, ierr))
65  PetscCallA(PetscFinalize(ierr))
66end program ex47f90
67
68!/*TEST
69!
70!  test:
71!    suffix: 0
72!    args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/quads-q2.msh
73!
74!TEST*/
75