THE API¶
This chapter contains the documentations for the main methods of the python-sscha code. It can be used both from advanced users, that wants to exploit python-sscha as a library, or developers, willing to add new features to the code (or adapt existing ones for their own purpouses).
The API is divided in Modules.
The Ensemble Module¶
This module deals with ensembles of configurations. It is used to generate random configurations from the dynamical matrix, to compute observables on the ensemble used in the SSCHA optimization. These includes the average force on atoms, the gradient of the SSCHA minimization, the quantum-thermal stress tensor, as well as properties of the ensemble, like reweighting.
The SchaMinimizer Module¶
This module is the main SSCHA minimizer. It allows to setup a single (one population) minimization. In this module the minimization algorithm is introduced, as well as stopping conditions and all the parameters usually located in the &inputscha namelist are read.
-
sscha.SchaMinimizer.
ApplyFCPrecond
(current_dyn, matrix, T=0)[source]¶ This function perform the precondition on a given matrix, by applying the inverse of the lambda function
- current_dyn3*nat x 3*nat
The current force-constant matrix to compute the Lambda tensor
- matrix3*nat x 3*nat
The matrix on which you want to apply the Lambda tensor.
- Tfloat
The temperature
- new_matrix3*nat x 3*nat
The matrix after the Lambda application
-
sscha.SchaMinimizer.
ApplyLambdaTensor
(current_dyn, matrix, T=0)[source]¶ This function perform the inverse of the precondtioning: it applies the Hessian matrix to the preconditioned gradient to obtain the real one. This is a test function.
- current_dyn3*nat x 3*nat
The current force-constant matrix to compute the Lambda tensor
- matrix3*nat x 3*nat
The matrix on which you want to apply the Lambda tensor.
- Tfloat
The temperature
- new_matrix3*nat x 3*nat
The matrix after the Lambda application
-
sscha.SchaMinimizer.
GetBestWykoffStep
(current_dyn)[source]¶ This is an alternative way to the preconditioning, in which the best wyckoff step is choosen rescaled on the current dynamical matrix.
NOTE: It works with real space matrices.
\[STEP = \frac{1}{\max \lambda(\Phi)}\]Where \(\lambda(\Phi)\) is the generic eigenvalue of the force constant matrix. This is because \(\Phi\) is correct Hessian of the free energy respect to the structure in the minimum.
The best step is returned in [Angstrom^2 / Ry].
- current_dynndarray 3n_at x 3n_at
The force constant matrix \(\Phi\). It should be in Ry/bohr^2.
-
sscha.SchaMinimizer.
GetStructPrecond
(current_dyn)[source]¶ NOTE: the Phi is in Ry/bohr^2 while the forces are in Ry/A NOTE: the output preconditioner (that must be interfaced with forces) is in A^2/Ry
The preconditioner of the structure minimization is computed directly from the dynamical matrix. It is the fake inverse (projected out the translations).
\[\Phi_{\alpha\beta}^{-1} = \frac{1}{\sqrt{M_\alpha M_\beta}}\sum_\mu \frac{e_\mu^\alpha e_\mu^\beta}{\omega_\mu^2}\]Where the sum is restricted to the non translational modes.
- current_dynPhonons()
The current dynamical matrix
- preconditionerndarray 3nat x 3nat
The inverse of the force constant matrix, it can be used as a preconditioner.
-
sscha.SchaMinimizer.
PerformRootStep
(dyn_q, grad_q, step_size=1, root_representation='sqrt', minimization_algorithm='sdes')[source]¶ As for the [Monacelli, Errea, Calandra, Mauri, PRB 2017], the nonlinear change of variable is used to perform the step.
It works as follows:
\[ \begin{align}\begin{aligned}\Phi \rightarrow \sqrt{x}{\Phi}\\\frac{\partial F}{\partial \Phi} \rightarrow \frac{\partial F}{\partial \sqrt{x}{\Phi}}\\\sqrt{x}{\Phi^{(n)}} \stackrel{\frac{\partial F}{\partial \sqrt{x}{\Phi}}}{\longrightarrow} \sqrt{x}{\Phi^{(n+1)}}\\\Phi^{(n+1)} = \left(\sqrt{x}{\Phi^{(n+1)}})^{x}\end{aligned}\end{align} \]Where the specific update step is determined by the minimization_algorithm, while the \(x\) order of the root representation is determined by the root_representation argument.
- dyn_qndarray( NQ x 3nat x 3nat )
The dynamical matrix in q space. The Nq are the total number of q.
- grad_qndarray( NQ x 3nat x 3nat )
The gradient of the dynamical matrix.
- step_sizefloat
The step size for the minimization
- root_representationstring
choice between “normal”, “sqrt” and “root4”. The value of \(x\) will be, respectively, 1, 2, 4.
- minimization_algorithmstring
The minimization algorithm to be used for the update.
- new_dyn_qndarray( NQ x 3nat x 3nat )
The updated dynamical matrix in q space
The Relax Module¶
This module deals with relaxations that are iterated over more populations. It includes the variable cell optimization algorithm. Here the parameters readed in the &relax namelist are read and setup.
-
sscha.Relax.
GetStaticBulkModulus
(structure, ase_calculator, eps=0.001)[source]¶ This method uses finite differences on the cell to compute the static bulk modulus. The cell is strained into several volumes, and the stress tensor is computed in orther to obtain the bulk modulus. Only the symmmetry relevant terms are computed.
- structureCC.Structure.Structure()
The structure on which you want to compute the static bulk modulus
- ase_calculatorase.calculators.calculator.Calculator()
One of the ase calculators to get the stress tensor in several strained cells.
- epsfloat
The strain module
- bk_modndarray (9x9)
The bulk modulus as a 9x9 matrix, expressed in eV / A^3
The Cluster Module¶
The Cluster module provide the interface between python-sscha and remote servers to witch you submit the energy and forces calculations. The input in &cluster namespace is interpreted in this module
-
sscha.Cluster.
units
= {'A': 63541.72207603944, 'AUT': 0.002375996331368385, 'Ang': 1.0, 'Angstrom': 1.0, 'Bohr': 0.5291772083535413, 'C': 6.241509647120417e+18, 'Debye': 0.20819435181122592, 'GPa': 0.006241509647120417, 'Ha': 27.21138386556469, 'Hartree': 27.21138386556469, 'J': 6.241509647120418e+18, 'Pascal': 6.241509647120417e-12, 'Ry': 13.605691932782346, 'Rydberg': 13.605691932782346, '_Grav': 6.67428e-11, '_Nav': 6.02214179e+23, '_amu': 1.660538782e-27, '_auf': 8.238722061327264e-08, '_aup': 29421010848651.156, '_aut': 2.418884325704007e-17, '_auv': 2187691.2538987417, '_c': 299792458.0, '_e': 1.602176487e-19, '_eps0': 8.85418781762039e-12, '_hbar': 1.054571628251774e-34, '_hplanck': 6.62606896e-34, '_k': 1.3806504e-23, '_me': 9.10938215e-31, '_mp': 1.672621637e-27, '_mu0': 1.2566370614359173e-06, 'alpha': 0.007297352536796447, 'eV': 1.0, 'fs': 0.09822695141392761, 'invcm': 0.00012398418754199978, 'kB': 8.617342790900664e-05, 'kJ': 6.241509647120418e+21, 'kcal': 2.611447636355183e+22, 'kg': 6.022141794216764e+26, 'm': 10000000000.0, 'mol': 6.02214179e+23, 'nm': 10.0, 's': 98226951413927.6, 'second': 98226951413927.6}¶ This is an untility script that is able to manage the submission into a cluster of an ensemble