hvs_vtk_open Subroutine

public subroutine hvs_vtk_open(vtk_file, use_iter, proc, time)

Open the output files in VTK format.

This will open VTU files and if multiple processes are used a PVTU file. We always write unstructured meshes, so we also write the header for the unstructured mesh here already. The actual mesh data is then to be written by hvs_vtk_write_meshdata.

Arguments

Type IntentOptional Attributes Name
type(hvs_vtk_file_type), intent(inout) :: vtk_file

The file description to open.

logical, intent(in) :: use_iter

User specified settings for the output Whether to use iteration as part of filename

type(tem_comm_env_type), intent(in) :: proc

Parallel environment to use for the output.

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.


Calls

proc~~hvs_vtk_open~~CallsGraph proc~hvs_vtk_open hvs_vtk_open proc~tem_open tem_open 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~tem_time_sim_stamp tem_time_sim_stamp proc~hvs_vtk_open->proc~tem_time_sim_stamp proc~newunit newunit proc~tem_open->proc~newunit proc~tem_abort tem_abort proc~tem_open->proc~tem_abort proc~upper_to_lower upper_to_lower proc~tem_open->proc~upper_to_lower mpi_abort mpi_abort proc~tem_abort->mpi_abort

Called by

proc~~hvs_vtk_open~~CalledByGraph proc~hvs_vtk_open hvs_vtk_open proc~hvs_dump_debug_array hvs_dump_debug_array proc~hvs_dump_debug_array->proc~hvs_vtk_open proc~hvs_output_open hvs_output_open proc~hvs_output_open->proc~hvs_vtk_open proc~tem_tracker tem_tracker proc~tem_tracker->proc~hvs_output_open

Source Code

  subroutine hvs_vtk_open(vtk_file, use_iter, proc, time)
    !> The file description to open.
    type(hvs_vtk_file_type), intent(inout) :: vtk_file

    !> User specified settings for the output
    ! type(hvs_vtk_config_type), intent(in) :: vtk_config
    !> Whether to use iteration as part of filename
    logical, intent(in) :: use_iter

    !> Parallel environment to use for  the output.
    type(tem_comm_env_type), intent(in) :: proc

    !> 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
    ! ----------------------------------------------------------------------!
    character(len=PathLen) :: filename
    character(len=PathLen) :: headerline
    character :: linebreak
    integer :: pos
    character(len=labelLen) :: byte_order
    ! ----------------------------------------------------------------------!

    if ( isLittleEndian ) then
      byte_order = 'LittleEndian'
    else
      byte_order = 'BigEndian'
    end if

    if (proc%comm_size > 1) then
      write(filename,'(a,i6.6)') trim(vtk_file%basename) // '_p', proc%rank
    else
      write(filename,'(a)') trim(vtk_file%basename)
    end if

    vtk_file%timestamp = ''
    if (present(time)) then
      if ( use_iter ) then
        write(vtk_file%timestamp, '(a)') &
          & '_t' // trim(tem_time_iter_stamp(time))
      else
        write(vtk_file%timestamp, '(a)') &
          & '_t' // trim(tem_time_sim_stamp(time))
      end if
    end if

    write(filename,'(a)') trim(filename) // trim(vtk_file%timestamp) // '.vtu'

    if (.not.vtk_file%write_pvtu) then
      write(logunit(3),*) 'Opening VTU file: ' // trim(filename)
      vtk_file%last_opened_file = trim(filename)
    end if

    ! Open the file (always unformatted stream, for ascii output numbers will
    ! be converted to strings before writing).
    call tem_open( newunit = vtk_file%outunit, &
      &            file    = trim(filename),   &
      &            action  = 'write',          &
      &            status  = 'replace',        &
      &            form    = 'unformatted',    &
      &            access  = 'stream'          )

    linebreak = new_line('x')

    ! Writing the header to the VTU file:
    write(headerline,'(a)') '<?xml version="1.0"?>'
    write(vtk_file%outunit) trim(headerline)//linebreak

    write(headerline,'(a)') '<VTKFile type="UnstructuredGrid" version="0.1"' &
      &                     //  ' byte_order="'&
      &                     //trim(byte_order)//'">'
    write(vtk_file%outunit) trim(headerline)//linebreak

    write(headerline,'(a)') ' <UnstructuredGrid>'
    write(vtk_file%outunit) trim(headerline)//linebreak

    ! Add point in time information to the VTU file.
    if (present(time)) then
      write(headerline,'(a)') '<FieldData>'
      write(vtk_file%outunit) trim(headerline)//linebreak

      write(headerline,'(a)') '<DataArray type="Float64" Name="TIME" '// &
        &                     'NumberOfTuples="1" format="ascii">'
      write(vtk_file%outunit) trim(headerline)//linebreak

      write(headerline,*) time%sim
      write(vtk_file%outunit) trim(headerline)//linebreak

      write(headerline,'(a)') '</DataArray>'
      write(vtk_file%outunit) trim(headerline)//linebreak

      write(headerline,'(a)') '<DataArray type="Int32" Name="CYCLE" '// &
        &                     'NumberOfTuples="1" format="ascii">'
      write(vtk_file%outunit) trim(headerline)//linebreak

      write(headerline,*) time%iter
      write(vtk_file%outunit) trim(headerline)//linebreak

      write(headerline,'(a)') '</DataArray>'
      write(vtk_file%outunit) trim(headerline)//linebreak

      write(headerline,'(a)') '</FieldData>'
      write(vtk_file%outunit) trim(headerline)//linebreak
    end if

    ! Open the PVTU file if necessary
    if (vtk_file%write_pvtu) then
      write(filename,'(a)') trim(vtk_file%basename)         &
        &                 //trim(vtk_file%timestamp)//'.pvtu'
      write(logunit(3),*) 'Opening PVTU file: ' // trim(filename)
      vtk_file%last_opened_file = trim(filename)

      call tem_open( newunit = vtk_file%punit, &
        &            file    = trim(filename), &
        &            action  = 'write',        &
        &            status  = 'replace',      &
        &            form    = 'unformatted',  &
        &            access  = 'stream'        )

      write(headerline,'(a)') '<?xml version="1.0"?>'
      write(vtk_file%punit) trim(headerline)//linebreak

      write(headerline,'(a)') '<VTKFile type="PUnstructuredGrid" ' &
        &                     // 'version="0.1" byte_order="'&
        &                     //trim(byte_order)//'">'
      write(vtk_file%punit) trim(headerline)//linebreak

      write(headerline,'(a)') ' <PUnstructuredGrid GhostLevel="0">'
      write(vtk_file%punit) trim(headerline)//linebreak
      write(vtk_file%punit) linebreak
    end if

    ! Append current filename in pvd file.
    ! if pvtu is present, write pvtu filename else vtu filename
    if ( vtk_file%write_pvd ) then
      pos = INDEX(trim(filename), pathSep, .true.)
      write(headerline,'(a)') '  <DataSet timestep="' &
        &                     //trim(tem_time_sim_stamp(time))//'" file="' &
        &                     //trim(filename(pos+1:))//'"/>'
      write(vtk_file%pvdunit) trim(headerline)//linebreak
      flush(vtk_file%pvdunit)
    end if
  end subroutine hvs_vtk_open