xref: /petsc/src/dm/tests/ex54f.F90 (revision 9b88ac225e01f016352a5f4cd90e158abe5f5675)
162212064STapashree Pradhan! test verifies DMShellSetCreateFieldDecomposition interface in Fortran
262212064STapashree Pradhan#include "petsc/finclude/petsc.h"
3*c5e229c2SMartin Diehlprogram main
462212064STapashree Pradhan  use petsc
562212064STapashree Pradhan  implicit none
662212064STapashree Pradhan  type(tDM)          :: dm
762212064STapashree Pradhan  PetscErrorCode     :: ierr
862212064STapashree Pradhan  interface
962212064STapashree Pradhan    subroutine myFieldDecomp(dm, nfields, fieldNames, isFields, subDms, ierr)
1062212064STapashree Pradhan      use petsc
1162212064STapashree Pradhan      implicit none
1262212064STapashree Pradhan      type(tDM), intent(in) :: dm
1362212064STapashree Pradhan      PetscInt, intent(out) :: nfields
1462212064STapashree Pradhan      character(len=30), allocatable, intent(out) :: fieldNames(:)
1562212064STapashree Pradhan      type(tIS), allocatable, intent(out) :: isFields(:)
1662212064STapashree Pradhan      type(tDM), allocatable, intent(out) :: subDms(:)
1762212064STapashree Pradhan      PetscErrorCode, intent(out) :: ierr
1862212064STapashree Pradhan    end subroutine myFieldDecomp
1962212064STapashree Pradhan  end interface
2062212064STapashree Pradhan  ! initializing PETSc
2162212064STapashree Pradhan  PetscCallA(PetscInitialize(PETSC_NULL_CHARACTER, ierr))
2262212064STapashree Pradhan  ! creating a DMShell object
2362212064STapashree Pradhan  PetscCallA(DMShellCreate(PETSC_COMM_WORLD, dm, ierr))
2462212064STapashree Pradhan  ! registering the Fortran field decomposition callback
2562212064STapashree Pradhan  PetscCallA(DMShellSetCreateFieldDecomposition(dm, myFieldDecomp, ierr))
2662212064STapashree Pradhan  ! for this minimal test, we simply print a success message to the console
2762212064STapashree Pradhan  print *, 'DMShellSetCreateFieldDecomposition set successfully.'
2862212064STapashree Pradhan  ! cleanup
2962212064STapashree Pradhan  PetscCallA(DMDestroy(dm, ierr))
3062212064STapashree Pradhan  PetscCallA(PetscFinalize(ierr))
3162212064STapashree Pradhanend program main
3262212064STapashree Pradhan
3362212064STapashree Pradhan! a simple Fortran callback for field decomposition.
3462212064STapashree Pradhansubroutine myFieldDecomp(dm, nfields, fieldNames, isFields, subDms, ierr)
3562212064STapashree Pradhan  use petsc
3662212064STapashree Pradhan  implicit none
3762212064STapashree Pradhan  type(tDM), intent(in) :: dm
3862212064STapashree Pradhan  PetscInt, intent(out) :: nfields
3962212064STapashree Pradhan  character(len=30), allocatable, intent(out) :: fieldNames(:)
4062212064STapashree Pradhan  type(tIS), allocatable, intent(out) :: isFields(:)
4162212064STapashree Pradhan  type(tDM), allocatable, intent(out) :: subDms(:)
4262212064STapashree Pradhan  PetscErrorCode, intent(out) :: ierr
4362212064STapashree Pradhan  PetscInt :: i
4462212064STapashree Pradhan  ! defining a simple decomposition with two fields
4562212064STapashree Pradhan  nfields = 2
4662212064STapashree Pradhan  allocate (fieldNames(nfields))
4762212064STapashree Pradhan  allocate (isFields(nfields))
4862212064STapashree Pradhan  allocate (subDms(nfields))
4962212064STapashree Pradhan  fieldNames(1) = 'field1'
5062212064STapashree Pradhan  fieldNames(2) = 'field2'
5162212064STapashree Pradhan  ! set the pointer arrays to NULL (using pointer assignment)
5262212064STapashree Pradhan  do i = 1, nfields
5362212064STapashree Pradhan    isFields(i) = PETSC_NULL_IS
5462212064STapashree Pradhan    subDms(i) = PETSC_NULL_DM
5562212064STapashree Pradhan  end do
5662212064STapashree Pradhan  ierr = 0
5762212064STapashree Pradhan  print *, 'myFieldDecomp callback invoked.'
5862212064STapashree Pradhanend subroutine myFieldDecomp
5962212064STapashree Pradhan!/*TEST
6062212064STapashree Pradhan!
6162212064STapashree Pradhan!   test:
6262212064STapashree Pradhan!TEST*/
63