kdotpy
Loading...
Searching...
No Matches
kdotpy.diagonalization.diagdata.DiagDataPoint Class Reference
Inheritance diagram for kdotpy.diagonalization.diagdata.DiagDataPoint:
Collaboration diagram for kdotpy.diagonalization.diagdata.DiagDataPoint:

Public Member Functions

 __init__ (self, k, eival=None, eivec=None, paramval=None, grid_index=None, opts=None)
 
 obsids (self)
 
 __str__ (self)
 
 hash_id (self, length=6, precision='%.12e')
 
 file_id (self)
 
 stitch_with (self, k, eival, eivec, targetenergy_old, targetenergy_new, inplace=False, accuracy=0.01)
 
 update (self, new_ddp)
 
 extend_by (self, k, eival, eivec, paramval=None, obsvals=None, obsids=None, char=None, llindex=None, bindex=None, accuracy=1e-6)
 
 extend (self, *args, **kwds)
 
 set_observables (self, obsvals, obsids=None)
 
 calculate_observables (self, params, obs, obs_prop=None, overlap_eivec=None, magn=None)
 
 add_observable (self, obsvals=None, obsid=None)
 
 delete_eivec (self)
 
 build_tuple_index_cache (self)
 
 get_index (self, val)
 
 get_index_with_llindex (self, val, llindex)
 
 get_ubindex (self)
 
 get_eival (self, val)
 
 get_eival0 (self)
 
 get_char (self, val)
 
 get_all_char (self)
 
 get_observable (self, obs, val=None)
 
 set_observable_value (self, obs, bandval, obsval)
 
 subset (self, sel)
 
 subset_inplace (self, sel)
 
 select_llindex (self, ll)
 
 select_bindex (self, b)
 
 select_obs (self, obs, val, accuracy=None)
 
 select_eival (self, val)
 
 select_char (self, which, inplace=False)
 
 sort_by_eival (self, inplace=False, reverse=False)
 
 sort_by_obs (self, obs, inplace=False)
 
 set_eivec_phase (self, accuracy=1e-6, inplace=False)
 
 get_eivec_coeff (self, norbitals, accuracy=1e-6, ll_full=False, ny=None)
 
 set_char (self, chardata, eival=None, llindex=None, eival_accuracy=1e-6)
 
 set_bindex (self, bindexdata, eival=None, llindex=None, aligned_with_e0=False)
 
 set_llindex (self, llindex)
 
 set_eivec (self, eivec, val=None)
 
 filter_transitions (self, ee, broadening=None, ampmin=100, inplace=False)
 
 to_binary_file (self, filename)
 
 __getitem__ (self, i)
 Compatibility / legacy functions; remove later.
 
 __len__ (self)
 
- Public Member Functions inherited from kdotpy.types.DiagDataPoint

Public Attributes

int k = k
 
list eival = eival
 
int neig = len(eival)
 
 eivec = None
 
 dim = None
 
 obsvals = None
 
 bindex = None
 
 llindex = None
 
bool aligned_with_e0 = False
 
list char = None
 
 transitions = None
 
 wffigure = None
 
int paramval = None
 
 current_step = None
 
 ham = None
 
 grid_index = grid_index
 
dict tuple_index = None
 
 opts = opts
 
int char = 1:
 
int tuple_index = 0:
 
int bindex = 0
 

Protected Attributes

list _obsids = None
 

Additional Inherited Members

- Static Public Attributes inherited from kdotpy.types.DiagDataPoint
Optional eival [np.ndarray]
 
Optional eivec [np.ndarray]
 
Optional dim [int]
 
Optional obsvals [np.ndarray]
 
Optional bindex [np.ndarray]
 
Optional llindex [np.ndarray]
 
Union char [np.ndarray, list, None]
 
Optional current_step [int]
 
Optional tuple_index [dict]
 
- Static Protected Attributes inherited from kdotpy.types.DiagDataPoint
Optional _obsids [list[str]]
 

Detailed Description

Container class for eigenvalue and eigenstate properties for a single k or B point.

    Attributes:
    k             Float or Vector instance. Momentum value.
    paramval      None, float or Vector instance. 'Parameter' value, currently
                  used only as magnetic field value.
    neig          Integer. Number of eigenvalues stored at this point.
    dim           Integer. Dimensionality, i.e., the size of the eigenvectors.
    eival         Numpy array or length neig. The eigenvalues in meV.
    eivec         Numpy array of shape (dim, neig). The eigenvectors belonging
                  to the eigenvectors. May be set to None to save memory
                  consumption.
    obsids        List of strings. These contain the observable ids for the
                  observables that are stored in obsvals.
    obsvals       Numpy array of shape (nobs, neig), where nobs = len(obsids).
                  The array contains the observable values for the eigenstates.
    bindex        Numpy array of length neig with integer elements. The band
                  indices belonging to the eigenstates.
    llindex       Numpy array of length neig with integer elements, or None. The
                  Landau-level indices belonging to the eigenstates.
    char          Numpy array of length neig with string elements, or None. The
                  band characters of the eigenstates.
    transitions   TransitionsData instance or None. Stores data for optical
                  transitions. See transitions.py.
    wffigure      Integer, string, or matplotlib Figure object. Identifier for
                  a matplotlib figure. None is used for absence of a figure. The
                  value is a list if separate figures are made for each state.
    current_step  Integer. Progress for this instance in the calculation model.
    ham           Sparse matrix (or tuple of matrices) for Hamiltonian(s)
                  evaluated for this instance's parameters.
    grid_index    Integer. Position of this instance in the flattened
                  VectorGrid. Used for priority ordering during diagonalization.
    tuple_index   Dict instance. Stores a mapping of band indices to array
                  indices.

Constructor & Destructor Documentation

◆ __init__()

kdotpy.diagonalization.diagdata.DiagDataPoint.__init__ ( self,
k,
eival = None,
eivec = None,
paramval = None,
grid_index = None,
opts = None )

Member Function Documentation

◆ __getitem__()

kdotpy.diagonalization.diagdata.DiagDataPoint.__getitem__ ( self,
i )

Compatibility / legacy functions; remove later.

◆ __len__()

kdotpy.diagonalization.diagdata.DiagDataPoint.__len__ ( self)

◆ __str__()

kdotpy.diagonalization.diagdata.DiagDataPoint.__str__ ( self)
Return human readable description.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ add_observable()

kdotpy.diagonalization.diagdata.DiagDataPoint.add_observable ( self,
obsvals = None,
obsid = None )
Add the values and id of an observable.

        Arguments:
        obsvals  Numpy array or None. If set, the observable values that will be
                 added. This array must have length neig. If None, add "NaN"
                 values.
        obsid    String or None. If set, add this observable id. If None, add an
                 empty string as observable id for the new observable.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ build_tuple_index_cache()

kdotpy.diagonalization.diagdata.DiagDataPoint.build_tuple_index_cache ( self)
Build and store a dict instance which maps tuple indices to array indices.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ calculate_observables()

kdotpy.diagonalization.diagdata.DiagDataPoint.calculate_observables ( self,
params,
obs,
obs_prop = None,
overlap_eivec = None,
magn = None )
Calculate observables.

        Arguments:
        params    PhysParams instance. Needed to calculate the observables, see
                  observables.py.
        obs       List of strings. The observables that will be calculated.
        obs_prop  ObservableList instance containing all observable properties.
        overlap_eivec  None or a dict instance, whose keys are the band labels
                       and values are the eigenvectors (numpy arrays). If set,
                       calculate overlap observables. If None, no overlap
                       observables are calculated.
        magn      Float, Vector instance, or None. If not None, the magnetic
                  field strength

        Note:
        Eigenvectors are required, i.e., self.eivec must not be None.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ delete_eivec()

kdotpy.diagonalization.diagdata.DiagDataPoint.delete_eivec ( self)
Delete the eigenvector data

Reimplemented from kdotpy.types.DiagDataPoint.

◆ extend()

kdotpy.diagonalization.diagdata.DiagDataPoint.extend ( self,
* args,
** kwds )
Extend data point; deal with either a DiagDataPoint or separate arguments.

        Argument:
        *args      Either a DiagDataPoint or an argument list that is passed to
                   self.extend_by().
        **kwds     Keyword arguments passed to self.extend_by().

        Note:
        If the first argument is a DiagDataPoint, all following arguments and
        keyword arguments are ignored.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ extend_by()

kdotpy.diagonalization.diagdata.DiagDataPoint.extend_by ( self,
k,
eival,
eivec,
paramval = None,
obsvals = None,
obsids = None,
char = None,
llindex = None,
bindex = None,
accuracy = 1e-6 )
Extend DiagDataPoint with additional states; prevent duplicates.

        Arguments:
        k, ...   See DiagDataPoint class info. The data that is added.

        Note:
        Arguments k and paramval serve as a check that the momentum and
        parameter value of the added data match that of the existing data. If
        not, an error is raised.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ file_id()

kdotpy.diagonalization.diagdata.DiagDataPoint.file_id ( self)
Provides a human readable id derived from k and paramval.

        Use the default str() method (no optional arguments) for the type (being
        either Vector or float).

Reimplemented from kdotpy.types.DiagDataPoint.

◆ filter_transitions()

kdotpy.diagonalization.diagdata.DiagDataPoint.filter_transitions ( self,
ee,
broadening = None,
ampmin = 100,
inplace = False )
Filter transitions by energy

        See DiagData.filter_transitions() and TransitionsData.at_energy() for
        more information.

◆ get_all_char()

kdotpy.diagonalization.diagdata.DiagDataPoint.get_all_char ( self)
Get all band characters.

        Returns:
        A dict, whose keys are the character labels and whose values are the
        (energy) eigenvalues.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ get_char()

kdotpy.diagonalization.diagdata.DiagDataPoint.get_char ( self,
val )
Look for state and return band character

Reimplemented from kdotpy.types.DiagDataPoint.

◆ get_eival()

kdotpy.diagonalization.diagdata.DiagDataPoint.get_eival ( self,
val )
Look for state and return eigenvalue.

        Argument:
        val    Any value that self.get_index can handle. Specifically, if val is
               a float, then return the eigenvalue closest to that value.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ get_eival0()

kdotpy.diagonalization.diagdata.DiagDataPoint.get_eival0 ( self)
Get energy of charge neutrality (using universal band indices)

Reimplemented from kdotpy.types.DiagDataPoint.

◆ get_eivec_coeff()

kdotpy.diagonalization.diagdata.DiagDataPoint.get_eivec_coeff ( self,
norbitals,
accuracy = 1e-6,
ll_full = False,
ny = None )
Get complex coefficients for each orbital, for each eigenvector
        The coefficients are extracted for each orbital as the eigenvector
        component where the absolute value is maximal. If this happens at
        multiple locations, then choose the value at the largest index
        (equivalent to largest z value).

        Arguments:
        norbitals   6 or 8. The number of orbitals.
        accuracy    Float. The 'fuzziness' of determining which values are
                    considered maximal. This is a relative number in terms of
                    the maximal absolute value.
        ll_full     True or False. If True, take a section of the eigenvector
                    corresponding to the Landau level with the largest weight.
                    If False (default), use the full eigenvector.
        ny          None or integer. The size in the 'y direction'; for LL mode,
                    this value serves as number of LLs in the basis. Required to
                    be set if ll_full is True, otherwise it is ignored.

        Returns:
        coeff       Numpy array of shape (neig, norbitals) and type complex.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ get_index()

kdotpy.diagonalization.diagdata.DiagDataPoint.get_index ( self,
val )
Get index (position of eigenstate) in the data arrays.

        Argument:
        val   If an integer, return this value. If a float, return the index of
              the nearest eigenvalue. If a string, return the index of the state
              with this character label. If a 1-tuple, return the index of the
              state with this band index. If a 2-tuple, return the index of the
              state with this LL index and band index.

        Returns:
        Integer index (from 0 to neig-1) or None if there is no match.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ get_index_with_llindex()

kdotpy.diagonalization.diagdata.DiagDataPoint.get_index_with_llindex ( self,
val,
llindex )
Get index of state near energy with a specific LL index.

        Arguments:
        val      Float. Energy value.
        llindex  Integer. The LL index to which the search is restricted.

        Returns:
        Integer index or None.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ get_observable()

kdotpy.diagonalization.diagdata.DiagDataPoint.get_observable ( self,
obs,
val = None )
Get observable values

        Arguments:
        obs   Integer, string, or None. If integer, take the n-th observable. If
              a string, take the observable with that obsid. If None, take all
              observables.
        val   None or a value that self.get_index() can handle. If set, then
              return the observable value(s) for that state. If None, return
              values for all states.

        Returns:
        A float (if both obs and val are None) or an array of floats (1- or
        2-dimensional (as approriate for the inputs). The value None may be
        returned on error, i.e., if obs is not a valid observable and/or if val
        does not refer to a valid state.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ get_ubindex()

kdotpy.diagonalization.diagdata.DiagDataPoint.get_ubindex ( self)
Get universal band index, i.e., an array of integers != 0 increasing in energy.
        In absence of llindex, return the bindex. With llindex, take into
        account the electrons and holes for all Landau levels.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ hash_id()

kdotpy.diagonalization.diagdata.DiagDataPoint.hash_id ( self,
length = 6,
precision = '%.12e' )
Provides a stable hash value derived from k and paramval.

        Arguments:
        length      Length of the returned hexadecimal string (default: 6).
        precision   Format specifier defining the precision of input string
                    values for k and paramval (default: '%.12e')

        Returns:
        Hexadecimal hash string

Reimplemented from kdotpy.types.DiagDataPoint.

◆ obsids()

kdotpy.diagonalization.diagdata.DiagDataPoint.obsids ( self)

Reimplemented from kdotpy.types.DiagDataPoint.

◆ select_bindex()

kdotpy.diagonalization.diagdata.DiagDataPoint.select_bindex ( self,
b )
Select states with a specific band index.

        Argument:
        b     Integer. The band index.

        Returns:
        A new DiagDataPoint instance.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ select_char()

kdotpy.diagonalization.diagdata.DiagDataPoint.select_char ( self,
which,
inplace = False )
Select states by band character.

        Arguments:
        which  String or list. If a string, look for band characters that start
               with this string. If a list, match any string in the list.

        Returns:
        A new DiagDataPoint instance (inplace = False) or the present instance
        (inplace = True).

Reimplemented from kdotpy.types.DiagDataPoint.

◆ select_eival()

kdotpy.diagonalization.diagdata.DiagDataPoint.select_eival ( self,
val )
Select states by eigenvalue.

        Arguments:
        val   Number, 2-tuple, or list. If a number, match the value exactly. If
              a 2-tuple treat the two values (numeric or None) as lower and
              upper bound for a search interval. If a list, match any value in
              the list.

        Returns:
        A new DiagDataPoint instance.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ select_llindex()

kdotpy.diagonalization.diagdata.DiagDataPoint.select_llindex ( self,
ll )
Select states with a specific LL index.

        Argument:
        ll    Integer. The LL index.

        Returns:
        A new DiagDataPoint instance.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ select_obs()

kdotpy.diagonalization.diagdata.DiagDataPoint.select_obs ( self,
obs,
val,
accuracy = None )
Select states by observable value.

        Arguments:
        obs       String. The observable id.
        val       Number, 2-tuple, or list. If a number, match the value
                  exactly or approximately. If a 2-tuple treat the two values
                  (numeric or None) as lower and upper bound for a search
                  interval. If a list, match any value in the list.
        accuracy  None or positive float. Test equality with this accuracy. If
                  None, match exactly. This only applies to testing equalities,
                  i.e., if val is a number.

        Returns:
        A new DiagDataPoint instance.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ set_bindex()

kdotpy.diagonalization.diagdata.DiagDataPoint.set_bindex ( self,
bindexdata,
eival = None,
llindex = None,
aligned_with_e0 = False )
Set band indices.

        Arguments:
        bindexdata  List or array of integers. The band indices.
        eival       List/array or None. If None, set bindexdata to self.bindex
                    as is. If a list or array of numbers, then match these
                    values to the eigenvalues (self.eival).
        llindex     Integer or None. If set, match only states with this LL
                    index. Only works if eival is not None.
        aligned_with_e0  True or False. Whether the band indices were aligned
                         with the zero energy. This should be set to True if the
                         band indices were set directly from e0, or is the band
                         indices are obtained from a BandAlignPoint with
                         aligned_with_e0 set to True.

        Returns:
        The present DiagDataPoint instance.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ set_char()

kdotpy.diagonalization.diagdata.DiagDataPoint.set_char ( self,
chardata,
eival = None,
llindex = None,
eival_accuracy = 1e-6 )
Set band characters.

        Arguments:
        chardata   List or array of strings. The character data. If a DiagDataPoint
                   is given, extract all arguments in a recursive call.
        eival      List/array or None. If None, set chardata to self.char as is.
                   If a list or array of numbers, then match these values to the
                   eigenvalues (self.eival).
        llindex    Integer or None. If set, match only states with this LL
                   index. Only works if eival is not None.

        Returns:
        The present DiagDataPoint instance.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ set_eivec()

kdotpy.diagonalization.diagdata.DiagDataPoint.set_eivec ( self,
eivec,
val = None )
Set eigenvectors.

        Arguments:
        eivec   Numpy array or DiagDataPoint instance. If an array, this is the
                eigenvector data that will be set to self.eivec. If a
                DiagDataPoint instance, copy the eigenvector data from there.
        val     None or list/array of values that match using self.get_index().
                If set, the specified input data is applied to the matching
                states. If none, then the data is applied as is.

        Returns:
        The present DiagDataPoint instance.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ set_eivec_phase()

kdotpy.diagonalization.diagdata.DiagDataPoint.set_eivec_phase ( self,
accuracy = 1e-6,
inplace = False )
Multiply each eigenvector by a phase factor to fix the arbitrary phase.

        For each eigenvector, look for the largest absolute component |psi_i|
        and divide by the phase psi_i / |psi_i|. The result is that the
        resulting eigenvector will have Im(psi_i) = 0 and Re(psi_i) > 0. If
        there are multiple values psi_i of almost the same size, choose the
        largest i.

        Arguments:
        accuracy  Float. Fuzziness of determining which psi_i are considered
                  maximal. The value is relative to the maximum |psi_i|.
        inplace   True or False. Whether to return a new instance (False) or the
                  present one (True).

        Returns:
        New or present DiagDataPoint instance.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ set_llindex()

kdotpy.diagonalization.diagdata.DiagDataPoint.set_llindex ( self,
llindex )
Set band indices.

        Arguments:
        llindex  List or array of integers. The LL indices. These are set to
                 self.llindex as is.

        Returns:
        The present DiagDataPoint instance.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ set_observable_value()

kdotpy.diagonalization.diagdata.DiagDataPoint.set_observable_value ( self,
obs,
bandval,
obsval )
Set observable values to specific states.

        Arguments:
        obs      Integer or string. Observable index or id, respectively.
        bandval  Float or integer number or a list or array. If numeric, look
                 for state using self.getindex(). If a list or array, look for
                 multiple states using self.getindex().
        obsval   Float or array. The observable value(s). If an array, the shape
                 must be set appropriately.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ set_observables()

kdotpy.diagonalization.diagdata.DiagDataPoint.set_observables ( self,
obsvals,
obsids = None )
Set observable values.

        Argument:
        obsvals   List or array. Observable values.
        obsids    None, list or array. Set self._obsids to this value (also
                  valid for None).

Reimplemented from kdotpy.types.DiagDataPoint.

◆ sort_by_eival()

kdotpy.diagonalization.diagdata.DiagDataPoint.sort_by_eival ( self,
inplace = False,
reverse = False )
Sort by eigenvalues.

        Arguments:
        inplace   True or False. Whether to return a new instance (False) or the
                  present one (True).
        reverse   True or False. Reverse or standard sorting direction.

        Returns:
        New or present DiagDataPoint instance.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ sort_by_obs()

kdotpy.diagonalization.diagdata.DiagDataPoint.sort_by_obs ( self,
obs,
inplace = False )
Sort by eigenvalues.

        Arguments:
        obs       String. Observable id.
        inplace   True or False. Whether to return a new instance (False) or the
                  present one (True).

        Returns:
        New or present DiagDataPoint instance.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ stitch_with()

kdotpy.diagonalization.diagdata.DiagDataPoint.stitch_with ( self,
k,
eival,
eivec,
targetenergy_old,
targetenergy_new,
inplace = False,
accuracy = 0.01 )
Stitch together multiple diagonalization solutions.
        Overlapping duplicate eigenvalues are recalculated by weighted mean.
        From duplicate eigenvectors, the one with the eigenvalue closer to its
        its target value (higher weight) is chosen.

        Arguments:
        k, ...          See DiagDataPoint class info. The data that is added.
        targetenergies  Target energies that have been used to calculate
                        solutions. Used to calculate weights.
        inplace         Replace eigenvectors of current DiagDataPoint instance
                        if True, otherwise return a new instance (default).
accuracy        Estimate of solver precision. Used to determine
                degeneracy of states.

        Note:
        Currently only supports bare diagonalization results without
        observables, ll indices, etc. as those quantities are usually not yet
        calculated. This DiagDataPoint is expected to be sorted by eival.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ subset()

kdotpy.diagonalization.diagdata.DiagDataPoint.subset ( self,
sel )
Take subset; can also be used for reordering

        Argument:
        sel   Integer or array. Anything that can be used as index to a numpy
              array.

        Returns:
        A new DiagDataPoint instance.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ subset_inplace()

kdotpy.diagonalization.diagdata.DiagDataPoint.subset_inplace ( self,
sel )
Take subset and discard other states.

        Argument:
        sel   Integer or array. Anything that can be used as index to a numpy
              array.

        Returns:
        The present DiagDataPoint instance with only the selected states.

Reimplemented from kdotpy.types.DiagDataPoint.

◆ to_binary_file()

kdotpy.diagonalization.diagdata.DiagDataPoint.to_binary_file ( self,
filename )
Save data to a binary file (Numpy npz or HDF5) file.
        This function saves all fields (member variables) specified in global
        variable binfile_ddp_fields as well as the x values (momentum and/or
        parameter value).

        For Numpy format: The file is a compressed npz file with a collection of
        numpy arrays. For more information on the file format, consult:
        https://numpy.org/doc/stable/reference/generated/numpy.lib.format.html

        For HDF5 format: The file is a HDF5 container and the data is saved in a
        separate group for each DiagDataPoint. The values for k and b are stored
        as attributes. We do not use compression because it would reduce the
        file size only minimally. See also: https://docs.h5py.org

        Argument:
        filename   String. The file name. The output type is extracted from the
                   file name extension.

        No return value

Reimplemented from kdotpy.types.DiagDataPoint.

◆ update()

kdotpy.diagonalization.diagdata.DiagDataPoint.update ( self,
new_ddp )
Update the current DiagDataPoint instance from another instance.
        This is useful if the current instance is already linked to a DiagData
        instance. Keeps the attributes 'grid_index' and 'current_step' from the
        current instance if not set in the new instance.

Reimplemented from kdotpy.types.DiagDataPoint.

Member Data Documentation

◆ _obsids

list kdotpy.diagonalization.diagdata.DiagDataPoint._obsids = None
protected

◆ aligned_with_e0

bool kdotpy.diagonalization.diagdata.DiagDataPoint.aligned_with_e0 = False

◆ bindex [1/2]

kdotpy.diagonalization.diagdata.DiagDataPoint.bindex = None

◆ bindex [2/2]

int kdotpy.diagonalization.diagdata.DiagDataPoint.bindex = 0

◆ char [1/2]

list kdotpy.diagonalization.diagdata.DiagDataPoint.char = None

◆ char [2/2]

int kdotpy.diagonalization.diagdata.DiagDataPoint.char = 1:

◆ current_step

kdotpy.diagonalization.diagdata.DiagDataPoint.current_step = None

◆ dim

kdotpy.diagonalization.diagdata.DiagDataPoint.dim = None

◆ eival

list kdotpy.diagonalization.diagdata.DiagDataPoint.eival = eival

◆ eivec

kdotpy.diagonalization.diagdata.DiagDataPoint.eivec = None

◆ grid_index

kdotpy.diagonalization.diagdata.DiagDataPoint.grid_index = grid_index

◆ ham

kdotpy.diagonalization.diagdata.DiagDataPoint.ham = None

◆ k

int kdotpy.diagonalization.diagdata.DiagDataPoint.k = k

◆ llindex

kdotpy.diagonalization.diagdata.DiagDataPoint.llindex = None

◆ neig

int kdotpy.diagonalization.diagdata.DiagDataPoint.neig = len(eival)

◆ obsvals

kdotpy.diagonalization.diagdata.DiagDataPoint.obsvals = None

◆ opts

kdotpy.diagonalization.diagdata.DiagDataPoint.opts = opts

◆ paramval

int kdotpy.diagonalization.diagdata.DiagDataPoint.paramval = None

◆ transitions

kdotpy.diagonalization.diagdata.DiagDataPoint.transitions = None

◆ tuple_index [1/2]

dict kdotpy.diagonalization.diagdata.DiagDataPoint.tuple_index = None

◆ tuple_index [2/2]

int kdotpy.diagonalization.diagdata.DiagDataPoint.tuple_index = 0:

◆ wffigure

kdotpy.diagonalization.diagdata.DiagDataPoint.wffigure = None

The documentation for this class was generated from the following file: