62    class(*), 
target, 
intent(in) :: data
 
   63    real(kind=
rp), 
intent(in), 
optional :: t
 
   64    type(
mesh_t), 
pointer :: msh => null()
 
   65    type(
field_t), 
pointer :: fld => null()
 
   67    type(
dofmap_t), 
pointer :: dm => null()
 
   70    character(len=10) :: id_str
 
   93       write(id_str,
'(i10.10)') 
pe_rank 
   95       open(unit=9, 
file=trim(this%fname(1:suffix_pos-1))//id_str//
'.vtk')
 
   97       open(unit=9, 
file=trim(this%fname))
 
  101    write(9, fmt=
'(A)') 
'# vtk DataFile Version 2.0' 
  102    write(9, fmt=
'(A)') 
'Neko' 
  103    write(9, fmt=
'(A)') 
'ASCII' 
  105    if (
associated(msh)) 
then 
  106       write(9, fmt=
'(A)') 
'DATASET UNSTRUCTURED_GRID' 
  110       if (
associated(mfld)) 
then 
  112       else if (
associated(fld)) 
then 
  115    else if (
associated(dm)) 
then 
  116       write(9, fmt=
'(A)') 
'DATASET POLYDATA' 
  121    else if (
associated(tet_msh)) 
then 
  122       write(9, fmt=
'(A)') 
'DATASET UNSTRUCTURED_GRID' 
  124    else if (
associated(tri_msh)) 
then 
  125       write(9, fmt=
'(A)') 
'DATASET UNSTRUCTURED_GRID' 
 
  136    class(*), 
target, 
intent(inout) :: data
 
  138    call neko_error(
'VTK file read not implemented')
 
 
  144    type(
mesh_t), 
intent(inout) :: msh
 
  145    integer :: i, j, vtk_type
 
  146    integer,  
dimension(8), 
parameter :: vcyc_to_sym = (/1, 2, 4, 3, &
 
  149    write(unit, fmt=
'(A,I8,A)') 
'POINTS', msh%mpts,
' double' 
  151       write(unit, fmt=
'(F15.8,F15.8,F15.8)') 
real(msh%points(i)%x,
dp)
 
  155    write(unit, fmt=
'(A,I8,I8)')  
'CELLS', msh%nelv, msh%nelv*(msh%npts+1)
 
  158       write(unit, fmt=
'(I8,8I8)') msh%npts, &
 
  159            (msh%get_local(msh%elements(i)%e%pts(vcyc_to_sym(j))%p) - 1, &
 
  164    write(unit, fmt=
'(A,I8)') 
'CELL_TYPES', msh%nelv
 
  166    if (msh%gdim .eq. 3) vtk_type = 12
 
  168       write(unit, fmt=
'(I2)') vtk_type
 
 
  179    write(unit, fmt=
'(A,I8)') 
'CELL_DATA', mfld%msh%nelv
 
  180    write(unit, fmt=
'(A,A,A,I8)') 
'SCALARS ', trim(mfld%name), 
' int', 1
 
  181    write(unit, fmt=
'(A)') 
'LOOKUP_TABLE default' 
  183    do i = 1, mfld%msh%nelv
 
  184       write(unit, fmt=
'(I8)') mfld%data(i)
 
 
  194    type(
field_t), 
intent(inout) :: fld
 
  195    real(kind=
dp), 
allocatable :: point_data(:)
 
  196    integer :: i, j, lx, ly, lz, id(8)
 
  198    if ( (fld%Xh%lx - 1 .gt. 1) .or. &
 
  199         (fld%Xh%ly - 1 .gt. 1) .or. &
 
  200         (fld%Xh%lz - 1 .gt. 1)) 
then 
  201       call neko_log%warning(
"Interpolate high-order data onto a low-order mesh")
 
  204    write(unit, fmt=
'(A,I8)') 
'POINT_DATA', fld%msh%mpts
 
  205    write(unit, fmt=
'(A,A,A,I8)') 
'SCALARS ', trim(fld%name), 
' double', 1
 
  206    write(unit, fmt=
'(A)') 
'LOOKUP_TABLE default' 
  211    allocate(point_data(fld%msh%mpts))
 
  213    do i = 1, fld%msh%nelv
 
  214       do j = 1, fld%msh%npts
 
  215          id(j) = fld%msh%get_local(fld%msh%elements(i)%e%pts(j)%p)
 
  218       point_data(id(1)) = 
real(fld%x(1,1,1,i),
dp)
 
  219       point_data(id(2)) = 
real(fld%x(lx,1,1,i),
dp)
 
  220       point_data(id(3)) = 
real(fld%x(1,ly,1,i),
dp)
 
  221       point_data(id(4)) = 
real(fld%x(lx,ly,1,i),
dp)
 
  223       if (fld%msh%gdim .eq. 3) 
then 
  224          point_data(id(5)) = 
real(fld%x(1,1,lz,i),
dp)
 
  225          point_data(id(6)) = 
real(fld%x(lx,1,lz,i),
dp)
 
  226          point_data(id(7)) = 
real(fld%x(1,ly,lz,i),
dp)
 
  227          point_data(id(8)) = 
real(fld%x(lx,ly,lz,i),
dp)
 
  232    write(unit, *) point_data
 
  234    deallocate(point_data)
 
 
  244    write(unit, fmt=
'(A,I8,A)') 
'POINTS', 
size(dm%x),
' double' 
  246    do i = 1, dm%msh%nelv
 
  250                write(unit, fmt=
'(F15.8,F15.8,F15.8)') &
 
  251                     real(dm%x(j,k,l,i),dp),&
 
  252                     real(dm%y(j,k,l,i),dp),&
 
  253                     real(dm%z(j,k,l,i),dp)
 
  259    write(unit, fmt=
'(A,I8,I8)') 
'VERTICES', 
size(dm%x), 2*
size(dm%x)
 
  261       write(unit, fmt=
'(I8,I8)') 1,i-1
 
 
  271    integer :: i, j, k, l
 
  273    write(unit, fmt=
'(A,I8)') 
'POINT_DATA', 
size(dm%dof)
 
  274    write(unit, fmt=
'(A,A,A,I8)') 
'SCALARS ', 
'dof_id', 
' integer', 1
 
  275    write(unit, fmt=
'(A)') 
'LOOKUP_TABLE default' 
  277    do i = 1, dm%msh%nelv
 
  281                write(unit, fmt=
'(I8)') 
real(dm%dof(j,k,l,i),
dp)
 
  287    write(unit, fmt=
'(A,A,A,I8)') 
'SCALARS ', 
'shared_dof', 
' integer', 1
 
  288    write(unit, fmt=
'(A)') 
'LOOKUP_TABLE default' 
  290    do i = 1, dm%msh%nelv
 
  294                if (dm%shared_dof(j,k,l,i)) 
then 
  295                   write(unit, fmt=
'(I8)') 1
 
  297                   write(unit, fmt=
'(I8)') 0
 
 
  310    integer, 
parameter :: npts = 4
 
  311    integer :: i, j, vtk_type
 
  314    write(unit, fmt=
'(A,I8,A)') 
'POINTS', tet_msh%msh%mpts,
' double' 
  315    do i = 1, tet_msh%msh%mpts
 
  316       write(unit, fmt=
'(F15.8,F15.8,F15.8)') 
real(tet_msh%msh%points(i)%x,
dp)
 
  320    write(unit, fmt=
'(A,I8,I8)')  
'CELLS', tet_msh%nelv, tet_msh%nelv*(npts+1)
 
  322    do i = 1, tet_msh%nelv
 
  323       write(unit, fmt=
'(I8,8I8)') npts, &
 
  324            (tet_msh%msh%get_local(tet_msh%el(i)%pts(j)%p) - 1, &
 
  329    write(unit, fmt=
'(A,I8)') 
'CELL_TYPES', tet_msh%nelv
 
  331    do i = 1, tet_msh%nelv
 
  332       write(unit, fmt=
'(I2)') vtk_type
 
 
  341    integer, 
parameter :: npts = 3
 
  342    integer :: i, j, vtk_type
 
  345    write(unit, fmt=
'(A,I8,A)') 
'POINTS', tri_msh%mpts,
' double' 
  346    do i = 1, tri_msh%mpts
 
  347       write(unit, fmt=
'(F15.8,F15.8,F15.8)') 
real(tri_msh%points(i)%x,
dp)
 
  351    write(unit, fmt=
'(A,I8,I8)')  
'CELLS', tri_msh%nelv, tri_msh%nelv*(npts+1)
 
  353    do i = 1, tri_msh%nelv
 
  354       write(unit, fmt=
'(I8,8I8)') npts, &
 
  355            (tri_msh%el(i)%pts(j)%p%id() - 1, j=1, npts)
 
  359    write(unit, fmt=
'(A,I8)') 
'CELL_TYPES', tri_msh%nelv
 
  361    do i = 1, tri_msh%nelv
 
  362       write(unit, fmt=
'(I2)') vtk_type
 
 
integer pe_size
MPI size of communicator.
 
Defines a mapping of the degrees of freedom.
 
Module for file I/O operations.
 
type(log_t), public neko_log
Global log stream.
 
integer, parameter, public dp
 
integer, parameter, public rp
Global precision used in computations.
 
Defines a tetrahedral mesh.
 
Defines a triangular surface mesh.
 
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
 
subroutine vtk_file_write(this, data, t)
Write data in legacy VTK.
 
subroutine vtk_file_write_tet_mesh(unit, tet_msh)
Write a tetrahedral mesh in legacy VTK format.
 
subroutine vtk_file_read(this, data)
 
subroutine vtk_file_write_dofmap_data(unit, dm)
Write a dofmap dm data as point data.
 
subroutine vtk_file_write_tri_mesh(unit, tri_msh)
Write a triangular mesh in legacy VTK format.
 
subroutine vtk_file_write_dofmap_coordinates(unit, dm)
Write xyz-coordinates of a dofmap dm as points.
 
subroutine vtk_file_write_cell_data(unit, mfld)
Write a mesh field mfld as cell data.
 
subroutine vtk_file_write_point_data(unit, fld)
Write a field fld as point data.
 
subroutine vtk_file_write_mesh(unit, msh)
Write a mesh in legacy VTK format.
 
Interface for legacy VTK files.