xref: /petsc/src/dm/impls/plex/tutorials/ex4f90.F90 (revision 9b88ac225e01f016352a5f4cd90e158abe5f5675)
1c4762a1bSJed Brown! setting up DMPlex for finite elements
2c4762a1bSJed Brown! Contributed by Pratheek Shanthraj <p.shanthraj@mpie.de>
3ce78bad3SBarry Smith#include <petsc/finclude/petsc.h>
4*c5e229c2SMartin Diehlprogram main
5ce78bad3SBarry Smith  use petsc
6c4762a1bSJed Brown  implicit none
7c4762a1bSJed Brown  DM :: dm
8ce78bad3SBarry Smith  PetscDS :: ds
9ce78bad3SBarry Smith  PetscInt :: dim = 3, zero = 0
10c4762a1bSJed Brown  PetscBool :: simplex = PETSC_TRUE
11c4762a1bSJed Brown  PetscBool :: interpolate = PETSC_TRUE
12c4762a1bSJed Brown  PetscReal :: refinementLimit = 0.0
13c4762a1bSJed Brown  PetscErrorCode :: ierr
14ce78bad3SBarry Smith  PetscTabulation, pointer :: tab(:)
15ce78bad3SBarry Smith  PetscFE fe, rfe
16ce78bad3SBarry Smith  PetscObject obj
17ce78bad3SBarry Smith  PetscInt :: one = 1, mone = -1
18c4762a1bSJed Brown
19d8606c27SBarry Smith  PetscCallA(PetscInitialize(PETSC_NULL_CHARACTER, ierr))
20ce78bad3SBarry Smith  PetscCallA(DMPlexCreateDoublet(PETSC_COMM_WORLD, dim, simplex, interpolate, refinementLimit, dm, ierr))
21ccfd86f1SBarry Smith  PetscCallA(PetscFECreateDefault(PETSC_COMM_WORLD, dim, one, simplex, 'name', mone, fe, ierr))
22ccfd86f1SBarry Smith  PetscCallA(PetscObjectSetName(fe, 'name', ierr))
23ccfd86f1SBarry Smith  PetscCallA(DMSetField(dm, zero, PETSC_NULL_DMLABEL, PetscObjectCast(fe), ierr))
24ccfd86f1SBarry Smith  PetscCallA(DMSetField(dm, one, PETSC_NULL_DMLABEL, PetscObjectCast(fe), ierr))
25c4762a1bSJed Brown
26ce78bad3SBarry Smith  PetscCallA(DMSetUp(dm, ierr))
27ce78bad3SBarry Smith  PetscCallA(DMCreateDS(dm, ierr))
28ce78bad3SBarry Smith  PetscCallA(DMGetDS(dm, ds, ierr))
29ce78bad3SBarry Smith  PetscCallA(PetscDSGetTabulation(ds, tab, ierr))
30ce78bad3SBarry Smith  print *, tab(1)%ptr%T(1)%ptr
31ce78bad3SBarry Smith  print *, tab(1)%ptr%T(2)%ptr
32ce78bad3SBarry Smith  print *, tab(2)%ptr%T(1)%ptr
33ce78bad3SBarry Smith  print *, tab(2)%ptr%T(2)%ptr
34ce78bad3SBarry Smith  PetscCallA(PetscDSRestoreTabulation(ds, tab, ierr))
35ce78bad3SBarry Smith
36ce78bad3SBarry Smith  PetscCallA(PetscDSGetDiscretization(ds, zero, obj, ierr))
37ce78bad3SBarry Smith  PetscObjectSpecificCast(rfe, obj)
38ccfd86f1SBarry Smith  PetscCallA(PetscFEDestroy(fe, ierr))
39d8606c27SBarry Smith  PetscCallA(DMDestroy(dm, ierr))
40d8606c27SBarry Smith  PetscCallA(PetscFinalize(ierr))
41c4762a1bSJed Brownend program main
42ce78bad3SBarry Smith!/*TEST
43ce78bad3SBarry Smith!
44ce78bad3SBarry Smith!  test:
45ce78bad3SBarry Smith!    nsize: 1
46ce78bad3SBarry Smith!
47ce78bad3SBarry Smith!TEST*/
48