xref: /petsc/src/dm/impls/plex/tutorials/ex3f90.F90 (revision 9b88ac225e01f016352a5f4cd90e158abe5f5675)
1c4762a1bSJed Brown! setting up 3-D DMPlex using DMPlexCreateFromDAG(), DMPlexInterpolate() and
2c4762a1bSJed Brown! DMPlexComputeCellGeometryFVM()
3c4762a1bSJed Brown! Contributed by Adrian Croucher <a.croucher@auckland.ac.nz>
4c4762a1bSJed Brown#include <petsc/finclude/petscdmplex.h>
5*c5e229c2SMartin Diehlprogram main
6ce78bad3SBarry Smith  use petscdmplex
7c4762a1bSJed Brown  implicit none
8c4762a1bSJed Brown  DM :: dm, dmi
9c4762a1bSJed Brown  PetscFV :: fvm
10c4762a1bSJed Brown  PetscInt, parameter :: dim = 3
11c4762a1bSJed Brown
12c4762a1bSJed Brown  PetscInt :: depth = 1
13c4762a1bSJed Brown  PetscErrorCode :: ierr
14c4762a1bSJed Brown  PetscInt, dimension(2) :: numPoints
15c4762a1bSJed Brown  PetscInt, dimension(14) :: coneSize
16c4762a1bSJed Brown  PetscInt, dimension(16) :: cones
17c4762a1bSJed Brown  PetscInt, dimension(16) :: coneOrientations
18c4762a1bSJed Brown  PetscScalar, dimension(36) :: vertexCoords
19c4762a1bSJed Brown  PetscReal ::  vol = 0.
20c4762a1bSJed Brown  PetscReal, target, dimension(dim)  :: centroid
21c4762a1bSJed Brown  PetscReal, target, dimension(dim)  :: normal
22c4762a1bSJed Brown  PetscReal, pointer :: pcentroid(:)
23c4762a1bSJed Brown  PetscReal, pointer :: pnormal(:)
24c4762a1bSJed Brown
25c4762a1bSJed Brown  PetscReal, target, dimension(dim)  :: v0
26c4762a1bSJed Brown  PetscReal, target, dimension(dim*dim)  :: J
27c4762a1bSJed Brown  PetscReal, target, dimension(dim*dim)  :: invJ
28c4762a1bSJed Brown  PetscReal, pointer :: pv0(:)
29c4762a1bSJed Brown  PetscReal, pointer :: pJ(:)
30c4762a1bSJed Brown  PetscReal, pointer :: pinvJ(:)
31c4762a1bSJed Brown  PetscReal :: detJ
32c4762a1bSJed Brown
33c4762a1bSJed Brown  PetscInt :: i
34c4762a1bSJed Brown
35c4762a1bSJed Brown  numPoints = [12, 2]
36c4762a1bSJed Brown  coneSize = [8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0]
37c4762a1bSJed Brown  cones = [2, 5, 4, 3, 6, 7, 8, 9, 3, 4, 11, 10, 7, 12, 13, 8]
38c4762a1bSJed Brown  coneOrientations = [0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0]
39c4762a1bSJed Brown  vertexCoords = [                                                  &
40c4762a1bSJed Brown &  -0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, -0.5, 1.0, 0.0,           &
41c4762a1bSJed Brown &  -0.5, 0.0, 1.0, 0.0, 0.0, 1.0, 0.0, 1.0, 1.0, -0.5, 1.0, 1.0,           &
42c4762a1bSJed Brown &   0.5, 0.0, 0.0, 0.5, 1.0, 0.0, 0.5, 0.0, 1.0, 0.5, 1.0, 1.0]
43c4762a1bSJed Brown  pcentroid => centroid
44c4762a1bSJed Brown  pnormal => normal
45c4762a1bSJed Brown
46c4762a1bSJed Brown  pv0 => v0
47c4762a1bSJed Brown  pJ => J
48c4762a1bSJed Brown  pinvJ => invJ
49c4762a1bSJed Brown
50d8606c27SBarry Smith  PetscCallA(PetscInitialize(ierr))
51d8606c27SBarry Smith  PetscCallA(DMPlexCreate(PETSC_COMM_WORLD, dm, ierr))
52d8606c27SBarry Smith  PetscCallA(PetscObjectSetName(dm, 'testplex', ierr))
53d8606c27SBarry Smith  PetscCallA(DMSetDimension(dm, dim, ierr))
54c4762a1bSJed Brown
55d8606c27SBarry Smith  PetscCallA(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords, ierr))
56c4762a1bSJed Brown
57d8606c27SBarry Smith  PetscCallA(DMPlexInterpolate(dm, dmi, ierr))
58d8606c27SBarry Smith  PetscCallA(DMPlexCopyCoordinates(dm, dmi, ierr))
59d8606c27SBarry Smith  PetscCallA(DMDestroy(dm, ierr))
60c4762a1bSJed Brown  dm = dmi
61c4762a1bSJed Brown
62d8606c27SBarry Smith  PetscCallA(DMView(dm, PETSC_VIEWER_STDOUT_WORLD, ierr))
63c4762a1bSJed Brown
64c4762a1bSJed Brown  do i = 0, 1
65d8606c27SBarry Smith    PetscCallA(DMPlexComputeCellGeometryFVM(dm, i, vol, pcentroid, pnormal, ierr))
66c4762a1bSJed Brown    write (*, '(a, i2, a, f8.4, a, 3(f8.4, 1x))') 'cell: ', i, ' volume: ', vol, ' centroid: ', pcentroid(1), pcentroid(2), pcentroid(3)
67d8606c27SBarry Smith    PetscCallA(DMPlexComputeCellGeometryAffineFEM(dm, i, pv0, pJ, pinvJ, detJ, ierr))
68c4762a1bSJed Brown  end do
69c4762a1bSJed Brown
70d8606c27SBarry Smith  PetscCallA(PetscFVCreate(PETSC_COMM_WORLD, fvm, ierr))
71d8606c27SBarry Smith  PetscCallA(PetscFVSetUp(fvm, ierr))
72d8606c27SBarry Smith  PetscCallA(PetscFVDestroy(fvm, ierr))
73c4762a1bSJed Brown
74d8606c27SBarry Smith  PetscCallA(DMDestroy(dm, ierr))
75d8606c27SBarry Smith  PetscCallA(PetscFinalize(ierr))
76c4762a1bSJed Brownend program main
77c4762a1bSJed Brown
78c4762a1bSJed Brown!/*TEST
79c4762a1bSJed Brown!
80c4762a1bSJed Brown!   test:
81c4762a1bSJed Brown!      suffix: 0
82c4762a1bSJed Brown!
83c4762a1bSJed Brown!TEST*/
84