Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
vtk_file.f90
Go to the documentation of this file.
1! Copyright (c) 2019-2025, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
36 use num_types, only : rp, dp
39 use mesh, only : mesh_t
40 use field, only : field_t
41 use dofmap, only : dofmap_t
42 use mesh_field, only : mesh_fld_t
43 use tet_mesh, only : tet_mesh_t
44 use tri_mesh, only : tri_mesh_t
45 use logger, only : neko_log
46 use comm, only : pe_size, pe_rank
47 implicit none
48 private
49
51 type, public, extends(generic_file_t) :: vtk_file_t
52 contains
53 procedure :: read => vtk_file_read
54 procedure :: write => vtk_file_write
55 end type vtk_file_t
56
57contains
58
60 subroutine vtk_file_write(this, data, t)
61 class(vtk_file_t), intent(inout) :: this
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()
66 type(mesh_fld_t), pointer :: mfld => null()
67 type(dofmap_t), pointer :: dm => null()
68 type(tet_mesh_t), pointer :: tet_msh => null()
69 type(tri_mesh_t), pointer :: tri_msh => null()
70 character(len=10) :: id_str
71 integer:: suffix_pos, file_unit
72 character(len=1024) :: fname
73
74 select type (data)
75 type is (mesh_t)
76 msh => data
77 type is (field_t)
78 msh => data%msh
79 fld => data
80 type is (mesh_fld_t)
81 msh => data%msh
82 mfld => data
83 type is (dofmap_t)
84 dm => data
85 type is (tet_mesh_t)
86 tet_msh => data
87 type is (tri_mesh_t)
88 tri_msh => data
89 class default
90 call neko_log%error('Invalid data')
91 end select
92
93 fname = trim(this%get_base_fname())
94 if (pe_size .gt. 1) then
95 write(id_str, '(i10.10)') pe_rank
96 suffix_pos = filename_suffix_pos(fname)
97 open(newunit = file_unit, &
98 file = trim(fname(1:suffix_pos-1))//id_str//'.vtk')
99 else
100 open(newunit = file_unit, file = trim(fname))
101 end if
102
103 ! Write legacy header
104 write(file_unit, fmt = '(A)') '# vtk DataFile Version 2.0'
105 write(file_unit, fmt = '(A)') 'Neko'
106 write(file_unit, fmt = '(A)') 'ASCII'
107
108 if (associated(msh)) then
109 write(file_unit, fmt = '(A)') 'DATASET UNSTRUCTURED_GRID'
110
111 call vtk_file_write_mesh(9, msh)
112
113 if (associated(mfld)) then
114 call vtk_file_write_cell_data(9, mfld)
115 else if (associated(fld)) then
116 call vtk_file_write_point_data(9, fld)
117 end if
118 else if (associated(dm)) then
119 write(file_unit, fmt = '(A)') 'DATASET POLYDATA'
120
122
124 else if (associated(tet_msh)) then
125 write(file_unit, fmt = '(A)') 'DATASET UNSTRUCTURED_GRID'
126 call vtk_file_write_tet_mesh(9, tet_msh)
127 else if (associated(tri_msh)) then
128 write(file_unit, fmt = '(A)') 'DATASET UNSTRUCTURED_GRID'
129 call vtk_file_write_tri_mesh(9, tri_msh)
130 else
131 call neko_error('Invalid data')
132 end if
133
134 close(file_unit)
135 end subroutine vtk_file_write
136
137 subroutine vtk_file_read(this, data)
138 class(vtk_file_t) :: this
139 class(*), target, intent(inout) :: data
140
141 call neko_error('VTK file read not implemented')
142 end subroutine vtk_file_read
143
145 subroutine vtk_file_write_mesh(unit, msh)
146 integer :: unit
147 type(mesh_t), intent(inout) :: msh
148 integer :: i, j, vtk_type
149 integer, dimension(8), parameter :: vcyc_to_sym = [1, 2, 4, 3, &
150 5, 6, 8, 7]
151 ! Dump coordinates
152 write(unit, fmt = '(A,I8,A)') 'POINTS', msh%mpts, ' double'
153 do i = 1, msh%mpts
154 write(unit, fmt = '(F15.8,F15.8,F15.8)') real(msh%points(i)%x, dp)
155 end do
156
157 ! Dump cells
158 write(unit, fmt = '(A,I8,I8)') 'CELLS', msh%nelv, msh%nelv*(msh%npts+1)
159 j = 0
160 do i = 1, msh%nelv
161 write(unit, fmt = '(I8,8I8)') msh%npts, &
162 (msh%get_local(msh%elements(i)%e%pts(vcyc_to_sym(j))%p) - 1, &
163 j = 1, msh%npts)
164 end do
165
166 ! Dump cell type for each element
167 write(unit, fmt = '(A,I8)') 'CELL_TYPES', msh%nelv
168 vtk_type = 9
169 if (msh%gdim .eq. 3) vtk_type = 12
170 do i = 1, msh%nelv
171 write(unit, fmt = '(I2)') vtk_type
172 end do
173
174 end subroutine vtk_file_write_mesh
175
177 subroutine vtk_file_write_cell_data(unit, mfld)
178 integer :: unit
179 type(mesh_fld_t), intent(in) :: mfld
180 integer :: i
181
182 write(unit, fmt = '(A,I8)') 'CELL_DATA', mfld%msh%nelv
183 write(unit, fmt = '(A,A,A,I8)') 'SCALARS ', trim(mfld%name), ' int', 1
184 write(unit, fmt = '(A)') 'LOOKUP_TABLE default'
185
186 do i = 1, mfld%msh%nelv
187 write(unit, fmt = '(I8)') mfld%data(i)
188 end do
189
190 end subroutine vtk_file_write_cell_data
191
195 subroutine vtk_file_write_point_data(unit, fld)
196 integer :: unit
197 type(field_t), intent(inout) :: fld
198 real(kind=dp), allocatable :: point_data(:)
199 integer :: i, j, lx, ly, lz, id(8)
200
201 if ( (fld%Xh%lx - 1 .gt. 1) .or. &
202 (fld%Xh%ly - 1 .gt. 1) .or. &
203 (fld%Xh%lz - 1 .gt. 1)) then
204 call neko_log%warning("Interpolate high-order data onto a " // &
205 "low-order mesh")
206 end if
207
208 write(unit, fmt = '(A,I8)') 'POINT_DATA', fld%msh%mpts
209 write(unit, fmt = '(A,A,A,I8)') 'SCALARS ', trim(fld%name), ' double', 1
210 write(unit, fmt = '(A)') 'LOOKUP_TABLE default'
211
212 lx = fld%Xh%lx
213 ly = fld%Xh%ly
214 lz = fld%Xh%lz
215 allocate(point_data(fld%msh%mpts))
216
217 do i = 1, fld%msh%nelv
218 do j = 1, fld%msh%npts
219 id(j) = fld%msh%get_local(fld%msh%elements(i)%e%pts(j)%p)
220 end do
221
222 point_data(id(1)) = real(fld%x(1, 1, 1, i), dp)
223 point_data(id(2)) = real(fld%x(lx, 1, 1, i), dp)
224 point_data(id(3)) = real(fld%x(1, ly, 1, i), dp)
225 point_data(id(4)) = real(fld%x(lx, ly, 1, i), dp)
226
227 if (fld%msh%gdim .eq. 3) then
228 point_data(id(5)) = real(fld%x(1, 1, lz, i), dp)
229 point_data(id(6)) = real(fld%x(lx, 1, lz, i), dp)
230 point_data(id(7)) = real(fld%x(1, ly, lz, i), dp)
231 point_data(id(8)) = real(fld%x(lx, ly, lz, i), dp)
232 end if
233
234 end do
235
236 write(unit, *) point_data
237
238 deallocate(point_data)
239
240 end subroutine vtk_file_write_point_data
241
244 integer :: unit
245 type(dofmap_t), intent(inout) :: dm
246 integer :: i, j, k, l
247
248 write(unit, fmt = '(A,I8,A)') 'POINTS', size(dm%x), ' double'
249
250 do i = 1, dm%msh%nelv
251 do l = 1, dm%Xh%lz
252 do k = 1, dm%Xh%ly
253 do j = 1, dm%Xh%lx
254 write(unit, fmt = '(F15.8,F15.8,F15.8)') &
255 real(dm%x(j,k,l,i), dp), &
256 real(dm%y(j,k,l,i), dp), &
257 real(dm%z(j,k,l,i), dp)
258 end do
259 end do
260 end do
261 end do
262
263 write(unit, fmt = '(A,I8,I8)') 'VERTICES', size(dm%x), 2*size(dm%x)
264 do i = 1, size(dm%x)
265 write(unit, fmt = '(I8,I8)') 1, i-1
266 end do
267
268
270
272 subroutine vtk_file_write_dofmap_data(unit, dm)
273 integer :: unit
274 type(dofmap_t), intent(inout) :: dm
275 integer :: i, j, k, l
276
277 write(unit, fmt = '(A,I8)') 'POINT_DATA', size(dm%dof)
278 write(unit, fmt = '(A,A,A,I8)') 'SCALARS ', 'dof_id', ' integer', 1
279 write(unit, fmt = '(A)') 'LOOKUP_TABLE default'
280
281 do i = 1, dm%msh%nelv
282 do l = 1, dm%Xh%lz
283 do k = 1, dm%Xh%ly
284 do j = 1, dm%Xh%lx
285 write(unit, fmt = '(I8)') real(dm%dof(j,k,l,i), dp)
286 end do
287 end do
288 end do
289 end do
290
291 write(unit, fmt = '(A,A,A,I8)') 'SCALARS ', 'shared_dof', ' integer', 1
292 write(unit, fmt = '(A)') 'LOOKUP_TABLE default'
293
294 do i = 1, dm%msh%nelv
295 do l = 1, dm%Xh%lz
296 do k = 1, dm%Xh%ly
297 do j = 1, dm%Xh%lx
298 if (dm%shared_dof(j,k,l,i)) then
299 write(unit, fmt = '(I8)') 1
300 else
301 write(unit, fmt = '(I8)') 0
302 end if
303 end do
304 end do
305 end do
306 end do
307
308 end subroutine vtk_file_write_dofmap_data
309
311 subroutine vtk_file_write_tet_mesh(unit, tet_msh)
312 integer :: unit
313 type(tet_mesh_t), intent(inout) :: tet_msh
314 integer, parameter :: npts = 4
315 integer :: i, j, vtk_type
316
317 ! Dump coordinates
318 write(unit, fmt = '(A,I8,A)') 'POINTS', tet_msh%msh%mpts, ' double'
319 do i = 1, tet_msh%msh%mpts
320 write(unit, fmt = '(F15.8,F15.8,F15.8)') &
321 real(tet_msh%msh%points(i)%x, dp)
322 end do
323
324 ! Dump cells
325 write(unit, fmt = '(A,I8,I8)') 'CELLS', tet_msh%nelv, tet_msh%nelv*(npts+1)
326 j = 0
327 do i = 1, tet_msh%nelv
328 write(unit, fmt = '(I8,8I8)') npts, &
329 (tet_msh%msh%get_local(tet_msh%el(i)%pts(j)%p) - 1, &
330 j = 1, npts)
331 end do
332
333 ! Dump cell type for each element
334 write(unit, fmt = '(A,I8)') 'CELL_TYPES', tet_msh%nelv
335 vtk_type = 10
336 do i = 1, tet_msh%nelv
337 write(unit, fmt = '(I2)') vtk_type
338 end do
339
340 end subroutine vtk_file_write_tet_mesh
341
343 subroutine vtk_file_write_tri_mesh(unit, tri_msh)
344 integer :: unit
345 type(tri_mesh_t), intent(inout) :: tri_msh
346 integer, parameter :: npts = 3
347 integer :: i, j, vtk_type
348
349 ! Dump coordinates
350 write(unit, fmt = '(A,I8,A)') 'POINTS', tri_msh%mpts, ' double'
351 do i = 1, tri_msh%mpts
352 write(unit, fmt = '(F15.8,F15.8,F15.8)') real(tri_msh%points(i)%x, dp)
353 end do
354
355 ! Dump cells
356 write(unit, fmt = '(A,I8,I8)') 'CELLS', tri_msh%nelv, tri_msh%nelv*(npts+1)
357 j = 0
358 do i = 1, tri_msh%nelv
359 write(unit, fmt = '(I8,8I8)') npts, &
360 (tri_msh%el(i)%pts(j)%p%id() - 1, j = 1, npts)
361 end do
362
363 ! Dump cell type for each element
364 write(unit, fmt = '(A,I8)') 'CELL_TYPES', tri_msh%nelv
365 vtk_type = 5
366 do i = 1, tri_msh%nelv
367 write(unit, fmt = '(I2)') vtk_type
368 end do
369
370 end subroutine vtk_file_write_tri_mesh
371
372end module vtk_file
double real
Definition comm.F90:1
integer, public pe_size
MPI size of communicator.
Definition comm.F90:59
integer, public pe_rank
MPI rank.
Definition comm.F90:56
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a field.
Definition field.f90:34
Module for file I/O operations.
Definition file.f90:34
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:77
Defines a mesh field.
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines a tetrahedral mesh.
Definition tet_mesh.f90:35
Defines a triangular surface mesh.
Definition tri_mesh.f90:35
Utilities.
Definition utils.f90:35
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
Definition utils.f90:58
Legacy VTK file format.
Definition vtk_file.f90:35
subroutine vtk_file_write(this, data, t)
Write data in legacy VTK.
Definition vtk_file.f90:61
subroutine vtk_file_write_tet_mesh(unit, tet_msh)
Write a tetrahedral mesh in legacy VTK format.
Definition vtk_file.f90:312
subroutine vtk_file_read(this, data)
Definition vtk_file.f90:138
subroutine vtk_file_write_dofmap_data(unit, dm)
Write a dofmap dm data as point data.
Definition vtk_file.f90:273
subroutine vtk_file_write_tri_mesh(unit, tri_msh)
Write a triangular mesh in legacy VTK format.
Definition vtk_file.f90:344
subroutine vtk_file_write_dofmap_coordinates(unit, dm)
Write xyz-coordinates of a dofmap dm as points.
Definition vtk_file.f90:244
subroutine vtk_file_write_cell_data(unit, mfld)
Write a mesh field mfld as cell data.
Definition vtk_file.f90:178
subroutine vtk_file_write_point_data(unit, fld)
Write a field fld as point data.
Definition vtk_file.f90:196
subroutine vtk_file_write_mesh(unit, msh)
Write a mesh in legacy VTK format.
Definition vtk_file.f90:146
A generic file handler.
Interface for legacy VTK files.
Definition vtk_file.f90:51