xref: /petsc/src/dm/tests/ex54f.F90 (revision c5e229c2f66f66995aed5443a26600af2aec4a3f)
1! test verifies DMShellSetCreateFieldDecomposition interface in Fortran
2#include "petsc/finclude/petsc.h"
3#include "petsc/finclude/petscdmshell.h"
4program main
5  use petsc
6  use petscdmshell
7  implicit none
8  type(tDM)          :: dm
9  PetscErrorCode     :: ierr
10  interface
11    subroutine myFieldDecomp(dm, nfields, fieldNames, isFields, subDms, ierr)
12      use petsc
13      implicit none
14      type(tDM), intent(in) :: dm
15      PetscInt, intent(out) :: nfields
16      character(len=30), allocatable, intent(out) :: fieldNames(:)
17      type(tIS), allocatable, intent(out) :: isFields(:)
18      type(tDM), allocatable, intent(out) :: subDms(:)
19      PetscErrorCode, intent(out) :: ierr
20    end subroutine myFieldDecomp
21  end interface
22  ! initializing PETSc
23  PetscCallA(PetscInitialize(PETSC_NULL_CHARACTER, ierr))
24  ! creating a DMShell object
25  PetscCallA(DMShellCreate(PETSC_COMM_WORLD, dm, ierr))
26  ! registering the Fortran field decomposition callback
27  PetscCallA(DMShellSetCreateFieldDecomposition(dm, myFieldDecomp, ierr))
28  ! for this minimal test, we simply print a success message to the console
29  print *, 'DMShellSetCreateFieldDecomposition set successfully.'
30  ! cleanup
31  PetscCallA(DMDestroy(dm, ierr))
32  PetscCallA(PetscFinalize(ierr))
33end program main
34
35! a simple Fortran callback for field decomposition.
36subroutine myFieldDecomp(dm, nfields, fieldNames, isFields, subDms, ierr)
37  use petsc
38  implicit none
39  type(tDM), intent(in) :: dm
40  PetscInt, intent(out) :: nfields
41  character(len=30), allocatable, intent(out) :: fieldNames(:)
42  type(tIS), allocatable, intent(out) :: isFields(:)
43  type(tDM), allocatable, intent(out) :: subDms(:)
44  PetscErrorCode, intent(out) :: ierr
45  PetscInt :: i
46  ! defining a simple decomposition with two fields
47  nfields = 2
48  allocate (fieldNames(nfields))
49  allocate (isFields(nfields))
50  allocate (subDms(nfields))
51  fieldNames(1) = 'field1'
52  fieldNames(2) = 'field2'
53  ! set the pointer arrays to NULL (using pointer assignment)
54  do i = 1, nfields
55    isFields(i) = PETSC_NULL_IS
56    subDms(i) = PETSC_NULL_DM
57  end do
58  ierr = 0
59  print *, 'myFieldDecomp callback invoked.'
60end subroutine myFieldDecomp
61!/*TEST
62!
63!   test:
64!TEST*/
65