ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
hamilt::HContainer< T > Class Template Reference

#include <hcontainer.h>

Collaboration diagram for hamilt::HContainer< T >:

Public Member Functions

 ~HContainer ()
 
 HContainer ()
 
 HContainer (const HContainer< T > &HR_in, T *data_array=nullptr)
 copy constructor when data_array is not nullptr, new HContainer will be wrapper for data_array data of HR_in will not be copied, please call add() after this constructor to copy data.
 
HContaineroperator= (const HContainer< T > &HR_in)=delete
 
 HContainer (HContainer< T > &&HR_in) noexcept
 
HContaineroperator= (HContainer< T > &&HR_in) noexcept
 
 HContainer (int natom)
 
 HContainer (const UnitCell &ucell_, const Parallel_Orbitals *paraV=nullptr)
 
 HContainer (const Parallel_Orbitals *paraV, T *data_pointer=nullptr, const std::vector< int > *ijr_info=nullptr)
 use 2d-block-recycle parallel case to initialize atom_pairs, mainly used now. pass a data pointer to HContainer, which means HContainer is a wrapper it will not allocate memory for atom_pairs this case will forbit inserting empty atom_pair
 
void set_paraV (const Parallel_Orbitals *paraV_in)
 set parallel orbital pointer to check parallel information
 
const Parallel_Orbitalsget_paraV () const
 get parallel orbital pointer to check parallel information
 
void allocate (T *data_array=nullptr, bool if_zero=false)
 allocate memory for all <IJR> matrix if data_array is not nullptr, use memory after data_array for each BaseMatrix; if BaseMatrix has memory allocated before, it will be freed first. if data_array is nullptr, allocate memory for each BaseMatrix
 
void set_zero ()
 set values of all <IJR> matrix to zero
 
void insert_pair (const AtomPair< T > &atom_ij)
 a AtomPair object can be inserted into HContainer, two steps: 1, find AtomPair with atom index atom_i and atom_j 2.1, if found, add to exist AtomPair, 2.2, if not found, insert new one and sort.
 
AtomPair< T > * find_pair (int i, int j) const
 find AtomPair with atom index atom_i and atom_j This interface can be used to find AtomPair, if found, return pointer will be the exist one, if not found, return pointer will be nullptr.
 
BaseMatrix< T > * find_matrix (int i, int j, int rx, int ry, int rz)
 find BaseMatrix with atom index atom_i and atom_j and R index (rx, ry, rz) This interface can be used to find BaseMatrix in AtomPair, if found, return pointer will be the exist one, if not found, return pointer will be nullptr.
 
const BaseMatrix< T > * find_matrix (int i, int j, int rx, int ry, int rz) const
 
BaseMatrix< T > * find_matrix (int i, int j, const ModuleBase::Vector3< int > &R_index)
 
const BaseMatrix< T > * find_matrix (int i, int j, const ModuleBase::Vector3< int > &R_index) const
 
int find_matrix_offset (int i, int j, int rx, int ry, int rz) const
 find the offset of BaseMatrix with atom index atom_i and atom_j and R index (rx, ry, rz) if found, return this->find_matrix(i, j, rx, ry, rz)->get_pointer() - this->get_wrapper(); if not found, return -1
 
int find_matrix_offset (int i, int j, const ModuleBase::Vector3< int > &R_index) const
 
AtomPair< T > & get_atom_pair (int i, int j) const
 return a reference of AtomPair with index of atom I and J in atom_pairs
 
AtomPair< T > & get_atom_pair (int index) const
 return a reference of AtomPair with index in atom_pairs (R is not fixed) tmp_atom_pairs (R is fixed)
 
void add (const HContainer< T > &other)
 operator() for accessing value with indexes this interface is not implemented now, because it is too expensive to access data
 
bool fix_R (int rx_in, int ry_in, int rz_in) const
 save atom-pair pointers into this->tmp_atom_pairs for selected R index
 
void unfix_R () const
 set current_R to -1, which means R is not fixed clear this->tmp_atom_pairs
 
void fix_gamma ()
 restrict R indexes of all atom-pair to 0, 0, 0 add BaseMatrix<T> with non-zero R index to this->atom_pairs[i].values[0] set gamma_only = true in this mode:
 
void loop_R (const size_t &index, int &rx, int &ry, int &rz) const
 interface for call a R loop for HContainer it can return a new R-index with (rx,ry,rz) for each loop if index==0, a new loop of R will be initialized
 
size_t size_R_loop () const
 calculate number of R index which has counted AtomPairs
 
int find_R (const int &rx_in, const int &ry_in, const int &rz_in) const
 find index of R in tmp_R_index, used when current_R is fixed
 
int find_R (const ModuleBase::Vector3< int > &R_in) const
 
size_t size_atom_pairs () const
 calculate number of AtomPairs for current R index
 
Tdata (int i, int j) const
 get data pointer of AtomPair with index of I, J
 
Tdata (int i, int j, int *R) const
 get data pointer of AtomPair with index of I, J, Rx, Ry, Rz
 
int get_current_R () const
 get current_R
 
bool is_gamma_only () const
 judge if gamma_only
 
size_t get_memory_size () const
 get total memory bites of HContainer
 
size_t get_nnr () const
 calculate total size of data in HContainer, named nnr inherited from history all AtomPairs and BaseMatrixs are counted
 
std::vector< int > get_ijr_info () const
 get infomation of IJR pairs in HContainer the return vector format is {size_I, I1, size_J, J1, size_R, R1x, R1y, R1z, ..., J2, ...}
 
void insert_ijrs (const std::vector< int > *ijrs)
 use infomation of IJ pairs to expand HContainer the input vector format is {size_IJ_pairs, I1, J1, size_R, R1x, R1y, R1z, ..., I2, J2, ...} HContainer has not been allocated after this function, user should call allocate(...) to allocate memory.
 
void insert_ijrs (const std::vector< int > *ijrs, const UnitCell &ucell, const int npol=1)
 use infomation of IJ pairs to expand HContainer the number of wavefunctions are stored in UnitCell. HContainer has not been allocated after this function, user should call allocate(...) to allocate memory.
 
Tget_wrapper () const
 return the wrapper_pointer
 
void shape_synchron (const HContainer< T > &other)
 synchronization of atom-pairs for read-in HContainer new <IJR> pair from read-in HContainer will be inserted into this->atom-pairs
 
const std::vector< std::vector< int > > & get_sparse_ap () const
 get sparse_ap
 
const std::vector< std::vector< int > > & get_sparse_ap_index () const
 get sparse_ap_index
 
int get_nbasis () const
 get number of basis in each H matrix
 

Private Attributes

std::vector< AtomPair< T > > atom_pairs
 
std::vector< std::vector< int > > sparse_ap
 
std::vector< std::vector< int > > sparse_ap_index
 
std::vector< const AtomPair< T > * > tmp_atom_pairs
 temporary atom-pair lists to loop selected R index
 
std::vector< ModuleBase::Vector3< int > > tmp_R_index
 
int current_R = -1
 
bool gamma_only = false
 
Twrapper_pointer = nullptr
 if wrapper_pointer is not nullptr, this HContainer is a wrapper there is only one case that "wrapper_pointer == nullptr": HContainer is constructed by allocated AtomPairs by insert_pair() in other cases, wrapper_pointer is not nullptr and is memory pool of AtomPairs
 
bool allocated = false
 
const Parallel_OrbitalsparaV = nullptr
 pointer of Parallel_Orbitals, which is used to get atom-pair information
 

Detailed Description

template<typename T>
class hamilt::HContainer< T >

class HContainer

used to store a matrix for atom-pair local Hamiltonian with specific R-index

<Psi_{mu_I,0}|H|Psi_{nu_J,R}>

template T can be double or std::complex<double>

examples for using this class:

  1. initialize a HContainer a. use unitcell to initialize atom_pairs
    // ucell is a UnitCell object
    // in this method, all empty atom-pairs will be initialized in HR
    // atom-pairs are sorted by matrix of (i, j)
    HContainer<double> HR(ucell, nullptr); // serial case
    HContainer<double> HR(ucell, paraV); // 2d-block-recycle parallel case
    Definition hcontainer.h:144
    const Parallel_Orbitals * paraV
    pointer of Parallel_Orbitals, which is used to get atom-pair information
    Definition hcontainer.h:502
    b. use insert_pair() to insert atom-pair
    AtomPair<double> atom_ij...
    HR.insert_pair(atom_ij);
    // allocate is nessasary if atom-pairs are inserted
    HR.allocate(nullptr, 0); // arrange memory by HContainer
    HR.allocate(data_array, 0); // use data_array as memory pool
    Definition atom_pair.h:42
    void allocate(T *data_array=nullptr, bool if_zero=false)
    allocate memory for all <IJR> matrix if data_array is not nullptr, use memory after data_array for ea...
    Definition hcontainer.cpp:175
    c. use Parallel_Orbital to initialize atom_pairs and HContainer
    // paraV is a Parallel_Orbital object, 2d-block-recycle parallel case
    HCotainer<double> HR(paraV);
    // initialize a new AtomPair with atom paraV
    AtomPair<double> atom_ij(0, 1, paraV);
    // insert atom_ij into HR
    HR.insert_pair(atom_ij);
    @icode
    2. get target AtomPair with index of atom I and J, or with index in atom_pairs
    a. use interface find_pair() to get pointer of target AtomPair
    std::vector< AtomPair< T > > atom_pairs
    Definition hcontainer.h:473
    void insert_pair(const AtomPair< T > &atom_ij)
    a AtomPair object can be inserted into HContainer, two steps: 1, find AtomPair with atom index atom_i...
    Definition hcontainer.cpp:624
    AtomPair< T > * find_pair(int i, int j) const
    find AtomPair with atom index atom_i and atom_j This interface can be used to find AtomPair,...
    Definition hcontainer.cpp:219
    // HR is a HContainer object AtomPair<double>* atom_ij = HR.find_pair(0, 1); // check if atom_ij is found if (atom_ij != nullptr) { // do something }
    b. use interface get_atom_pair() to get reference of target AtomPair
    AtomPair< T > & get_atom_pair(int i, int j) const
    return a reference of AtomPair with index of atom I and J in atom_pairs
    Definition hcontainer.cpp:353
    // HR is a HContainer object AtomPair<double>& atom_ij_1 = HR.get_atom_pair(0, 1); AtomPair<double>& atom_ij_2 = HR.get_atom_pair(1);//suppose 0,1 is the second atom-pair in atom_pairs
    3. get data pointer of target local matrix <Psi_{mu_I,R}|H|Psi_{nu_J,0}>
    a. use interface data() with atom_i and atom_j and R index
    T * data(int i, int j) const
    get data pointer of AtomPair with index of I, J
    Definition hcontainer.cpp:592
    // HR is a HContainer object // suppose atom_i = 0, atom_j = 1, int[3] target_R = {0, 0, 0} double* target_data = HR.data(0, 1, target_R);
    b. fix_R and use data() with atom_i and atom_j without R index
    bool fix_R(int rx_in, int ry_in, int rz_in) const
    save atom-pair pointers into this->tmp_atom_pairs for selected R index
    Definition hcontainer.cpp:417
    HR.fix_R(0, 0, 0); double* target_data = HR.data(0, 1); HR.unfix_R();
    4. use for-loop to do something with atom-pairs with specific R index
    a. loop R-index first and then loop atom-pairs
    // HR is a const HContainer object, which has been initialized int rx, ry, rz; // call size_R_loop() to get number of different R indexes, // be careful, it will cost some time to loop all atom-pairs to gather R indexes int size_for_loop_R = HR.size_R_loop(); for (int iR = 0; iR < size_for_loop_R ; iR++) { // call loop_R() to get R coordinates (rx, ry, rz) HR.loop_R(iR, rx, ry, rz); // call fix_R() to save atom-pairs with R index (rx, ry, rz) into tmp_atom_pairs HR.fix_R(rx, ry, rz); // loop fixed atom-pairs for (int i = 0; i < HR.size_atom_pairs(); i++) { // get pointer of target atom-pair double* data_pointer = HR.data(i); // or get reference of target atom-pair AtomPair<double>& atom_ijR = HR.get_atom_pair(i); // do something with atom_ijR or data_pointer ... } } // call unfix_R() to clear tmp_atom_pairs HR.unfix_R();
    b. loop atom-pairs first and then loop R-index
    // HR is a const HContainer object, which has been initialized // loop atom-pairs for (int i = 0; i < HR.size_atom_pairs(); i++) { // get reference of target atom-pair AtomPair<double>& atom_ij = HR.get_atom_pair(i); // loop R-index for (int iR = 0; iR < atom_ij.size_R(); iR++) { const ModuleBase::Vector3<int> r_index = atom_ij.get_R_index(iR); auto tmp_matrix = atom_ij.get_HR_values(r_index.x, r_index.y, r_index.z); // do something with tmp_matrix ... } }
    c. loop atom-pairs with gamma_only case
    bool gamma_only
    Definition hcontainer.h:488
    ... HR.fix_gamma(); // HR is a const HContainer object, which has been initialized and fixed to gamma // loop atom-pairs directly without R index for (int i = 0; i < HR.size_atom_pairs(); i++) { // get data pointer of target atom-pair double* data_pointer = HR.get_pointer(i); // do something with data_pointer ... }

Constructor & Destructor Documentation

◆ ~HContainer()

template<typename T >
hamilt::HContainer< T >::~HContainer ( )

◆ HContainer() [1/6]

template<typename T >
hamilt::HContainer< T >::HContainer ( )

◆ HContainer() [2/6]

template<typename T >
hamilt::HContainer< T >::HContainer ( const HContainer< T > &  HR_in,
T data_array = nullptr 
)

copy constructor when data_array is not nullptr, new HContainer will be wrapper for data_array data of HR_in will not be copied, please call add() after this constructor to copy data.

◆ HContainer() [3/6]

template<typename T >
hamilt::HContainer< T >::HContainer ( HContainer< T > &&  HR_in)
noexcept

◆ HContainer() [4/6]

template<typename T >
hamilt::HContainer< T >::HContainer ( int  natom)

◆ HContainer() [5/6]

template<typename T >
hamilt::HContainer< T >::HContainer ( const UnitCell ucell_,
const Parallel_Orbitals paraV = nullptr 
)
Here is the call graph for this function:

◆ HContainer() [6/6]

template<typename T >
hamilt::HContainer< T >::HContainer ( const Parallel_Orbitals paraV,
T data_pointer = nullptr,
const std::vector< int > *  ijr_info = nullptr 
)

use 2d-block-recycle parallel case to initialize atom_pairs, mainly used now. pass a data pointer to HContainer, which means HContainer is a wrapper it will not allocate memory for atom_pairs this case will forbit inserting empty atom_pair

Member Function Documentation

◆ add()

template<typename T >
void hamilt::HContainer< T >::add ( const HContainer< T > &  other)

operator() for accessing value with indexes this interface is not implemented now, because it is too expensive to access data

Parameters
atom_iindex of atom i
atom_jindex of atom j
rx_inindex of R in x direction
ry_inindex of R in y direction
rz_inindex of R in z direction
muindex of orbital mu
nuindex of orbital nu
Returns
T&

add another HContainer to this HContainer

Here is the call graph for this function:
Here is the caller graph for this function:

◆ allocate()

template<typename T >
void hamilt::HContainer< T >::allocate ( T data_array = nullptr,
bool  if_zero = false 
)

allocate memory for all <IJR> matrix if data_array is not nullptr, use memory after data_array for each BaseMatrix; if BaseMatrix has memory allocated before, it will be freed first. if data_array is nullptr, allocate memory for each BaseMatrix

Parameters
data_arraypointer of data array
if_zeroif true, set all values to zero
Here is the caller graph for this function:

◆ data() [1/2]

template<typename T >
T * hamilt::HContainer< T >::data ( int  i,
int  j 
) const

get data pointer of AtomPair with index of I, J

Parameters
iindex of atom i
jindex of atom j
Returns
T* pointer of data
Here is the call graph for this function:
Here is the caller graph for this function:

◆ data() [2/2]

template<typename T >
T * hamilt::HContainer< T >::data ( int  i,
int  j,
int *  R 
) const

get data pointer of AtomPair with index of I, J, Rx, Ry, Rz

Parameters
iindex of atom i
jindex of atom j
Rint[3] of R index
Returns
T* pointer of data
Here is the call graph for this function:

◆ find_matrix() [1/4]

template<typename T >
BaseMatrix< T > * hamilt::HContainer< T >::find_matrix ( int  i,
int  j,
const ModuleBase::Vector3< int > &  R_index 
)

◆ find_matrix() [2/4]

template<typename T >
const BaseMatrix< T > * hamilt::HContainer< T >::find_matrix ( int  i,
int  j,
const ModuleBase::Vector3< int > &  R_index 
) const
Here is the call graph for this function:

◆ find_matrix() [3/4]

template<typename T >
BaseMatrix< T > * hamilt::HContainer< T >::find_matrix ( int  i,
int  j,
int  rx,
int  ry,
int  rz 
)

find BaseMatrix with atom index atom_i and atom_j and R index (rx, ry, rz) This interface can be used to find BaseMatrix in AtomPair, if found, return pointer will be the exist one, if not found, return pointer will be nullptr.

Here is the caller graph for this function:

◆ find_matrix() [4/4]

template<typename T >
const BaseMatrix< T > * hamilt::HContainer< T >::find_matrix ( int  i,
int  j,
int  rx,
int  ry,
int  rz 
) const

◆ find_matrix_offset() [1/2]

template<typename T >
int hamilt::HContainer< T >::find_matrix_offset ( int  i,
int  j,
const ModuleBase::Vector3< int > &  R_index 
) const

◆ find_matrix_offset() [2/2]

template<typename T >
int hamilt::HContainer< T >::find_matrix_offset ( int  i,
int  j,
int  rx,
int  ry,
int  rz 
) const

find the offset of BaseMatrix with atom index atom_i and atom_j and R index (rx, ry, rz) if found, return this->find_matrix(i, j, rx, ry, rz)->get_pointer() - this->get_wrapper(); if not found, return -1

Here is the call graph for this function:
Here is the caller graph for this function:

◆ find_pair()

template<typename T >
AtomPair< T > * hamilt::HContainer< T >::find_pair ( int  i,
int  j 
) const

find AtomPair with atom index atom_i and atom_j This interface can be used to find AtomPair, if found, return pointer will be the exist one, if not found, return pointer will be nullptr.

Parameters
iatom index of atom i
jatom index of atom j
Returns
AtomPair<T>*
Here is the call graph for this function:
Here is the caller graph for this function:

◆ find_R() [1/2]

template<typename T >
int hamilt::HContainer< T >::find_R ( const int &  rx_in,
const int &  ry_in,
const int &  rz_in 
) const

find index of R in tmp_R_index, used when current_R is fixed

Parameters
rx_inindex of R in x direction
ry_inindex of R in y direction
rz_inindex of R in z direction
Returns
int index of R in tmp_R_index

◆ find_R() [2/2]

template<typename T >
int hamilt::HContainer< T >::find_R ( const ModuleBase::Vector3< int > &  R_in) const

◆ fix_gamma()

template<typename T >
void hamilt::HContainer< T >::fix_gamma ( )

restrict R indexes of all atom-pair to 0, 0, 0 add BaseMatrix<T> with non-zero R index to this->atom_pairs[i].values[0] set gamma_only = true in this mode:

  1. fix_R() can not be used
  2. there is no interface to set gamma_only = false, user should create a new HContainer if needed
  3. get_size_for_loop_R() and loop_R() can not be used
  4. get_AP_size() can be used
  5. data(i, j) can be used to get pointer of target atom-pair with R = 0, 0, 0 , data(i,j,R) can not be used
  6. insert_pair() can be safely used, but the R index will be ignored
  7. find_matrix() can be safely used, but the R index will be ignored
  8. operator() can be used, but the R index will be ignored
  9. get_atom_pair(), find_atom_pair() can be used, be careful that AtomPair::get_HR_values() is dangerous in this mode.
Here is the caller graph for this function:

◆ fix_R()

template<typename T >
bool hamilt::HContainer< T >::fix_R ( int  rx_in,
int  ry_in,
int  rz_in 
) const

save atom-pair pointers into this->tmp_atom_pairs for selected R index

Parameters
rx_inindex of R in x direction
ry_inindex of R in y direction
rz_inindex of R in z direction
Returns
true if success
Here is the caller graph for this function:

◆ get_atom_pair() [1/2]

template<typename T >
AtomPair< T > & hamilt::HContainer< T >::get_atom_pair ( int  i,
int  j 
) const

return a reference of AtomPair with index of atom I and J in atom_pairs

Parameters
iindex of atom i
jindex of atom j
Here is the caller graph for this function:

◆ get_atom_pair() [2/2]

template<typename T >
AtomPair< T > & hamilt::HContainer< T >::get_atom_pair ( int  index) const

return a reference of AtomPair with index in atom_pairs (R is not fixed) tmp_atom_pairs (R is fixed)

Parameters
indexindex of atom-pair

◆ get_current_R()

template<typename T >
int hamilt::HContainer< T >::get_current_R ( ) const

get current_R

◆ get_ijr_info()

template<typename T >
std::vector< int > hamilt::HContainer< T >::get_ijr_info ( ) const

get infomation of IJR pairs in HContainer the return vector format is {size_I, I1, size_J, J1, size_R, R1x, R1y, R1z, ..., J2, ...}

Here is the caller graph for this function:

◆ get_memory_size()

template<typename T >
size_t hamilt::HContainer< T >::get_memory_size ( ) const

get total memory bites of HContainer

Here is the caller graph for this function:

◆ get_nbasis()

template<typename T >
int hamilt::HContainer< T >::get_nbasis ( ) const
inline

get number of basis in each H matrix

Returns
int
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_nnr()

template<typename T >
size_t hamilt::HContainer< T >::get_nnr ( ) const
inline

calculate total size of data in HContainer, named nnr inherited from history all AtomPairs and BaseMatrixs are counted

Here is the caller graph for this function:

◆ get_paraV()

template<typename T >
const Parallel_Orbitals * hamilt::HContainer< T >::get_paraV ( ) const
inline

get parallel orbital pointer to check parallel information

Returns
const Parallel_Orbitals* , if return is nullptr, it means HContainer is not in parallel mode
Here is the caller graph for this function:

◆ get_sparse_ap()

template<typename T >
const std::vector< std::vector< int > > & hamilt::HContainer< T >::get_sparse_ap ( ) const
inline

get sparse_ap

Returns
std::vector<std::vector<int>>&

◆ get_sparse_ap_index()

template<typename T >
const std::vector< std::vector< int > > & hamilt::HContainer< T >::get_sparse_ap_index ( ) const
inline

get sparse_ap_index

Returns
std::vector<std::vector<int>>&

◆ get_wrapper()

template<typename T >
T * hamilt::HContainer< T >::get_wrapper ( ) const
inline

return the wrapper_pointer

Here is the caller graph for this function:

◆ insert_ijrs() [1/2]

template<typename T >
void hamilt::HContainer< T >::insert_ijrs ( const std::vector< int > *  ijrs)

use infomation of IJ pairs to expand HContainer the input vector format is {size_IJ_pairs, I1, J1, size_R, R1x, R1y, R1z, ..., I2, J2, ...} HContainer has not been allocated after this function, user should call allocate(...) to allocate memory.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ insert_ijrs() [2/2]

template<typename T >
void hamilt::HContainer< T >::insert_ijrs ( const std::vector< int > *  ijrs,
const UnitCell ucell,
const int  npol = 1 
)

use infomation of IJ pairs to expand HContainer the number of wavefunctions are stored in UnitCell. HContainer has not been allocated after this function, user should call allocate(...) to allocate memory.

◆ insert_pair()

template<typename T >
void hamilt::HContainer< T >::insert_pair ( const AtomPair< T > &  atom_ij)

a AtomPair object can be inserted into HContainer, two steps: 1, find AtomPair with atom index atom_i and atom_j 2.1, if found, add to exist AtomPair, 2.2, if not found, insert new one and sort.

Parameters
atom_ijAtomPair object
Here is the call graph for this function:
Here is the caller graph for this function:

◆ is_gamma_only()

template<typename T >
bool hamilt::HContainer< T >::is_gamma_only ( ) const

judge if gamma_only

◆ loop_R()

template<typename T >
void hamilt::HContainer< T >::loop_R ( const size_t &  index,
int &  rx,
int &  ry,
int &  rz 
) const

interface for call a R loop for HContainer it can return a new R-index with (rx,ry,rz) for each loop if index==0, a new loop of R will be initialized

Parameters
indexindex of R loop
rxindex of R in x direction, would be set in the function
ryindex of R in y direction, would be set in the function
rzindex of R in z direction, would be set in the function
Here is the caller graph for this function:

◆ operator=() [1/2]

template<typename T >
HContainer & hamilt::HContainer< T >::operator= ( const HContainer< T > &  HR_in)
delete

◆ operator=() [2/2]

template<typename T >
HContainer< T > & hamilt::HContainer< T >::operator= ( HContainer< T > &&  HR_in)
noexcept

◆ set_paraV()

template<typename T >
void hamilt::HContainer< T >::set_paraV ( const Parallel_Orbitals paraV_in)
inline

set parallel orbital pointer to check parallel information

◆ set_zero()

template<typename T >
void hamilt::HContainer< T >::set_zero ( )

set values of all <IJR> matrix to zero

Here is the caller graph for this function:

◆ shape_synchron()

template<typename T >
void hamilt::HContainer< T >::shape_synchron ( const HContainer< T > &  other)

synchronization of atom-pairs for read-in HContainer new <IJR> pair from read-in HContainer will be inserted into this->atom-pairs

Here is the call graph for this function:

◆ size_atom_pairs()

template<typename T >
size_t hamilt::HContainer< T >::size_atom_pairs ( ) const

calculate number of AtomPairs for current R index

Here is the caller graph for this function:

◆ size_R_loop()

template<typename T >
size_t hamilt::HContainer< T >::size_R_loop ( ) const

calculate number of R index which has counted AtomPairs

start a new iteration of loop_R there is different R-index in this->atom_pairs[i].R_values search them one by one and store them in this->tmp_R_index

search (rx, ry, rz) with (it->R_values[i].x, it->R_values[i].y, it->R_values[i].z) if (rx, ry, rz) not found in this->tmp_R_index, insert the (rx, ry, rz) into end of this->tmp_R_index no need to sort this->tmp_R_index, using find_R() to find the (rx, ry, rz) -> int in tmp_R_index

Here is the caller graph for this function:

◆ unfix_R()

template<typename T >
void hamilt::HContainer< T >::unfix_R ( ) const

set current_R to -1, which means R is not fixed clear this->tmp_atom_pairs

Here is the caller graph for this function:

Member Data Documentation

◆ allocated

template<typename T >
bool hamilt::HContainer< T >::allocated = false
private

◆ atom_pairs

template<typename T >
std::vector<AtomPair<T> > hamilt::HContainer< T >::atom_pairs
private

◆ current_R

template<typename T >
int hamilt::HContainer< T >::current_R = -1
mutableprivate

◆ gamma_only

template<typename T >
bool hamilt::HContainer< T >::gamma_only = false
private

◆ paraV

template<typename T >
const Parallel_Orbitals* hamilt::HContainer< T >::paraV = nullptr
private

pointer of Parallel_Orbitals, which is used to get atom-pair information

◆ sparse_ap

template<typename T >
std::vector<std::vector<int> > hamilt::HContainer< T >::sparse_ap
private

◆ sparse_ap_index

template<typename T >
std::vector<std::vector<int> > hamilt::HContainer< T >::sparse_ap_index
private

◆ tmp_atom_pairs

template<typename T >
std::vector<const AtomPair<T>*> hamilt::HContainer< T >::tmp_atom_pairs
mutableprivate

temporary atom-pair lists to loop selected R index

◆ tmp_R_index

template<typename T >
std::vector<ModuleBase::Vector3<int> > hamilt::HContainer< T >::tmp_R_index
mutableprivate

◆ wrapper_pointer

template<typename T >
T* hamilt::HContainer< T >::wrapper_pointer = nullptr
private

if wrapper_pointer is not nullptr, this HContainer is a wrapper there is only one case that "wrapper_pointer == nullptr": HContainer is constructed by allocated AtomPairs by insert_pair() in other cases, wrapper_pointer is not nullptr and is memory pool of AtomPairs


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