grogupy package
Submodules
grogupy.core module
- grogupy.core.build_hh_ss(dh)[source]
It builds the Hamiltonian and Overlap matrix from the sisl.dh class.
It restructures the data in the SPIN BOX representation, where NS is the number of supercells and NO is the number of orbitals.
- Parameters:
dh – sisl.physics.Hamiltonian Hamiltonian read in by sisl
- Returns:
- (NS, NO, NO) np.array_like
Hamiltonian in SPIN BOX representation
- ss(NS, NO, NO) np.array_like
Overlap matrix in SPIN BOX representation
- Return type:
hh
- grogupy.core.calc_Vu(H, Tu)[source]
Calculates the local perturbation in case of a spin rotation.
- Parameters:
H – (NO, NO) np.array_like Hamiltonian
Tu – (NO, NO) array_like Rotation around u
- Returns:
- (NO, NO) np.array_like
First order perturbed matrix
- Vu2(NO, NO) np.array_like
Second order perturbed matrix
- Return type:
Vu1
- grogupy.core.onsite_projection(matrix, idx1, idx2)[source]
It produces the slices of a matrix for the on site projection.
The slicing is along the last two axes as these contains the orbital indexing.
- Parameters:
matrix – (…, :, :) np.array_like Some matrix
idx – np.array_like The indexes of the orbitals
- Returns:
- np.array_like
Reduced matrix based on the projection
- grogupy.core.parallel_Gk(HK, SK, eran, eset)[source]
Calculates the Greens function by inversion.
It calculates the Greens function on all the energy levels at the same time.
- Parameters:
HK – (NO, NO), np.array_like Hamiltonian at a given k point
SK – (NO, NO), np.array_like Overlap Matrix at a given k point
eran – (eset) np.array_like Energy sample along the contour
eset – int Number of energy samples along the contour
- Returns:
- (eset, NO, NO), np.array_like
Green’s function at a given k point
- Return type:
Gk
- grogupy.core.remove_clutter_for_save(pairs, magnetic_entities)[source]
Removes unimportant data from the dictionaries.
It is used before saving to throw away data that is not needed for post processing.
- Parameters:
pairs – dict Contains all the pair information
magnetic_entities – dict Contains all the magnetic entity information
- Returns:
- dict
Contains all the reduced pair information
- magnetic_entitiesdict
Contains all the reduced magnetic entity information
- Return type:
pairs
- grogupy.core.sequential_GK(HK, SK, eran, eset)[source]
Calculates the Greens function by inversion.
It calculates sequentially over the energy levels.
- Parameters:
HK – (NO, NO), np.array_like Hamiltonian at a given k point
SK – (NO, NO), np.array_like Overlap Matrix at a given k point
eran – (eset) np.array_like Energy sample along the contour
eset – int Number of energy samples along the contour
- Returns:
- (eset, NO, NO), np.array_like
Green’s function at a given k point
- Return type:
Gk
- grogupy.core.setup_pairs_and_magnetic_entities(magnetic_entities, pairs, dh, simulation_parameters)[source]
It creates the complete structure of the dictionaries and fills some basic data.
It creates orbital indexes, spin box indexes, coordinates and tags for magnetic entities. Furthermore it creates the structures for the energies, the perturbed potentials and the Greens function calculation. It dose the same for the pairs.
- Parameters:
pairs – dict Contains the initial pair information
magnetic_entities – dict Contains the initial magnetic entity information
dh – sisl.physics.Hamiltonian Hamiltonian read in by sisl
simulation_parameters – dict A set of parameters from the simulation
- Returns:
- dict
Contains the initial information and the complete structure
- magnetic_entitiesdict
Contains the initial information and the complete structure
- Return type:
pairs
grogupy.grogu module
grogupy.io module
- grogupy.io.load_pickle(infile)[source]
Loads the data from the infile with pickle.
- Parameters:
infile – str Path to infile
- Returns:
- dict
A dictionary of data
- Return type:
data
- grogupy.io.print_atoms_and_pairs(magnetic_entities, pairs)[source]
It prints the pair and magnetic entity information for the grogu out.
- Parameters:
magnetic_entities – dict It contains the data on the magnetic entities
pairs – dict It contains the data on the pairs
- grogupy.io.print_job_description(simulation_parameters)[source]
It prints the parameters and the description of the job.
- Parameters:
simulation_parameters – dict It contains the simulations parameters
- grogupy.io.print_parameters(simulation_parameters)[source]
It prints the simulation parameters for the grogu out.
- Parameters:
simulation_parameters – dict It contains the simulations parameters
grogupy.magnetism module
- grogupy.magnetism.blow_up_orbindx(orb_indices)[source]
Function to blow up orbital indices to make SPIN BOX indices.
- Parameters:
orb_indices – np.array_like These are the indices in ORBITAL BOX
- Returns:
- np.array_like
These are the indices in SPIN BOX
- Return type:
orb_indices
- grogupy.magnetism.calculate_anisotropy_tensor(mag_ent)[source]
Calculates the renormalized anisotropy tensor from the energies.
It uses the grogu convention for output.
- Parameters:
mag_ent – dict An element from the magnetic entities
- Returns:
- np.array_like
elements of the anisotropy tensor
- Return type:
K
- grogupy.magnetism.calculate_exchange_tensor(pair)[source]
Calculates the exchange tensor from the energies.
It produces the isotropic exchange, the relevant elements from the Dzyaloshinskii-Morilla (Dm) tensor, the symmetric-anisotropy and the complete exchange tensor.
- Parameters:
pair – dict An element from the pairs
- Returns:
- float
Isotropic exchange (Tr[J] / 3)
- J_Snp.array_like
Symmetric-anisotropy (J_S = J - J_iso * I ––> Jxx, Jyy, Jxy, Jxz, Jyz)
- Dnp.array_like
DM elements (Dx, Dy, Dz)
- Jnp.array_like
Complete exchange tensor flattened (Jxx, Jxy, Jxz, Jyx, Jyy, Jyz, Jzx, Jzy, Jzz)
- Return type:
J_iso
- grogupy.magnetism.parse_magnetic_entity(dh, atom=None, l=None, **kwargs)[source]
Function to define orbital indexes of a given magnetic entity.
- Parameters:
dh – sisl.physics.Hamiltonian Hamiltonian from sisl
atom – integer or list of integers, optional Defining atom (or atoms) in the unit cell forming the magnetic entity. Defaults to None
l – integer, optional Defining the angular momentum channel. Defaults to None
- Returns:
- list
The orbital indexes of the given magnetic entity
- grogupy.magnetism.spin_tracer(M)[source]
Spin tracer utility.
This takes an operator with the orbital-spin sequence: orbital 1 up, orbital 1 down, orbital 2 up, orbital 2 down, that is in the SPIN-BOX representation, and extracts orbital dependent Pauli traces.
- Parameters:
M – np.array_like Traceable matrix
- Returns:
- dict
It contains the traced matrix with “x”, “y”, “z” and “c”
grogupy.utilities module
- grogupy.utilities.RotM(theta, u, eps=1e-10)[source]
Definition of rotation matrix with angle theta around direction u.
- Parameters:
theta – float The angle of rotation
u – np.array_like The rotation axis
eps – float, optional Cutoff for small elements in the resulting matrix. Defaults to 1e-10
- Returns:
- np.array_like
The rotation matrix
- grogupy.utilities.RotMa2b(a, b, eps=1e-10)[source]
Definition of rotation matrix rotating unit vector a to unit vector b.
Function returns array R such that R @ a = b holds.
- Parameters:
a – np.array_like First vector
b – np.array_like Second vector
eps – float, optional Cutoff for small elements in the resulting matrix. Defaults to 1e-10
- Returns:
- np.array_like
The rotation matrix with the above property
- grogupy.utilities.commutator(a, b)[source]
Shorthand for commutator.
Commutator of two matrices in the mathematical sense.
- Parameters:
a – np.array_like The first matrix
b – np.array_like The second matrix
- Returns:
- np.array_like
The commutator of a and b
- grogupy.utilities.crossM(u)[source]
Definition for the cross-product matrix.
It acts as a cross product with vector u.
- Parameters:
u – list or np.array_like The second vector in the cross product
- Returns:
- np.array_like
The matrix that represents teh cross product with a vector
- grogupy.utilities.hsk(H, ss, sc_off, k=(0, 0, 0))[source]
Speed up Hk and Sk generation.
Calculates the Hamiltonian and the Overlap matrix at a given k point. It is faster that the sisl version.
- Parameters:
H – np.array_like Hamiltonian in spin box form
ss – np.array_like Overlap matrix in spin box form
sc_off – list supercell indexes of the Hamiltonian
k – tuple, optional The k point where the matrices are set up. Defaults to (0, 0, 0)
- Returns:
- np.array_like
Hamiltonian at the given k point
- np.array_like
Overlap matrix at the given k point
- grogupy.utilities.int_de_ke(traced, we)[source]
It numerically integrates the traced matrix.
It is a wrapper from numpy.trapz and it contains the relevant constants to calculate the energy integral from equation 93 or 96.
- Parameters:
traced – np.array_like The trace of a matrix or a matrix product
we – float The weight of a point on the contour
- Returns:
- float
The energy calculated from the integral formula
- grogupy.utilities.make_contour(emin=-20, emax=0.0, enum=42, p=150)[source]
A more sophisticated contour generator.
Calculates the parameters for the complex contour integral. It uses the Legendre-Gauss quadrature method. It returns a class that contains the information for the contour integral.
- Parameters:
emin – int, optional Energy minimum of the contour. Defaults to -20
emax – float, optional Energy maximum of the contour. Defaults to 0.0, so the Fermi level
enum – int, optional Number of sample points along the contour. Defaults to 42
p – int, optional Shape parameter that describes the distribution of the sample points. Defaults to 150
- Returns:
- ccont
Contains all the information for the contour integral
- grogupy.utilities.make_kset(dirs='xyz', NUMK=20)[source]
Simple k-grid generator to sample the Brillouin zone.
Depending on the value of the dirs argument k sampling in 1,2 or 3 dimensions is generated. If dirs argument does not contain either of x, y or z a kset of a single k-pont at the origin is returned. The total number of k points is the NUMK**(dimensions)
- Parameters:
dirs – str, optional Directions of the k points in the Brillouin zone. They are the three lattice vectors. Defaults to “xyz”
NUMK – int, optional The number of k points in a direction. Defaults to 20
- Returns:
- np.array_like
An array of k points that uniformly sample the Brillouin zone in the given directions