SPEED
WRITE_VTK_MESH.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine write_vtk_mesh (nn_loc, loc_n_num, nm, tag_mat, prop_mat, sdeg, xx_loc, yy_loc, zz_loc, cs_nnz_loc, cs_loc, mpi_id)
 ...Writing VTK file to visualise in Paraview
 

Function/Subroutine Documentation

◆ write_vtk_mesh()

subroutine write_vtk_mesh ( integer*4  nn_loc,
integer*4, dimension(nn_loc)  loc_n_num,
integer*4  nm,
integer*4, dimension(nm)  tag_mat,
real*8, dimension(nm,4)  prop_mat,
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,
integer*4  mpi_id 
)

...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_MESH.f90.

36
38
39 implicit none
40
41 include 'SPEED.MPI'
42
43 integer*4 :: nm, nel_loc, nn_loc, cs_nnz_loc, mpi_id
44 integer*4 :: ie, inode, nn, im
45 integer*4 :: unit_mpi
46 integer*4 :: ic1, ic2, ic3, ic4, ic5, ic6, ic7, ic8
47
48 real*8 :: lambda, mu, qs, qp, dum, gamma
49
50 integer*4, dimension(nm) :: sdeg, tag_mat
51 integer*4, dimension(nn_loc) :: loc_n_num, mat_id
52 integer*4, dimension(0:cs_nnz_loc) :: cs_loc
53
54 real*8, dimension(nm,4) :: prop_mat
55 real*8, dimension(nn_loc) :: xx_loc, yy_loc, zz_loc
56 real*8, dimension(nn_loc) :: rho, vs, vp, poisn
57
58 character*70 :: file_nhe_proc, file_nhe_new, mpi_file
59
60 rho = 0.d0
61 lambda = 0.d0
62 mu = 0.d0
63 qs = 0.d0
64 qp = 0.d0
65
66 vs = 0.d0
67 vp = 0.d0
68
69 nel_loc = cs_loc(0) - 1
70 !Tagging Each node to blovkID in *.mate file
71 do ie=1,nel_loc
72 im = cs_loc(cs_loc(ie -1) + 0 )
73 mat_id(cs_loc(cs_loc(ie-1) + 1 : cs_loc(ie)-1)) = im
74 enddo
75
76 do inode = 1,nn_loc
77
78 rho(inode) = prop_mat(mat_id(inode),1);
79 lambda = prop_mat(mat_id(inode),2);
80 mu = prop_mat(mat_id(inode),3);
81 gamma = prop_mat(mat_id(inode),4)
82
83 vs(inode) = sqrt(mu/rho(inode))
84 vp(inode) = sqrt((lambda/rho(inode)) + (2.d0*vs(inode)**2.d0))
85 dum = (vp(inode)/vs(inode))**2.d0
86 poisn(inode) = 0.5d0*(dum-2)/(dum-1)
87 enddo
88
89
90 if (mpi_id.eq.0) write(*,'(A)')
91 if (mpi_id.eq.0) write(*,'(A)')'------Writing VTK file - Mechanical Properties----------'
92
93 mpi_file = 'MONITORS'
94 unit_mpi = 2500 + mpi_id
95 ! write(file_nhe_proc,'(A6,I6.6,A4)') 'NHMECH', mpi_id, '.vtk'
96 write(file_nhe_proc,'(A10,I5.5,A4)') 'DIS_X_PROC', mpi_id, '.vtk'
97
98 if(len_trim(mpi_file) .ne. 70) then
99 file_nhe_new = mpi_file(1:len_trim(mpi_file)) // '/' // file_nhe_proc
100 else
101 file_nhe_new = file_nhe_proc
102 endif
103
104 !----------------------------------------------------------------------
105 open(unit_mpi,file=file_nhe_new,status='replace')
106 write(unit_mpi,'(a)') '# vtk DataFile Version 3.1'
107 write(unit_mpi,'(a)') 'material model VTK file'
108 write(unit_mpi,'(a)') 'ASCII'
109 write(unit_mpi,'(a)') 'DATASET UNSTRUCTURED_GRID'
110 write(unit_mpi, '(a,i12,a)') 'POINTS ', nn_loc, ' float'
111
112 ! Node Coordinates
113 do inode=1,nn_loc
114 write(unit_mpi,'(3e20.12)') xx_loc(inode),yy_loc(inode),zz_loc(inode)
115 enddo
116 write(unit_mpi,*) ''
117
118 ! Connectivity (note: node indices for vtk start at 0)
119 write(unit_mpi,'(a,i12,i12)') "CELLS ",nel_loc,nel_loc*9
120 do ie=1,nel_loc
121 im = cs_loc(cs_loc(ie -1) + 0 )
122 nn = sdeg(im) +1
123
124 ic1 = cs_loc(cs_loc(ie -1) + nn*nn*(1 -1) + nn*(1 -1) + 1) - 1
125 ic2 = cs_loc(cs_loc(ie -1) + nn*nn*(1 -1) + nn*(1 -1) + nn) - 1
126 ic3 = cs_loc(cs_loc(ie -1) + nn*nn*(1 -1) + nn*(nn -1) + nn) - 1
127 ic4 = cs_loc(cs_loc(ie -1) + nn*nn*(1 -1) + nn*(nn -1) + 1) - 1
128 ic5 = cs_loc(cs_loc(ie -1) + nn*nn*(nn -1) + nn*(1 -1) + 1) - 1
129 ic6 = cs_loc(cs_loc(ie -1) + nn*nn*(nn -1) + nn*(1 -1) + nn) - 1
130 ic7 = cs_loc(cs_loc(ie -1) + nn*nn*(nn -1) + nn*(nn -1) + nn) - 1
131 ic8 = cs_loc(cs_loc(ie -1) + nn*nn*(nn -1) + nn*(nn -1) + 1) - 1
132
133 write(unit_mpi,'(9i12)') 8, ic1, ic2, ic3, ic4, ic5, ic6, ic7, ic8
134 enddo
135 write(unit_mpi,*) ''
136
137 ! vtkCellType: hexahedrons (ID = 12)
138 write(unit_mpi,'(a,i12)') "CELL_TYPES ",nel_loc
139 write(unit_mpi,'(6i12)') (12,ie=1,nel_loc)
140 write(unit_mpi,*) ''
141
142 ! ! Writing Scalar data
143 ! write(unit_mpi,'(a,i12)') "POINT_DATA ",nn_loc
144
145 ! ! RHO
146 ! write(unit_mpi,'(a)') "SCALARS Rho float"
147 ! write(unit_mpi,'(a)') "LOOKUP_TABLE default"
148 ! do inode = 1,nn_loc
149 ! write(unit_mpi,*) rho(inode)
150 ! enddo
151 ! write(unit_mpi,*) ''
152
153 ! ! VS
154 ! write(unit_mpi,'(a)') "SCALARS VS float"
155 ! write(unit_mpi,'(a)') "LOOKUP_TABLE default"
156 ! do inode = 1,nn_loc
157 ! write(unit_mpi,*) VS(inode)
158 ! enddo
159 ! write(unit_mpi,*) ''
160
161 ! ! VP
162 ! ! write(unit_mpi,'(a)') "SCALARS VP float"
163 ! write(unit_mpi,'(a)') "SCALARS Thick float"
164 ! write(unit_mpi,'(a)') "LOOKUP_TABLE default"
165 ! do inode = 1,nn_loc
166 ! ! write(unit_mpi,*) VP(inode)
167 ! write(unit_mpi,*) thick_nodes(inode)
168 ! enddo
169 ! write(unit_mpi,*) ''
170
171 ! ! Poissons Ratio
172 ! ! write(unit_mpi,'(a)') "SCALARS Poisson float"
173 ! write(unit_mpi,'(a)') "SCALARS depth float"
174 ! write(unit_mpi,'(a)') "LOOKUP_TABLE default"
175 ! do inode = 1,nn_loc
176 ! ! write(unit_mpi,*) poisn(inode)
177 ! write(unit_mpi,*) zs_elev(inode)
178 ! enddo
179 ! write(unit_mpi,*) ''
180
181 close(unit_mpi)
182 !------------------------------------------------------------------
183
184 if (mpi_id.eq.0) write(*,'(A)')'Completed.'
185
SPEED exit codes.
Definition MODULES.f90:25