forked from pseudospectators/mpi2vis
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathconvert_vtk.f90
321 lines (247 loc) · 10.8 KB
/
convert_vtk.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
!=========================================================
! converts *.mpiio files to paraview readable *.vtk format
!
! can use either scalar or vector mode
!
! SCALAR:
! ./convert_vtk.out nx ny nz -scalar mask_00010.mpiio
! this creates the mask_00010.vtk file
!
! VECTOR:
! ./convert_vtk.out nx ny nz -vector ux_00010.mpiio uy_00010.mpiio uz_00010.mpiio
! this creates ONE file u_00010.vtk
!
! PROGRAM ALSO WORKS WITH *.BINARY FILES!
!
!=========================================================
! THOMAS ENGELS, Aix-Marseille Université and Technische Universität Berlin
! APRIL 2013
!=========================================================
program convert_mpiio
use mpi
implicit none
integer, parameter :: pr_in = 8, pr_out = 4 ! precision for input and output (INPUT PRECISION MUST BE KNOWN!!!)
integer :: nx, ny, nz, mpicode, nunit
integer :: ix, iy, iz, iSaveAscii, iArg, ii, j1,j2
real (kind=pr_out), dimension (:,:,:), allocatable :: field_out
real (kind=pr_out), dimension (:,:,:), allocatable :: ux,uy,uz
real (kind=pr_out), dimension (:,:), allocatable :: tmp
character (len=128) :: fname, nx_str, ny_str, nz_str, fname_out, modus
character (len=128) :: fname_x, fname_y, fname_z
character (len=128) :: cbuffer, frmt
!--------------------------------------------------
! read in command line arguments
!--------------------------------------------------
call get_command_argument(1, nx_str)
if (len_trim(nx_str).ne.0) then
read (nx_str, *) nx
else
call help
write(*,*) "nx not specified!"
stop
endif
call get_command_argument(2, ny_str)
if (len_trim(ny_str).ne.0) then
read (ny_str, *) ny
else
call help
write (*,*) "ny not specified!"
stop
endif
call get_command_argument(3, nz_str)
if (len_trim(nz_str).ne.0) then
read (nz_str, *) nz
else
call help
write(*,*) "nz not specified!"
stop
endif
call get_command_argument(4, modus)
if (len_trim(modus)==0) then
call help
write (*,*) "you forgot to set the modus, either -scalar or -vector?"
stop
endif
call MPI_INIT (mpicode)
!-----------------------------------------
! scalar modus
!-----------------------------------------
if ( modus=='-scalar' ) then
!-------------------------------------------
! get filename (only one in this case)
!-------------------------------------------
call get_command_argument( 5, fname )
if (len_trim(fname)==0) then
write(*,*) ("no file given")
stop
endif
call check(fname)
! fetch variable name "mask_00010.mpiio"
j1 = index ( fname,"_" ) ! name is fname(1:j1-1) "mask"
j2 = index ( fname,"." ) ! filename is fname(1:j2-1) "mask_00010"
!---------------------------------------------------------------
! read in the files and write the VTK file
!---------------------------------------------------------------
nunit = 11
open (nunit, file=fname(1:j2-1)//".vtk", status='replace', access='STREAM',convert='big_endian')
!----------------------------
! write VTK file
!----------------------------
write(nunit) "# vtk DataFile Version 3.0"//char(10) ! standard header
write(nunit) "Fluid Field"//char(10) ! description
write(nunit) "BINARY"//char(10) ! it's a binary file
write(nunit) "DATASET STRUCTURED_POINTS"//char(10) ! we use uniform grids
write(cbuffer,fmt="('DIMENSIONS ',3(i5,1x))") nx,ny,nz ! this is the grid dimensions
write(nunit) cbuffer//char(10)
write(cbuffer,fmt="('ORIGIN ',3(i1,1x))") 0, 0, 0 ! our origin is 0,0,0
write(nunit) cbuffer//char(10)
write(cbuffer,fmt="('SPACING ',3(i1,1x))") 1, 1, 1 ! our spacing is 1,1,1 (note it's possible to use the real spacing with the domain size)
write(nunit) cbuffer//char(10)
write(cbuffer,fmt="('POINT_DATA',I10)") nx*ny*nz ! tell paraview to expect point data, not cell data
write(nunit) cbuffer//char(10)
write(cbuffer,fmt="('SCALARS ',A10,' float 1')") fname(1:j1-1) ! save ONE scalar field, named as the variable name (the last "1" is the dim of scalar field)
write(nunit) cbuffer//char(10)
write(cbuffer,fmt="('LOOKUP_TABLE default')") ! standard lookup table (don't touch)
write(nunit) cbuffer//char(10)
allocate ( field_out(0:nx-1, 0:ny-1, 0:nz-1) ) ! alloc output field
call ReadMPIIO (nx, ny, nz, mpicode, fname, field_out) ! read in file (field_out is single precision)
write(nunit) field_out ! dump to vtk file
close(unit=nunit) ! close vtk file
deallocate ( field_out )
!-----------------------------------------
! vector modus
!-----------------------------------------
elseif ( modus=='-vector') then
!-------------------------------------------
! get filenames (three files in this case)
!-------------------------------------------
call get_command_argument( 5, fname_x )
if (len_trim(fname_x)==0) then
write (*,*) "no x-file given"
stop
endif
call check(fname_x)
call get_command_argument( 6, fname_y )
if (len_trim(fname_y)==0) then
write (*,*) "no y-file given"
stop
endif
call check(fname_x)
call get_command_argument( 7, fname_z )
if (len_trim(fname_z)==0) then
write (*,*) "no z-file given"
stop
endif
call check(fname_x)
! fetch variable name "vorx_00010.mpiio" "vory_00010.mpiio" "vorz_00010.mpiio"
j1 = index ( fname_x,"x" ) ! name is fname(1:j1-1) "vor"
ii = index ( fname_x,"_" )
j2 = index ( fname_x,"." ) ! filename is fname(1:j2-1) "vorx_00010"
fname_out = fname_x(1:j1-1)//fname_x(ii:j2-1)//".vtk" ! out name is "vor_00010.vtk"
!---------------------------------------------------------------
! read in the files and write the VTK file
!---------------------------------------------------------------
nunit = 11
open (nunit, file=fname_out, status='replace', access='STREAM',convert='big_endian')
!----------------------------
! write VTK file
!----------------------------
write(nunit) "# vtk DataFile Version 3.0"//char(10) ! standard header
write(nunit) "Fluid Field"//char(10) ! description
write(nunit) "BINARY"//char(10) ! it's a binary file
write(nunit) "DATASET STRUCTURED_POINTS"//char(10) ! we use uniform grids
write(cbuffer,fmt="('DIMENSIONS ',3(i5,1x))") nx,ny,nz ! this is the grid dimensions
write(nunit) cbuffer//char(10)
write(cbuffer,fmt="('ORIGIN ',3(i1,1x))") 0, 0, 0 ! our origin is 0,0,0
write(nunit) cbuffer//char(10)
write(cbuffer,fmt="('SPACING ',3(i1,1x))") 1, 1, 1 ! our spacing is 1,1,1 (note it's possible to use the real spacing with the domain size)
write(nunit) cbuffer//char(10)
write(cbuffer,fmt="('POINT_DATA',I10)") nx*ny*nz ! tell paraview to expect point data, not cell data
write(nunit) cbuffer//char(10)
write(cbuffer,fmt="('VECTORS ',A10,' float')") fname_x(1:j1-1) ! save ONE scalar field, named as the variable name
write(nunit) cbuffer//char(10)
allocate ( ux(0:nx-1, 0:ny-1, 0:nz-1) )
allocate ( uy(0:nx-1, 0:ny-1, 0:nz-1) )
allocate ( uz(0:nx-1, 0:ny-1, 0:nz-1) )
allocate ( tmp(3,nx*ny*nz) )
call ReadMPIIO (nx, ny, nz, mpicode, fname_x, ux) ! read in X-file (single precision)
call ReadMPIIO (nx, ny, nz, mpicode, fname_y, uy) ! read in Y-file (single precision)
call ReadMPIIO (nx, ny, nz, mpicode, fname_z, uz) ! read in Z-file (single precision)
ii = 0 ! counter
do iz=0,nz-1
do iy=0,ny-1
do ix=0,nx-1
ii = ii + 1
tmp(1,ii) = ux(ix,iy,iz)
tmp(2,ii) = uy(ix,iy,iz)
tmp(3,ii) = uz(ix,iy,iz)
end do
end do
end do
write (nunit) tmp ! dump to disk
close(unit=nunit) ! close vtk file
deallocate (ux,uy,uz,tmp)
else
write (*,*) "wrong modus, choose -scalar or -vector"
stop
endif
call MPI_FINALIZE (mpicode)
end program convert_mpiio
subroutine help
write (*,*) "--------------------------------------------------------------"
write (*,*) " converter "
write (*,*) " *.mpiio -> *.vtk (for PARAVIEW)"
write (*,*) "--------------------------------------------------------------"
write (*,*) " usage: ./convert_vtk.out nx ny nz -modus file1 [file2] [file3]"
write (*,*) ""
write (*,*) " for example (SCALAR):"
write (*,*) " ./convert_vtk.out 64 128 64 -scalar mask_00010.mpiio"
write (*,*) ""
write (*,*) " for example (VECTOR):"
write (*,*) " ./convert_vtk.out 64 128 64 -vector vorx_00010.mpiio vory_00010.mpiio vorz_00010.mpiio"
write (*,*) "--------------------------------------------------------------"
end subroutine
subroutine ReadMPIIO (nx, ny, nz, mpicode, fname, field_out)
use mpi
implicit none
character(len=128), intent (in) :: fname
integer, intent(in) :: nx,ny,nz
integer, intent(inout) :: mpicode
integer, dimension(MPI_STATUS_SIZE) :: mpistatus
integer :: filedesc, j1,j2
real(kind=4), intent(out) :: field_out (0:nx-1, 0:ny-1, 0:nz-1)
real(kind=8), dimension(:,:,:), allocatable :: field_in
integer, parameter :: mpireal = MPI_DOUBLE_PRECISION ! double precision array for input
j1=index(fname,".")
j2=len_trim(fname)
!--------------------------------------------------
! read MPI data or BINARY data
!--------------------------------------------------
if ( fname(j1:j2) == ".mpiio") then
allocate ( field_in(0:nx-1, 0:ny-1, 0:nz-1) )
call MPI_FILE_OPEN (MPI_COMM_WORLD,trim(fname),MPI_MODE_RDONLY,MPI_INFO_NULL,filedesc,mpicode)
call MPI_FILE_READ_ORDERED (filedesc,field_in,nx*ny*nz,mpireal,mpistatus,mpicode)
call MPI_FILE_CLOSE (filedesc,mpicode)
field_out = real(field_in, kind=4) ! convert field to output precision
deallocate (field_in)
elseif ( fname(j1:j2) == ".binary") then
open (14, file=trim(fname), status='old', action='read', form='unformatted')
read (14) field_out
close(14)
endif
write(*,'("Reading ",A,"... Min:Max=",es12.4,":",es12.4,1x,"nx:ny:nz=",i3,":",i3,":",i3)') &
trim(fname),minval (field_out), maxval (field_out),nx,ny,nz
end subroutine
!-------------------------------------------
! check if file exists (returns true if so)
!-------------------------------------------
subroutine check(fname)
implicit none
character(len=128), intent (in) :: fname
logical :: file_exists
inquire( file=trim(fname),exist=file_exists )
if (file_exists .eqv. .false.) then
write (*,'("ERROR!! File ",A," NOT found.")') trim(fname)
stop
endif
end subroutine