ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
LRI_CV.h
Go to the documentation of this file.
1//=======================
2// AUTHOR : Peize Lin
3// DATE : 2022-08-17
4//=======================
5
6#ifndef LRI_CV_H
7#define LRI_CV_H
8
9#include "Matrix_Orbs11.h"
10#include "Matrix_Orbs21.h"
14
15#include <RI/global/Tensor.h>
16#include <RI/global/Global_Func-2.h>
17
18#include <vector>
19#include <map>
20#include <functional>
21#include <pthread.h>
22
23template<typename Tdata>
24class LRI_CV
25{
26private:
27 using TA = int;
28 using TC = std::array<int,3>;
29 using TAC = std::pair<TA,TC>;
30 using Tdata_real = RI::Global_Func::To_Real_t<Tdata>;
31
32public:
33 LRI_CV();
34 ~LRI_CV();
35
36 void set_orbitals(
37 const UnitCell &ucell,
38 const LCAO_Orbitals& orb,
39 const std::vector<std::vector<std::vector<Numerical_Orbital_Lm>>> &lcaos_in,
40 const std::vector<std::vector<std::vector<Numerical_Orbital_Lm>>> &abfs_in,
41 const std::vector<std::vector<std::vector<Numerical_Orbital_Lm>>> &abfs_ccp_in,
42 const double &kmesh_times,
43 ORB_gaunt_table& MGT,
44 const bool& init_MGT,
45 const bool& init_C);
46 inline std::map<TA,std::map<TAC,RI::Tensor<Tdata>>>
47 cal_Vs(
48 const UnitCell &ucell,
49 const std::vector<TA> &list_A0,
50 const std::vector<TAC> &list_A1,
51 const std::map<std::string,bool> &flags); // "writable_Vws"
52 inline std::map<TA, std::map<TAC, std::array<RI::Tensor<Tdata>, 3>>>
53 cal_dVs(
54 const UnitCell &ucell,
55 const std::vector<TA> &list_A0,
56 const std::vector<TAC> &list_A1,
57 const std::map<std::string,bool> &flags); // "writable_dVws"
58 std::pair<std::map<TA, std::map<TAC, RI::Tensor<Tdata>>>,
59 std::map<TA, std::map<TAC, std::array<RI::Tensor<Tdata>, 3>>>>
61 const UnitCell &ucell,
62 const std::vector<TA> &list_A0,
63 const std::vector<TAC> &list_A1,
64 const std::map<std::string,bool> &flags); // "cal_dC", "writable_Cws", "writable_dCws", "writable_Vws", "writable_dVws"
65
66 size_t get_index_abfs_size(const size_t &iat){return this->index_abfs[iat].count_size; }
67
68private:
69 std::vector<std::vector<std::vector<Numerical_Orbital_Lm>>> lcaos;
70 std::vector<std::vector<std::vector<Numerical_Orbital_Lm>>> abfs;
71 std::vector<std::vector<std::vector<Numerical_Orbital_Lm>>> abfs_ccp;
74 std::vector<double> lcaos_rcut;
75 std::vector<double> abfs_ccp_rcut;
76
77public:
78 std::map<int,std::map<int,std::map<Abfs::Vector3_Order<double>,RI::Tensor<Tdata>>>> Vws;
79 std::map<int,std::map<int,std::map<Abfs::Vector3_Order<double>,RI::Tensor<Tdata>>>> Cws;
80 std::map<int,std::map<int,std::map<Abfs::Vector3_Order<double>,std::array<RI::Tensor<Tdata>,3>>>> dVws;
81 std::map<int,std::map<int,std::map<Abfs::Vector3_Order<double>,std::array<RI::Tensor<Tdata>,3>>>> dCws;
82private:
83 pthread_rwlock_t rwlock_Vw;
84 pthread_rwlock_t rwlock_Cw;
85 pthread_rwlock_t rwlock_dVw;
86 pthread_rwlock_t rwlock_dCw;
87
90
91 template<typename Tresult>
92 using T_func_DPcal_data = std::function<Tresult(
93 const int it0,
94 const int it1,
96 const std::map<std::string,bool> &flags)>;
97 using T_func_cal_Rcut = std::function<double(const int it0, const int it1)>;
98 template<typename Tresult>
99 std::map<TA,std::map<TAC,Tresult>>
101 const UnitCell &ucell,
102 const std::vector<TA>& list_A0,
103 const std::vector<TAC>& list_A1,
104 const std::map<std::string, bool>& flags,
105 const T_func_cal_Rcut& func_cal_Rcut,
106 const T_func_DPcal_data<Tresult>& func_DPcal_data);
107
108 inline double cal_V_Rcut(const int it0, const int it1);
109 inline double cal_C_Rcut(const int it0, const int it1);
110
111 inline RI::Tensor<Tdata>
112 DPcal_V(
113 const int it0,
114 const int it1,
116 const std::map<std::string,bool> &flags); // "writable_Vws"
117 inline std::array<RI::Tensor<Tdata>,3>
118 DPcal_dV(
119 const int it0,
120 const int it1,
122 const std::map<std::string,bool> &flags); // "writable_dVws"
123 std::pair<RI::Tensor<Tdata>, std::array<RI::Tensor<Tdata>,3>>
125 const int it0,
126 const int it1,
128 const std::map<std::string,bool> &flags); // "cal_dC", "writable_Cws", "writable_dCws", "writable_Vws", "writable_dVws"
129
130 template<typename To11, typename Tfunc>
131 To11 DPcal_o11(
132 const int it0,
133 const int it1,
135 const bool &flag_writable_o11ws,
136 pthread_rwlock_t &rwlock_o11,
137 std::map<int,std::map<int,std::map<Abfs::Vector3_Order<double>,To11>>> &o11ws,
138 const Tfunc &func_cal_o11);
139};
140
141#include "LRI_CV.hpp"
142
143#endif
Definition abfs-vector3_order.h:16
Definition ORB_read.h:19
Definition LRI_CV.h:25
RI::Global_Func::To_Real_t< Tdata > Tdata_real
Definition LRI_CV.h:30
std::array< RI::Tensor< Tdata >, 3 > DPcal_dV(const int it0, const int it1, const Abfs::Vector3_Order< double > &R, const std::map< std::string, bool > &flags)
Definition LRI_CV.hpp:306
std::pair< RI::Tensor< Tdata >, std::array< RI::Tensor< Tdata >, 3 > > DPcal_C_dC(const int it0, const int it1, const Abfs::Vector3_Order< double > &R, const std::map< std::string, bool > &flags)
Definition LRI_CV.hpp:336
RI::Tensor< Tdata > DPcal_V(const int it0, const int it1, const Abfs::Vector3_Order< double > &R, const std::map< std::string, bool > &flags)
Definition LRI_CV.hpp:291
std::vector< double > abfs_ccp_rcut
Definition LRI_CV.h:75
std::map< int, std::map< int, std::map< Abfs::Vector3_Order< double >, RI::Tensor< Tdata > > > > Vws
Definition LRI_CV.h:78
std::map< TA, std::map< TAC, RI::Tensor< Tdata > > > cal_Vs(const UnitCell &ucell, const std::vector< TA > &list_A0, const std::vector< TAC > &list_A1, const std::map< std::string, bool > &flags)
Definition LRI_CV.hpp:157
std::map< int, std::map< int, std::map< Abfs::Vector3_Order< double >, std::array< RI::Tensor< Tdata >, 3 > > > > dVws
Definition LRI_CV.h:80
std::function< double(const int it0, const int it1)> T_func_cal_Rcut
Definition LRI_CV.h:97
~LRI_CV()
Definition LRI_CV.hpp:30
std::vector< std::vector< std::vector< Numerical_Orbital_Lm > > > lcaos
Definition LRI_CV.h:69
std::map< TA, std::map< TAC, Tresult > > cal_datas(const UnitCell &ucell, const std::vector< TA > &list_A0, const std::vector< TAC > &list_A1, const std::map< std::string, bool > &flags, const T_func_cal_Rcut &func_cal_Rcut, const T_func_DPcal_data< Tresult > &func_DPcal_data)
std::function< Tresult(const int it0, const int it1, const Abfs::Vector3_Order< double > &R, const std::map< std::string, bool > &flags)> T_func_DPcal_data
Definition LRI_CV.h:96
LRI_CV()
Definition LRI_CV.hpp:21
int TA
Definition LRI_CV.h:27
pthread_rwlock_t rwlock_dVw
Definition LRI_CV.h:85
std::vector< double > lcaos_rcut
Definition LRI_CV.h:74
std::pair< std::map< TA, std::map< TAC, RI::Tensor< Tdata > > >, std::map< TA, std::map< TAC, std::array< RI::Tensor< Tdata >, 3 > > > > cal_Cs_dCs(const UnitCell &ucell, const std::vector< TA > &list_A0, const std::vector< TAC > &list_A1, const std::map< std::string, bool > &flags)
Definition LRI_CV.hpp:200
ModuleBase::Element_Basis_Index::IndexLNM index_lcaos
Definition LRI_CV.h:72
Matrix_Orbs21 m_abfslcaos_lcaos
Definition LRI_CV.h:89
Matrix_Orbs11 m_abfs_abfs
Definition LRI_CV.h:88
std::map< TA, std::map< TAC, std::array< RI::Tensor< Tdata >, 3 > > > cal_dVs(const UnitCell &ucell, const std::vector< TA > &list_A0, const std::vector< TAC > &list_A1, const std::map< std::string, bool > &flags)
Definition LRI_CV.hpp:178
pthread_rwlock_t rwlock_Cw
Definition LRI_CV.h:84
std::vector< std::vector< std::vector< Numerical_Orbital_Lm > > > abfs
Definition LRI_CV.h:70
std::pair< TA, TC > TAC
Definition LRI_CV.h:29
ModuleBase::Element_Basis_Index::IndexLNM index_abfs
Definition LRI_CV.h:73
size_t get_index_abfs_size(const size_t &iat)
Definition LRI_CV.h:66
std::array< int, 3 > TC
Definition LRI_CV.h:28
pthread_rwlock_t rwlock_dCw
Definition LRI_CV.h:86
double cal_C_Rcut(const int it0, const int it1)
Definition LRI_CV.hpp:103
std::map< int, std::map< int, std::map< Abfs::Vector3_Order< double >, std::array< RI::Tensor< Tdata >, 3 > > > > dCws
Definition LRI_CV.h:81
std::map< int, std::map< int, std::map< Abfs::Vector3_Order< double >, RI::Tensor< Tdata > > > > Cws
Definition LRI_CV.h:79
void set_orbitals(const UnitCell &ucell, const LCAO_Orbitals &orb, const std::vector< std::vector< std::vector< Numerical_Orbital_Lm > > > &lcaos_in, const std::vector< std::vector< std::vector< Numerical_Orbital_Lm > > > &abfs_in, const std::vector< std::vector< std::vector< Numerical_Orbital_Lm > > > &abfs_ccp_in, const double &kmesh_times, ORB_gaunt_table &MGT, const bool &init_MGT, const bool &init_C)
Definition LRI_CV.hpp:40
double cal_V_Rcut(const int it0, const int it1)
Definition LRI_CV.hpp:98
pthread_rwlock_t rwlock_Vw
Definition LRI_CV.h:83
To11 DPcal_o11(const int it0, const int it1, const Abfs::Vector3_Order< double > &R, const bool &flag_writable_o11ws, pthread_rwlock_t &rwlock_o11, std::map< int, std::map< int, std::map< Abfs::Vector3_Order< double >, To11 > > > &o11ws, const Tfunc &func_cal_o11)
Definition LRI_CV.hpp:237
std::vector< std::vector< std::vector< Numerical_Orbital_Lm > > > abfs_ccp
Definition LRI_CV.h:71
Definition Matrix_Orbs11.h:22
Definition Matrix_Orbs21.h:21
std::vector< Index_T > IndexLNM
Definition element_basis_index.h:41
Definition ORB_gaunt_table.h:9
Definition unitcell.h:16