matscipy.pressurecoupling

Classes to be used with ASE in order to perform pressure relaxation and/or sliding simulations under pressure.

A usage example can be found in docs/applications/tribology.ipynb.

Some parts are based on L. Pastewka, S. Moser, and M. Moseler, Tribol. Lett. 39, 49 (2010) as indicated again below.

Classes

AutoDamping(C11[, p_c])

Automatic damping.

FixedDamping(gamma[, M_factor])

Damping with fixed damping constant and fixed mass.

FixedMassCriticalDamping(C11[, M_factor])

Damping with fixed mass and critical damping constant.

SlideLog(handle)

Reader for logs written with SlideLogger instance.

SlideLogger(handle, atoms, slider, integrator)

Logger to be attached to an ASE integrator.

SlideWithNormalPressureCuboidCell(top_mask, ...)

ASE constraint used for sliding with pressure coupling.

class matscipy.pressurecoupling.AutoDamping(C11, p_c=0.01)

Bases: object

Automatic damping.

Following L. Pastewka, S. Moser, and M. Moseler, Tribol. Lett. 39, 49 (2010).

Parameters:
  • C11 (float) – Elastic material constant.

  • p_c (float) – Empirical cut-off parameter.

Methods

get_M_gamma(slider, atoms)

Calculate mass M and dissipation constant gamma.

__init__(C11, p_c=0.01)
get_M_gamma(slider, atoms)

Calculate mass M and dissipation constant gamma.

Parameters:
Returns:

  • M (float) – Mass parameter.

  • gamma (float) – Dissipation constant parameter.

class matscipy.pressurecoupling.FixedDamping(gamma, M_factor=1.0)

Bases: object

Damping with fixed damping constant and fixed mass.

Parameters:
  • gamma (float) – Damping constant.

  • M_factor (float) – Multiplicative factor to increase actual mass of upper rigid atoms.

Methods

get_M_gamma(slider, atoms)

Calculate mass M and damping constant gamma.

__init__(gamma, M_factor=1.0)
get_M_gamma(slider, atoms)

Calculate mass M and damping constant gamma.

Parameters:
Returns:

  • M (float) – Mass parameter.

  • gamma (float) – Damping parameter.

class matscipy.pressurecoupling.FixedMassCriticalDamping(C11, M_factor=1.0)

Bases: object

Damping with fixed mass and critical damping constant.

Useful for fast pressure equilibration with small lid mass.

Parameters:
  • C11 (float) – Elastic material constant.

  • M_factor (float) – Multiplicative factor to increase actual mass of upper rigid atoms.

Methods

get_M_gamma(slider, atoms)

Calculate mass M and damping constant gamma.

__init__(C11, M_factor=1.0)
get_M_gamma(slider, atoms)

Calculate mass M and damping constant gamma.

Parameters:
Returns:

  • M (float) – Mass parameter.

  • gamma (float) – Damping parameter.

class matscipy.pressurecoupling.SlideWithNormalPressureCuboidCell(top_mask, bottom_mask, Pdir, P, vdir, v, damping)

Bases: object

ASE constraint used for sliding with pressure coupling.

Following L. Pastewka, S. Moser, and M. Moseler, Tribol. Lett. 39, 49 (2010).

Parameters:
  • top_mask (boolean numpy array) – Array a with a[i] == True for each index i of the constraint top atoms (the atoms which slide with constant speed).

  • bottom_mask (boolean numpy array) – same as top_mask but for completely fixed bottom atoms.

  • Pdir (int) – Index of cell axis (0, 1, 2) along which normal pressure is applied.

  • P (int) – Normal pressure in ASE units (e.g. 10.0 * ase.units.GPa).

  • vdir (int) – Index of cell axis (0, 1, 2) along which to slide.

  • v (float) – Constant sliding speed in ASE units (e.g. 100.0 * ase.units.m / ase.units.s).

  • damping – Damping object (e.g. matscipy.pressurecoupling.AutoDamping instance).

Attributes:
Tdir

Direction used for thermostatting.

middle_mask

Mask of free atoms.

Methods

adjust_forces(atoms, forces)

Adjust forces of upper and lower rigid atoms.

adjust_momenta(atoms, momenta)

Adjust momenta of upper and lower rigid atoms.

adjust_positions(atoms, positions)

Do not adjust positions.

adjust_potential_energy(atoms)

Do not adjust energy.

get_A(atoms)

Calculate cell area normal to applied load.

__init__(top_mask, bottom_mask, Pdir, P, vdir, v, damping)
property Tdir

Direction used for thermostatting.

Thermostat direction is normal to the sliding direction and the direction of the applied load direction.

Returns:

Direction used for thermostatting.

Return type:

int

property middle_mask

Mask of free atoms.

Returns:

Array a with a[i] == True for each index i of the atoms not being part of lower or upper rigid group.

Return type:

numpy boolean array

adjust_positions(atoms, positions)

Do not adjust positions.

get_A(atoms)

Calculate cell area normal to applied load.

Returns:

Cell area normal to applied load.

Return type:

float

Raises:

NotImplementedError – If atoms.get_cell() is non-orthogonal, SlideWithNormalPressureCuboidCell only works for orthogonal cells.

adjust_forces(atoms, forces)

Adjust forces of upper and lower rigid atoms.

Raises:

NotImplementedError – If atoms.get_cell() is non-orthogonal, SlideWithNormalPressureCuboidCell only works for orthogonal cells.

adjust_momenta(atoms, momenta)

Adjust momenta of upper and lower rigid atoms.

Raises:

NotImplementedError – If atoms.get_cell() is non-orthogonal, SlideWithNormalPressureCuboidCell only works for orthogonal cells.

adjust_potential_energy(atoms)

Do not adjust energy.

class matscipy.pressurecoupling.SlideLogger(handle, atoms, slider, integrator, step_offset=0)

Bases: object

Logger to be attached to an ASE integrator.

For new files (not restart jobs), the write_header method should be called once in order to write the header to the file.

Parameters:
  • handle (filehandle) – Filehandle e.g. pointing to a file opened in w or a mode.

  • atoms (ase.Atoms) – Atomic configuration.

  • slider (slider object) – Instance of SlideWithNormalPressureCuboidCell.

  • integrator (ASE integrator object,) – Instance of ASE integrator e.g. ase.md.langevin.Langevin.

  • step_offset (int) – Last step already written to log file, useful for restarts.

Examples

  1. For new runs:

    log_handle = open(logfn, ‘w’, 1) # line buffered

    logger = SlideLogger(log_handle, …)

    logger.write_header()

    integrator.attach(logger)

    integrator.run(steps_integrate)

    log_handle.close()

  2. For restarts:

    with open(logfn, ‘r’) as log_handle:

    step_offset = SlideLog(log_handle).step[-1]

    log_handle = open(logfn, ‘a’, 1) # line buffered append

    logger = SlideLogger(log_handle, …, step_offset=step_offset)

    integrator.attach(logger)

    integrator.run(steps_integrate)

    log_handle.close()

Methods

__call__()

Write current status (time, T, P, ...) to log-file.

write_header()

Write header of log-file.

__init__(handle, atoms, slider, integrator, step_offset=0)
write_header()

Write header of log-file.

class matscipy.pressurecoupling.SlideLog(handle)

Bases: object

Reader for logs written with SlideLogger instance.

Parameters:

instance (handle matscipy.pressurecoupling.SlideLogger) – Handle or filename pointing to log file.

step

Step indices 0, 1, 2, ….

Type:

ndarray

time

Simulation time in fs at step.

Type:

ndarray

T_thermostat

Instantantaneous temperature in K from thermostat-region only from degrees of freedom along thermalized direction.

Type:

ndarray

P_top

Normal pressure on lid in GPa.

Type:

ndarray

P_bottom

Normal pressure on base in GPa.

Type:

ndarray

h

Separation of lid and base in Ang.

Type:

ndarray

v

Normal speed of lid in Ang / fs.

Type:

ndarray

a

Normal acceleration of lid in Ang / fs ** 2.

Type:

ndarray

tau_top

Shear stress on lid in GPa.

Type:

ndarray

tau_bottom

Shear stress on base in GPa.

Type:

ndarray

rows

All data in a 2d array with axis 0 step and axis 1 the values ordered as above.

Type:

ndarray

__init__(handle)