Next: Error Codes, Previous: Summary of Fortran 77 Interface, Up: Top [Index]
character(len=80) function nf90mpi_inq_libvers()
character(len=80) function nf90mpi_strerror(ncerr)
integer, intent(in) :: ncerr
logical function nf90mpi_issyserr(ncerr)
integer, intent(in) :: ncerr
character(len=80) function nf90mpi_strerrno(ncerr)
integer, intent(in) :: ncerr
integer function nf90mpi_create(mpi_comm, path, cmode, mpi_info, ncid)
character(len=*), intent( in) :: path
integer, intent( in) :: mpi_comm
integer, intent( in) :: cmode
integer, intent( in) :: mpi_info
integer, intent(out) :: ncid
integer function nf90mpi_open(mpi_comm, path, omode, mpi_info, ncid)
character(len=*), intent( in) :: path
integer, intent( in) :: mpi_comm
integer, intent( in) :: omode
integer, intent( in) :: mpi_info
integer, intent(out) :: ncid
integer function nf90mpi_enddef(ncid, h_minfree, v_align, v_minfree, r_align)
integer, intent(in) :: ncid
integer(kind=MPI_OFFSET_KIND), optional, intent(in) :: h_minfree
integer(kind=MPI_OFFSET_KIND), optional, intent(in) :: v_align
integer(kind=MPI_OFFSET_KIND), optional, intent(in) :: v_minfree
integer(kind=MPI_OFFSET_KIND), optional, intent(in) :: r_align
integer function nf90mpi_redef(ncid)
integer, intent(in) :: ncid
integer function nf90mpi_flush(ncid)
integer, intent(in) :: ncid
integer function nf90mpi_sync(ncid)
integer, intent(in) :: ncid
integer function nf90mpi_sync_numrecs(ncid)
integer, intent(in) :: ncid
integer function nf90mpi_close(ncid)
integer, intent(in) :: ncid
integer function nf90mpi_set_fill(ncid, fillmode, old_mode)
integer, intent( in) :: ncid
integer, intent( in) :: fillmode
integer, intent(out) :: old_mode
integer function nf90mpi_set_default_format(new_format, old_format)
integer, intent( in) :: new_format
integer, intent(out) :: old_format
integer function nf90mpi_abort(ncid)
integer, intent(in) :: ncid
integer function nf90mpi_delete(name, mpi_info)
character(len=*), intent(in) :: name
integer, intent(in) :: mpi_info
integer function nf90mpi_begin_indep_data(ncid)
integer, intent(in) :: ncid
integer function nf90mpi_end_indep_data(ncid)
integer, intent(in) :: ncid
integer function nf90mpi_wait(ncid, num, req, st)
integer, intent(in) :: ncid
integer, intent(in) :: num
integer, dimension(:), intent(inout) :: req
integer, dimension(:), intent(out) :: st
integer function nf90mpi_wait_all(ncid, num, req, st)
integer, intent(in) :: ncid
integer, intent(in) :: num
integer, dimension(:), intent(inout) :: req
integer, dimension(:), intent(out) :: st
integer function nf90mpi_cancel(ncid, num, req, st)
integer, intent(in) :: ncid
integer, intent(in) :: num
integer, dimension(:), intent(inout) :: req
integer, dimension(:), intent(out) :: st
integer function nf90mpi_inquire(ncid, nDimensions, nVariables, nAttributes, unlimitedDimId, format)
integer, intent( in) :: ncid
integer, optional, intent(out) :: nDimensions
integer, optional, intent(out) :: nVariables
integer, optional, intent(out) :: nAttributes
integer, optional, intent(out) :: unlimitedDimId
integer, optional, intent(out) :: format
integer function nf90mpi_inq_file_format(path, format)
character(len=*), intent( in) :: path
integer, intent(out) :: format
integer function nf90mpi_inq_default_format(default_format)
integer, intent(out) :: default_format
integer function nf90mpi_inq_header_size(ncid, h_size)
integer, intent( in) :: ncid
integer(kind=MPI_OFFSET_KIND), intent(out) :: h_size
integer function nf90mpi_inq_header_extent(ncid, h_extent)
integer, intent( in) :: ncid
integer(kind=MPI_OFFSET_KIND), intent(out) :: h_extent
integer function nf90mpi_inq_files_opened(nfiles, ncids)
integer, intent(out) :: nfiles
integer, dimension(:), intent(out) :: ncids
integer function nf90mpi_inq_path(ncid, pathlen, path)
integer, intent( in) :: ncid
integer, intent(out) :: pathlen
character(len=*), intent(out) :: path
integer function nf90mpi_inq_num_rec_vars(ncid, nvars)
integer, intent( in) :: ncid
integer, intent(out) :: nvars
integer function nf90mpi_inq_num_fix_vars(ncid, nvars)
integer, intent( in) :: ncid
integer, intent(out) :: nvars
integer function nf90mpi_inq_recsize(ncid, recsize)
integer, intent( in) :: ncid
integer(kind=MPI_OFFSET_KIND), intent(out) :: recsize
integer function nf90mpi_inq_file_info(ncid, mpi_info)
integer, intent( in) :: ncid
integer, intent(out) :: mpi_info
integer function nf90mpi_get_file_info(ncid, mpi_info)
integer, intent( in) :: ncid
integer, intent(out) :: mpi_info
integer function nf90mpi_inq_striping(ncid, striping_size, striping_count)
integer, intent( in) :: ncid
integer, intent(out) :: striping_size
integer, intent(out) :: striping_count
integer function nf90mpi_inq_put_size(ncid, put_size)
integer, intent( in) :: ncid
integer(kind=MPI_OFFSET_KIND), intent(out) :: put_size
integer function nf90mpi_inq_get_size(ncid, get_size)
integer, intent( in) :: ncid
integer(kind=MPI_OFFSET_KIND), intent(out) :: get_size
integer function nf90mpi_inq_nreqs(ncid, nreqs)
integer, intent( in) :: ncid
integer, intent(out) :: nreqs
integer function nf90mpi_inq_malloc_size(size)
integer(kind=MPI_OFFSET_KIND), intent(out) :: size
integer function nf90mpi_inq_malloc_max_size(size)
integer(kind=MPI_OFFSET_KIND), intent(out) :: size
integer function nf90mpi_inq_malloc_list()
integer :: nf90mpi_inq_malloc_list
integer function nf90mpi_def_dim(ncid, name, len, dimid)
integer, intent( in) :: ncid
character(len=*), intent( in) :: name
integer(kind=MPI_OFFSET_KIND), intent( in) :: len
integer, intent(out) :: dimid
integer function nf90mpi_inq_dimid(ncid, name, dimid)
integer, intent( in) :: ncid
character(len=*), intent( in) :: name
integer, intent(out) :: dimid
integer function nf90mpi_inquire_dimension(ncid, dimid, name, len)
integer, intent( in) :: ncid
integer, intent( in) :: dimid
character(len=*), optional, intent(out) :: name
integer(kind=MPI_OFFSET_KIND), optional, intent(out) :: len
integer function nf90mpi_rename_dim(ncid, dimid, name)
integer, intent( in) :: ncid
character(len=*), intent( in) :: name
integer, intent( in) :: dimid
integer function nf90mpi_inquire_attribute(ncid, varid, name, xtype, len, attnum)
integer, intent( in) :: ncid
integer, intent( in) :: varid
character(len = *), intent( in) :: name
integer, optional, intent(out) :: xtype
integer, optional, intent(out) :: attnum
integer(kind=MPI_OFFSET_KIND), optional, intent(out) :: len
integer function nf90mpi_inq_attname(ncid, varid, attnum, name)
integer, intent( in) :: ncid
integer, intent( in) :: varid
integer, intent( in) :: attnum
character(len=*), intent(out) :: name
integer function nf90mpi_copy_att(ncid_in, varid_in, name, ncid_out, varid_out)
integer, intent(in) :: ncid_in
integer, intent(in) :: varid_in
character(len=*), intent(in) :: name
integer, intent(in) :: ncid_out
integer, intent(in) :: varid_out
integer function nf90mpi_rename_att(ncid, varid, curname, newname)
integer, intent(in) :: ncid
integer, intent(in) :: varid
character(len=*), intent(in) :: curname
character(len=*), intent(in) :: newname
integer function nf90mpi_del_att(ncid, varid, name)
integer, intent(in) :: ncid
integer, intent(in) :: varid
character(len=*), intent(in) :: name
integer function nf90mpi_put_att(ncid, varid, name, values)
integer, intent(in) :: ncid
integer, intent(in) :: varid
character(len=*), intent(in) :: name
<any valid type, scalar or array of rank 1>, intent(in) :: values
integer function nf90mpi_get_att(ncid, varid, name, values)
integer, intent( in) :: ncid
integer, intent( in) :: varid
character(len=*), intent( in) :: name
<any valid type, scalar or array of rank 1>, intent(out) :: values
integer function nf90mpi_def_var(ncid, name, xtype, dimids, varid)
integer, intent( in) :: ncid
character(len = *), intent( in) :: name
integer, intent( in) :: xtype
integer, dimension(:), intent( in) :: dimids ! omitted for scalar
integer function nf90mpi_inq_varid(ncid, name, varid)
integer, intent( in) :: ncid
character(len=*), intent( in) :: name
integer, intent(out) :: varid
integer function nf90mpi_inquire_variable(ncid, varid, name, xtype, ndims, dimids, nAtts)
integer, intent( in) :: ncid
integer, intent( in) :: varid
character(len=*), optional, intent(out) :: name
integer, optional, intent(out) :: xtype
integer, optional, intent(out) :: ndims
integer, dimension(*), optional, intent(out) :: dimids
integer, optional, intent(out) :: nAtts
integer function nf90mpi_inq_varoffset(ncid, varid, offset)
integer, intent( in) :: ncid, varid
integer(kind=MPI_OFFSET_KIND), intent(out) :: offset
integer function nf90mpi_rename_var(ncid, varid, newname)
integer, intent( in) :: ncid
integer, intent( in) :: varid
character(len=*), intent( in) :: newname
integer function nf90mpi_def_var_fill(ncid, varid, no_fill, fill_value)
integer, intent( in) :: ncid
integer, intent( in) :: varid
integer, intent( in) :: no_fill
<any valid type>, intent( in) :: fill_value
integer function nf90mpi_inq_var_fill(ncid, varid, no_fill, fill_value)
integer, intent( in) :: ncid
integer, intent( in) :: varid
integer, intent(out) :: no_fill
<any valid type>, intent(out) :: fill_value
integer function nf90mpi_fill_var_rec(ncid, varid, recno)
integer, intent( in) :: ncid
integer, intent( in) :: varid
integer(kind=MPI_OFFSET_KIND), intent( in) :: recno
integer function nf90mpi_put_var(ncid, varid, values, start, stride, map, bufcount, buftype)
integer, intent(in ) :: ncid
integer, intent(in ) :: varid
<any valid type, scalar or array of any rank>, intent(inout) :: values
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: start
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: count
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: stride
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: map
integer (kind=MPI_OFFSET_KIND), optional, intent(in ) :: bufcount
integer, optional, intent(in ) :: buftype
integer function nf90mpi_put_var_all(ncid, varid, values, start, stride, map, bufcount, buftype)
integer, intent(in ) :: ncid
integer, intent(in ) :: varid
<any valid type, scalar or array of any rank>, intent(inout) :: values
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: start
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: count
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: stride
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: map
integer (kind=MPI_OFFSET_KIND), optional, intent(in ) :: bufcount
integer, optional, intent(in ) :: buftype
integer function nf90mpi_get_var(ncid, varid, values, start, stride, map, bufcount, buftype)
integer, intent( in) :: ncid
integer, intent( in) :: varid
<any valid type, scalar or array of any rank>, intent(out) :: values
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent( in) :: start
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent( in) :: count
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent( in) :: stride
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent( in) :: map
integer (kind=MPI_OFFSET_KIND), optional, intent( in) :: bufcount
integer, optional, intent( in) :: buftype
integer function nf90mpi_get_var_all(ncid, varid, values, start, stride, map, bufcount, buftype)
integer, intent( in) :: ncid
integer, intent( in) :: varid
<any valid type, scalar or array of any rank>, intent(out) :: values
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent( in) :: start
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent( in) :: count
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent( in) :: stride
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent( in) :: map
integer (kind=MPI_OFFSET_KIND), optional, intent( in) :: bufcount
integer, optional, intent( in) :: buftype
integer function nf90mpi_iput_var(ncid, varid, values, req, start, stride, map, bufcount, buftype)
integer, intent(in ) :: ncid
integer, intent(in ) :: varid
<any valid type, scalar or array of any rank>, intent(inout) :: values
integer, intent( out) :: req
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: start
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: count
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: stride
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: map
integer (kind=MPI_OFFSET_KIND), optional, intent(in ) :: bufcount
integer, optional, intent(in ) :: buftype
integer function nf90mpi_iget_var(ncid, varid, values, req, start, stride, map, bufcount, buftype)
integer, intent( in) :: ncid
integer, intent( in) :: varid
<any valid type, scalar or array of any rank>, intent(out) :: values
integer, intent(out) :: req
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent( in) :: start
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent( in) :: count
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent( in) :: stride
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent( in) :: map
integer (kind=MPI_OFFSET_KIND), optional, intent( in) :: bufcount
integer, optional, intent( in) :: buftype
integer function nf90mpi_bput_var(ncid, varid, values, req, start, stride, map, bufcount, buftype)
integer, intent(in ) :: ncid
integer, intent(in ) :: varid
<any valid type, scalar or array of any rank>, intent(inout) :: values
integer, intent( out) :: req
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: start
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: count
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: stride
integer (kind=MPI_OFFSET_KIND), dimension(:), optional, intent(in ) :: map
integer (kind=MPI_OFFSET_KIND), optional, intent(in ) :: bufcount
integer, optional, intent(in ) :: buftype
integer function nf90mpi_buffer_attach(ncid, bufsize)
integer, intent( in) :: ncid
integer (kind=MPI_OFFSET_KIND), intent( in) :: bufsize
integer function nf90mpi_buffer_detach(ncid)
integer, intent( in) :: ncid
integer function nf90mpi_inq_buffer_size(ncid, buf_size)
integer, intent( in) :: ncid
integer (kind=MPI_OFFSET_KIND), intent(out) :: buf_size
integer function nf90mpi_inq_buffer_usage(ncid, usage)
integer, intent( in) :: ncid
integer (kind=MPI_OFFSET_KIND), intent(out) :: usage
integer function nf90mpi_put_vard(ncid, varid, filetype, buf, bufcount, buftype)
integer, intent(in ) :: ncid
integer, intent(in ) :: varid
integer, intent(in ) :: filetype
<any valid type, scalar or array of any rank>, intent(inout) :: buf
integer (kind=MPI_OFFSET_KIND), intent(in ) :: bufcount
integer, intent(in ) :: buftype
integer function nf90mpi_put_vard_all(ncid, varid, filetype, buf, bufcount, buftype)
integer, intent(in ) :: ncid
integer, intent(in ) :: varid
integer, intent(in ) :: filetype
<any valid type, scalar or array of any rank>, intent(inout) :: buf
integer (kind=MPI_OFFSET_KIND), intent(in ) :: bufcount
integer, intent(in ) :: buftype
integer function nf90mpi_get_vard(ncid, varid, filetype, buf, bufcount, buftype)
integer, intent( in) :: ncid
integer, intent( in) :: varid
integer, intent( in) :: filetype
<any valid type, scalar or array of any rank>, intent(out) :: buf
integer (kind=MPI_OFFSET_KIND), intent( in) :: bufcount
integer, intent( in) :: buftype
integer function nf90mpi_get_vard_all(ncid, varid, filetype, buf, bufcount, buftype)
integer, intent( in) :: ncid
integer, intent( in) :: varid
integer, intent( in) :: filetype
<any valid type, scalar or array of any rank>, intent(out) :: buf
integer (kind=MPI_OFFSET_KIND), intent( in) :: bufcount
integer, intent( in) :: buftype
integer function nf90mpi_put_varn(ncid, varid, values, num, start, count)
integer, intent(in ) :: ncid
integer, intent(in ) :: varid
<any valid type, scalar or array of any rank>, intent(inout) :: values
integer, intent(in ) :: num
integer (kind=MPI_OFFSET_KIND), dimension(:,:), intent(in ) :: start
integer (kind=MPI_OFFSET_KIND), dimension(:,:), optional, intent(in ) :: count
integer function nf90mpi_put_varn_all(ncid, varid, values, num, start, count)
integer, intent(in ) :: ncid
integer, intent(in ) :: varid
<any valid type, scalar or array of any rank>, intent(inout) :: values
integer, intent(in ) :: num
integer (kind=MPI_OFFSET_KIND), dimension(:,:), intent(in ) :: start
integer (kind=MPI_OFFSET_KIND), dimension(:,:), optional, intent(in ) :: count
integer function nf90mpi_get_varn(ncid, varid, values, num, start, count)
integer, intent( in) :: ncid
integer, intent( in) :: varid
<any valid type, scalar or array of any rank>, intent(out) :: values
integer, intent( in) :: num
integer (kind=MPI_OFFSET_KIND), dimension(:,:), intent( in) :: start
integer (kind=MPI_OFFSET_KIND), dimension(:,:), optional, intent( in) :: count
integer function nf90mpi_get_varn_all(ncid, varid, values, num, start, count)
integer, intent( in) :: ncid
integer, intent( in) :: varid
<any valid type, scalar or array of any rank>, intent(out) :: values
integer, intent( in) :: num
integer (kind=MPI_OFFSET_KIND), dimension(:,:), intent( in) :: start
integer (kind=MPI_OFFSET_KIND), dimension(:,:), optional, intent( in) :: count
integer function nf90mpi_iput_varn(ncid, varid, values, req, num, start, count)
integer, intent(in ) :: ncid
integer, intent(in ) :: varid
<any valid type, scalar or array of any rank>, intent(inout) :: values
integer, intent(out ) :: req
integer, intent(in ) :: num
integer (kind=MPI_OFFSET_KIND), dimension(:,:), intent(in ) :: start
integer (kind=MPI_OFFSET_KIND), dimension(:,:), optional, intent(in ) :: count
integer function nf90mpi_iget_varn(ncid, varid, values, req, num, start, count)
integer, intent( in) :: ncid
integer, intent( in) :: varid
<any valid type, scalar or array of any rank>, intent(out) :: values
integer, intent(out) :: req
integer, intent( in) :: num
integer (kind=MPI_OFFSET_KIND), dimension(:,:), intent( in) :: start
integer (kind=MPI_OFFSET_KIND), dimension(:,:), optional, intent( in) :: count
integer function nf90mpi_bput_varn(ncid, varid, values, req, num, start, count)
integer, intent(in ) :: ncid
integer, intent(in ) :: varid
<any valid type, scalar or array of any rank>, intent(inout) :: values
integer, intent( out) :: req
integer, intent(in ) :: num
integer (kind=MPI_OFFSET_KIND), dimension(:,:), intent(in ) :: start
integer (kind=MPI_OFFSET_KIND), dimension(:,:), optional, intent(in ) :: count
Next: Error Codes, Previous: Summary of Fortran 77 Interface, Up: Top [Index]