Reading a VMAP file

The most common use case is to get the element nodal stress tensor for a certain geometry 1 and a certain load state STATE-2 out of the vmap file. The vmap interface provides you the nodal geometry (node coordinates), the mesh connectivity index and the field variables.

You can retrieve a DataFrame of the mesh with the desired variables in just one statement.

>>> (pylife.vmap.VMAPImport('demos/plate_with_hole.vmap')
     .make_mesh('1', 'STATE-2')
     .join_coordinates()
     .join_variable('STRESS_CAUCHY')
     .join_variable('E')
     .to_frame())
                            x         y    z        S11       S22  S33        S12  S13  S23       E11       E22  E33       E12  E13  E23
element_id node_id
1          1734     14.897208  5.269875  0.0  27.080811  6.927080  0.0 -13.687358  0.0  0.0  0.000119 -0.000006  0.0 -0.000169  0.0  0.0
           1582     14.555333  5.355806  0.0  28.319006  1.178649  0.0 -10.732705  0.0  0.0  0.000133 -0.000035  0.0 -0.000133  0.0  0.0
           1596     14.630658  4.908741  0.0  47.701195  5.512213  0.0 -17.866833  0.0  0.0  0.000219 -0.000042  0.0 -0.000221  0.0  0.0
           4923     14.726271  5.312840  0.0  27.699907  4.052865  0.0 -12.210032  0.0  0.0  0.000126 -0.000020  0.0 -0.000151  0.0  0.0
           4924     14.592996  5.132274  0.0  38.010101  3.345431  0.0 -14.299768  0.0  0.0  0.000176 -0.000038  0.0 -0.000177  0.0  0.0
...                       ...       ...  ...        ...       ...  ...        ...  ...  ...       ...       ...  ...       ...  ...  ...
4770       3812    -13.189782 -5.691876  0.0  36.527439  2.470588  0.0 -14.706686  0.0  0.0  0.000170 -0.000040  0.0 -0.000182  0.0  0.0
           12418   -13.560289 -5.278386  0.0  32.868889  3.320898  0.0 -14.260107  0.0  0.0  0.000152 -0.000031  0.0 -0.000177  0.0  0.0
           14446   -13.673285 -5.569107  0.0  34.291058  3.642457  0.0 -13.836027  0.0  0.0  0.000158 -0.000032  0.0 -0.000171  0.0  0.0
           14614   -13.389065 -5.709927  0.0  36.063541  2.828889  0.0 -13.774759  0.0  0.0  0.000168 -0.000038  0.0 -0.000171  0.0  0.0
           14534   -13.276068 -5.419206  0.0  33.804211  2.829817  0.0 -14.580153  0.0  0.0  0.000157 -0.000035  0.0 -0.000181  0.0  0.0

[37884 rows x 15 columns]

Supported features

So far the following data can be read from a vmap file

Geometry

  • node positions

  • node element index

Field variables

Any field variables can be read and joined to the node element index from the following locations:

  • element

  • node

  • element nodal

In particular, field variables at integration point location cannot cannot be read, as that would require extrapolating them to the node positions. This functionality is not available in pyLife.

The VMAPImport Class

class pylife.vmap.VMAPImport(filename)[source]

The interface class to import a vmap file

Parameters:

filename (string) – The path to the vmap file to be read

Raises:

Exception – if the file cannot be read an exception is raised. So far any exception from the h5py module is passed through.

element_sets(geometry)[source]

Returns a list of the element_sets present in the vmap file

filter_element_set(element_set)[source]

Filters a node set out of the current mesh

Parameters:

element_set (string, optional) – The element set defined in the vmap file as geometry set

Return type:

self

Raises:

APIUseError – if the mesh has not been initialized using make_mesh()

filter_node_set(node_set)[source]

Filters a node set out of the current mesh

Parameters:

node_set (string) – The node set defined in the vmap file as geometry set

Return type:

self

Raises:

APIUseError – if the mesh has not been initialized using make_mesh()

geometries()[source]

Returns a list of geometry strings of geometries present in the vmap data

join_coordinates()[source]

Join the coordinates of the predefined geometry in the mesh

Return type:

self

Raises:

APIUseError – if the mesh has not been initialized using make_mesh()

Examples

Receive the mesh with the node coordinates

>>> pylife.vmap.VMAPImport('demos/plate_with_hole.vmap').make_mesh('1').join_coordinates().to_frame()
                            x         y    z
element_id node_id
1          1734     14.897208  5.269875  0.0
           1582     14.555333  5.355806  0.0
           1596     14.630658  4.908741  0.0
           4923     14.726271  5.312840  0.0
           4924     14.592996  5.132274  0.0
...                       ...       ...  ...
4770       3812    -13.189782 -5.691876  0.0
           12418   -13.560289 -5.278386  0.0
           14446   -13.673285 -5.569107  0.0
           14614   -13.389065 -5.709927  0.0
           14534   -13.276068 -5.419206  0.0

[37884 rows x 3 columns]

join_variable(var_name, state=None, column_names=None)[source]

Joins a field output variable to the mesh

Parameters:
  • var_name (string) – The name of the field variables

  • state (string, opional) – The load state of which the field variable is to be read If not given, the last defined state, either defined in make_mesh() or defeined in join_variable() is used.

  • column_names (list of string, optional) – The names of the columns names to be used in the DataFrame If not provided, it will be chosen according to the list shown below. The length of the list must match the dimension of the variable.

Return type:

self

Raises:
  • APIUseError – if the mesh has not been initialized using make_mesh()

  • KeyError – if the geometry, state or varname is not found of if the vmap file is corrupted

  • KeyError – if there are no column names given and known for the variable.

  • ValueError – if the length of the column_names does not match the dimension of the variable

Notes

The mesh must be initialized with make_mesh(). The final DataFrame can be retrieved with to_frame().

If the column_names argument is not provided the following column names are chosen

  • ‘DISPLACEMENT’: ['dx', 'dy', 'dz']

  • ‘STRESS_CAUCHY’: ['S11', 'S22', 'S33', 'S12', 'S13', 'S23']

  • ‘E’: ['E11', 'E22', 'E33', 'E12', 'E13', 'E23']

If that fails a KeyError exception is risen.

Examples

Receiving the ‘DISPLACEMENT’ of ‘STATE-1’ , the stress and strain tensors of ‘STATE-2’

>>> (pylife.vmap.VMAPImport('demos/plate_with_hole.vmap')
     .make_mesh('1')
     .join_variable('DISPLACEMENT', 'STATE-1')
     .join_variable('STRESS_CAUCHY', 'STATE-2')
     .join_variable('E').to_frame())
                     dx   dy   dz        S11       S22  S33        S12  S13  S23       E11       E22  E33       E12  E13  E23
element_id node_id
1          1734     0.0  0.0  0.0  27.080811  6.927080  0.0 -13.687358  0.0  0.0  0.000119 -0.000006  0.0 -0.000169  0.0  0.0
           1582     0.0  0.0  0.0  28.319006  1.178649  0.0 -10.732705  0.0  0.0  0.000133 -0.000035  0.0 -0.000133  0.0  0.0
           1596     0.0  0.0  0.0  47.701195  5.512213  0.0 -17.866833  0.0  0.0  0.000219 -0.000042  0.0 -0.000221  0.0  0.0
           4923     0.0  0.0  0.0  27.699907  4.052865  0.0 -12.210032  0.0  0.0  0.000126 -0.000020  0.0 -0.000151  0.0  0.0
           4924     0.0  0.0  0.0  38.010101  3.345431  0.0 -14.299768  0.0  0.0  0.000176 -0.000038  0.0 -0.000177  0.0  0.0
...                 ...  ...  ...        ...       ...  ...        ...  ...  ...       ...       ...  ...       ...  ...  ...
4770       3812     0.0  0.0  0.0  36.527439  2.470588  0.0 -14.706686  0.0  0.0  0.000170 -0.000040  0.0 -0.000182  0.0  0.0
           12418    0.0  0.0  0.0  32.868889  3.320898  0.0 -14.260107  0.0  0.0  0.000152 -0.000031  0.0 -0.000177  0.0  0.0
           14446    0.0  0.0  0.0  34.291058  3.642457  0.0 -13.836027  0.0  0.0  0.000158 -0.000032  0.0 -0.000171  0.0  0.0
           14614    0.0  0.0  0.0  36.063541  2.828889  0.0 -13.774759  0.0  0.0  0.000168 -0.000038  0.0 -0.000171  0.0  0.0
           14534    0.0  0.0  0.0  33.804211  2.829817  0.0 -14.580153  0.0  0.0  0.000157 -0.000035  0.0 -0.000181  0.0  0.0

[37884 rows x 15 columns]

Todo

Write a more central document about pyLife’s column names.

make_mesh(geometry, state=None)[source]

Makes the initial mesh

Parameters:
  • geometry (string) – The geometry defined in the vmap file

  • state (string, optional) – The load state of which the field variable is to be read. If not given, the state must be defined in join_variable().

Return type:

self

Raises:
  • KeyError – if the geometry is not found of if the vmap file is corrupted

  • KeyError – if the node_set or element_set is not found in the geometry.

  • APIUseError – if both, a node_set and an element_set are given

Notes

This methods defines the initial mesh to which coordinate data can be joined by join_coordinates() and field variables can be joined by join_variable()

Examples

Get the mesh data with the coordinates of geometry ‘1’ and the stress tensor of ‘STATE-2’

>>> (pylife.vmap.VMAPImport('demos/plate_with_hole.vmap')
     .make_mesh('1', 'STATE-2')
     .join_coordinates()
     .join_variable('STRESS_CAUCHY')
     .to_frame()
                            x         y    z        S11       S22  S33        S12  S13  S23
element_id node_id
1          1734     14.897208  5.269875  0.0  27.080811  6.927080  0.0 -13.687358  0.0  0.0
           1582     14.555333  5.355806  0.0  28.319006  1.178649  0.0 -10.732705  0.0  0.0
           1596     14.630658  4.908741  0.0  47.701195  5.512213  0.0 -17.866833  0.0  0.0
           4923     14.726271  5.312840  0.0  27.699907  4.052865  0.0 -12.210032  0.0  0.0
           4924     14.592996  5.132274  0.0  38.010101  3.345431  0.0 -14.299768  0.0  0.0
...                       ...       ...  ...        ...       ...  ...        ...  ...  ...
4770       3812    -13.189782 -5.691876  0.0  36.527439  2.470588  0.0 -14.706686  0.0  0.0
           12418   -13.560289 -5.278386  0.0  32.868889  3.320898  0.0 -14.260107  0.0  0.0
           14446   -13.673285 -5.569107  0.0  34.291058  3.642457  0.0 -13.836027  0.0  0.0
           14614   -13.389065 -5.709927  0.0  36.063541  2.828889  0.0 -13.774759  0.0  0.0
           14534   -13.276068 -5.419206  0.0  33.804211  2.829817  0.0 -14.580153  0.0  0.0
node_sets(geometry)[source]

Returns a list of the node_sets present in the vmap file

nodes(geometry)[source]

Retrieves the node positions

Parameters:

geometry (string) – The geometry defined in the vmap file

Returns:

node_positions – a DataFrame with the node numbers as index and the columns ‘x’, ‘y’ and ‘z’ for the node coordinates.

Return type:

DataFrame

Raises:

KeyError – if the geometry is not found of if the vmap file is corrupted

states()[source]

Returns a list of state strings of states present in the vmap data

to_frame()[source]

Returns the mesh and resets the mesh

Returns:

mesh – The mesh data joined so far

Return type:

DataFrame

Raises:

APIUseError – if there is no mesh present, i.e. make_mesh() has not been called yet or the mesh has been reset in the meantime.

Notes

This method resets the mesh, i.e. make_mesh() must be called again in order to fetch more mesh data in another mesh.

try_get_geometry_set(geometry_name, geometry_set_name)[source]
try_get_vmap_object(group_full_path)[source]
variables(geometry, state)[source]

Ask for available variables for a certain geometry and state.

Parameters:
  • geometry (string) – Name of the geometry

  • state (string) – Name of the state

Returns:

variables – List of available variable names for the geometry state combination

Return type:

list

Raises:

KeyError – if the geometry state combination is not available.