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