Density modules#

The density modules of MAICoS are tools for computing density, temperature, and chemical potential profiles from molecular simulation trajectory files. Profiles can be extracted either in Cartesian or cylindrical coordinate systems. Units for the density are the same as GROMACS, i.e. mass, number or charge. See the gmx density manual for details.

From the command line

You can extract a density profile from molecular dynamics trajectory files directly from the terminal. For this example, we use the airwater data file of MAICoS. First, go to the directory

cd tests/data/airwater/

then type:

maicos DensityPlanar -s conf.gro -traj traj.trr

Here conf.gro and traj.trr are GROMACS configuration and trajectory files, respectively. The density profile appears in a .dat file. You can visualise all the options of the module DensityPlanar by typing

maicos DensityPlanar -h

From the Python interpreter

In order to calculate the density using MAICoS in a Python environment, first import MAICoS and MDAnalysis:

import MDAnalysis as mda
import maicos

Then create a MDAnalysis universe:

u = mda.Universe('conf.gro', 'traj.trr')
group_H2O = u.select_atoms('type O or type H')

And run MAICoS’ DensityPlanar module:

dplan = maicos.DensityPlanar(group_H2O)
dplan.run()

Results can be accessed from dplan.results. More details are given in the Tutorials.

Density planar#

Compute partial densities/temperature profiles in the Cartesian systems.

Parameters#

atomgroupslist[AtomGroup]

a list of AtomGroup for which the densities are calculated

densstr

Density: mass, number, charge, temperature.

dimint

Dimension for binning (x=0, y=1, z=2)

binwidthfloat

binwidth (nanometer)

comgroupAtomGroup

Perform the binning relative to the center of mass of the selected group.

centerbool

Perform the binning relative to the center of the (changing) box.

mubool

Calculate the chemical potential (requires dens=’number’)

muoutstr

Prefix for output filename for chemical potential

temperaturefloat

temperature (K) for chemical potential

massfloat

Mass (u) for the chemical potential. By default taken from topology.

zposfloat

position (nm) at which the chemical potential will be computed. By default average over box.

outputstr

Output filename

concfreqint

Default number of frames after which results are calculated and files refreshed. If 0 results are only calculated at the end of the analysis and not saved by default.

verbosebool

Turn on more logging and debugging

Attributes#

results.zlist

bins

results.dens_meannp.ndarray

calculated densities

results.dens_stdnp.ndarray

density standard deviation

results.dens_errnp.ndarray

density error

results.mufloat

chemical potential (only if mu=True)

results.dmufloat

error of chemical potential (only if mu=True)

Density cylinder#

Compute partial densities across a cylinder.

Parameters#

atomgroupslist[AtomGroup]

A list of AtomGroup for which the densities are calculated.

densstr

Density: mass, number, charge, temp

dimint

Dimension for binning (x=0, y=1, z=2)

centerstr

Perform the binning relative to the center of this selection string of teh first AtomGroup. If None center of box is used.

radiusfloat

Radius of the cylinder (nm). If None smallest box extension is taken.

lengthfloat

Length of the cylinder (nm). If None length of box in the binning dimension is taken.

binwidthfloat

binwidth (nanometer)

outputstr

Output filename

concfreqint

Default number of frames after which results are calculated and files refreshed. If 0 results are only calculated at the end of the analysis and not saved by default.

verbosebool

Turn on more logging and debugging

Attributes#

results.rnp.ndarray

bins

results.dens_meannp.ndarray

calculated densities

results.dens_mean_sqnp.ndarray

squared calculated density

results.dens_stdnp.ndarray

density standard deviation

results.dens_errnp.ndarray

density error