-
Notifications
You must be signed in to change notification settings - Fork 12
3.1 Input data format
We use a standardised and human-readable ASCII format as input for CrackPy. This has the advantage that output from almost any data source can be easily converted. Therefore, it is advantageous for academic or research use. Each input set consists of two files:
- nodemap.txt
- connections.txt (optional)
The nodemap file stores the nodal or facet data. This file is crucial as all processes are based on it. Metadata or comments can be stored in the header, identified by the symbol #
. The data section includes individual columns containing nodal coordinates, displacement vector, and strain tensor defined for plane strain conditions.
# Process data:
# specimen : Dummy2
####################################################################################################
# SIGNALS:
# force : value_element : 14985.6119156
# displacement : value_element : -9.17953395844
#
# Kraft : analog_input : 15102.5390625
# Verschiebung : analog_input : -7.6751708984375
####################################################################################################
# ID; x_undef [mm]; y_undef [mm]; z_undef [mm]; u [mm]; v [mm]; w [mm]; epsx [%]; epsy [%]; epsxy [1]; epseqv [%]
1; -80.0564567626; 84.9779203814; 0.1509636527; 0.003655774286017; 0.059203740209341; -0.158300071954727; -0.185958549380302; -0.112336091697216; -0.000613267708104; 0.310294896364212
2; -79.2904753387; 85.0053386055; 0.1975326408; 0.001983210910112; 0.060109287500381; -0.155325025320053; -0.209498479962349; -0.005866464227438; -0.000913561263587; 0.267783910036087
Where:
-
ID
is the node number. -
x_undef
,y_undef
,z_undef
are the nodal cooridnates in the undeformed reference coordinate system, given in milimeters. -
u
,v
,w
are the displacements in x, y and z direction, given in milimeters. -
epsx
,epsy
,epsxy
are the componentens of the Cauchy's strain tensor. Note that$\varepsilon_{\mathrm{x}}$ and$\varepsilon_{\mathrm{y}}$ are given in %, while$\varepsilon_{\mathrm{xy}}$ is given without a ratio! -
epseqv
is the equivalent von Mises strain and is optional as it can be computed in CrackPy. - The separator is a semicolon
;
The connection file is optional and is mainly used for plotting. It stores how the nodes within the elements of the mesh are connected to each other. This file is based in the vtk
format. Within CrackPy we use pyvista.UnstructuredGrid for precessing the data.
Type; Element #; Node 1; Node 2; Node 3
3; 1; 108; 1; 0
3; 2; 109; 2; 1
3; 3; 110; 3; 2
For export of data from Zeiss GOM ARAMIS see: ARAMIS Exporter
We developed the following export script for Ansys Mechanical using PyMAPDL. See fe_model.py for the full code.
Caution! The shear strains (epto_xy
by 0.5.
def export_nodemap(self, filename) -> None:
"""Saves the current Ansys node data into a txt-file called 'Nodemap'.
Args:
"filename" (str): Name of the nodemap file.
"""
# Change to post-processing
self.mapdl.post1()
self.mapdl.set('last')
self.mapdl.allsel()
# Get coordinates, displacements, strains, and stresses
node_list = self.mapdl.get_array(entity='NODE', item1='NLIST')
node_loc_x = self.mapdl.get_array(entity='NODE', item1='LOC', it1num='X')
node_loc_y = self.mapdl.get_array(entity='NODE', item1='LOC', it1num='Y')
node_loc_z = self.mapdl.get_array(entity='NODE', item1='LOC', it1num='Z')
u_x = self.mapdl.get_array(entity='NODE', item1='U', it1num='X')
u_y = self.mapdl.get_array(entity='NODE', item1='U', it1num='Y')
u_z = self.mapdl.get_array(entity='NODE', item1='U', it1num='Z')
epto_xx = self.mapdl.get_array(entity='NODE', item1='EPTO', it1num='X') * 100.0
epto_yy = self.mapdl.get_array(entity='NODE', item1='EPTO', it1num='Y') * 100.0
epto_xy = self.mapdl.get_array(entity='NODE', item1='EPTO', it1num='XY') * 0.5
epto_eqv = self.mapdl.get_array(entity='NODE', item1='EPTO', it1num='EQV') * 100.0
s_xx = self.mapdl.get_array(entity='NODE', item1='S', it1num='X')
s_yy = self.mapdl.get_array(entity='NODE', item1='S', it1num='Y')
s_xy = self.mapdl.get_array(entity='NODE', item1='S', it1num='XY')
s_eqv = self.mapdl.get_array(entity='NODE', item1='S', it1num='EQV')
export_data = np.asarray([node_list, node_loc_x, node_loc_y, node_loc_z, u_x, u_y, u_z,
epto_xx, epto_yy, epto_xy, epto_eqv, s_xx, s_yy, s_xy, s_eqv]).T
with open(os.path.join(self.OUTPUT_PATH, f'{filename}_nodemap.txt'), 'w') as file:
file.write('# Header for nodemap file\n')
current_time = datetime.now().strftime("%d/%m/%Y %H:%M:%S")
file.write(f'# {"date":>16}: {current_time}\n')
file.write(f'# \n')
file.write(f'# {"index":>10}; '
f'{"x_undef":>12s}; {"y_undef":>12s}; {"z_undef":>12s}; '
f'{"ux":>12s}; {"uy":>12s}; {"uz":>12s}; '
f'{"eps_x":>12s}; {"eps_y":>12s}; {"eps_xy":>12s}; {"eps_eqv":>12s}; '
f'{"s_x":>12s}; {"s_y":>12s}; {"s_xy":>12s}; {"s_eqv":>12s}\n')
for i, line in enumerate(export_data):
file.write(f'{line[0]:12.1f}; {line[1]:12.6f}; {line[2]:12.6f}; {line[3]:12.6f}; '
f'{line[4]:12.6f}; {line[5]:12.6f}; {line[6]:12.6f}; '
f'{line[7]:12.6f}; {line[8]:12.6f}; {line[9]:12.6f}; {line[10]:12.6f}; '
f'{line[11]:12.6f}; {line[12]:12.6f}; {line[13]:12.6f}; {line[14]:12.6f}\n')
def export_vtk(self, filename):
"""Saves the current Ansys results as vtk file.
Args:
"filename" (str): Name of the vtk file.
"""
self.mapdl.post1()
self.mapdl.set("LAST")
self.mapdl.allsel()
# Define PyVista unstructured grid object and fill it with data
grid = self.mapdl.mesh.grid
grid.point_data['node_list'] = self.mapdl.get_array(entity='NODE', item1='NLIST')
grid.point_data['u_x'] = self.mapdl.get_array(entity='NODE', item1='U', it1num='X')
grid.point_data['u_y'] = self.mapdl.get_array(entity='NODE', item1='U', it1num='Y')
grid.point_data['u_z'] = self.mapdl.get_array(entity='NODE', item1='U', it1num='Z')
grid.point_data['epto_xx'] = self.mapdl.get_array(entity='NODE', item1='EPTO', it1num='X') * 100.0
grid.point_data['epto_yy'] = self.mapdl.get_array(entity='NODE', item1='EPTO', it1num='Y') * 100.0
grid.point_data['epto_xy'] = self.mapdl.get_array(entity='NODE', item1='EPTO', it1num='XY') * 0.5
grid.point_data['epto_eqv'] = self.mapdl.get_array(entity='NODE', item1='EPTO', it1num='EQV') * 100.0
grid.point_data['s_xx'] = self.mapdl.get_array(entity='NODE', item1='S', it1num='X')
grid.point_data['s_yy'] = self.mapdl.get_array(entity='NODE', item1='S', it1num='Y')
grid.point_data['s_xy'] = self.mapdl.get_array(entity='NODE', item1='S', it1num='XY')
grid.point_data['s_eqv'] = self.mapdl.get_array(entity='NODE', item1='S', it1num='EQV')
grid.save(os.path.join(self.OUTPUT_PATH, f'{filename}_mesh.vtk'), binary=False)
return grid
For a full example see:
- PyAnsys: 01_finite_element_model.py
- Ansys Parametric Design Language (APDL): 04_nodedata_export.inp