load_datafile Subroutine

private subroutine load_datafile(me)

This subroutine reads the data from disc and stores it in the tem_from_file_temporal_type.

Arguments

Type IntentOptional Attributes Name
type(tem_from_file_temporal_type), intent(inout) :: me

temporal predefined from_file type


Calls

proc~~load_datafile~~CallsGraph proc~load_datafile load_datafile interface~append~24 append proc~load_datafile->interface~append~24 interface~init~24 init proc~load_datafile->interface~init~24 proc~tem_abort tem_abort proc~load_datafile->proc~tem_abort proc~tem_open tem_open proc~load_datafile->proc~tem_open proc~append_arrayga2d_real append_arrayga2d_real interface~append~24->proc~append_arrayga2d_real proc~append_singlega2d_real append_singlega2d_real interface~append~24->proc~append_singlega2d_real proc~init_ga2d_real init_ga2d_real interface~init~24->proc~init_ga2d_real mpi_abort mpi_abort proc~tem_abort->mpi_abort proc~tem_open->proc~tem_abort proc~newunit newunit proc~tem_open->proc~newunit proc~upper_to_lower upper_to_lower proc~tem_open->proc~upper_to_lower interface~expand~22 expand proc~append_arrayga2d_real->interface~expand~22 proc~append_singlega2d_real->interface~expand~22 proc~expand_ga2d_real expand_ga2d_real interface~expand~22->proc~expand_ga2d_real

Called by

proc~~load_datafile~~CalledByGraph proc~load_datafile load_datafile proc~load_temporal_from_file load_temporal_from_file proc~load_temporal_from_file->proc~load_datafile proc~tem_load_temporal tem_load_temporal proc~tem_load_temporal->proc~load_temporal_from_file proc~load_spacetime_predefined load_spacetime_predefined proc~load_spacetime_predefined->proc~tem_load_temporal proc~tem_load_spacetime_single tem_load_spacetime_single proc~tem_load_spacetime_single->proc~load_spacetime_predefined

Source Code

  subroutine load_datafile( me )
    ! --------------------------------------------------------------------------
    !> temporal predefined from_file type
    type(tem_from_file_temporal_type),intent(inout) :: me
    ! --------------------------------------------------------------------------
    ! unit of the datafile
    integer :: datafile_unit
    ! status of IO
    integer :: stat = 0
    ! tmp time value
    real(kind=rk) :: tmp_time
    ! tmp vel value
    real(kind=rk) :: tmp_vel
    ! --------------------------------------------------------------------------

    ! initialize the data array
    call init( me = me%signal, width = 2 )

    ! open the datafile
    call tem_open( newunit = datafile_unit,     &
      &            file    = trim(me%datafile), &
      &            action  = 'read',            &
      &            status  = 'old',             &
      &            form    = 'formatted',       &
      &            access  = 'sequential'       )

    ! as long as there is no io error (end of file is not reached)
    do
      ! ... read the data
      read( unit = datafile_unit, fmt='(2e15.8)', iostat=stat ) tmp_time, &
        &                                                       tmp_vel
      ! ... check wether the file has reached its end to prevent storing
      !     the last item twice
      if ( stat .ne. 0) then
        exit
      end if
      ! ... and append them to the growing array of tuples
      call append( me = me%signal, val = (/ tmp_time, tmp_vel/) )
    end do

    ! multiply the data with the provided factor
    me%signal%val(2,:)=me%fac*me%signal%val(2,:)

    if( stat .gt. 0 ) then
      write(logUnit(1),*)'Error when reading the inlet data from file!'
      write(logUnit(1),*)'STOPPING'
      call tem_abort()
    end if

    ! close the file
    close( unit = datafile_unit )

  end subroutine load_datafile