hvs_vtk_dump_data Subroutine

public subroutine hvs_vtk_dump_data(vtk_file, varpos, varsys, mesh, time, subtree)

Dump the given data (input) with the given name in the given format (vtu) to the given unit.

Arguments

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

VTK file to write data to.

integer, intent(in) :: varpos

Position of the variable to write

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

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

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

Mesh to write the data on.

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

Point in time to use for this data.

Can be important for space-time function evaluations.

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

Optional restriction of the elements to output.


Calls

proc~~hvs_vtk_dump_data~2~~CallsGraph proc~hvs_vtk_dump_data~2 hvs_vtk_dump_data interface~convert_to_base64 convert_to_Base64 proc~hvs_vtk_dump_data~2->interface~convert_to_base64 proc~int8_to_base64 int8_to_base64 interface~convert_to_base64->proc~int8_to_base64 proc~int64_to_base64 int64_to_base64 interface~convert_to_base64->proc~int64_to_base64 proc~real64_to_base64 real64_to_base64 interface~convert_to_base64->proc~real64_to_base64 proc~char_to_base64 char_to_base64 interface~convert_to_base64->proc~char_to_base64 proc~real32_to_base64 real32_to_base64 interface~convert_to_base64->proc~real32_to_base64 proc~int32_to_base64 int32_to_base64 interface~convert_to_base64->proc~int32_to_base64 interface~convert_to_base64_single convert_to_Base64_single proc~int8_to_base64->interface~convert_to_base64_single interface~encodeindex EncodeIndex proc~int8_to_base64->interface~encodeindex interface~encodeblock EncodeBlock proc~int8_to_base64->interface~encodeblock proc~int64_to_base64->interface~convert_to_base64_single proc~int64_to_base64->interface~encodeindex proc~int64_to_base64->interface~encodeblock proc~real64_to_base64->interface~convert_to_base64_single proc~real64_to_base64->interface~encodeindex proc~real64_to_base64->interface~encodeblock proc~char_to_base64->interface~convert_to_base64_single proc~char_to_base64->interface~encodeindex proc~char_to_base64->interface~encodeblock proc~real32_to_base64->interface~convert_to_base64_single proc~real32_to_base64->interface~encodeindex proc~real32_to_base64->interface~encodeblock proc~int32_to_base64->interface~convert_to_base64_single proc~int32_to_base64->interface~encodeindex proc~int32_to_base64->interface~encodeblock proc~real64_to_base64_single real64_to_base64_single interface~convert_to_base64_single->proc~real64_to_base64_single proc~int8_to_base64_single int8_to_base64_single interface~convert_to_base64_single->proc~int8_to_base64_single proc~int64_to_base64_single int64_to_base64_single interface~convert_to_base64_single->proc~int64_to_base64_single proc~real32_to_base64_single real32_to_base64_single interface~convert_to_base64_single->proc~real32_to_base64_single proc~char_to_base64_single char_to_base64_single interface~convert_to_base64_single->proc~char_to_base64_single proc~int32_to_base64_single int32_to_base64_single interface~convert_to_base64_single->proc~int32_to_base64_single

Contents

Source Code


Source Code

  subroutine hvs_vtk_dump_data( vtk_file, varpos, varSys, mesh, time, subtree )
    ! --------------------------------------------------------------------------!
    !> VTK file to write data to.
    type(hvs_vtk_file_type), intent(in) :: vtk_file

    !> Position of the variable to write
    integer, intent(in) :: varpos

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

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

    !> Point in time to use for this data.
    !!
    !! Can be important for space-time function evaluations.
    type(tem_time_type), intent(in) :: time

    !> Optional restriction of the elements to output.
    type(tem_subtree_type), optional, intent(in) :: subtree
    ! --------------------------------------------------------------------------!
    integer :: nElems
    integer :: iElem
    integer :: nComp
    integer :: epos(1)
    character(len=PathLen) :: headerline
    character :: linebreak
    character(len=LabelLen) :: dataname
    real(kind=c_double), target, allocatable :: indata(:)
    integer :: spatialComp
    logical :: use_ascii
    ! --------------------------------------------------------------------------!

    nComp = varsys%method%val(varPos)%nComponents
    dataname = varsys%varname%val(varPos)
    if (present(subtree)) then
      nElems = subtree%nElems
    else
      nElems = mesh%nElems
    end if

    linebreak = new_line('x')
    use_ascii = (trim(vtk_file%dataform) == 'ascii')

    spatialComp = nComp

    ! Promote 2D data to 3D, as we always create 3D meshes.
    ! 3rd component will be set to 0.
    if (nComp == 2) then
      spatialComp = 3
    end if

    ! Open the data array.
    write(headerline,'(a,i0,a)') '    <DataArray type="Float64" Name="'        &
      &                          // trim(dataname)                             &
      &                          // '" NumberOfComponents="', spatialComp,     &
      &                          '" format="' // trim(vtk_file%dataform) // '">'
    write(vtk_file%outunit) trim(headerline)//linebreak

    if (vtk_file%write_pvtu) then
      write(headerline,'(a,i0,a)') '    <PDataArray type="Float64" Name="'   &
        &                          // trim(dataname)                         &
        &                          // '" NumberOfComponents="', spatialComp, &
        &                          '" format="' // trim(vtk_file%dataform)   &
        &                          //'"/>'
      write(vtk_file%punit) trim(headerline)//linebreak
    end if

    ! Dump the actual values into the data array.
    if (use_ascii) then

      allocate(indata(spatialComp))
      indata = 0.0_rk
      do iElem = 1, nElems
        if (present(subtree)) then
          epos = subtree%map2global(iElem)
        else
          epos = iElem
        end if

        call varsys%method%val(varpos)                    &
          &        %get_element( varsys = varsys,         &
          &                      elempos = epos,          &
          &                      time    = time,          &
          &                      tree    = mesh,          &
          &                      nElems  = 1,             &
          &                      nDofs   = 1,             &
          &                      res     = indata(:nComp) )
        write(headerline,*) indata
      end do

    else

      allocate(indata(nElems*spatialComp))
      if (nComp == 2) then
        do iElem = 1, nElems
          if (present(subtree)) then
            epos = subtree%map2global(iElem)
          else
            epos = iElem
          end if
          call varsys%method%val(varpos)                          &
            &        %get_element( varsys  = varsys,              &
            &                      elempos = epos,                &
            &                      time    = time,                &
            &                      tree    = mesh,                &
            &                      nElems  = 1,                   &
            &                      nDofs   = 1,                   &
            &                      res     = indata(1+(iElem-1)*3 &
            &                                       :iElem*3-1)   )
          indata(iElem*3) = 0.0_c_double
        end do
      else
        do iElem = 1, nElems
          if (present(subtree)) then
            epos = subtree%map2global(iElem)
          else
            epos = iElem
          end if
          call varsys%method%val(varpos)                              &
            &        %get_element( varsys  = varsys,                  &
            &                      elempos = epos,                    &
            &                      time    = time,                    &
            &                      tree    = mesh,                    &
            &                      nElems  = 1,                       &
            &                      nDofs   = 1,                       &
            &                      res     = indata(1+(iElem-1)*nComp &
            &                                       :iElem*nComp)     )
        end do
      end if
      call convert_to_base64( indata, nElems*spatialComp, vtk_file%outunit )
      write(vtk_file%outunit) linebreak
      deallocate(indata)

    end if

    ! Close the data array again.
    write(headerline,'(a)') '    </DataArray>'
    write(vtk_file%outunit) trim(headerline)//linebreak

  end subroutine hvs_vtk_dump_data