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: |
|
---|---|
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: |
|
---|---|
Returns: |
|
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
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: |
|
---|---|
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: |
|
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: |
|
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: |
|
---|---|
Returns: | err – Error. |
Return type: | float |
Note
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: |
|
---|---|
Returns: |
|
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: |
|
---|
amfe.tools.
resulting_force
(mechanical_system, force_vec, ref_point=None)[source]¶Compute the resulting force and moment of a given force vector.
Parameters: |
|
---|---|
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) |