replace_short_sequences_global

mdtools.dynamics.replace_short_sequences_global(list_of_arrays, min_len, val=0, inplace=True, verbose=False, debug=False)[source]

Replace consecutive occurrences of the same value that are shorter than a given minimum length by a given value at any location of the arrays.

Parameters:
  • list_of_arrays (tuple or list or numpy.ndarray) – Sequence of arrays (usually one array per trajectory frame). All arrays in the sequence must have the same shape and all arrays must be NumPy arrays.

  • min_len (int) – All consecutive occurrences shorter than min_len will be replaced by val. If min_len is less than two, list_of_arrays will be returned without any changes.

  • val (scalar, optional) – The value used to replace the values in consecutive occurrences shorter than min_len. Be aware that the arrays in list_of_arrays must have a suitable dtype to be able to hold val.

  • inplace (bool, optional) – If True, modify list_of_arrays inplace instead of creating a copy of list_of_arrays.

  • verbose (bool, optional) – If True, print progress information to standard output.

  • debug (bool, optional) – If True, check the input arguments.

    Deprecated since version 0.0.0.dev0: This argument is without use and will be removed in a future release.

See also

mdtools.dynamics.correct_intermittency()

Correct for intermittent behavior of discrete variables stored in a sequence of arrays

mdtools.dynamics.correct_intermittency_1d()

Correct for intermittent behavior of discrete variables stored in a 1-dimensional array

mdtools.dynamics.replace_short_sequences()

Replace consecutive occurrences of the same value that are shorter than a given minimum length by a given value

Notes

This function is very similar to mdtools.dynamics.replace_short_sequences(), but here consecutive occurrences of the same value are not checked between elements at the same location of the arrays in list_of_arrays but all elements of the arrays in list_of_arrays are taken into account.

min_len must not exceed 65535 (np.uint16(-1)). If you really need a higher value for min_len, change the dtype of seen_nframes_ago and seen_nframes_ahead in the source code of this function. Be aware that high values of min_len increase the computational cost.

Examples

>>> a = np.array([[1, 2, 3],
...               [1, 2, 3]])
>>> b = np.array([[2, 2, 3],
...               [2, 2, 3]])
>>> c = np.array([[3, 3, 3],
...               [3, 3, 3]])
>>> x = (a, b, c)
>>> y = mdt.dyn.replace_short_sequences_global(x,
...                                            min_len=2,
...                                            val=0,
...                                            inplace=False,
...                                            verbose=False)
>>> for arr in y:
...     print(arr)
[[0 2 3]
 [0 2 3]]
[[2 2 3]
 [2 2 3]]
[[3 3 3]
 [3 3 3]]
>>> y = mdt.dyn.replace_short_sequences_global(x,
...                                            min_len=3,
...                                            val=0,
...                                            inplace=False,
...                                            verbose=False)
>>> for arr in y:
...     print(arr)
[[0 0 3]
 [0 0 3]]
[[0 0 3]
 [0 0 3]]
[[3 3 3]
 [3 3 3]]

If min_len is greater than len(list_of_arrays), all elements of all arrays are replaced by val:

>>> y = mdt.dyn.replace_short_sequences_global(x,
...                                            min_len=len(x)+1,
...                                            val=0,
...                                            inplace=False,
...                                            verbose=False)
>>> for arr in y:
...     print(arr)
[[0 0 0]
 [0 0 0]]
[[0 0 0]
 [0 0 0]]
[[0 0 0]
 [0 0 0]]

If min_len is less than two, the input list_of_arrays is returned without any changes:

>>> y = mdt.dyn.replace_short_sequences_global(x,
...                                            min_len=1,
...                                            val=0,
...                                            inplace=False,
...                                            verbose=False)
>>> all(np.array_equal(x[i], arr) for i, arr in enumerate(y))
True

Inplace modification:

>>> y = mdt.dyn.replace_short_sequences_global(x,
...                                            min_len=3,
...                                            val=0,
...                                            inplace=True,
...                                            verbose=False)
>>> x is y
True