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
71 integer:: suffix_pos, file_unit
72 character(len=1024) :: fname
93 fname = trim(this%get_base_fname())
95 write(id_str,
'(i10.10)')
pe_rank
97 open(newunit=file_unit,
file=trim(fname(1:suffix_pos-1))//id_str//
'.vtk')
99 open(newunit=file_unit,
file=trim(fname))
103 write(file_unit, fmt=
'(A)')
'# vtk DataFile Version 2.0'
104 write(file_unit, fmt=
'(A)')
'Neko'
105 write(file_unit, fmt=
'(A)')
'ASCII'
107 if (
associated(msh))
then
108 write(file_unit, fmt=
'(A)')
'DATASET UNSTRUCTURED_GRID'
112 if (
associated(mfld))
then
114 else if (
associated(fld))
then
117 else if (
associated(dm))
then
118 write(file_unit, fmt=
'(A)')
'DATASET POLYDATA'
123 else if (
associated(tet_msh))
then
124 write(file_unit, fmt=
'(A)')
'DATASET UNSTRUCTURED_GRID'
126 else if (
associated(tri_msh))
then
127 write(file_unit, fmt=
'(A)')
'DATASET UNSTRUCTURED_GRID'
138 class(*),
target,
intent(inout) :: data
140 call neko_error(
'VTK file read not implemented')
146 type(
mesh_t),
intent(inout) :: msh
147 integer :: i, j, vtk_type
148 integer,
dimension(8),
parameter :: vcyc_to_sym = (/1, 2, 4, 3, &
151 write(unit, fmt=
'(A,I8,A)')
'POINTS', msh%mpts,
' double'
153 write(unit, fmt=
'(F15.8,F15.8,F15.8)')
real(msh%points(i)%x,
dp)
157 write(unit, fmt=
'(A,I8,I8)')
'CELLS', msh%nelv, msh%nelv*(msh%npts+1)
160 write(unit, fmt=
'(I8,8I8)') msh%npts, &
161 (msh%get_local(msh%elements(i)%e%pts(vcyc_to_sym(j))%p) - 1, &
166 write(unit, fmt=
'(A,I8)')
'CELL_TYPES', msh%nelv
168 if (msh%gdim .eq. 3) vtk_type = 12
170 write(unit, fmt=
'(I2)') vtk_type
181 write(unit, fmt=
'(A,I8)')
'CELL_DATA', mfld%msh%nelv
182 write(unit, fmt=
'(A,A,A,I8)')
'SCALARS ', trim(mfld%name),
' int', 1
183 write(unit, fmt=
'(A)')
'LOOKUP_TABLE default'
185 do i = 1, mfld%msh%nelv
186 write(unit, fmt=
'(I8)') mfld%data(i)
196 type(
field_t),
intent(inout) :: fld
197 real(kind=
dp),
allocatable :: point_data(:)
198 integer :: i, j, lx, ly, lz, id(8)
200 if ( (fld%Xh%lx - 1 .gt. 1) .or. &
201 (fld%Xh%ly - 1 .gt. 1) .or. &
202 (fld%Xh%lz - 1 .gt. 1))
then
203 call neko_log%warning(
"Interpolate high-order data onto a low-order mesh")
206 write(unit, fmt=
'(A,I8)')
'POINT_DATA', fld%msh%mpts
207 write(unit, fmt=
'(A,A,A,I8)')
'SCALARS ', trim(fld%name),
' double', 1
208 write(unit, fmt=
'(A)')
'LOOKUP_TABLE default'
213 allocate(point_data(fld%msh%mpts))
215 do i = 1, fld%msh%nelv
216 do j = 1, fld%msh%npts
217 id(j) = fld%msh%get_local(fld%msh%elements(i)%e%pts(j)%p)
220 point_data(id(1)) =
real(fld%x(1,1,1,i),
dp)
221 point_data(id(2)) =
real(fld%x(lx,1,1,i),
dp)
222 point_data(id(3)) =
real(fld%x(1,ly,1,i),
dp)
223 point_data(id(4)) =
real(fld%x(lx,ly,1,i),
dp)
225 if (fld%msh%gdim .eq. 3)
then
226 point_data(id(5)) =
real(fld%x(1,1,lz,i),
dp)
227 point_data(id(6)) =
real(fld%x(lx,1,lz,i),
dp)
228 point_data(id(7)) =
real(fld%x(1,ly,lz,i),
dp)
229 point_data(id(8)) =
real(fld%x(lx,ly,lz,i),
dp)
234 write(unit, *) point_data
236 deallocate(point_data)
246 write(unit, fmt=
'(A,I8,A)')
'POINTS',
size(dm%x),
' double'
248 do i = 1, dm%msh%nelv
252 write(unit, fmt=
'(F15.8,F15.8,F15.8)') &
253 real(dm%x(j,k,l,i),dp),&
254 real(dm%y(j,k,l,i),dp),&
255 real(dm%z(j,k,l,i),dp)
261 write(unit, fmt=
'(A,I8,I8)')
'VERTICES',
size(dm%x), 2*
size(dm%x)
263 write(unit, fmt=
'(I8,I8)') 1,i-1
273 integer :: i, j, k, l
275 write(unit, fmt=
'(A,I8)')
'POINT_DATA',
size(dm%dof)
276 write(unit, fmt=
'(A,A,A,I8)')
'SCALARS ',
'dof_id',
' integer', 1
277 write(unit, fmt=
'(A)')
'LOOKUP_TABLE default'
279 do i = 1, dm%msh%nelv
283 write(unit, fmt=
'(I8)')
real(dm%dof(j,k,l,i),
dp)
289 write(unit, fmt=
'(A,A,A,I8)')
'SCALARS ',
'shared_dof',
' integer', 1
290 write(unit, fmt=
'(A)')
'LOOKUP_TABLE default'
292 do i = 1, dm%msh%nelv
296 if (dm%shared_dof(j,k,l,i))
then
297 write(unit, fmt=
'(I8)') 1
299 write(unit, fmt=
'(I8)') 0
312 integer,
parameter :: npts = 4
313 integer :: i, j, vtk_type
316 write(unit, fmt=
'(A,I8,A)')
'POINTS', tet_msh%msh%mpts,
' double'
317 do i = 1, tet_msh%msh%mpts
318 write(unit, fmt=
'(F15.8,F15.8,F15.8)')
real(tet_msh%msh%points(i)%x,
dp)
322 write(unit, fmt=
'(A,I8,I8)')
'CELLS', tet_msh%nelv, tet_msh%nelv*(npts+1)
324 do i = 1, tet_msh%nelv
325 write(unit, fmt=
'(I8,8I8)') npts, &
326 (tet_msh%msh%get_local(tet_msh%el(i)%pts(j)%p) - 1, &
331 write(unit, fmt=
'(A,I8)')
'CELL_TYPES', tet_msh%nelv
333 do i = 1, tet_msh%nelv
334 write(unit, fmt=
'(I2)') vtk_type
343 integer,
parameter :: npts = 3
344 integer :: i, j, vtk_type
347 write(unit, fmt=
'(A,I8,A)')
'POINTS', tri_msh%mpts,
' double'
348 do i = 1, tri_msh%mpts
349 write(unit, fmt=
'(F15.8,F15.8,F15.8)')
real(tri_msh%points(i)%x,
dp)
353 write(unit, fmt=
'(A,I8,I8)')
'CELLS', tri_msh%nelv, tri_msh%nelv*(npts+1)
355 do i = 1, tri_msh%nelv
356 write(unit, fmt=
'(I8,8I8)') npts, &
357 (tri_msh%el(i)%pts(j)%p%id() - 1, j=1, npts)
361 write(unit, fmt=
'(A,I8)')
'CELL_TYPES', tri_msh%nelv
363 do i = 1, tri_msh%nelv
364 write(unit, fmt=
'(I2)') vtk_type
integer, public pe_size
MPI size of communicator.
integer, public pe_rank
MPI rank.
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.