load_temporal_from_file Subroutine

private subroutine load_temporal_from_file(me, conf, thandle)

This subroutine loads the information needed to read data from a file.

Example: Read time periodic velocity signals from file 'data.dat' as inlet condition without interpolation between the velocity signals.

boundary_condition =
        {
          { label = 'inlet',
            kind = 'inlet_ubb',
            velocityX = { kind = 'combined',
                          temporal = { predefined = 'datafile',
                                        filename   = 'data.dat',
                                        intp       = 'none',
                                        periodic   = true
                                      }
                        },
            velocityY = 0.0
            velocityZ = 0.0
          }
        }

The options for interpolating the values are \li \a 'none': No interpolation is used. The signals will be chosen by: [t1,t2[ = v1, [t2,t3[ = v2, ... \li \a 'linear': A linear interpolation is done between neighboring values according to the following equation: \f[ \frac{v_{i+1}-v_{i}}{t_{i+1}-t_{i}} \cdot (t_{sim} \% t_{period} - t_{i}) + v_{i} \f]

Arguments

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

temporal predefined from_file type

type(flu_State) :: conf

lua state type

integer, intent(in) :: thandle

aotus parent handle


Calls

proc~~load_temporal_from_file~~CallsGraph proc~load_temporal_from_file load_temporal_from_file aot_get_val aot_get_val proc~load_temporal_from_file->aot_get_val aot_table_open aot_table_open proc~load_temporal_from_file->aot_table_open proc~load_datafile load_datafile proc~load_temporal_from_file->proc~load_datafile proc~tem_abort tem_abort proc~load_temporal_from_file->proc~tem_abort proc~load_datafile->proc~tem_abort interface~append~24 append proc~load_datafile->interface~append~24 interface~init~24 init proc~load_datafile->interface~init~24 proc~tem_open tem_open proc~load_datafile->proc~tem_open mpi_abort mpi_abort proc~tem_abort->mpi_abort 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 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

Called by

proc~~load_temporal_from_file~~CalledByGraph proc~load_temporal_from_file load_temporal_from_file 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 proc~tem_load_spacetime_single->proc~tem_load_spacetime_single interface~tem_load_spacetime tem_load_spacetime interface~tem_load_spacetime->proc~tem_load_spacetime_single proc~tem_load_spacetime_table tem_load_spacetime_table proc~tem_load_spacetime_table->proc~tem_load_spacetime_single

Source Code

  subroutine load_temporal_from_file( me, conf, thandle )
    ! ---------------------------------------------------------------------------
    !> temporal predefined from_file type
    type(tem_from_file_temporal_type),intent(inout) :: me
    !> lua state type
    type(flu_State) :: conf
    !> aotus parent handle
    integer, intent(in) :: thandle
    ! ---------------------------------------------------------------------------
    ! local variables
    integer :: iError
    ! aotus handle for ramping table
    integer :: rampHandle
    ! ---------------------------------------------------------------------------

    ! read the filename of the datafile
    call aot_get_val( L       = conf,                                          &
      &               thandle = thandle,                                       &
      &               key     = 'filename',                                    &
      &               val     = me%datafile,                                   &
      &               ErrCode = iError )
    if (btest(iError, aoterr_Fatal)) then
      write(logUnit(1),*)'FATAL Error occured, while retrieving datafile:'
      if (btest(iError, aoterr_NonExistent))                                   &
        & write(logUnit(1),*)'Variable not existent!'
      if (btest(iError, aoterr_WrongType))                                     &
        & write(logUnit(1),*)'Variable has wrong type!'
      write(logUnit(1),*)'STOPPING'
      call tem_abort()
    end if

    ! read the interpolation kind for the data
    call aot_get_val( L       = conf,                                          &
      &               thandle = thandle,                                       &
      &               key     = 'intp',                                        &
      &               val     = me%intp,                                       &
      &               ErrCode = iError,                                        &
      &               default = 'none' )
    if (btest(iError, aoterr_Fatal)) then
      write(logUnit(1),*)'FATAL Error occured, while retrieving interpolation:'
      if (btest(iError, aoterr_NonExistent))                                   &
        & write(logUnit(1),*)'Variable not existent!'
      if (btest(iError, aoterr_WrongType))                                     &
        & write(logUnit(1),*)'Variable has wrong type!'
      write(logUnit(1),*)'STOPPING'
      call tem_abort()
    end if

    ! open the ramping table
    call aot_table_open( L       = conf,                                       &
      &                  thandle = rampHandle,                                 &
      &                  parent  = thandle,                                    &
      &                  key     = 'ramping')
    if (rampHandle /= 0) then
      me%ramp = .true.
      ! read the ramping time
      call aot_get_val( L       = conf,                                        &
        &               thandle = rampHandle,                                  &
        &               key     = 'rampVal',                                   &
        &               val     = me%rampVal,                                  &
        &               ErrCode = iError)
      if (btest(iError, aoterr_Fatal)) then
        write(logUnit(1),*)'FATAL Error occured, while retrieving ramping '//  &
          &                'value:'
        if (btest(iError, aoterr_NonExistent))                                 &
          & write(logUnit(1),*)'Variable not existent!'
        if (btest(iError, aoterr_WrongType))                                   &
          & write(logUnit(1),*)'Variable has wrong type!'
        write(logUnit(1),*)'STOPPING'
        call tem_abort()
      end if

      ! read the ramping time
      call aot_get_val( L       = conf,                                        &
        &               thandle = rampHandle,                                  &
        &               key     = 'rampT',                                     &
        &               val     = me%rampT,                                    &
        &               ErrCode = iError)
      if (btest(iError, aoterr_Fatal)) then
        write(logUnit(1),*)'FATAL Error occured, while retrieving ramping time:'
        if (btest(iError, aoterr_NonExistent))                                 &
          & write(logUnit(1),*)'Variable not existent!'
        if (btest(iError, aoterr_WrongType))                                   &
          & write(logUnit(1),*)'Variable has wrong type!'
        write(logUnit(1),*)'STOPPING'
        call tem_abort()
      end if
    end if

    ! read the factor
    call aot_get_val( L       = conf,                                          &
      &               thandle = thandle,                                       &
      &               key     = 'fac',                                         &
      &               val     = me%fac,                                        &
      &               ErrCode = iError,                                        &
      &               default = 1.0_rk )
    if (btest(iError, aoterr_Fatal)) then
      write(logUnit(1),*)'FATAL Error occured, while retrieving factor: '
      if (btest(iError, aoterr_NonExistent))                                   &
        & write(logUnit(1),*)'Variable not existent!'
      if (btest(iError, aoterr_WrongType))                                     &
        & write(logUnit(1),*)'Variable has wrong type!'
      write(logUnit(1),*)'STOPPING'
      call tem_abort()
    end if

    ! read in wether the data should be treated periodic
    call aot_get_val( L       = conf,                                          &
      &               thandle = thandle,                                       &
      &               key     = 'periodic',                                    &
      &               val     = me%periodic,                                   &
      &               ErrCode = iError,                                        &
      &               default = .false. )
    if (btest(iError, aoterr_Fatal)) then
      write(logUnit(1),*)'FATAL Error occured, while retrieving periodic:'
      if (btest(iError, aoterr_NonExistent))                                   &
        & write(logUnit(1),*)'Variable not existent!'
      if (btest(iError, aoterr_WrongType))                                     &
        & write(logUnit(1),*)'Variable has wrong type!'
      write(logUnit(1),*)'STOPPING'
      call tem_abort()
    end if

    write(logUnit(3),*)"Using data from file: "//trim(me%datafile)
    write(logUnit(3),*)"  interpolation:      "//trim(me%intp)
    if (rampHandle /= 0) then
      write(logUnit(3),*)"  ramping time:       ",me%rampT
      write(logUnit(3),*)"  ramping value:      ",me%rampVal
    end if
    write(logUnit(3),*)"  periodic:           ",me%periodic

    ! load the data from file
    call load_datafile( me )

  end subroutine load_temporal_from_file