1
1
openmpi/ompi/mpi/f90/scripts/mpi_sizeof.f90.sh
Jeff Squyres faf63c68f8 Merge over from the /tmp/fortran-stuff branch
- split mpif.h into mpif.h and mpif-common.h[.in]
- mpif-common.h is included by various f90 things and contains output
  from configure
- mpif.h defines some f77-specific stuff and then includes
  mpif-common.h 

This commit was SVN r9997.
2006-05-20 02:15:49 +00:00

175 строки
5.2 KiB
Bash
Исполняемый файл

#! /bin/sh
#
# Copyright (c) 2004-2006 The Trustees of Indiana University and Indiana
# University Research and Technology
# Corporation. All rights reserved.
# Copyright (c) 2004-2005 The University of Tennessee and The University
# of Tennessee Research Foundation. All rights
# reserved.
# Copyright (c) 2004-2005 High Performance Computing Center Stuttgart,
# University of Stuttgart. All rights reserved.
# Copyright (c) 2004-2005 The Regents of the University of California.
# All rights reserved.
# Copyright (c) 2006 Cisco Systems, Inc. All rights reserved.
# $COPYRIGHT$
#
# Additional copyrights may follow
#
# $HEADER$
#
. "$1/fortran_kinds.sh"
procedure='MPI_Sizeof'
rank=0
for kind in $lkinds
do
proc="${procedure}${rank}DL${kind}"
echo "subroutine ${proc}(x, size, ierr)"
echo " implicit none"
echo " include 'mpif-common.h'"
echo " include 'fortran_sizes.h'"
echo " logical*${kind}, intent(in) :: x"
echo " integer, intent(out) :: size"
echo " integer, intent(out) :: ierr"
echo " size = OMPI_SIZEOF_F90_LOGICAL${kind}"
echo " ierr = 0"
echo "end subroutine ${proc}"
echo
done
for kind in $ikinds
do
proc="${procedure}${rank}DI${kind}"
echo "subroutine ${proc}(x, size, ierr)"
echo " implicit none"
echo " include 'mpif-common.h'"
echo " include 'fortran_sizes.h'"
echo " integer*${kind}, intent(in) :: x"
echo " integer, intent(out) :: size"
echo " integer, intent(out) :: ierr"
echo " size = OMPI_SIZEOF_F90_INT${kind}"
echo " ierr = 0"
echo "end subroutine ${proc}"
echo
done
for kind in $rkinds
do
proc="${procedure}${rank}DR${kind}"
echo "subroutine ${proc}(x, size, ierr)"
echo " implicit none"
echo " include 'mpif-common.h'"
echo " include 'fortran_sizes.h'"
echo " real*${kind}, intent(in) :: x"
echo " integer, intent(out) :: size"
echo " integer, intent(out) :: ierr"
echo " size = OMPI_SIZEOF_F90_REAL${kind}"
echo " ierr = 0"
echo "end subroutine ${proc}"
echo
done
for kind in $ckinds
do
case "$kind" in 4) size_kind='8' ; esac
case "$kind" in 8) size_kind='16' ; esac
case "$kind" in 16) size_kind='32' ; esac
proc="${procedure}${rank}DC${kind}"
echo "subroutine ${proc}(x, size, ierr)"
echo " implicit none"
echo " include 'mpif-common.h'"
echo " include 'fortran_sizes.h'"
echo " complex*${kind}, intent(in) :: x"
echo " integer, intent(out) :: size"
echo " integer, intent(out) :: ierr"
echo " size = OMPI_SIZEOF_F90_COMPLEX${size_kind}"
echo " ierr = 0"
echo "end subroutine ${proc}"
echo
done
for rank in $ranks
do
case "$rank" in 1) dim=':' ; esac
case "$rank" in 2) dim=':,:' ; esac
case "$rank" in 3) dim=':,:,:' ; esac
case "$rank" in 4) dim=':,:,:,:' ; esac
case "$rank" in 5) dim=':,:,:,:,:' ; esac
case "$rank" in 6) dim=':,:,:,:,:,:' ; esac
case "$rank" in 7) dim=':,:,:,:,:,:,:' ; esac
for kind in $lkinds
do
proc="${procedure}${rank}DL${kind}"
echo "subroutine ${proc}(x, size, ierr)"
echo " implicit none"
echo " include 'mpif-common.h'"
echo " include 'fortran_sizes.h'"
echo " logical*${kind}, dimension(${dim}), intent(in) :: x"
echo " integer, intent(out) :: size"
echo " integer, intent(out) :: ierr"
echo " size = OMPI_SIZEOF_F90_LOGICAL${kind}"
echo " ierr = 0"
echo "end subroutine ${proc}"
echo
done
for kind in $ikinds
do
proc="${procedure}${rank}DI${kind}"
echo "subroutine ${proc}(x, size, ierr)"
echo " implicit none"
echo " include 'mpif-common.h'"
echo " include 'fortran_sizes.h'"
echo " integer*${kind}, dimension(${dim}), intent(in) :: x"
echo " integer, intent(out) :: size"
echo " integer, intent(out) :: ierr"
echo " size = OMPI_SIZEOF_F90_INT${kind}"
echo " ierr = 0"
echo "end subroutine ${proc}"
echo
done
for kind in $rkinds
do
proc="${procedure}${rank}DR${kind}"
echo "subroutine ${proc}(x, size, ierr)"
echo " implicit none"
echo " include 'mpif-common.h'"
echo " include 'fortran_sizes.h'"
echo " real*${kind}, dimension(${dim}), intent(in) :: x"
echo " integer, intent(out) :: size"
echo " integer, intent(out) :: ierr"
echo " size = OMPI_SIZEOF_F90_REAL${kind}"
echo " ierr = 0"
echo "end subroutine ${proc}"
echo
done
for kind in $ckinds
do
case "$kind" in 4) size_kind='8' ; esac
case "$kind" in 8) size_kind='16' ; esac
case "$kind" in 16) size_kind='32' ; esac
proc="${procedure}${rank}DC${kind}"
echo "subroutine ${proc}(x, size, ierr)"
echo " implicit none"
echo " include 'mpif-common.h'"
echo " include 'fortran_sizes.h'"
echo " complex*${kind}, dimension(${dim}), intent(in) :: x"
echo " integer, intent(out) :: size"
echo " integer, intent(out) :: ierr"
echo " size = OMPI_SIZEOF_F90_COMPLEX${size_kind}"
echo " ierr = 0"
echo "end subroutine ${proc}"
echo
done
echo
done
echo
echo
echo