amfe.mechanical_system module

Module handling the whole mechanical system, no matter if it’s a finite element system, defined by certain parameters or a multibody system.

class amfe.mechanical_system.MechanicalSystem(stress_recovery=False)[source]

Bases: object

Master class for mechanical systems with the goal to black-box the routines of assembly and element selection.

mesh_class

instance of Mesh() – Class handling the mesh.

assembly_class

instance of Assembly() – Class handling the assembly.

dirichlet_class

instance of DirichletBoundary – Class handling the Dirichlet boundary conditions.

T_output

list of floats – List of timesteps saved.

u_output

list of ndarrays – List of unconstrained displacement arrays corresponding to the timesteps in T_output.

S_output

list of ndarrays – List of stress arrays corresponding to the timesteps in T_output.

E_output

list of ndarrays – List of strain arrays corresponding to the timesteps in T_output.

stress

ndarray – Array of nodal stress of the last assembly run. Shape is (no_of_nodes, 6).

strain

ndarray – Array of nodal strain of the last assembly run. Shape is (no_of_nodes, 6).

stress_recovery

bool – Flag for option stress_recovery.

iteration_info

ndarray – array containing the information of an iterative solution procedure. iteration_info[:,0] is the time information, iteration_info[:,1] is the number of iteations, iteration_info[:,3] is the residual.

D(u=None, t=0)[source]

Return the damping matrix of the mechanical system

Parameters:
  • u (ndarray, optional) – Displacement field in voigt notation
  • t (float, optional) – Time
Returns:

K – Stiffness matrix with applied constraints in sparse csr-format

Return type:

sp.sparse.sparse_matrix

K(u=None, t=0)[source]

Compute the stiffness matrix of the mechanical system

Parameters:
  • u (ndarray, optional) – Displacement field in voigt notation
  • t (float, optional) – Time
Returns:

K – Stiffness matrix with applied constraints in sparse csr-format

Return type:

sp.sparse.sparse_matrix

K_and_f(u=None, t=0)[source]

Compute tangential stiffness matrix and nonlinear force vector in one assembly run.

M(u=None, t=0)[source]

Compute the Mass matrix of the dynamical system.

Parameters:
  • u (ndarray, optional) – array of the displacement
  • t (float) – time
Returns:

M – Mass matrix with applied constraints in sparse csr-format

Return type:

sp.sparse.sparse_matrix

S_and_res(u, du, ddu, dt, t, beta, gamma)[source]

Compute jacobian and residual for implicit time integration.

Parameters:
  • u (ndarray) – displacement; dimension (ndof,)
  • du (ndarray) – velocity; dimension (ndof,)
  • ddu (ndarray) – acceleration; dimension (ndof,)
  • dt (float) – time step width
  • t (float) – time of current time step (for time dependent loads)
  • beta (float) – weighting factor for position in generalized-\(\alpha\) scheme
  • gamma (float) – weighting factor for velocity in generalized-\(\alpha\) scheme
Returns:

  • S (ndarray) – jacobian matrix of residual; dimension (ndof, ndof)
  • res (ndarray) – residual; dimension (ndof,)

Notes

Time integration scheme: The iteration matrix is composed using the generalized-\(\alpha\) scheme:

\[\mathbf S = \frac{1}{h^2\beta}\mathbf{M} + \frac{\gamma}{h\beta} \mathbf D + \mathbf K\]

which bases on the time discretization of the velocity and the displacement:

\[\mathbf{\dot{q}}_{n+1} & = \mathbf{\dot{q}}_{n} + (1-\gamma)h \mathbf{\ddot{q}}_{n} + \gamma h \mathbf{\ddot{q}}_{n+1}\]
\[\mathbf{q}_{n+1} & = \mathbf{q}_n + h \mathbf{\dot{q}}_n + \left(\frac{1}{2} - \beta\right)h^2\mathbf{\ddot{q}}_n + h^2\beta\mathbf{\ddot{q}}_{n+1}\]

This method is using the variables/methods

  • self.M()
  • self.M_constr
  • self.K_and_f()
  • self.f_ext()

If these methods are implemented correctly in a daughter class, the time integration interface should work properly.

apply_dirichlet_boundaries(key, coord, mesh_prop='phys_group')[source]

Apply dirichlet-boundaries to the system.

Parameters:
  • key (int) – Key for mesh property which is to be chosen. Matches the group given in the gmsh file. For help, the function mesh_information or boundary_information gives the groups
  • coord (str {'x', 'y', 'z', 'xy', 'xz', 'yz', 'xyz'}) – coordinates which should be fixed
  • mesh_prop (str {'phys_group', 'geom_entity', 'el_type'}, optional) – label of which the element should be chosen from. Default is ‘phys_group’.
Returns:

Return type:

None

apply_neumann_boundaries(key, val, direct, time_func=None, shadow_area=False, mesh_prop='phys_group')[source]

Apply neumann boundaries to the system via skin elements.

Parameters:
  • key (int) – Key of the physical domain to be chosen for the neumann bc
  • val (float) – value for the pressure/traction onto the element
  • direct (ndarray or str 'normal') – Direction, in which force should act at. If
  • time_func (function object) –

    Function object returning a value between -1 and 1 given the input t:

    >>> val = time_func(t)
    
  • shadow_area (bool, optional) – flag, if force should be proportional to shadow area of surface with respect to direction. Default: False.
  • mesh_prop (str {'phys_group', 'geom_entity', 'el_type'}, optional) – label of which the element should be chosen from. Default is phys_group.
Returns:

Return type:

None

apply_rayleigh_damping(alpha, beta)[source]

Apply Rayleigh damping to the system.

The damping matrix D is defined as

D = alpha*M + beta*K(0)

Thus, it is Rayleigh Damping applied to the linearized system around zero deformation.

Parameters:
  • alpha (float) – damping coefficient for the mass matrix
  • beta (float) – damping coefficient for the stiffness matrix
clear_timesteps()[source]

Clear the timesteps gathered internally

deflate_mesh()[source]

Remove free floating nodes not connected to a selected element from the mesh.

Parameters:None
Returns:
Return type:None
export_paraview(filename, field_list=None)[source]

Export the system with the given information to paraview.

Parameters:
  • filename (str) – filename to which the xdmf file and the hdf5 file will be saved.
  • field_list (list, optional) –

    list of tuples containing a field to be exported as well as a dictionary with the attribute information of the hdf5 file. Example:

    >>> # example field list with reduced displacement not to export
    >>> # ParaView and strain epsilon to be exported to ParaView
    >>> field_list = [(q_red, {'ParaView':False, 'Name':'q_red'}),
                      (eps, {'ParaView':True,
                             'Name':'epsilon',
                             'AttributeType':'Tensor6',
                             'Center':'Node',
                             'NoOfComponents':6})]
    
Returns:

Return type:

None

f_ext(u, du, t)[source]

Return the nonlinear external force of the right hand side of the equation, i.e. the excitation.

f_int(u, t=0)[source]

Return the elastic restoring force of the system

gen_alpha(q, dq, ddq, q_old, dq_old, ddq_old, f_ext_old, dt, t, alpha_m, alpha_f, beta, gamma)[source]

Computation of Jacobian and residual for generalized-alpha time integration scheme.

TODO

load_mesh_from_csv(node_list_csv, element_list_csv, no_of_dofs_per_node=2, explicit_node_numbering=False, ele_type=False)[source]

Loads the mesh from two csv-files containing the node and the element list.

Parameters:
  • node_list_csv (str) – filename of the csv-file containing the coordinates of the nodes (x, y, z)
  • element_list_csv (str) – filename of the csv-file containing the nodes which belong to one element
  • no_of_dofs_per_node (int, optional) – degree of freedom per node as saved in the csv-file
  • explicit_node_numbering (bool, optional) – flag stating, if the node numbers are explcitly numbered in the csv file, i.e. if the first column gives the numbers of the nodes.
  • ele_type (str) – Spezifiy elements type of the mesh (e.g. for a Tri-Mesh different elements types as Tri3, Tri4, Tri6 can be used) If not spezified value is set to ‘False’
Returns:

Return type:

None

Examples

todo

load_mesh_from_gmsh(msh_file, phys_group, material, scale_factor=1)[source]

Load the mesh from a msh-file generated by gmsh.

Parameters:
  • msh_file (str) – file name to an existing .msh file
  • phys_group (int) – integer key of the physical group which is considered as the mesh part
  • material (amfe.Material) – Material associated with the physical group to be computed
  • scale_factor (float, optional) – scale factor for the mesh to adjust the units. The default value is 1, i.e. no scaling is done.
Returns:

Return type:

None

tie_mesh(master_key, slave_key, master_prop='phys_group', slave_prop='phys_group', tying_type='fixed', verbose=False, conform_slave_mesh=False, fix_mesh_dist=0.001)[source]

Tie nonconforming meshes for a given master and slave side.

Parameters:
  • master_key (int or string) – mesh key of the master face mesh. The master face mesh has to be at least the size of the slave mesh. It is better, when the master mesh is larger than the slave mesh.
  • slave_key (int or string) – mesh key of the slave face mesh or point cloud
  • master_prop (string, optional) – mesh property for which master_key is specified. Default value: ‘phys_group’
  • slave_prop (string, optional) – mesh property for which slave_key is specified. Default value: ‘phys_group’
  • tying_type (string {'fixed', 'slide'}) – Mesh tying type. ‘fixed’ glues the meshes together while ‘slide’ allows for a sliding motion between the meshes.
Returns:

Return type:

None

Notes

The master mesh has to embrace the full slave mesh. If this is not the case, the routine will fail, a slave point outside the master mesh cannot be addressed to a specific element.

write_timestep(t, u)[source]

write the timestep to the mechanical_system class

class amfe.mechanical_system.ReducedSystem(V_basis=None, assembly='indirect', **kwargs)[source]

Bases: amfe.mechanical_system.MechanicalSystem

Class for reduced systems. It is directly inherited from MechanicalSystem. Provides the interface for an integration scheme and so on where a basis vector is to be chosen...

Notes

The Basis V is a Matrix with x = V*q mapping the reduced set of coordinates q onto the physical coordinates x. The coordinates x are constrained, i.e. the x denotes the full system in the sense of the problem set and not of the pure finite element set.

The system runs without providing a V_basis when constructing the method only for the unreduced routines.

Examples

TODO

K(u=None, t=0)[source]
K_and_f(u=None, t=0)[source]
K_unreduced(u=None, t=0)[source]

Unreduced Stiffness Matrix.

Parameters:
  • u (ndarray, optional) – Displacement of constrained system. Default is zero vector.
  • t (float, optionial) – Time. Default is 0.
Returns:

K – Stiffness matrix

Return type:

sparse csr matrix

M(u=None, t=0)[source]
M_unreduced()[source]

Unreduced mass matrix.

clear_timesteps()[source]
export_paraview(filename, field_list=None)[source]

Export the produced results to ParaView via XDMF format.

f_ext(u, du, t)[source]
f_int(u, t=0)[source]
f_int_unreduced(u, t=0)[source]

Internal nonlinear force of the unreduced system.

Parameters:
  • u (ndarray) – displacement of unreduces system.
  • t (float, optional) – time, default value: 0.
Returns:

f_nl – nonlinear force of unreduced system.

Return type:

ndarray

write_timestep(t, u)[source]
class amfe.mechanical_system.ExternalForce(force_basis, force_series, t_series)[source]

Bases: object

Class for mimicking the external forces based on a force basis and time values. The force values are linearly interpolated.

f_ext(u, du, t)[source]

Mimicked external force for the given force time series