read_sparse_real Subroutine

public subroutine read_sparse_real(filename, npts, ids_loc, values_loc, nskip)

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: filename
integer, intent(in) :: npts
integer, intent(in) :: ids_loc(:)
real, intent(out), allocatable :: values_loc(:)
integer, intent(in), optional :: nskip

Calls

proc~~read_sparse_real~~CallsGraph proc~read_sparse_real read_sparse_real mpi_bcast mpi_bcast proc~read_sparse_real->mpi_bcast

Called by

proc~~read_sparse_real~~CalledByGraph proc~read_sparse_real read_sparse_real proc~init_vegetation init_vegetation proc~init_vegetation->proc~read_sparse_real program~udales uDALES program~udales->proc~init_vegetation

Source Code

  subroutine read_sparse_real(filename, npts, ids_loc, values_loc, nskip)
    implicit none

    character(len=*), intent(in)      :: filename
    integer, intent(in)               :: npts
    integer, intent(in)               :: ids_loc(:)
    real, allocatable, intent(out)    :: values_loc(:)
    integer, intent(in), optional     :: nskip

    real, allocatable :: values_glob(:)
    integer :: n, nh, ierr
    character(80) :: chmess

    nh = 1
    if (present(nskip)) nh = nskip

    allocate(values_glob(npts))
    values_glob = 0.

    if (myid == 0) then
      open(ifinput, file=filename, status='old', iostat=ierr)
      if (ierr /= 0) then
        write(*, '(A,A)') 'ERROR: Cannot open file: ', trim(filename)
        write(*, '(A,I0)') '  iostat error code: ', ierr
        stop 1
      end if

      do n = 1, nh
        read(ifinput, '(a80)', iostat=ierr) chmess
        if (ierr /= 0) then
          write(*, '(A,I0,A,A)') 'ERROR: Cannot read header line ', n, ' from ', trim(filename)
          stop 1
        end if
      end do

      do n = 1, npts
        read(ifinput, *, iostat=ierr) values_glob(n)
        if (ierr /= 0) then
          write(*, '(A,I0,A,A)') 'ERROR: Cannot read data line ', n, ' from ', trim(filename)
          write(*, '(A,I0)') '  iostat error code: ', ierr
          stop 1
        end if
      end do

      close(ifinput)
    end if

    call MPI_BCAST(values_glob, npts, MY_REAL, 0, comm3d, mpierr)

    allocate(values_loc(size(ids_loc)))
    do n = 1, size(ids_loc)
      values_loc(n) = values_glob(ids_loc(n))
    end do

    deallocate(values_glob)
  end subroutine read_sparse_real