hvs_vtk_close Subroutine

public subroutine hvs_vtk_close(vtk_file, proc)

This routine finalizes the vtu file i.e closing cellData xml and creating pvtu file to combile all parallel vtu files

Arguments

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

The file descriptor to close again.

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

Communicator for the parallel environment.


Called by

proc~~hvs_vtk_close~2~~CalledByGraph proc~hvs_vtk_close~2 hvs_vtk_close proc~hvs_dump_debug_array hvs_dump_debug_array proc~hvs_dump_debug_array->proc~hvs_vtk_close~2

Contents

Source Code


Source Code

  subroutine hvs_vtk_close( vtk_file, proc )
    !> The file descriptor to close again.
    type(hvs_vtk_file_type), intent(in) :: vtk_file
    !> Communicator for the parallel environment.
    type(tem_comm_env_type), intent(in) :: proc
    ! ----------------------------------------------------------------------!
    character(len=PathLen) :: headerline
    character :: linebreak
    character(len=PathLen) :: filename
    integer :: irank
    integer :: pos
    ! ----------------------------------------------------------------------!

    linebreak = new_line('x')

    if (vtk_file%has_celldata) then
      write(headerline,'(a)') '   </CellData>'
      write(vtk_file%outunit) trim(headerline)//linebreak
      write(vtk_file%outunit) linebreak
      if (vtk_file%write_pvtu) then
        write(headerline,'(a)') '  </PCellData>'
        write(vtk_file%punit) trim(headerline)//linebreak
        write(vtk_file%punit) linebreak
      end if
    end if

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

    close(vtk_file%outunit)

    ! write pvtu file from root process
    if (vtk_file%write_pvtu) then

      pos = INDEX(trim(vtk_file%basename), pathSep, .true.)
      do irank = 0, proc%comm_size-1
        write(filename,'(a,i6.6)') trim(vtk_file%basename(pos+1:)) &
          &                        // '_p', irank
        write(filename,'(a)') trim(filename) // trim(vtk_file%timestamp) &
          &                   // '.vtu'
        write(headerline,'(a)') '<Piece Source="'//trim(filename)//'"/>'
        write(vtk_file%punit) trim(headerline)//linebreak
      end do

      write(headerline,'(a)') '</PUnstructuredGrid>'
      write(vtk_file%punit) trim(headerline)//linebreak
      write(headerline,'(a)') '</VTKFile>'
      write(vtk_file%punit) trim(headerline)//linebreak

      !close the unit
      close(vtk_file%punit)

    end if

  end subroutine hvs_vtk_close