SPEED
WRITE_VTK_VOLUME_TIMESERIES.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine write_vtk_volume_timeseries (its, mpi_id, nn_loc, loc_n_num, nm, sdeg, xx_loc, yy_loc, zz_loc, cs_nnz_loc, cs_loc, u2, stress)
 ...Writing VTK file to visualise in Paraview
 

Function/Subroutine Documentation

◆ write_vtk_volume_timeseries()

subroutine write_vtk_volume_timeseries ( integer*4  its,
integer*4  mpi_id,
integer*4  nn_loc,
integer*4, dimension(nn_loc)  loc_n_num,
integer*4  nm,
integer*4, dimension(nm)  sdeg,
real*8, dimension(nn_loc)  xx_loc,
real*8, dimension(nn_loc)  yy_loc,
real*8, dimension(nn_loc)  zz_loc,
integer*4  cs_nnz_loc,
integer*4, dimension(0:cs_nnz_loc)  cs_loc,
real*8, dimension(nn_loc*3)  u2,
real*8, dimension(nn_loc*6)  stress 
)

...Writing VTK file to visualise in Paraview

Author
Srihari Sangaraju
Date
July, 2021
Version
1.0
Parameters
[in]loc_n_num.Global node number of 'i'th local node is loc_n_num(i)
[in]nn_loc.No. of nodes in Local/Current Partition
[in]nmat_nheNo. of Blocks specified with NHE case
[in]nhe_matTag/Labels of Blocks where NHE has to be implemented
[in]xs_locx-coordinate of spectral nodes
[in]ys_locy-coordinate of spectral nodes
[in]zs_locz-coordinate of spectral nodes

Definition at line 33 of file WRITE_VTK_VOLUME_TIMESERIES.f90.

35
37
38 implicit none
39
40 include 'SPEED.MPI'
41
42 integer*4 :: nm, nel_loc, nn_loc, cs_nnz_loc, mpi_id, its
43 integer*4 :: ie, inode, nn, im
44 integer*4 :: unit_mpi
45 integer*4 :: ic1, ic2, ic3, ic4, ic5, ic6, ic7, ic8
46
47 integer*4, dimension(nm) :: sdeg
48 integer*4, dimension(nn_loc) :: loc_n_num
49 integer*4, dimension(0:cs_nnz_loc) :: cs_loc
50
51 real*8, dimension(nn_loc) :: xx_loc, yy_loc, zz_loc
52 real*8, dimension(nn_loc*3) :: u2
53 real*8, dimension(nn_loc*6) :: stress
54
55 character*70 :: file_nhe_proc, file_nhe_new, mpi_file
56
57 nel_loc = cs_loc(0) - 1
58
59 mpi_file = 'MONITORS'
60 unit_mpi = 2500 + mpi_id
61 write(file_nhe_proc,'(A10,I5.5,A1,I8.8,A4)') 'DIS_X_PROC', mpi_id, '_', its, '.vtk'
62
63 if(len_trim(mpi_file) .ne. 70) then
64 file_nhe_new = mpi_file(1:len_trim(mpi_file)) // '/' // file_nhe_proc
65 else
66 file_nhe_new = file_nhe_proc
67 endif
68
69 !----------------------------------------------------------------------
70 open(unit_mpi,file=file_nhe_new,status='replace')
71 write(unit_mpi,'(a)') '# vtk DataFile Version 3.1'
72 write(unit_mpi,'(a)') 'material model VTK file'
73 write(unit_mpi,'(a)') 'ASCII'
74 write(unit_mpi,'(a)') 'DATASET UNSTRUCTURED_GRID'
75 write(unit_mpi, '(a,i12,a)') 'POINTS ', nn_loc, ' float'
76
77 ! Node Coordinates
78 do inode=1,nn_loc
79 write(unit_mpi,'(3e20.12)') xx_loc(inode),yy_loc(inode),zz_loc(inode)
80 enddo
81 write(unit_mpi,*) ''
82
83 ! Connectivity (note: node indices for vtk start at 0)
84 write(unit_mpi,'(a,i12,i12)') "CELLS ",nel_loc,nel_loc*9
85 do ie=1,nel_loc
86 im = cs_loc(cs_loc(ie -1) + 0 )
87 nn = sdeg(im) +1
88
89 ic1 = cs_loc(cs_loc(ie -1) + nn*nn*(1 -1) + nn*(1 -1) + 1) - 1
90 ic2 = cs_loc(cs_loc(ie -1) + nn*nn*(1 -1) + nn*(1 -1) + nn) - 1
91 ic3 = cs_loc(cs_loc(ie -1) + nn*nn*(1 -1) + nn*(nn -1) + nn) - 1
92 ic4 = cs_loc(cs_loc(ie -1) + nn*nn*(1 -1) + nn*(nn -1) + 1) - 1
93 ic5 = cs_loc(cs_loc(ie -1) + nn*nn*(nn -1) + nn*(1 -1) + 1) - 1
94 ic6 = cs_loc(cs_loc(ie -1) + nn*nn*(nn -1) + nn*(1 -1) + nn) - 1
95 ic7 = cs_loc(cs_loc(ie -1) + nn*nn*(nn -1) + nn*(nn -1) + nn) - 1
96 ic8 = cs_loc(cs_loc(ie -1) + nn*nn*(nn -1) + nn*(nn -1) + 1) - 1
97
98 write(unit_mpi,'(9i12)') 8, ic1, ic2, ic3, ic4, ic5, ic6, ic7, ic8
99 enddo
100 write(unit_mpi,*) ''
101
102 ! vtkCellType: hexahedrons (ID = 12)
103 write(unit_mpi,'(a,i12)') "CELL_TYPES ",nel_loc
104 write(unit_mpi,'(6i12)') (12,ie=1,nel_loc)
105 write(unit_mpi,*) ''
106
107 ! Writing Scalar data
108 write(unit_mpi,'(a,i12)') "POINT_DATA ",nn_loc
109
110 ! UX
111 write(unit_mpi,'(a)') "SCALARS stress_xy float"
112 write(unit_mpi,'(a)') "LOOKUP_TABLE default"
113 do inode = 1,nn_loc
114 !write(unit_mpi,*) U2(3*(inode-1)+1)
115 write(unit_mpi,*) stress(6*(inode -1) +4)
116 enddo
117 write(unit_mpi,*) ''
118
119
120 close(unit_mpi)
121 !------------------------------------------------------------------
122
123 ! if (mpi_id.eq.0) write(*,'(A)')'Completed.'
124
SPEED exit codes.
Definition MODULES.f90:25