Neko 1.99.1
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, file=trim(fname(1:suffix_pos-1))//id_str//'.vtk')
98 else
99 open(newunit=file_unit, file=trim(fname))
100 end if
101
102 ! Write legacy header
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'
106
107 if (associated(msh)) then
108 write(file_unit, fmt='(A)') 'DATASET UNSTRUCTURED_GRID'
109
110 call vtk_file_write_mesh(9, msh)
111
112 if (associated(mfld)) then
113 call vtk_file_write_cell_data(9, mfld)
114 else if (associated(fld)) then
115 call vtk_file_write_point_data(9, fld)
116 end if
117 else if (associated(dm)) then
118 write(file_unit, fmt='(A)') 'DATASET POLYDATA'
119
121
123 else if (associated(tet_msh)) then
124 write(file_unit, fmt='(A)') 'DATASET UNSTRUCTURED_GRID'
125 call vtk_file_write_tet_mesh(9, tet_msh)
126 else if (associated(tri_msh)) then
127 write(file_unit, fmt='(A)') 'DATASET UNSTRUCTURED_GRID'
128 call vtk_file_write_tri_mesh(9, tri_msh)
129 else
130 call neko_error('Invalid data')
131 end if
132
133 close(file_unit)
134 end subroutine vtk_file_write
135
136 subroutine vtk_file_read(this, data)
137 class(vtk_file_t) :: this
138 class(*), target, intent(inout) :: data
139
140 call neko_error('VTK file read not implemented')
141 end subroutine vtk_file_read
142
144 subroutine vtk_file_write_mesh(unit, msh)
145 integer :: unit
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, &
149 5, 6, 8, 7/)
150 ! Dump coordinates
151 write(unit, fmt='(A,I8,A)') 'POINTS', msh%mpts,' double'
152 do i = 1, msh%mpts
153 write(unit, fmt='(F15.8,F15.8,F15.8)') real(msh%points(i)%x,dp)
154 end do
155
156 ! Dump cells
157 write(unit, fmt='(A,I8,I8)') 'CELLS', msh%nelv, msh%nelv*(msh%npts+1)
158 j = 0
159 do i = 1, msh%nelv
160 write(unit, fmt='(I8,8I8)') msh%npts, &
161 (msh%get_local(msh%elements(i)%e%pts(vcyc_to_sym(j))%p) - 1, &
162 j=1, msh%npts)
163 end do
164
165 ! Dump cell type for each element
166 write(unit, fmt='(A,I8)') 'CELL_TYPES', msh%nelv
167 vtk_type = 9
168 if (msh%gdim .eq. 3) vtk_type = 12
169 do i = 1, msh%nelv
170 write(unit, fmt='(I2)') vtk_type
171 end do
172
173 end subroutine vtk_file_write_mesh
174
176 subroutine vtk_file_write_cell_data(unit, mfld)
177 integer :: unit
178 type(mesh_fld_t), intent(in) :: mfld
179 integer :: i
180
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'
184
185 do i = 1, mfld%msh%nelv
186 write(unit, fmt='(I8)') mfld%data(i)
187 end do
188
189 end subroutine vtk_file_write_cell_data
190
194 subroutine vtk_file_write_point_data(unit, fld)
195 integer :: unit
196 type(field_t), intent(inout) :: fld
197 real(kind=dp), allocatable :: point_data(:)
198 integer :: i, j, lx, ly, lz, id(8)
199
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")
204 end if
205
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'
209
210 lx = fld%Xh%lx
211 ly = fld%Xh%ly
212 lz = fld%Xh%lz
213 allocate(point_data(fld%msh%mpts))
214
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)
218 end do
219
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)
224
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)
230 end if
231
232 end do
233
234 write(unit, *) point_data
235
236 deallocate(point_data)
237
238 end subroutine vtk_file_write_point_data
239
242 integer :: unit
243 type(dofmap_t), intent(inout) :: dm
244 integer :: i,j,k,l
245
246 write(unit, fmt='(A,I8,A)') 'POINTS', size(dm%x),' double'
247
248 do i = 1, dm%msh%nelv
249 do l = 1, dm%Xh%lz
250 do k = 1, dm%Xh%ly
251 do j = 1, dm%Xh%lx
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)
256 end do
257 end do
258 end do
259 end do
260
261 write(unit, fmt='(A,I8,I8)') 'VERTICES', size(dm%x), 2*size(dm%x)
262 do i = 1, size(dm%x)
263 write(unit, fmt='(I8,I8)') 1,i-1
264 end do
265
266
268
270 subroutine vtk_file_write_dofmap_data(unit, dm)
271 integer :: unit
272 type(dofmap_t), intent(inout) :: dm
273 integer :: i, j, k, l
274
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'
278
279 do i = 1, dm%msh%nelv
280 do l = 1, dm%Xh%lz
281 do k = 1, dm%Xh%ly
282 do j = 1, dm%Xh%lx
283 write(unit, fmt='(I8)') real(dm%dof(j,k,l,i),dp)
284 end do
285 end do
286 end do
287 end do
288
289 write(unit, fmt='(A,A,A,I8)') 'SCALARS ', 'shared_dof', ' integer', 1
290 write(unit, fmt='(A)') 'LOOKUP_TABLE default'
291
292 do i = 1, dm%msh%nelv
293 do l = 1, dm%Xh%lz
294 do k = 1, dm%Xh%ly
295 do j = 1, dm%Xh%lx
296 if (dm%shared_dof(j,k,l,i)) then
297 write(unit, fmt='(I8)') 1
298 else
299 write(unit, fmt='(I8)') 0
300 end if
301 end do
302 end do
303 end do
304 end do
305
306 end subroutine vtk_file_write_dofmap_data
307
309 subroutine vtk_file_write_tet_mesh(unit, tet_msh)
310 integer :: unit
311 type(tet_mesh_t), intent(inout) :: tet_msh
312 integer, parameter :: npts = 4
313 integer :: i, j, vtk_type
314
315 ! Dump coordinates
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)
319 end do
320
321 ! Dump cells
322 write(unit, fmt='(A,I8,I8)') 'CELLS', tet_msh%nelv, tet_msh%nelv*(npts+1)
323 j = 0
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, &
327 j=1, npts)
328 end do
329
330 ! Dump cell type for each element
331 write(unit, fmt='(A,I8)') 'CELL_TYPES', tet_msh%nelv
332 vtk_type = 10
333 do i = 1, tet_msh%nelv
334 write(unit, fmt='(I2)') vtk_type
335 end do
336
337 end subroutine vtk_file_write_tet_mesh
338
340 subroutine vtk_file_write_tri_mesh(unit, tri_msh)
341 integer :: unit
342 type(tri_mesh_t), intent(inout) :: tri_msh
343 integer, parameter :: npts = 3
344 integer :: i, j, vtk_type
345
346 ! Dump coordinates
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)
350 end do
351
352 ! Dump cells
353 write(unit, fmt='(A,I8,I8)') 'CELLS', tri_msh%nelv, tri_msh%nelv*(npts+1)
354 j = 0
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)
358 end do
359
360 ! Dump cell type for each element
361 write(unit, fmt='(A,I8)') 'CELL_TYPES', tri_msh%nelv
362 vtk_type = 5
363 do i = 1, tri_msh%nelv
364 write(unit, fmt='(I2)') vtk_type
365 end do
366
367 end subroutine vtk_file_write_tri_mesh
368
369end module vtk_file
double real
Definition comm.F90:1
integer, public pe_size
MPI size of communicator.
Definition comm.F90:58
integer, public pe_rank
MPI rank.
Definition comm.F90:55
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:70
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:310
subroutine vtk_file_read(this, data)
Definition vtk_file.f90:137
subroutine vtk_file_write_dofmap_data(unit, dm)
Write a dofmap dm data as point data.
Definition vtk_file.f90:271
subroutine vtk_file_write_tri_mesh(unit, tri_msh)
Write a triangular mesh in legacy VTK format.
Definition vtk_file.f90:341
subroutine vtk_file_write_dofmap_coordinates(unit, dm)
Write xyz-coordinates of a dofmap dm as points.
Definition vtk_file.f90:242
subroutine vtk_file_write_cell_data(unit, mfld)
Write a mesh field mfld as cell data.
Definition vtk_file.f90:177
subroutine vtk_file_write_point_data(unit, fld)
Write a field fld as point data.
Definition vtk_file.f90:195
subroutine vtk_file_write_mesh(unit, msh)
Write a mesh in legacy VTK format.
Definition vtk_file.f90:145
A generic file handler.
Interface for legacy VTK files.
Definition vtk_file.f90:51