load_variable_operation Subroutine

public subroutine load_variable_operation(me, conf, parent, iError)

Loading of a variable, that is defined in terms of an operation on other variables.

Arguments

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

The variable to read from the Lua script(conf) and fill

type(flu_State) :: conf

Lua handle connected to the script to read the table from

integer, intent(in) :: parent

A parent table handle in which to look the current operation table

integer, intent(out) :: iError

iError .ne. 0 if operation is not loaded successfully.


Calls

proc~~load_variable_operation~~CallsGraph proc~load_variable_operation load_variable_operation proc~aot_table_open aot_table_open proc~load_variable_operation->proc~aot_table_open proc~aot_table_close aot_table_close proc~load_variable_operation->proc~aot_table_close proc~aot_get_val aot_get_val proc~load_variable_operation->proc~aot_get_val proc~tem_abort tem_abort proc~load_variable_operation->proc~tem_abort proc~tem_reduction_transient_load tem_reduction_transient_load proc~load_variable_operation->proc~tem_reduction_transient_load mpi_abort mpi_abort proc~tem_abort->mpi_abort proc~tem_reduction_transient_load->proc~aot_table_open proc~tem_reduction_transient_load->proc~aot_table_close proc~tem_reduction_transient_load->proc~tem_abort proc~aot_get_val~2 aot_get_val proc~tem_reduction_transient_load->proc~aot_get_val~2

Called by

proc~~load_variable_operation~~CalledByGraph proc~load_variable_operation load_variable_operation proc~tem_variable_load_single tem_variable_load_single proc~tem_variable_load_single->proc~load_variable_operation proc~tem_variable_load_vector tem_variable_load_vector proc~tem_variable_load_vector->proc~tem_variable_load_single interface~tem_variable_load tem_variable_load interface~tem_variable_load->proc~tem_variable_load_single interface~tem_variable_load->proc~tem_variable_load_vector

Contents


Source Code

  subroutine load_variable_operation( me, conf, parent, iError)
    ! --------------------------------------------------------------------------!
    !> The variable to read from the Lua script(conf) and fill
    type(tem_variable_type), intent(inout) :: me

    !> Lua handle connected to the script to read the table from
    type(flu_state) :: conf

    !> A parent table handle in which to look the current operation table
    integer, intent(in) :: parent

    !> iError .ne. 0 if operation is not loaded successfully.
    integer, intent(out) :: iError
    ! --------------------------------------------------------------------------!
    integer :: oper_handle, iIn, nInvar
    integer, allocatable :: vError(:)
    ! --------------------------------------------------------------------------!

    call aot_table_open( L       = conf,          &
      &                  parent  = parent,        &
      &                  thandle = oper_handle,   &
      &                  key     = 'operation'    )
    if (oper_handle == 0) then
      iError = -1
      write(logUnit(1),*) 'operation table not defined.'
      return
    end if
    write(logUnit(5),*) '  Loading operation table '

    ! get operation kind
    call aot_get_val( L       = conf,        &
      &               thandle = oper_handle, &
      &               val     = me%operType, &
      &               ErrCode = iError,      &
      &               key     = 'kind'       )

    if (btest(iError, aoterr_Fatal)) then
      write(logUnit(1),*) 'FATAL Error occured, while retrieving "kind" '//    &
        &                 'from operation table :'
      if ( btest( iError, aotErr_NonExistent ))                                &
        & write(logUnit(1),*)'Variable not existent!'
      if (btest(iError, aoterr_WrongType))                                     &
        & write(logUnit(1),*)'Variable has wrong type!'
    end if
    write(logUnit(5),*) '    kind = '//trim(me%operType)

    ! Get the dependent variables
    call aot_get_val( val       = me%input_varName, &
      &               ErrCode   = vError,           &
      &               maxlength = 100,              &
      &               L         = conf,             &
      &               thandle   = oper_handle,      &
      &               key       = 'input_varname'   )

    if ( any(btest(vError, aoterr_Fatal)) ) then

      write(logUnit(1),*) 'FATAL Error occured, while retrieving '
      write(logUnit(1),*) '"input_varName" for operation' // trim(me%operType)
      call tem_abort()

    end if

    ! number of variables in input_varname
    nInVar = size(me%input_varName)

    write(logUnit(5),*) '    input_varname:'
    do iIn = 1, nInVar
      write(logUnit(5), *) '      ', iIN, trim(me%input_varName(iIn))
    end do

    deallocate(vError)

    ! Load operation specific infos
    select case(trim(me%operType))
    case ('reduction_transient')
      call tem_reduction_transient_load(me     = me%redTransConfig, &
        &                               conf   = conf,              &
        &                               parent = oper_handle        )

    case ('extract')
      allocate(me%input_varIndex(me%nComponents))
      allocate(vError(me%nComponents))
      ! Extraction of component index from input_varname is possible
      ! only when there is only one variable name in input_varname
      if ( nInVar == 1 ) then
        call aot_get_val( L         = conf,               &
          &               thandle   = oper_handle,        &
          &               val       = me%input_varIndex, &
          &               ErrCode   = vError,             &
          &               key       = 'input_varindex'   )

        if ( any(btest(vError, aoterr_Fatal)) .or. &
          &  any(me%input_varIndex <= 0) ) then

          write(logUnit(1),*) 'FATAL Error occured, while retrieving '
          write(logUnit(1),*) '"input_varindex" for extract operation'
          write(logUnit(1),*) 'input_varindex: ',me%input_varIndex
          call tem_abort()

        end if
      else
        write(logUnit(1),*) 'Error loading component index for operation'//&
          &                 ' kind extract'
        write(logUnit(1),*) 'Size of input_varname: ', nInVar
        write(logUnit(1),*) 'Extract component index from input_varname'//&
          &                 'works only when size(input_varname) == 1'
        call tem_abort()
      end if

      write(logUnit(5),*) '     Component index to extract:', me%input_varIndex
    end select
    call aot_table_close( L = conf, thandle = oper_handle)

  end subroutine load_variable_operation