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