grogupy.physics.MagneticEntity
- class grogupy.physics.MagneticEntity(infile: str | tuple[Hamiltonian, DensityMatrix | None], atom: None | int | list[int] = None, l: None | int | list[int] | list[list[int]] = None, orb: None | int | list[int] | list[list[int]] = None)[source]
This class contains the data and the methods related to the magnetic entities.
It sets up the instance based on the indexing of the Hamiltonian by the
atom,land orbital (orb) parameters.There are four possible input types: 1. Cluster: a list of atoms 2. AtomShell: one atom and a list of shells indexed in the atom or a list of atoms and a list of lists containing the shells 3. AtomOrbital: one atom and a list of orbitals indexed in the atom or a list of atoms and a list of lists containing the orbitals 4. Orbitals: a list of orbitals indexed in the Hamiltonian
Parameters
- infile: Union[str, tuple[Union[sisl.physics.Hamiltonian, Hamiltonian], Union[sisl.physics.DensityMatrix, None]]]
Either the path to the .fdf file or a tuple of sisl or grogupy hamiltonian and a sisl density matrix
- atom: Union[None, int, list[int]], optional
Defining atom (or atoms) in the unit cell forming the magnetic entity, by default None
- l: Union[None, int, list[int], list[list[int]]], optional
Defining the angular momentum channel, by default None
- orb: Union[None, int, list[int], list[list[int]]], optional
Defining the orbital index in the Hamiltonian or on the atom, by default None
Examples
Creating a magnetic entity can be done by giving the Hamiltonian from the DFT calculation and somehow specifying the corresponding atoms and orbitals.
The following examples show you how to create magnetic entities in the Fe3GeTe2 system. You can compare the tags of the
MagneticEntityto the input parameters, to understand how to build the magnetic entity, that suits your needs.>>> fdf_path = "/Users/danielpozsar/Downloads/Fe3GeTe2/Fe3GeTe2.fdf"
To define a Cluster of atoms use a dictionary that only contains atoms.
>>> magnetic_entity = MagneticEntity(fdf_path, atom=[3,4,5]) >>> print(magnetic_entity.tag) 3Fe(l:All)--4Fe(l:All)--5Fe(l:All)
To define a magnetic entity with a single atom, but with specific shells use both the
atomandlkey in the dictionary.>>> magnetic_entity = MagneticEntity(fdf_path, atom=5, l=[[1,2,3]]) >>> print(magnetic_entity.tag) 5Fe(l:1-2-3)
Or you can define multiple atoms with different shells:
>>> magnetic_entity = MagneticEntity(fdf_path, atom=[4,5], l=[[1],[1,2,3]]) >>> print(magnetic_entity.tag) 4Fe(l:1)--5Fe(l:1-2-3)
To define a magnetic entity with a single atom, but with specific orbitals use both the
atomandorbkey in the dictionary. Be aware that these orbitals are indexed inside the atom, not in the total Hamiltonian.>>> magnetic_entity = MagneticEntity(fdf_path, atom=5, orb=[[1,2,3,4,5,6,7,8,9,10]]) >>> print(magnetic_entity.tag) 5Fe(o:1-2-3-4-5-6-7-8-9-10)
Or you can define multiple atoms with different orbitals:
>>> magnetic_entity = MagneticEntity(fdf_path, atom=[4,5], orb=[[1],[1,2,3,4,5,6,7,8,9,10]]) >>> print(magnetic_entity.tag) 4Fe(o:1)--5Fe(o:1-2-3-4-5-6-7-8-9-10)
And finally you can use only the
orbkey to directly index the orbitals from the Hamiltonian.>>> magnetic_entity = MagneticEntity(fdf_path, orb=[1,10,30,40,50]) >>> print(magnetic_entity.tag) 0Te(o:1-10)--2Ge(o:4)--3Fe(o:1-11)
Methods
- calculate_energies(weights) :
Calculates the energies of the infinitesimal rotations.
- calculate_anisotropy() :
Calculates the anisotropy matrix and the consistency from the energies.
- copy() :
Return a copy of this MagneticEntity
Attributes
- _orbital_box_indicesNDArray
The ORBITAL BOX indices
- _total_orbital_box_indicesNDArray
The ORBITAL BOX indices for the whole atom
- _atomNDArray
The list of atoms in the magnetic entity
- _llist[list[Union[None, int]]]
The list of l in the magnetic entity, None if it is incomplete
- _spin_box_indicesNDArray
The SPIN BOX indices
- _total_mulliken: Union[None, NDArray]
Sisl Mulliken charges from the total atom
- _local_mulliken: Union[None, NDArray]
Sisl Mulliken charges from the given orbitals
- SBSint
Length of the SPIN BOX indices
- _xyzNDArray
The coordinates of the magnetic entity (it can consist of many atoms)
- _xyz_centerNDArray
The center of coordinates for the magnetic entity
- _tagstr
The description of the magnetic entity
- _Vu1list[list[float]]
The list of the first order rotations
- _Vu2list[list[float]]
The list of the second order rotations
- _Giilist[NDArray]
The list of the projected Greens functions
- energiesUnion[None, NDArray]
The calculated energies for each direction
- KUnion[None, NDArray]
The magnetic anisotropy, by default None
- K_consistencyUnion[None, float]
Consistency check on the diagonal K elements, by default None
- dh_ds_id: str
Hash of the underlying system so MagneticEntities cannot be mixed between Hamiltonians
- cell: NDArray
Cell size of the underlying system so MagneticEntities cannot be mixed between Hamiltonians
- __init__(infile: str | tuple[Hamiltonian, DensityMatrix | None], atom: None | int | list[int] = None, l: None | int | list[int] | list[list[int]] = None, orb: None | int | list[int] | list[list[int]] = None) None[source]
Initialize the magnetic entity.
Methods
__init__(infile[, atom, l, orb])Initialize the magnetic entity.
calculate_anisotropy()Calculates the anisotropy matrix and the consistency from the energies.
calculate_energies(weights[, append, ...])Calculates the energies of the infinitesimal rotations.
copy()Returns the deepcopy of the instance.
fit_anisotropy_tensor(ref_xcf)Fits the anisotropy tensor to the energies.
reset()Resets the simulation results.
Attributes
KThe anisotropy tensor, in meV.
K_JThe anisotropy tensor, but in J.
K_consistencyThe consistency check, in meV.
K_consistency_JThe consistency check, but in J.
K_consistency_mRyThe consistency check, but in mRy.
K_consistency_meVThe consistency check, but in meV.
K_mRyThe anisotropy tensor, but in mRy.
K_meVThe anisotropy tensor, but in meV.
SBIThe spin box indices of the magnetic entity
SBSThe spin box size of the magnetic entity
cellThe cell of the Hamiltonian.
dh_ds_idThe ID of the Hamiltonian and the Density Matrix.
energies_JThe energies, but in J.
energies_mRyThe energies, but in mRy.
energies_meVThe energies, but in meV.
local_QThe charge of the magnetic entity.
local_SSpin moment of the magnetic entity.
local_SxSx of the magnetic entity.
local_SySy of the magnetic entity.
local_SzSz of the magnetic entity.
number_of_entitiestagThe description of the magnetic entity
total_QThe total charge of the atom or the atoms of the magnetic entity.
total_STotal spin moment of the atom or the atoms of the magnetic entity.
total_SxTotal Sx of the atom or the atoms of the magnetic entity.
total_SyTotal Sy of the atom or the atoms of the magnetic entity.
total_SzTotal Sz of the atom or the atoms of the magnetic entity.
xyz_centerThe mean of the position of the atoms that are in the magnetic entity.