mdtraj.compute_nematic_order

mdtraj.compute_nematic_order(traj, indices='chains')

Compute the nematic order parameter of a group in every frame.

The nematic order parameter describes the orientational order of a system with a value between 0 and 1. A completely isotropic system, such as a liquid, has no preferred direction and a nematic order parameter of 0. An anisotropic system, such as many liquid crystals, monolayers or bilayers, have a preferred orientation and will have a positive order parameter where a value of 1 signifies perfect ordering.

Parameters:

traj : Trajectory

Trajectory to compute ordering in.

indices : {‘chains’, ‘residues’, list of lists}, optional, default=’chains’

The group to consider. Users can pass their own indices as a list of lists with the “shape” (n_compounds, len(each_compound)). Recognized string keywords are ‘chains’ and ‘residues’.

Returns:

S2 : np.ndarray, shape=(traj.n_frames,), dtype=float64

Nematic order parameter values in every frame.

References

[R6466]Allen, M. P.; Tildesley , D. J. (1987), “Computer Simulation of Liquids”, Ch. 11.5
[R6566]http://cmt.dur.ac.uk/sjc/thesis_dlc/node65.html
[R6666]http://cmt.dur.ac.uk/sjc/thesis_dlc/node19.html

Examples

Ordering of chains in an alkylsilane monolayer of C10H31-Si(OH)2-:

>>> import mdtraj as md
>>> from mdtraj.testing import get_fn
>>> traj = md.load(get_fn('monolayer.xtc'), top=get_fn('monolayer.pdb'))
>>> # Each of the 100 chains contains 36 atoms.
>>> chain_indices = [[n+x for x in range(36)] for n in range(0, 3600, 36)]
>>> S2 = md.compute_nematic_order(traj, indices=chain_indices)

The chains were attached to a flat, crystalline substrate and are thus highly ordered with a mean S2 of ~0.996.

>>> traj = md.load(get_fn('tip3p_300K_1ATM.xtc'), top=get_fn('tip3p_300K_1ATM.pdb'))
>>> water_indices = [[n+x for x in range(3)] for n in range(0, 3600, 3)]
>>> S2 = md.compute_nematic_order(traj, indices=water_indices)

This water box is essentially isotropic and has a mean S2 of ~0.042.