amfe.tools module

A collection of tools which to not fit to one topic of the other modules.

Some tools here might be experimental.

amfe.tools.node2total(node_index, coordinate_index, ndof_node=2)[source]

Converts the node index and the corresponding coordinate index to the index of the total dof.

Parameters:
  • node_index (int) – Index of the node as shown in tools like paraview
  • coordinate_index (int) – Index of the coordinate; 0 if it’s x, 1 if it’s y etc.
  • ndof_node (int, optional) – Number of degrees of freedom per node
Returns:

total_index – Index of the total dof

Return type:

int

amfe.tools.total2node(total_index, ndof_node=2)[source]

Converts the total index in the global dofs to the coordinate index and the index fo the coordinate.

Parameters:
  • total_index (int) – Index of the total dof
  • ndof_node (int, optional) – Number of degrees of freedom per node
Returns:

  • node_index (int) – Index of the node as shown in tools like paraview
  • coordinate_index (int) – Index of the coordinate; 0 if it’s x, 1 if it’s y etc.

amfe.tools.read_hbmat(filename)[source]

Reads the hbmat file and returns it in the csc format.

Parameters:filename (string) – string of the filename
Returns:matrix – matrix which is saved in harwell-boeing format
Return type:sp.sparse.csc_matrix

Notes

Information on the Harwell Boeing format: http://people.sc.fsu.edu/~jburkardt/data/hb/hb.html

When the hbmat file is exported as an ASCII-file, the truncation of the numerical values can cause issues, for example

  • eigenvalues change
  • zero eigenvalues vanish
  • stiffness matrix becomes indefinite
  • etc.

Thus do not trust matrices which are imported with this method. The method is correct, but the truncation error in the hbmat file of the floating point digits might cause some issues.

amfe.tools.append_interactively(filename)[source]

Open an input dialog for interactively appending a string to a filename.

This filename function should make it easy to save output files from numerical experiments containing a time stamp with an additonal tag requested at time of saving.

Parameters:filename (string) – filename path, e.g. for saving
Returns:filename – filename path with additional stuff, maybe added for convenience or better understanding
Return type:string
amfe.tools.matshow_3d(A, thickness=0.8, cmap=<matplotlib.colors.ListedColormap object>, alpha=1.0)[source]

Show a matrix as bar-plot using matplotlib.bar3d plotting tools similar to pyplot.matshow.

Parameters:
  • A (ndarray) – Array to be plotted
  • thickness (float, optional) – thickness of the bar. Default: 0.8
  • cmap (matplotlib.cm function, optional) – Colormap-function of matplotlib. Default. mpl.cm.jet
  • alpha (float) – alpha channel value (transparency): alpha=1.0 is not transparent at all, alpha=0.0 is full transparent and thus invisible.
Returns:

barplot

Return type:

instance of mpl_toolkits.mplot3d.art3d.Poly3DCollection

See also

matplotlib.pyplot.matshow()

amfe.tools.amfe_dir(filename='')[source]

Return the absolute path of the filename given relative to the amfe directory.

Parameters:filename (string, optional) – relative path to something inside the amfe directory.
Returns:dir – string of the filename inside the AMFE-directory. Default value is ‘’, so the AMFE-directory is returned.
Return type:string
amfe.tools.h5_read_u(h5filename)[source]

Extract the displacement field of a given hdf5-file.

Parameters:h5filename (str) – Full filename (with e.g. .hdf5 ending) of the hdf5 file.
Returns:
  • u_constr (ndarray) – Displacement time series of the dofs with constraints implied. Shape is (ndof_constr, no_of_timesteps), i.e. u_constr[:,0] is the first timestep.
  • u_unconstr (ndarray) – Displacement time series of the dofs without constraints. I.e. the dofs are as in the mesh file. Shape is (ndof_unconstr, no_of_timesteps).
  • T (ndarray) – Time. Shape is (no_of_timesteps,).
amfe.tools.test(*args, **kwargs)[source]

Run all tests for AMfe.

amfe.tools.reorder_sparse_matrix(A)[source]

Reorder the sparse matrix A such that the bandwidth of the matrix is minimized using the Cuthill–McKee (RCM) algorithm.

Parameters:A (CSR or CSC sprarse symmetric matrix) – Sparse and symmetric matrix
Returns:
  • A_new (CSR or CSC sparse symmetric matrix) – reordered sparse and symmetric matrix
  • perm (ndarray) – vector of row and column permutation

References

E. Cuthill and J. McKee, “Reducing the Bandwidth of Sparse Symmetric Matrices”, ACM ‘69 Proceedings of the 1969 24th national conference, (1969).

amfe.tools.eggtimer(fkt)[source]

Egg timer for functions which reminds via speech, when the function has terminated.

The intention of this function is, that the user gets reminded, when longer simulations are over.

Parameters:fkt (function) – any function
Returns:fkt – function decorated with eggtimer. It reminds you via speech, when the function has terminated.
Return type:function

Examples

Import eggtimer function:

>>> from amfe import eggtimer

working directly on function:

>>> def square(a):
...     return a**2
...
>>> timed_square = eggtimer(square)
>>> timed_square(6)
36

working as decorator:

>>> @eggtimer
... def square(a):
...     return a**2
...
>>> square(6)
36
amfe.tools.compute_relative_error(red_file, ref_file, M=None)[source]

Compute the farhat error of two given hdf5 files displacement sets.

Parameters:
  • red_file (str) – Filename of the reduced or to be investigated file.
  • ref_file (str) – Filename of the reference file.
  • M (array-like, optional) – mass matrix of the given system. If None, no mass matrix is used. Default: None.
Returns:

err – Error.

Return type:

float

Note

\[ER = \frac{\sqrt{\sum\limits_{t\in T} \Delta u(t)^T M \Delta u(t)}}{ \sqrt{\sum\limits_{t\in T} u_{ref}(t)^T M u_{ref}(t)}}\]
amfe.tools.principal_angles(V1, V2, cosine=True, principal_vectors=False)[source]

Return the cosine of the principal angles of the two bases V1 and V2.

Parameters:
  • V1 (ndarray) – array denoting n-dimensional subspace spanned by V1 (Mxn)
  • V2 (ndarray) – array denoting subspace 2. Dimension is (MxO)
  • cosine (bool, optional) – flag stating, if the cosine of the angles is to be used
  • principal_vectors (bool, optional) – Option flag for returning principal vectors. Default is False.
Returns:

  • sigma (ndarray) – cosine of subspace angles
  • F1 (ndarray) – array of principal vectors of subspace spanned by V1. The columns give the principal vectors, i.e. F1[:,0] is the first principal vector associated with theta[0] and so on. Only returned, if principal_vectors=True.
  • F2 (ndarray) – array of principal vectors of subspace spanned by V2. Only returned if principal_vectors=True.

Notes

Both matrices V1 and V2 have live in the same vector space, i.e. they have to have the same number of rows.

Examples

TODO

References

[1]G. H. Golub and C. F. Van Loan. Matrix computations, volume 3. JHU Press, 2012.
amfe.tools.query_yes_no(question, default='yes')[source]

Ask a yes/no question and return their answer.

Parameters:
  • question (String) – The question to be asked
  • default (String "yes" or "no") – The default answer
  • Returns
  • --------
  • answer (Boolean) – Answer: True if yes, False if no.
amfe.tools.resulting_force(mechanical_system, force_vec, ref_point=None)[source]

Compute the resulting force and moment of a given force vector.

Parameters:
  • mechanical_system – mechanical system ẃith FE mesh
  • force_vec (array) – constrained force vector
  • ref_point (array-like, shape: (ndim), optional) – reference point to which the resulting moment is computed. Default value is None, meaning that the reference point is at the origin
Returns:

resulting_force_and_moment – resulting force and moment vector with (F_x, F_y, F_z, M_x, M_y, M_z)

Return type:

array, shape(6)