hvs_output_open Subroutine

public subroutine hvs_output_open(out_file, use_iter, mesh, varsys, time, subtree)

Open the output for a given mesh.

This writes the mesh data to the output files. Subsequently it can then be enriched by restart information.

Arguments

Type IntentOptional Attributes Name
type(hvs_output_file_type), intent(inout) :: out_file
logical, intent(in) :: use_iter

The output configuration settings to use.

type(treelmesh_type), intent(in) :: mesh

Mesh of the data to visualize.

type(tem_varSys_type), intent(in) :: varsys

Description of the available variable system to get the given varnames from.

type(tem_time_type), intent(in), optional :: time

Time information.

If this is present, the filename will be built with a time stamp and the time point information is written into the vtu file.

type(tem_subTree_type), intent(inout), optional :: subtree

Optional restriction of the elements to output.


Calls

proc~~hvs_output_open~~CallsGraph proc~hvs_output_open hvs_output_open proc~hvs_asciispatial_open hvs_asciiSpatial_open proc~hvs_output_open->proc~hvs_asciispatial_open proc~hvs_vtk_open hvs_vtk_open proc~hvs_output_open->proc~hvs_vtk_open proc~hvs_vtk_write_meshdata hvs_vtk_write_meshdata proc~hvs_output_open->proc~hvs_vtk_write_meshdata proc~hvs_vtk_write_varsys hvs_vtk_write_varSys proc~hvs_output_open->proc~hvs_vtk_write_varsys proc~tem_restart_openwrite tem_restart_openWrite proc~hvs_output_open->proc~tem_restart_openwrite proc~tem_time_reset tem_time_reset proc~hvs_output_open->proc~tem_time_reset proc~tem_updatepropertybits tem_updatePropertyBits proc~hvs_output_open->proc~tem_updatepropertybits proc~getheader getHeader proc~hvs_asciispatial_open->proc~getheader proc~tem_open tem_open proc~hvs_asciispatial_open->proc~tem_open proc~tem_time_sim_stamp tem_time_sim_stamp proc~hvs_asciispatial_open->proc~tem_time_sim_stamp proc~hvs_vtk_open->proc~tem_open proc~tem_time_iter_stamp tem_time_iter_stamp proc~hvs_vtk_open->proc~tem_time_iter_stamp proc~hvs_vtk_open->proc~tem_time_sim_stamp interface~convert_to_base64 convert_to_Base64 proc~hvs_vtk_write_meshdata->interface~convert_to_base64 proc~hvs_vtk_write_cd_header hvs_vtk_write_cd_Header proc~hvs_vtk_write_varsys->proc~hvs_vtk_write_cd_header mpi_allreduce mpi_allreduce proc~tem_restart_openwrite->mpi_allreduce mpi_file_open mpi_file_open proc~tem_restart_openwrite->mpi_file_open mpi_file_set_view mpi_file_set_view proc~tem_restart_openwrite->mpi_file_set_view proc~check_mpi_error check_mpi_error proc~tem_restart_openwrite->proc~check_mpi_error proc~dump_treelmesh dump_treelmesh proc~tem_restart_openwrite->proc~dump_treelmesh proc~tem_create_endiansuffix tem_create_EndianSuffix proc~tem_restart_openwrite->proc~tem_create_endiansuffix proc~tem_dump_subtree tem_dump_subTree proc~tem_restart_openwrite->proc~tem_dump_subtree proc~tem_restart_writeheader tem_restart_writeHeader proc~tem_restart_openwrite->proc~tem_restart_writeheader proc~tem_restart_openwrite->proc~tem_time_sim_stamp mpi_wtime mpi_wtime proc~tem_time_reset->mpi_wtime

Called by

proc~~hvs_output_open~~CalledByGraph proc~hvs_output_open hvs_output_open proc~tem_tracker tem_tracker proc~tem_tracker->proc~hvs_output_open

Source Code

  subroutine hvs_output_open(out_file, use_iter, mesh, varsys, time, subtree )
    ! --------------------------------------------------------------------------!
    type(hvs_output_file_type), intent(inout) :: out_file

    !> The output configuration settings to use.
    ! type(hvs_output_config_type), intent(in) :: out_config
    ! Use iteration in filename?
    logical, intent(in) :: use_iter

    !> Mesh of the data to visualize.
    type(treelmesh_type), intent(in) :: mesh

    !> Description of the available variable system to get the given varnames
    !! from.
    type(tem_varSys_type), intent(in) :: varsys

    !> Optional restriction of the elements to output.
    type(tem_subtree_type), optional, intent(inout) :: subtree

    !> Time information.
    !!
    !! If this is present, the filename will be built with a time stamp and
    !! the time point information is written into the vtu file.
    type(tem_time_type), intent(in), optional :: time
    ! ----------------------------------------------------------------------!
    integer :: nElems
    ! ----------------------------------------------------------------------!
    if (present(subtree)) then
      nElems = subtree%nElems
    else
      nElems = mesh%nElems
    end if

    if (present(time)) then
      out_file%time = time
    else
      call tem_time_reset(out_file%time)
    end if

    select case(out_file%vis_kind)
    case(hvs_AsciiSpatial)
      call hvs_asciiSpatial_open( asciiSpatial = out_file%asciiSpatial, &
        &                         outProc      = out_file%proc,         &
        &                         time         = out_file%time,         &
        &                         varsys       = varsys,                &
        &                         varpos       = out_file%varpos,       &
        &                         nDofs        = out_file%nDofs         )
    case(hvs_Internal)
      ! prepare the header differently for using global tree and tracking
      ! subTree. For harvester format, prepare the header.
      if( present(subTree) ) then
        ! For harvester format, prepare the header
        ! Here we pass the tracking subTree and the global tree to the header
        ! In tem_restart_openWrite the subTree will be written to disk

        ! update the property bits according to the global mesh if this has
        ! changed
        call tem_updatePropertyBits( mesh, subTree )

        call tem_restart_openWrite(me      = out_file%restart,             &
          &                        tree    = mesh,                         &
          &                        timing  = out_file%time,                &
          &                        varSys  = varSys,                     &
          &                        subTree = subTree,                      &
          &                        label   = trim(out_file%basename)//'_'  )
      else
        call tem_restart_openWrite( me     = out_file%restart, &
          &                         tree   = mesh,             &
          &                         timing = out_file%time,    &
          &                         varSys = varSys          )
      end if ! present subTree
    case(hvs_VTK)
      call hvs_vtk_open( vtk_file = out_file%vtk,   &
        &                use_iter = use_iter,       &
        &                proc     = out_file%proc,  &
        &                time     = out_file%time   )

      call hvs_vtk_write_meshdata( vtk_file = out_file%vtk,  &
        &                          vrtx     = out_file%vrtx, &
        &                          nElems   = nElems         )

      call hvs_vtk_write_varSys( vtk_file = out_file%vtk,   &
        &                        varsys   = varsys,         &
        &                        varpos   = out_file%varpos )
    end select

  end subroutine hvs_output_open