hvs_output_load Subroutine

public subroutine hvs_output_load(me, conf, parent, isReduce)

Read the output configuration from a Lua script.

Arguments

Type IntentOptional Attributes Name
type(hvs_output_config_type), intent(out) :: me

The output configuration settings to fill.

type(flu_State) :: conf

Handle of the Lua script to load the configuration from.

integer, intent(in), optional :: parent

Table handle to the table providing the output settings.

logical, intent(in) :: isReduce

true if reduction is defined


Calls

proc~~hvs_output_load~~CallsGraph proc~hvs_output_load hvs_output_load proc~aot_table_open aot_table_open proc~hvs_output_load->proc~aot_table_open proc~upper_to_lower upper_to_lower proc~hvs_output_load->proc~upper_to_lower proc~hvs_vtk_config_load hvs_vtk_config_load proc~hvs_output_load->proc~hvs_vtk_config_load proc~aot_get_val~2 aot_get_val proc~hvs_output_load->proc~aot_get_val~2 proc~tem_abort tem_abort proc~hvs_output_load->proc~tem_abort proc~aot_table_close aot_table_close proc~hvs_output_load->proc~aot_table_close proc~hvs_vtk_config_load->proc~tem_abort mpi_abort mpi_abort proc~tem_abort->mpi_abort

Called by

proc~~hvs_output_load~~CalledByGraph proc~hvs_output_load hvs_output_load proc~tem_load_trackingconfig tem_load_trackingConfig proc~tem_load_trackingconfig->proc~hvs_output_load proc~tem_load_tracking tem_load_tracking proc~tem_load_tracking->proc~tem_load_trackingconfig

Contents

Source Code


Source Code

  subroutine hvs_output_load(me, conf, parent, isReduce)
    ! --------------------------------------------------------------------------!
    !> The output configuration settings to fill.
    type(hvs_output_config_type), intent(out) :: me

    !> Handle of the Lua script to load the configuration from.
    type(flu_state) :: conf

    !> Table handle to the table providing the output settings.
    integer, intent(in), optional :: parent

    !> true if reduction is defined
    logical, intent(in) :: isReduce
    ! ----------------------------------------------------------------------!
    character(len=labelLen) :: vkind
    integer :: thandle
    integer :: iError
    ! ----------------------------------------------------------------------!
    write(logUnit(3), *) '  Loading output table ...'

    call aot_table_open( L       = conf,    &
      &                  parent  = parent,  &
      &                  thandle = thandle, &
      &                  key     = 'output' )

    if (thandle > 0) then
      ! If reduction is active set output%vis_kind to hvs_asciiTransient
      ! else load vis_kind from output table
      if (isReduce) then
        write(logUnit(7),*) 'Spatial reduction is active, set output format to ascii'
        vkind = 'ascii'
      else
        call aot_get_val( L       = conf,     &
          &               thandle = thandle,  &
          &               key     = 'format', &
          &               val     = vkind,    &
          &               ErrCode = iError    )

        if ( btest(iError, aoterr_Fatal) ) then
          write(logUnit(0),*) 'Fatal Error: In reading format for output'
          if ( btest(iError, aoterr_NonExistent) ) then
            write(logUnit(0),*) 'NonExistent: No format for output provided.'
          end if
          if( btest(iError, aoterr_WrongType) ) then
            write(logUnit(0),*) 'WrongType: No format for output provided.'
          end if
          write(logUnit(0),*) 'STOPPING'
          call tem_abort()
        end if
      end if

      vkind = adjustl(vkind)
      vkind = upper_to_lower(vkind)
      select case(trim(vkind))
      case('vtk')
        ! Output data in paraview unstructured vtk format (.vtu)
        me%vis_kind = hvs_VTK
        call hvs_vtk_config_load( me      = me%vtk, &
          &                       conf    = conf,   &
          &                       thandle = thandle )
      case('ascii')
        me%vis_kind = hvs_AsciiTransient
      case('asciispatial')
        ! write an ascii file for each time step and write all three barycenters
        ! each line has corrdinates and variables of only ONE element
        ! #       coordX  coordY  coordZ    var1       var2        var3
        me%vis_kind = hvs_AsciiSpatial
      case('harvester')
        ! ... in the harvester format, meaning
        ! there is a mesh in treelm format and corresponding to it
        ! there are restart files which hold the state information for the
        ! elements in the mesh for the requested time step
        me%vis_kind = hvs_Internal
      case('precice')
        ! ... for preCICE
        ! there is a mesh in treelm format and corresponding to it
        ! there are restart files which hold the state information for the
        ! elements in the mesh for the requested time step
        me%vis_kind = hvs_PreciceSpatial
      case default
        write(logunit(0),*) 'ERROR in hvs_output_load: '     &
          &                 // 'unknown visualization kind ' &
          &                 // trim(vkind)
        write(logunit(0),*) 'format has to be one of:'
        write(logunit(0),*) '* vtk'
        write(logunit(0),*) '* ascii'
        write(logunit(0),*) '* asciiSpatial'
        write(logunit(0),*) '* harvester'
        write(logunit(0),*) '* precice'
        write(logunit(0),*)
        call tem_abort()
      end select
    else
      write(logUnit(0),*) 'FATAL Error: output table is not defined.'
      write(logunit(0),*) 'Please provide an output table with at least the'
      write(logunit(0),*) 'format defined in it!'
      write(logUnit(0),*) 'STOPPING'
      call tem_abort()
    end if
    write(logUnit(5),*) '  Output format: '//trim(vkind)

    ! To decide whether to use get_point or get_element
    call aot_get_val( L       = conf,           &
      &               thandle = thandle,        &
      &               val     = me%useGetPoint, &
      &               ErrCode = iError,         &
      &               default = .false.,        &
      &               key     = 'use_get_point' )
    write(logUnit(7),*) '  Use get_point: ', me%useGetPoint

    ! Get the number of Dofs to be written in the output
    ! The default is set to -1. If the dofs are not specified,
    ! all the dofs should be dumped
    call aot_get_val( L       = conf,     &
      &               thandle = thandle,  &
      &               val     = me%nDofs, &
      &               ErrCode = iError,   &
      &               default = -1,       &
      &               key     = 'ndofs'   )

    call aot_table_close( L       = conf,   &
      &                   thandle = thandle )

  end subroutine hvs_output_load