hvs_output_write Subroutine

public subroutine hvs_output_write(out_file, varsys, mesh, subtree)

Arguments

Type IntentOptional Attributes Name
type(hvs_output_file_type), intent(inout) :: out_file
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_subTree_type), intent(in), optional :: subtree

Optional restriction of the elements to output.


Source Code

  subroutine hvs_output_write(out_file, varsys, mesh, subtree)
    ! --------------------------------------------------------------------------!
    type(hvs_output_file_type), intent(inout) :: out_file

    !> 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

    !> Optional restriction of the elements to output.
    type(tem_subtree_type), optional, intent(in) :: subtree
    ! --------------------------------------------------------------------------!
    integer :: iVar
    ! --------------------------------------------------------------------------!

    select case(out_file%vis_kind)
    case(hvs_AsciiTransient)
      if (out_file%useGetPoint) then
        ! Call the get point routine and dump the data at the points
        call hvs_ascii_dump_point_data(                                      &
        &                       ascii           = out_file%ascii,            &
        &                       outProc         = out_file%proc,             &
        &                       varPos          = out_file%varPos,           &
        &                       varSys          = varSys,                    &
        &                       mesh            = mesh,                      &
        &                       time            = out_file%time,             &
        &                       subTree         = subtree )
        ! &                       transientReduce = out_file%transientReduce   )
      else
        ! call the get_element routine and dump all dofs
        call hvs_ascii_dump_elem_data(                                        &
          &                     ascii           = out_file%ascii,             &
          &                     outProc         = out_file%proc,              &
          &                     varPos          = out_file%varPos,            &
          &                     varSys          = varSys,                     &
          &                     mesh            = mesh,                       &
          &                     time            = out_file%time,              &
          &                     nDofs           = out_file%nDofs,             &
          &                     subTree         = subTree )
          ! &                     transientReduce = out_file%transientReduce    )

      endif
    case(hvs_AsciiSpatial)
      if (out_file%useGetPoint) then
        call hvs_asciiSpatial_dump_point_data(                            &
          &                   asciiSpatial    = out_file%asciiSpatial,    &
          &                   varPos          = out_file%varPos,          &
          &                   varSys          = varSys,                   &
          &                   mesh            = mesh,                     &
          &                   time            = out_file%time,            &
          &                   subTree         = subTree,                  &
          ! &                   transientReduce = out_file%transientReduce, &
          &                   bary            = out_file%bary             )
      else
        call hvs_asciiSpatial_dump_elem_data(                             &
          &                   asciiSpatial    = out_file%asciiSpatial,    &
          &                   varPos          = out_file%varPos,          &
          &                   varSys          = varSys,                   &
          &                   mesh            = mesh,                     &
          &                   time            = out_file%time,            &
          &                   nDofs           = out_file%nDofs,           &
          &                   subTree         = subTree,                  &
          ! &                   transientReduce = out_file%transientReduce, &
          &                   bary            = out_file%bary             )
      end if
    case(hvs_Internal)
       call tem_restart_dump_data( restart = out_file%restart,        &
         &                         varSys  = varSys,                  &
         &                         tree    = mesh,                    &
         &                         time    = out_file%time,           &
         &                         subTree = subTree )
         ! &                         transientReduce = out_file%transientReduce )
    case(hvs_VTK)
      do iVar=1,out_file%nVars
        call hvs_vtk_dump_data( vtk_file = out_file%vtk,          &
          &                     varpos   = out_file%varpos(iVar), &
          &                     varSys   = varSys,                &
          &                     mesh     = mesh,                  &
          &                     subtree  = subtree,               &
          &                     time     = out_file%time          )
      end do
    end select

  end subroutine hvs_output_write