support for non QUAD cell type
This commit is contained in:
parent
0cf27a8229
commit
bd4db68a12
|
@ -880,14 +880,33 @@ class DADF5():
|
|||
else:
|
||||
|
||||
nodes = vtk.vtkPoints()
|
||||
with h5py.File(self.fname) as f:
|
||||
with h5py.File(self.fname,'r') as f:
|
||||
nodes.SetData(numpy_support.numpy_to_vtk(f['/geometry/x_n'][()],deep=True))
|
||||
|
||||
vtk_geom = vtk.vtkUnstructuredGrid()
|
||||
vtk_geom.SetPoints(nodes)
|
||||
vtk_geom.Allocate(f['/geometry/T_c'].shape[0])
|
||||
|
||||
if self.version_major == 0 and self.version_minor <= 5:
|
||||
vtk_type = vtk.VTK_HEXAHEDRON
|
||||
n_nodes = 8
|
||||
else:
|
||||
if f['/geometry/T_c'].attrs['VTK_TYPE'] == b'TRIANGLE':
|
||||
vtk_type = vtk.VTK_TRIANGLE
|
||||
n_nodes = 3
|
||||
elif f['/geometry/T_c'].attrs['VTK_TYPE'] == b'QUAD':
|
||||
vtk_type = vtk.VTK_QUAD
|
||||
n_nodes = 4
|
||||
elif f['/geometry/T_c'].attrs['VTK_TYPE'] == b'TETRA': # not tested
|
||||
vtk_type = vtk.VTK_TETRA
|
||||
n_nodes = 4
|
||||
elif f['/geometry/T_c'].attrs['VTK_TYPE'] == b'HEXAHEDRON':
|
||||
vtk_type = vtk.VTK_HEXAHEDRON
|
||||
n_nodes = 8
|
||||
|
||||
for i in f['/geometry/T_c']:
|
||||
vtk_geom.InsertNextCell(vtk.VTK_HEXAHEDRON,8,i-1) # not for all elements!
|
||||
vtk_geom.InsertNextCell(n_nodes,vtk_type,i-1)
|
||||
|
||||
elif mode == 'Point':
|
||||
Points = vtk.vtkPoints()
|
||||
Vertices = vtk.vtkCellArray()
|
||||
|
|
|
@ -104,7 +104,7 @@ subroutine mesh_init(ip,el)
|
|||
ip_reshaped,&
|
||||
node0_cell)
|
||||
|
||||
call writeGeometry(0,connectivity_elem,&
|
||||
call writeGeometry(elem,connectivity_elem,&
|
||||
reshape(connectivity_cell,[elem%NcellNodesPerCell,elem%nIPs*nElems]),&
|
||||
node0_cell,ip_reshaped)
|
||||
|
||||
|
@ -122,11 +122,12 @@ end subroutine mesh_init
|
|||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief Write all information needed for the DADF5 geometry
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
subroutine writeGeometry(elemType, &
|
||||
subroutine writeGeometry(elem, &
|
||||
connectivity_elem,connectivity_cell, &
|
||||
coordinates_nodes,coordinates_points)
|
||||
|
||||
integer, intent(in) :: elemType
|
||||
type(tElement), intent(in) :: &
|
||||
elem
|
||||
integer, dimension(:,:), intent(in) :: &
|
||||
connectivity_elem, &
|
||||
connectivity_cell
|
||||
|
@ -149,6 +150,7 @@ subroutine writeGeometry(elemType, &
|
|||
connectivity_temp = connectivity_cell
|
||||
call results_writeDataset('geometry',connectivity_temp,'T_c', &
|
||||
'connectivity of the cells','-')
|
||||
call results_addAttribute('VTK_TYPE',elem%vtkType,'geometry/T_c')
|
||||
|
||||
coordinates_temp = coordinates_nodes
|
||||
call results_writeDataset('geometry',coordinates_temp,'x_n', &
|
||||
|
|
|
@ -70,7 +70,7 @@ subroutine results_init
|
|||
|
||||
resultsFile = HDF5_openFile(trim(getSolverJobName())//'.hdf5','w',.true.)
|
||||
call results_addAttribute('DADF5_version_major',0)
|
||||
call results_addAttribute('DADF5_version_minor',5)
|
||||
call results_addAttribute('DADF5_version_minor',6)
|
||||
call results_addAttribute('DAMASK_version',DAMASKVERSION)
|
||||
call get_command(commandLine)
|
||||
call results_addAttribute('call',trim(commandLine))
|
||||
|
|
Loading…
Reference in New Issue