All Classes Files Functions Variables Pages
module_vtk_io.f90
Go to the documentation of this file.
1 !=====================================================================================================================================
2 ! EFISPEC3D !
3 ! (Elements FInis SPECtraux 3D) !
4 ! !
5 ! http://efispec.free.fr !
6 ! !
7 ! !
8 ! This file is part of EFISPEC3D !
9 ! Please refer to http://efispec.free.fr if you use it or part of it !
10 ! !
11 ! !
12 !1 ---> French License: CeCILL V2 !
13 ! !
14 ! Copyright BRGM 2009 contributeurs : Florent DE MARTIN !
15 ! David MICHEA !
16 ! Philippe THIERRY !
17 ! !
18 ! Contact: f.demartin at brgm.fr !
19 ! !
20 ! Ce logiciel est un programme informatique servant a resoudre l'equation du !
21 ! mouvement en trois dimensions via une methode des elements finis spectraux. !
22 ! !
23 ! Ce logiciel est regi par la licence CeCILL soumise au droit francais et !
24 ! respectant les principes de diffusion des logiciels libres. Vous pouvez !
25 ! utiliser, modifier et/ou redistribuer ce programme sous les conditions de la !
26 ! licence CeCILL telle que diffusee par le CEA, le CNRS et l'INRIA sur le site !
27 ! "http://www.cecill.info". !
28 ! !
29 ! En contrepartie de l'accessibilite au code source et des droits de copie, de !
30 ! modification et de redistribution accordes par cette licence, il n'est offert !
31 ! aux utilisateurs qu'une garantie limitee. Pour les memes raisons, seule une !
32 ! responsabilite restreinte pese sur l'auteur du programme, le titulaire des !
33 ! droits patrimoniaux et les concedants successifs. !
34 ! !
35 ! A cet egard l'attention de l'utilisateur est attiree sur les risques associes !
36 ! au chargement, a l'utilisation, a la modification et/ou au developpement et a !
37 ! la reproduction du logiciel par l'utilisateur etant donne sa specificite de !
38 ! logiciel libre, qui peut le rendre complexe a manipuler et qui le reserve donc !
39 ! a des developpeurs et des professionnels avertis possedant des connaissances !
40 ! informatiques approfondies. Les utilisateurs sont donc invites a charger et !
41 ! tester l'adequation du logiciel a leurs besoins dans des conditions permettant !
42 ! d'assurer la securite de leurs systemes et ou de leurs donnees et, plus !
43 ! generalement, a l'utiliser et l'exploiter dans les memes conditions de !
44 ! securite. !
45 ! !
46 ! Le fait que vous puissiez acceder a cet en-tete signifie que vous avez pris !
47 ! connaissance de la licence CeCILL et que vous en avez accepte les termes. !
48 ! !
49 ! !
50 !2 ---> International license: GNU GPL V3 !
51 ! !
52 ! EFISPEC3D is a computer program that solves the three-dimensional equations of !
53 ! motion using a finite spectral-element method. !
54 ! !
55 ! Copyright (C) 2009 Florent DE MARTIN !
56 ! !
57 ! Contact: f.demartin at brgm.fr !
58 ! !
59 ! This program is free software: you can redistribute it and/or modify it under !
60 ! the terms of the GNU General Public License as published by the Free Software !
61 ! Foundation, either version 3 of the License, or (at your option) any later !
62 ! version. !
63 ! !
64 ! This program is distributed in the hope that it will be useful, but WITHOUT ANY !
65 ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A !
66 ! PARTICULAR PURPOSE. See the GNU General Public License for more details. !
67 ! !
68 ! You should have received a copy of the GNU General Public License along with !
69 ! this program. If not, see http://www.gnu.org/licenses/. !
70 ! !
71 ! !
72 !3 ---> Third party libraries !
73 ! !
74 ! EFISPEC3D uses the following of third party libraries: !
75 ! !
76 ! --> METIS 5.1.0 !
77 ! see http://glaros.dtc.umn.edu/gkhome/metis/metis/overview !
78 ! !
79 ! --> Lib_VTK_IO !
80 ! see S. Zaghi's website: https://github.com/szaghi/Lib_VTK_IO !
81 ! !
82 ! --> INTERP_LINEAR !
83 ! see J. Burkardt website: http://people.sc.fsu.edu/~jburkardt/ !
84 ! !
85 ! --> SAC !
86 ! http://ds.iris.edu/files/sac-manual/ !
87 ! !
88 ! --> EXODUS II !
89 ! http://sourceforge.net/projects/exodusii/ !
90 ! !
91 ! --> NETCDF !
92 ! http://www.unidata.ucar.edu/software/netcdf/ !
93 ! !
94 ! --> HDF5 !
95 ! http://www.hdfgroup.org/HDF5/ !
96 ! !
97 ! Some of these libraries are located in directory lib and pre-compiled !
98 ! with intel compiler for x86_64 platform. !
99 ! !
100 ! !
101 !4 ---> Related Articles !
102 ! !
103 ! De Martin, F., Matsushima, M., Kawase, H. (BSSA, 2013) !
104 ! Impact of geometric effects on near-surface Green's functions !
105 ! doi:10.1785/0120130039 !
106 ! !
107 ! Aochi, H., Ducellier, A., Dupros, F., Delatre, M., Ulrich, T., De Martin, F., !
108 ! Yoshimi, M., (Pure Appl. Geophys. 2013) !
109 ! Finite Difference Simulations of Seismic Wave Propagation for the 2007 Mw 6.6 !
110 ! Niigata-ken Chuetsu-Oki Earthquake: Validity of Models and Reliable !
111 ! Input Ground Motion in the Near-Field !
112 ! doi:10.1007/s00024-011-0429-5 !
113 ! !
114 ! De Martin, F. (BSSA, 2011) !
115 ! Verification of a Spectral-Element Method Code for the Southern California !
116 ! Earthquake Center LOH.3 Viscoelastic Case !
117 ! doi:10.1785/0120100305 !
118 ! !
119 !=====================================================================================================================================
120 
126 
131 
132  use mod_global_variables, only : get_newunit
133 
134  implicit none
135 
136  private
137 
138  !functions for VTK XML
139  public :: vtk_ini_xml
140  public :: vtk_geo_xml
141  public :: vtk_con_xml
142  public :: vtk_dat_xml
143  public :: vtk_var_xml
144  public :: vtk_end_xml
145 
146  !portable kind-precision
147  public :: r4p, fr4p
148  public :: r_p, fr_p
149  public :: i4p, fi4p
150  public :: i1p, fi1p
151  public :: i_p, fi_p
152 
153  !----------------------------------------------------------------------------------------------------------------------------------
154 
157  interface vtk_geo_xml
158  module procedure vtk_geo_xml_strg_r4, & ! real(R4P) StructuredGrid
159  vtk_geo_xml_rect_r4, & ! real(R4P) RectilinearGrid
160  vtk_geo_xml_unst_r4, & ! real(R4P) UnstructuredGrid
161  vtk_geo_xml_closep ! closing tag "Piece" function
162  endinterface
163 
166  interface vtk_var_xml
167  module procedure vtk_var_xml_scal_r4, & ! real (R4P) scalar
168  vtk_var_xml_scal_i4, & ! integer(I4P) scalar
169  vtk_var_xml_vect_r4, & ! real (R4P) vectorial
170  vtk_var_xml_vect_i4 ! integer(I4P) vectorial
171  endinterface
172 
173 
174 !********************************
175 !Real precision definitions
176 !********************************
177 
179  integer , parameter :: R4P = selected_real_kind(6,37)
180 
182  integer , parameter :: R_P = R4P
183 
184 
185 !********************************
186 !Integer precision definitions
187 !********************************
188 
190  integer , parameter :: I4P = selected_int_kind(9)
191 
193  integer , parameter :: I1P = selected_int_kind(2)
194 
196  integer , parameter :: I_P = I4P
197 
198 
199 !********************************
200 !Real output formats
201 !********************************
202 
204  character( 9), parameter :: FR4P = '(E14.6E2)'
205 
207  character(10), parameter :: FR_P = '(E23.15E3)'
208 
209 
210 !********************************
211 !Integer output formats
212 !********************************
213 
215  character(5) , parameter :: FI4P = '(I12)'
216 
218  character(4) , parameter :: FI1P = '(I5)'
219 
221  character(5) , parameter :: FI_P = '(I12)'
222 
223 
224 !********************************
225 !private variables
226 !********************************
227 
229  integer(I4P) , parameter :: maxlen = 500
230 
232  character(1) , parameter :: end_rec = char(10)
233 
235  integer(I4P) , parameter :: f_out_ascii = 0
236 
238  integer(I4P) , parameter :: f_out_binary = 1
239 
241  integer(I4P) :: f_out = f_out_ascii
242 
244  integer(I4P) :: Unit_VTK
245 
247  integer(I4P) :: Unit_VTK_Append
248 
250  integer(I4P) :: N_Byte
251 
253  integer(I4P) :: ioffset
254 
256  integer(I4P) :: indent
257 
259  integer(I4P) :: Tipo_I4
260 
262  integer(I1P) :: Tipo_I1
263 
265  real(R4P) :: Tipo_R4
266 
268  character(len=maxlen) :: topology
269 
270 !----------------------------------------------------------------------------------------------------------------------------------
271 
272  contains
273 
274 !
275 !
287 !***********************************************************************************************************************************************************************************
288  function vtk_ini_xml(output_format,filename,mesh_topology,nx1,nx2,ny1,ny2,nz1,nz2) result(E_IO)
289 !***********************************************************************************************************************************************************************************
290 
291  implicit none
292 
293  !--------------------------------------------------------------------------------------------------------------------------------
294  character(*), intent(IN) :: output_format ! output format: ASCII or BINARY
295  character(*), intent(IN) :: filename ! file name
296  character(*), intent(IN) :: mesh_topology ! mesh topology
297  integer(I4P), intent(IN), optional :: nx1,nx2 ! initial and final nodes of x axis
298  integer(I4P), intent(IN), optional :: ny1,ny2 ! initial and final nodes of y axis
299  integer(I4P), intent(IN), optional :: nz1,nz2 ! initial and final nodes of z axis
300 
301  integer(I4P) :: e_io ! Input/Output inquiring flag: $0$ if IO is done, $> 0$ if IO is not done
302  character(len=maxlen) :: s_buffer ! buffer string
303  !--------------------------------------------------------------------------------------------------------------------------------
304 
305  topology = trim(mesh_topology)
306 
307  unit_vtk = get_newunit()
308 
309  select case(trim(output_format))
310 
311  case('ASCII')
312  f_out = f_out_ascii
313  open(unit=unit_vtk,file=trim(filename),form='FORMATTED',access='SEQUENTIAL',action='WRITE',iostat=e_io)
314  ! writing header of file
315  write(unit=unit_vtk,fmt='(A)' ,iostat=e_io)'<?xml version="1.0"?>'
316  write(unit=unit_vtk,fmt='(A)' ,iostat=e_io)'<VTKFile type="'//trim(topology)//'" version="0.1" byte_order="BigEndian">'
317  indent = 2
318  select case(trim(topology))
319  case('RectilinearGrid','StructuredGrid')
320  write(unit=unit_vtk,fmt='(A,6'//fi4p//',A)',iostat=e_io)repeat(' ',indent)//'<'//trim(topology)//' WholeExtent="',nx1,nx2,ny1,ny2,nz1,nz2,'">'
321  case('UnstructuredGrid')
322  write(unit=unit_vtk,fmt='(A)' ,iostat=e_io)repeat(' ',indent)//'<'//trim(topology)//'>'
323  endselect
324  indent = indent + 2
325 
326  case('BINARY')
327  f_out = f_out_binary
328  open(unit=unit_vtk,file=trim(filename),form='UNFORMATTED',access='SEQUENTIAL',action='WRITE',convert='BIG_ENDIAN',buffered='YES',recordtype='STREAM',iostat=e_io)
329  ! writing header of file
330  write(unit=unit_vtk, iostat=e_io)'<?xml version="1.0"?>'//end_rec
331  write(unit=unit_vtk, iostat=e_io)'<VTKFile type="'//trim(topology)//'" version="0.1" byte_order="BigEndian">'//end_rec
332  indent = 2
333  select case(trim(topology))
334  case('RectilinearGrid','StructuredGrid')
335  write(s_buffer,fmt='(A,6'//fi4p//',A)',iostat=e_io)repeat(' ',indent)//'<'//trim(topology)//' WholeExtent="',nx1,nx2,ny1,ny2,nz1,nz2,'">'
336  case('UnstructuredGrid')
337  write(s_buffer,fmt='(A)', iostat=e_io)repeat(' ',indent)//'<'//trim(topology)//'>'
338  endselect
339  write(unit=unit_vtk, iostat=e_io)trim(s_buffer)//end_rec
340  indent = indent + 2
341 
342  !opening the SCRATCH file used for appending raw binary data
343  unit_vtk_append=get_newunit()
344  open(unit=unit_vtk_append,form='UNFORMATTED',access='SEQUENTIAL',action='WRITE',convert='BIG_ENDIAN',recordtype='STREAM',buffered='YES',status='SCRATCH',iostat=e_io)
345  ioffset = 0 ! initializing offset puntator
346 
347  endselect
348 
349  return
350 
351 !***********************************************************************************************************************************************************************************
352  endfunction vtk_ini_xml
353 !***********************************************************************************************************************************************************************************
354 
355 !
356 !
369 !***********************************************************************************************************************************************************************************
370  function vtk_geo_xml_strg_r4(nx1,nx2,ny1,ny2,nz1,nz2,NN,X,Y,Z) result(E_IO)
371 !***********************************************************************************************************************************************************************************
372 
373  implicit none
374 
375  !--------------------------------------------------------------------------------------------------------------------------------
376  integer(I4P), intent(IN) :: nx1,nx2 ! initial and final nodes of x axis
377  integer(I4P), intent(IN) :: ny1,ny2 ! initial and final nodes of y axis
378  integer(I4P), intent(IN) :: nz1,nz2 ! initial and final nodes of z axis
379  integer(I4P), intent(IN) :: nn ! number of all nodes
380  real (R4P), intent(IN) :: x(1:nn) ! x coordinates
381  real (R4P), intent(IN) :: y(1:nn) ! y coordinates
382  real (R4P), intent(IN) :: z(1:nn) ! z coordinates
383 
384  integer(I4P) :: e_io ! Input/Output inquiring flag: $0$ if IO is done, $> 0$ if IO is not done
385  integer(I4P) :: n1 ! counter
386 
387  character(len=maxlen) :: s_buffer ! buffer string
388  !--------------------------------------------------------------------------------------------------------------------------------
389 
390 
391  select case(f_out)
392 
393  case(f_out_ascii)
394  write(unit=unit_vtk,fmt='(A,6'//fi4p//',A)',iostat=e_io)repeat(' ',indent)//'<Piece Extent="',nx1,nx2,ny1,ny2,nz1,nz2,'">'
395  indent = indent + 2
396  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'<Points>'
397  indent = indent + 2
398  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'<DataArray type="Float32" NumberOfComponents="3" Name="Point" format="ascii">'
399  write(unit=unit_vtk,fmt='(3'//fr4p//')', iostat=e_io)(x(n1),y(n1),z(n1),n1=1,nn)
400  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'</DataArray>'
401  indent = indent - 2
402  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'</Points>'
403 
404  case(f_out_binary)
405  write(s_buffer,fmt='(A,6'//fi4p//',A)',iostat=e_io)repeat(' ',indent)//'<Piece Extent="',nx1,nx2,ny1,ny2,nz1,nz2,'">'
406  indent = indent + 2
407  write(unit=unit_vtk, iostat=e_io)trim(s_buffer)//end_rec
408  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<Points>'//end_rec
409  indent = indent + 2
410  write(s_buffer,fmt='(I10)', iostat=e_io)ioffset
411  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<DataArray type="Float32" NumberOfComponents="3" Name="Point" format="appended" offset="',trim(s_buffer),'">'//end_rec
412  n_byte = 3*nn*sizeof(tipo_r4)
413  ioffset = ioffset + sizeof(tipo_i4) + n_byte
414  write(unit=unit_vtk_append, iostat=e_io)n_byte,'R4',3*nn
415  write(unit=unit_vtk_append, iostat=e_io)(x(n1),y(n1),z(n1),n1=1,nn)
416  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</DataArray>'//end_rec
417  indent = indent - 2
418  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</Points>'//end_rec
419 
420  endselect
421 
422  return
423 !***********************************************************************************************************************************************************************************
424  endfunction vtk_geo_xml_strg_r4
425 !***********************************************************************************************************************************************************************************
426 
427 !
428 !
440 !***********************************************************************************************************************************************************************************
441  function vtk_geo_xml_rect_r4(nx1,nx2,ny1,ny2,nz1,nz2,X,Y,Z) result(E_IO)
442 !***********************************************************************************************************************************************************************************
443 
444  implicit none
445 
446  !--------------------------------------------------------------------------------------------------------------------------------
447  integer(I4P), intent(IN) :: nx1,nx2 ! initial and final nodes of x axis
448  integer(I4P), intent(IN) :: ny1,ny2 ! initial and final nodes of y axis
449  integer(I4P), intent(IN) :: nz1,nz2 ! initial and final nodes of z axis
450  real (R4P), intent(IN) :: x(nx1:nx2) ! x coordinates
451  real (R4P), intent(IN) :: y(ny1:ny2) ! y coordinates
452  real (R4P), intent(IN) :: z(nz1:nz2) ! z coordinates
453 
454  integer(I4P) :: e_io ! Input/Output inquiring flag: $0$ if IO is done, $> 0$ if IO is not done
455  integer(I4P) :: n1 ! counter
456 
457  character(len=maxlen) :: s_buffer ! buffer string
458  !--------------------------------------------------------------------------------------------------------------------------------
459 
460  select case(f_out)
461 
462  case(f_out_ascii)
463  write(unit=unit_vtk,fmt='(A,6'//fi4p//',A)',iostat=e_io)repeat(' ',indent)//'<Piece Extent="',nx1,nx2,ny1,ny2,nz1,nz2,'">'
464  indent = indent + 2
465  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'<Coordinates>'
466  indent = indent + 2
467  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'<DataArray type="Float32" Name="X" format="ascii">'
468  write(unit=unit_vtk,fmt=fr4p, iostat=e_io)(x(n1),n1=nx1,nx2)
469  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'</DataArray>'
470  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'<DataArray type="Float32" Name="Y" format="ascii">'
471  write(unit=unit_vtk,fmt=fr4p, iostat=e_io)(y(n1),n1=ny1,ny2)
472  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'</DataArray>'
473  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'<DataArray type="Float32" Name="Z" format="ascii">'
474  write(unit=unit_vtk,fmt=fr4p, iostat=e_io)(z(n1),n1=nz1,nz2)
475  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'</DataArray>'
476  indent = indent - 2
477  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'</Coordinates>'
478 
479  case(f_out_binary)
480  write(s_buffer,fmt='(A,6'//fi4p//',A)',iostat=e_io)repeat(' ',indent)//'<Piece Extent="',nx1,nx2,ny1,ny2,nz1,nz2,'">'
481  indent = indent + 2
482  write(unit=unit_vtk, iostat=e_io)trim(s_buffer)//end_rec
483  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<Coordinates>'//end_rec
484  indent = indent + 2
485  write(s_buffer,fmt='(I8)', iostat=e_io)ioffset
486  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<DataArray type="Float32" Name="X" format="appended" offset="',trim(s_buffer),'">'//end_rec
487  n_byte = (nx2-nx1+1)*sizeof(tipo_r4)
488  ioffset = ioffset + sizeof(tipo_i4) + n_byte
489  write(unit=unit_vtk_append, iostat=e_io)n_byte,'R4',nx2-nx1+1
490  write(unit=unit_vtk_append, iostat=e_io)(x(n1),n1=nx1,nx2)
491  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</DataArray>'//end_rec
492  write(s_buffer,fmt='(I8)', iostat=e_io)ioffset
493  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<DataArray type="Float32" Name="Y" format="appended" offset="',trim(s_buffer),'">'//end_rec
494  n_byte = (ny2-ny1+1)*sizeof(tipo_r4)
495  ioffset = ioffset + sizeof(tipo_i4) + n_byte
496  write(unit=unit_vtk_append, iostat=e_io)n_byte,'R4',ny2-ny1+1
497  write(unit=unit_vtk_append, iostat=e_io)(y(n1),n1=ny1,ny2)
498  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</DataArray>'//end_rec
499  write(s_buffer,fmt='(I8)', iostat=e_io)ioffset
500  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<DataArray type="Float32" Name="Z" format="appended" offset="',trim(s_buffer),'">'//end_rec
501  n_byte = (nz2-nz1+1)*sizeof(tipo_r4)
502  ioffset = ioffset + sizeof(tipo_i4) + n_byte
503  write(unit=unit_vtk_append, iostat=e_io)n_byte,'R4',nz2-nz1+1
504  write(unit=unit_vtk_append, iostat=e_io)(z(n1),n1=nz1,nz2)
505  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</DataArray>'//end_rec
506  indent = indent - 2
507  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</Coordinates>'//end_rec
508 
509  endselect
510 
511  return
512 !***********************************************************************************************************************************************************************************
513  endfunction vtk_geo_xml_rect_r4
514 !***********************************************************************************************************************************************************************************
515 
516 !
517 !
525 !***********************************************************************************************************************************************************************************
526  function vtk_geo_xml_unst_r4(NN,NC,X,Y,Z) result(E_IO)
527 !***********************************************************************************************************************************************************************************
528 
529  implicit none
530 
531  !--------------------------------------------------------------------------------------------------------------------------------
532  integer(I4P), intent(IN) :: nn ! number of nodes
533  integer(I4P), intent(IN) :: nc ! number of cells
534  real (R4P), intent(IN) :: x(1:nn) ! x coordinates
535  real (R4P), intent(IN) :: y(1:nn) ! y coordinates
536  real (R4P), intent(IN) :: z(1:nn) ! z coordinates
537 
538  integer(I4P) :: e_io ! Input/Output inquiring flag: $0$ if IO is done, $> 0$ if IO is not done
539  integer(I4P) :: n1 ! counter
540 
541  character(len=maxlen) :: s_buffer ! buffer string
542  !--------------------------------------------------------------------------------------------------------------------------------
543 
544  select case(f_out)
545 
546  case(f_out_ascii)
547  write(unit=unit_vtk,fmt='(A,'//fi4p//',A,'//fi4p//',A)',iostat=e_io)repeat(' ',indent)//'<Piece NumberOfPoints="',nn,'" NumberOfCells="',nc,'">'
548  indent = indent + 2
549  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'<Points>'
550  indent = indent + 2
551  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'<DataArray type="Float32" NumberOfComponents="3" Name="Point" format="ascii">'
552  write(unit=unit_vtk,fmt='(3'//fr4p//')', iostat=e_io)(x(n1),y(n1),z(n1),n1=1,nn)
553  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'</DataArray>'
554  indent = indent - 2
555  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'</Points>'
556 
557  case(f_out_binary)
558  write(s_buffer,fmt='(A,'//fi4p//',A,'//fi4p//',A)',iostat=e_io)repeat(' ',indent)//'<Piece NumberOfPoints="',nn,'" NumberOfCells="',nc,'">'
559  indent = indent + 2
560  write(unit=unit_vtk, iostat=e_io)trim(s_buffer)//end_rec
561  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<Points>'//end_rec
562  indent = indent + 2
563  write(s_buffer,fmt='(I8)', iostat=e_io)ioffset
564  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<DataArray type="Float32" NumberOfComponents="3" Name="Point" format="appended" offset="',trim(s_buffer),'">'//end_rec
565  n_byte = 3*nn*sizeof(tipo_r4)
566  ioffset = ioffset + sizeof(tipo_i4) + n_byte
567  write(unit=unit_vtk_append, iostat=e_io)n_byte,'R4',3*nn
568  write(unit=unit_vtk_append, iostat=e_io)(x(n1),y(n1),z(n1),n1=1,nn)
569  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</DataArray>'//end_rec
570  indent = indent - 2
571  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</Points>'//end_rec
572 
573  endselect
574 
575  return
576 !***********************************************************************************************************************************************************************************
577  endfunction vtk_geo_xml_unst_r4
578 !***********************************************************************************************************************************************************************************
579 
580 !
581 !
584 !***********************************************************************************************************************************************************************************
585  function vtk_geo_xml_closep() result(E_IO)
586 !***********************************************************************************************************************************************************************************
587 
588  implicit none
589 
590  !--------------------------------------------------------------------------------------------------------------------------------
591  integer(I4P) :: e_io ! Input/Output inquiring flag: $0$ if IO is done, $> 0$ if IO is not done
592  !--------------------------------------------------------------------------------------------------------------------------------
593 
594  indent = indent - 2
595 
596  select case(f_out)
597 
598  case(f_out_ascii)
599  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'</Piece>'
600 
601  case(f_out_binary)
602  write(unit=unit_vtk,iostat=e_io)repeat(' ',indent)//'</Piece>'//end_rec
603 
604  endselect
605 
606  return
607 !***********************************************************************************************************************************************************************************
608  endfunction vtk_geo_xml_closep
609 !***********************************************************************************************************************************************************************************
610 
611 !
612 !
619 !***********************************************************************************************************************************************************************************
620  function vtk_con_xml(NC,connect,offset,cell_type) result(E_IO)
621 !***********************************************************************************************************************************************************************************
622 
623  implicit none
624 
625  !--------------------------------------------------------------------------------------------------------------------------------
626  integer(I4P), intent(IN) :: nc ! number of cells
627  integer(I4P), intent(IN) :: connect(:) ! mesh connectivity
628  integer(I4P), intent(IN) :: offset(1:nc) ! cell offset
629  integer(I1P), intent(IN) :: cell_type(1:nc) ! VTK cell type
630 
631  integer(I4P) :: e_io ! Input/Output inquiring flag: $0$ if IO is done, $> 0$ if IO is not done
632  integer(I4P) :: n1 ! counter
633 
634  character(len=maxlen) :: s_buffer ! buffer string
635  !--------------------------------------------------------------------------------------------------------------------------------
636 
637  select case(f_out)
638 
639  case(f_out_ascii)
640  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'<Cells>'
641  indent = indent + 2
642  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'<DataArray type="Int32" Name="connectivity" format="ascii">'
643  write(unit=unit_vtk,fmt=fi4p, iostat=e_io)(connect(n1),n1=1,size(connect))
644  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'</DataArray>'
645  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'<DataArray type="Int32" Name="offsets" format="ascii">'
646  write(unit=unit_vtk,fmt=fi4p, iostat=e_io)(offset(n1),n1=1,nc)
647  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'</DataArray>'
648  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'<DataArray type="Int8" Name="types" format="ascii">'
649  write(unit=unit_vtk,fmt=fi1p, iostat=e_io)(cell_type(n1),n1=1,nc)
650  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'</DataArray>'
651  indent = indent - 2
652  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'</Cells>'
653 
654  case(f_out_binary)
655  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<Cells>'//end_rec
656  indent = indent + 2
657  write(s_buffer,fmt='(I10)',iostat=e_io)ioffset
658  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<DataArray type="Int32" Name="connectivity" format="appended" offset="',trim(s_buffer),'">'//end_rec
659  n_byte = size(connect)*sizeof(tipo_i4)
660  ioffset = ioffset + sizeof(tipo_i4) + n_byte
661  write(unit=unit_vtk_append,iostat=e_io)n_byte,'I4',size(connect)
662  write(unit=unit_vtk_append,iostat=e_io)(connect(n1),n1=1,size(connect))
663  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</DataArray>'//end_rec
664  write(s_buffer,fmt='(I10)',iostat=e_io)ioffset
665  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<DataArray type="Int32" Name="offsets" format="appended" offset="',trim(s_buffer),'">'//end_rec
666  n_byte = nc*sizeof(tipo_i4)
667  ioffset = ioffset + sizeof(tipo_i4) + n_byte
668  write(unit=unit_vtk_append,iostat=e_io)n_byte,'I4',nc
669  write(unit=unit_vtk_append,iostat=e_io)(offset(n1),n1=1,nc)
670  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</DataArray>'//end_rec
671  write(s_buffer,fmt='(I10)',iostat=e_io)ioffset
672  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<DataArray type="Int8" Name="types" format="appended" offset="',trim(s_buffer),'">'//end_rec
673  n_byte = nc*sizeof(tipo_i1)
674  ioffset = ioffset + sizeof(tipo_i4) + n_byte
675  write(unit=unit_vtk_append,iostat=e_io)n_byte,'I1',nc
676  write(unit=unit_vtk_append,iostat=e_io)(cell_type(n1),n1=1,nc)
677  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</DataArray>'//end_rec
678  indent = indent - 2
679  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</Cells>'//end_rec
680 
681  endselect
682 
683  return
684 !***********************************************************************************************************************************************************************************
685  endfunction vtk_con_xml
686 !***********************************************************************************************************************************************************************************
687 
688 !
689 !
694 !***********************************************************************************************************************************************************************************
695  function vtk_dat_xml(var_location,var_block_action) result(E_IO)
696 !***********************************************************************************************************************************************************************************
697 
698  implicit none
699 
700  !--------------------------------------------------------------------------------------------------------------------------------
701  character(*), intent(IN) :: var_location ! location of saving variables: CELL for cell-centered, NODE for node-centered
702  character(*), intent(IN) :: var_block_action ! variables block action: OPEN or CLOSE block
703  integer(I4P) :: e_io ! Input/Output inquiring flag: $0$ if IO is done, $> 0$ if IO is not done
704  !--------------------------------------------------------------------------------------------------------------------------------
705 
706  select case(f_out)
707 
708  case(f_out_ascii)
709  select case(trim(var_location))
710  case('CELL')
711  select case(trim(var_block_action))
712  case('OPEN')
713  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'<CellData>'
714  indent = indent + 2
715  case('CLOSE')
716  indent = indent - 2
717  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'</CellData>'
718  endselect
719  case('NODE')
720  select case(trim(var_block_action))
721  case('OPEN')
722  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'<PointData>'
723  indent = indent + 2
724  case('CLOSE')
725  indent = indent - 2
726  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'</PointData>'
727  endselect
728  endselect
729 
730  case(f_out_binary)
731  select case(trim(var_location))
732  case('CELL')
733  select case(trim(var_block_action))
734  case('OPEN')
735  write(unit=unit_vtk,iostat=e_io)repeat(' ',indent)//'<CellData>'//end_rec
736  indent = indent + 2
737  case('CLOSE')
738  indent = indent - 2
739  write(unit=unit_vtk,iostat=e_io)repeat(' ',indent)//'</CellData>'//end_rec
740  endselect
741  case('NODE')
742  select case(trim(var_block_action))
743  case('OPEN')
744  write(unit=unit_vtk,iostat=e_io)repeat(' ',indent)//'<PointData>'//end_rec
745  indent = indent + 2
746  case('CLOSE')
747  indent = indent - 2
748  write(unit=unit_vtk,iostat=e_io)repeat(' ',indent)//'</PointData>'//end_rec
749  endselect
750  endselect
751 
752  endselect
753 
754  return
755 !***********************************************************************************************************************************************************************************
756  endfunction vtk_dat_xml
757 !***********************************************************************************************************************************************************************************
758 
759 !
760 !
766 !***********************************************************************************************************************************************************************************
767  function vtk_var_xml_scal_r4(NC_NN,varname,var) result(E_IO)
768 !***********************************************************************************************************************************************************************************
769 
770  implicit none
771 
772  !--------------------------------------------------------------------------------------------------------------------------------
773  integer(I4P), intent(IN) :: nc_nn ! number of cells or nodes
774  character(*), intent(IN) :: varname ! variable name
775  real (R4P), intent(IN) :: var(1:nc_nn) ! variable to be saved
776 
777  integer(I4P) :: e_io ! Input/Output inquiring flag: $0$ if IO is done, $> 0$ if IO is not done
778  integer(I4P) :: n1 ! counter
779 
780  character(len=maxlen) :: s_buffer ! buffer string
781  !--------------------------------------------------------------------------------------------------------------------------------
782 
783  select case(f_out)
784 
785  case(f_out_ascii)
786  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'<DataArray type="Float32" Name="'//trim(varname)//'" NumberOfComponents="1" format="ascii">'
787  write(unit=unit_vtk,fmt=fr4p, iostat=e_io)var
788  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'</DataArray>'
789 
790  case(f_out_binary)
791  write(s_buffer,fmt='(I10)',iostat=e_io)ioffset
792  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<DataArray type="Float32" Name="'//trim(varname)//'" NumberOfComponents="1" format="appended" offset="',trim(s_buffer),'">'//end_rec
793  n_byte = nc_nn*sizeof(tipo_r4)
794  ioffset = ioffset + sizeof(tipo_i4) + n_byte
795  write(unit=unit_vtk_append,iostat=e_io)n_byte,'R4',nc_nn
796  write(unit=unit_vtk_append,iostat=e_io)(var(n1),n1=1,nc_nn)
797  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</DataArray>'//end_rec
798 
799  endselect
800 
801  return
802 !***********************************************************************************************************************************************************************************
803  endfunction vtk_var_xml_scal_r4
804 !***********************************************************************************************************************************************************************************
805 
806 !
807 !
813 !***********************************************************************************************************************************************************************************
814  function vtk_var_xml_scal_i4(NC_NN,varname,var) result(E_IO)
815 !***********************************************************************************************************************************************************************************
816 
817  implicit none
818 
819  !--------------------------------------------------------------------------------------------------------------------------------
820  integer(I4P), intent(IN) :: nc_nn ! number of cells or nodes
821  character(*), intent(IN) :: varname ! variable name
822  integer(I4P), intent(IN) :: var(1:nc_nn) ! variable to be saved
823 
824  integer(I4P) :: e_io ! Input/Output inquiring flag: $0$ if IO is done, $> 0$ if IO is not done
825  integer(I4P) :: n1 ! counter
826 
827  character(len=maxlen) :: s_buffer ! buffer string
828  !--------------------------------------------------------------------------------------------------------------------------------
829 
830  select case(f_out)
831 
832  case(f_out_ascii)
833  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'<DataArray type="Int32" Name="'//trim(varname)//'" NumberOfComponents="1" format="ascii">'
834  write(unit=unit_vtk,fmt=fi4p, iostat=e_io)var
835  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'</DataArray>'
836 
837  case(f_out_binary)
838  write(s_buffer,fmt='(I8)', iostat=e_io)ioffset
839  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<DataArray type="Int32" Name="'//trim(varname)//'" NumberOfComponents="1" format="appended" offset="',trim(s_buffer),'">'//end_rec
840  n_byte = nc_nn*sizeof(tipo_i4)
841  ioffset = ioffset + sizeof(tipo_i4) + n_byte
842  write(unit=unit_vtk_append,iostat=e_io)n_byte,'I4',nc_nn
843  write(unit=unit_vtk_append,iostat=e_io)(var(n1),n1=1,nc_nn)
844  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</DataArray>'//end_rec
845 
846  endselect
847 
848  return
849 !***********************************************************************************************************************************************************************************
850  endfunction vtk_var_xml_scal_i4
851 !***********************************************************************************************************************************************************************************
852 
853 !
854 !
862 !***********************************************************************************************************************************************************************************
863  function vtk_var_xml_vect_r4(NC_NN,varname,varX,varY,varZ) result(E_IO)
864 !***********************************************************************************************************************************************************************************
865 
866  implicit none
867 
868  !--------------------------------------------------------------------------------------------------------------------------------
869  integer(I4P), intent(IN) :: nc_nn ! number of cells or nodes
870  character(*), intent(IN) :: varname ! variable name
871  real (R4P), intent(IN) :: varx(1:nc_nn) ! x component
872  real (R4P), intent(IN) :: vary(1:nc_nn) ! y component
873  real (R4P), intent(IN) :: varz(1:nc_nn) ! z component
874 
875  integer(I4P) :: e_io ! Input/Output inquiring flag: $0$ if IO is done, $> 0$ if IO is not done
876  integer(I4P) :: n1 ! counter
877 
878  character(len=maxlen) :: s_buffer ! buffer string
879  !--------------------------------------------------------------------------------------------------------------------------------
880 
881  select case(f_out)
882 
883  case(f_out_ascii)
884  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'<DataArray type="Float32" Name="'//trim(varname)//'" NumberOfComponents="3" format="ascii">'
885  write(unit=unit_vtk,fmt='(3'//fr4p//')',iostat=e_io)(varx(n1),vary(n1),varz(n1),n1=1,nc_nn)
886  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'</DataArray>'
887 
888  case(f_out_binary)
889  write(s_buffer,fmt='(I8)', iostat=e_io)ioffset
890  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<DataArray type="Float32" Name="'//trim(varname)//'" NumberOfComponents="3" format="appended" offset="',trim(s_buffer),'">'//end_rec
891  n_byte = 3*nc_nn*sizeof(tipo_r4)
892  ioffset = ioffset + sizeof(tipo_i4) + n_byte
893  write(unit=unit_vtk_append,iostat=e_io)n_byte,'R4',3*nc_nn
894  write(unit=unit_vtk_append,iostat=e_io)(varx(n1),vary(n1),varz(n1),n1=1,nc_nn)
895  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</DataArray>'//end_rec
896 
897  endselect
898 
899  return
900 !***********************************************************************************************************************************************************************************
901  endfunction vtk_var_xml_vect_r4
902 !***********************************************************************************************************************************************************************************
903 
904 !
905 !
913 !***********************************************************************************************************************************************************************************
914  function vtk_var_xml_vect_i4(NC_NN,varname,varX,varY,varZ) result(E_IO)
915 !***********************************************************************************************************************************************************************************
916 
917  implicit none
918 
919  !--------------------------------------------------------------------------------------------------------------------------------
920  integer(I4P), intent(IN) :: nc_nn ! number of cells or nodes
921  character(*), intent(IN) :: varname ! variable name
922  integer(I4P), intent(IN) :: varx(1:nc_nn) ! x component
923  integer(I4P), intent(IN) :: vary(1:nc_nn) ! y component
924  integer(I4P), intent(IN) :: varz(1:nc_nn) ! z component
925 
926  integer(I4P) :: e_io ! Input/Output inquiring flag: $0$ if IO is done, $> 0$ if IO is not done
927  integer(I4P) :: n1 ! counter
928 
929  character(len=maxlen) :: s_buffer ! buffer string
930  !--------------------------------------------------------------------------------------------------------------------------------
931 
932  select case(f_out)
933 
934  case(f_out_ascii)
935  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'<DataArray type="Int32" Name="'//trim(varname)//'" NumberOfComponents="3" format="ascii">'
936  write(unit=unit_vtk,fmt='(3'//fi4p//')',iostat=e_io)(varx(n1),vary(n1),varz(n1),n1=1,nc_nn)
937  write(unit=unit_vtk,fmt='(A)', iostat=e_io)repeat(' ',indent)//'</DataArray>'
938 
939  case(f_out_binary)
940  write(s_buffer,fmt='(I8)', iostat=e_io)ioffset
941  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'<DataArray type="Int32" Name="'//trim(varname)//'" NumberOfComponents="3" format="appended" offset="',trim(s_buffer),'">'//end_rec
942  n_byte = 3*nc_nn*sizeof(tipo_i4)
943  ioffset = ioffset + sizeof(tipo_i4) + n_byte
944  write(unit=unit_vtk_append,iostat=e_io)n_byte,'I4',3*nc_nn
945  write(unit=unit_vtk_append,iostat=e_io)(varx(n1),vary(n1),varz(n1),n1=1,nc_nn)
946  write(unit=unit_vtk, iostat=e_io)repeat(' ',indent)//'</DataArray>'//end_rec
947 
948  endselect
949 
950  return
951 !***********************************************************************************************************************************************************************************
952  endfunction vtk_var_xml_vect_i4
953 !***********************************************************************************************************************************************************************************
954 
955 !
956 !
959 !***********************************************************************************************************************************************************************************
960  function vtk_end_xml() result(E_IO)
961 !***********************************************************************************************************************************************************************************
962 
963  implicit none
964 
965  !--------------------------------------------------------------------------------------------------------------------------------
966  integer(I4P) :: e_io ! Input/Output inquiring flag: $0$ if IO is done, $> 0$ if IO is not done
967  character(2) :: var_type ! var\_type = R4,I4
968  real (R4P), allocatable :: v_r4(:) ! R4 vector for IO in AppendData
969  integer(I4P), allocatable :: v_i4(:) ! I4 vector for IO in AppendData
970  integer(I1P), allocatable :: v_i1(:) ! I4 vector for IO in AppendData
971  integer(I4P) :: n_v ! vector dimension
972  integer(I4P) :: n1 ! counter
973  !--------------------------------------------------------------------------------------------------------------------------------
974 
975  select case(f_out)
976 
977  case(f_out_ascii)
978 
979  indent = indent - 2
980  write(unit=unit_vtk,fmt='(A)',iostat=e_io)repeat(' ',indent)//'</'//trim(topology)//'>'
981  write(unit=unit_vtk,fmt='(A)',iostat=e_io)'</VTKFile>'
982 
983  case(f_out_binary)
984 
985  indent = indent - 2
986  write(unit =unit_vtk, iostat=e_io)repeat(' ',indent)//'</'//trim(topology)//'>'//end_rec
987  write(unit =unit_vtk, iostat=e_io)repeat(' ',indent)//'<AppendedData encoding="raw">'//end_rec
988  write(unit =unit_vtk, iostat=e_io)'_'
989  endfile(unit=unit_vtk_append,iostat=e_io)
990  rewind(unit =unit_vtk_append,iostat=e_io)
991 
992  do
993 
994  read(unit=unit_vtk_append,iostat=e_io,end=100)n_byte,var_type,n_v
995 
996  select case(var_type)
997 
998  case('R4')
999  allocate(v_r4(1:n_v))
1000  read (unit=unit_vtk_append,iostat=e_io)(v_r4(n1),n1=1,n_v)
1001  write(unit=unit_vtk, iostat=e_io)n_byte,(v_r4(n1),n1=1,n_v)
1002  deallocate(v_r4)
1003 
1004  case('I4')
1005  allocate(v_i4(1:n_v))
1006  read (unit=unit_vtk_append,iostat=e_io)(v_i4(n1),n1=1,n_v)
1007  write(unit=unit_vtk, iostat=e_io)n_byte,(v_i4(n1),n1=1,n_v)
1008  deallocate(v_i4)
1009 
1010  case('I1')
1011  allocate(v_i1(1:n_v))
1012  read (unit=unit_vtk_append,iostat=e_io)(v_i1(n1),n1=1,n_v)
1013  write(unit=unit_vtk, iostat=e_io)n_byte,(v_i1(n1),n1=1,n_v)
1014  deallocate(v_i1)
1015 
1016  endselect
1017 
1018  enddo
1019 
1020  100 continue
1021  write(unit=unit_vtk,iostat=e_io)end_rec
1022  write(unit=unit_vtk,iostat=e_io)repeat(' ',indent)//'</AppendedData>'//end_rec
1023  write(unit=unit_vtk,iostat=e_io)'</VTKFile>'//end_rec
1024 
1025  !->closing AppendData file
1026  close(unit=unit_vtk_append,iostat=e_io)
1027 
1028  endselect
1029 
1030  close(unit=unit_vtk,iostat=e_io)
1031 
1032  return
1033 !***********************************************************************************************************************************************************************************
1034  endfunction vtk_end_xml
1035 !***********************************************************************************************************************************************************************************
1036 
1037 end module mod_vtk_io
integer(i4p) function, public vtk_end_xml()
The VTK_END_XML function finalizes opened files.
integer(i4p) function vtk_geo_xml_unst_r4(NN, NC, X, Y, Z)
The VTK_GEO_XML_UNST_R4 function is used for saving mesh with a VTK UnstructuredGrid R4P topology...
integer(i4p) function, public vtk_con_xml(NC, connect, offset, cell_type)
The VTK_CON_XML function must be called when unstructured grid topology is used. It saves the connect...
integer(i4p) function, public vtk_ini_xml(output_format, filename, mesh_topology, nx1, nx2, ny1, ny2, nz1, nz2)
The VTK_INI_XML function is used for initializing file. This function must be the first to be called...
integer(i4p) function vtk_var_xml_scal_i4(NC_NN, varname, var)
The VTK_VAR_XML_SCAL_I4 function saves I4P scalar variables.
overloading of VTK_GEO_XML
This module defines all global variables of EFISPEC3D. Scalar variables are initialized directly in t...
integer(i4p) function vtk_var_xml_scal_r4(NC_NN, varname, var)
The VTK_VAR_XML_SCAL_R4 function saves R4P scalar variables.
integer(i4p) function, public vtk_dat_xml(var_location, var_block_action)
The VTK_DAT_XML function opens or closes CellData/PointData structures.
integer(i4p) function vtk_geo_xml_rect_r4(nx1, nx2, ny1, ny2, nz1, nz2, X, Y, Z)
The VTK_GEO_XML_RECT_R4 function is used for saving mesh with a VTK RectilinearGrid R4P topology...
This module programmed by S. Zaghi (https://github.com/szaghi/Lib_VTK_IO) contains functions to write...
integer(i4p) function vtk_var_xml_vect_r4(NC_NN, varname, varX, varY, varZ)
The VTK_VAR_XML_VECT_R4 function saves R4P vectorial variables.
overloading of VTK_VAR_XML
integer(i4p) function vtk_geo_xml_closep()
The VTK_GEO_XML_CLOSE function is used for closing mesh block data.
integer(i4p) function vtk_geo_xml_strg_r4(nx1, nx2, ny1, ny2, nz1, nz2, NN, X, Y, Z)
The VTK_GEO_XML_STRG_R4 function is used for saving mesh with a VTK StructuredGrid R4P topology...
integer(i4p) function vtk_var_xml_vect_i4(NC_NN, varname, varX, varY, varZ)
The VTK_VAR_XML_VECT_I4 function saves I4P vectorial variables.