forked from sanshar/Dice
-
Notifications
You must be signed in to change notification settings - Fork 0
/
integral.h
188 lines (155 loc) · 6.21 KB
/
integral.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
/*
Developed by Sandeep Sharma with contributions from James E. T. Smith and Adam A. Holmes, 2017
Copyright (c) 2017, Sandeep Sharma
This file is part of DICE.
This program is free software: you can redistribute it and/or modify it under the terms
of the GNU General Public License as published by the Free Software Foundation,
either version 3 of the License, or (at your option) any later version.
This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with this program.
If not, see <http://www.gnu.org/licenses/>.
*/
#ifndef INTEGRAL_HEADER_H
#define INTEGRAL_HEADER_H
#include <vector>
#include <string>
#include <iostream>
#include <Eigen/Dense>
#include <map>
#include <utility>
#include "iowrapper.h"
#include <boost/interprocess/managed_shared_memory.hpp>
#include <boost/serialization/complex.hpp>
#include <boost/serialization/vector.hpp>
using namespace std;
using namespace Eigen;
class schedule;
bool myfn(double i, double j);
class compAbs {
public:
bool operator()(const float& a, const float& b) const { return fabs(a) < fabs(b); }
bool operator()(const complex<double>& a, const complex<double>& b) const { return std::abs(a) < std::abs(b); }
};
class oneInt {
private:
friend class boost::serialization::access;
template<class Archive>
void serialize(Archive & ar, const unsigned int version){ ar & store & norbs; }
public:
std::vector<CItype> store;
int norbs;
//I explicitly store all elements of the matrix
//so for normal operator if i and j dont have the same spin
//then it will just return zero. If we have SOC and
// i and j have different spin then it can be a complex number.
inline CItype& operator()(int i, int j) { return store.at(i*norbs+j); }
};
class twoInt {
private:
friend class boost::serialization::access;
template<class Archive>
void serialize(Archive & ar, const unsigned int version){
ar & maxEntry \
& Direct \
& Exchange \
& zero \
& norbs \
& ksym;
}
public:
double* store;
double maxEntry;
MatrixXd Direct, Exchange;
double zero ;
size_t norbs;
bool ksym;
twoInt() :zero(0.0),maxEntry(100.) {}
inline double& operator()(int i, int j, int k, int l) {
zero = 0.0;
if (!((i%2 == j%2) && (k%2 == l%2))) return zero;
int I=i/2;int J=j/2;int K=k/2;int L=l/2;
if(!ksym) {
int IJ = max(I,J)*(max(I,J)+1)/2 + min(I,J);
int KL = max(K,L)*(max(K,L)+1)/2 + min(K,L);
int A = max(IJ,KL), B = min(IJ,KL);
return store[A*(A+1)/2+B];
} else {
int IJ = I*norbs+J, KL = K*norbs+L;
int A = max(IJ,KL), B = min(IJ,KL);
return store[A*(A+1)/2+B];
}
}
};
class twoIntHeatBath {
public:
//i,j,a,b are spatial orbitals
//first pair is i,j (i>j)
//the map contains a list of integral which are equal to (ia|jb) where a,b are the second pair
//if the integrals are for the same spin you have to separately store (ia|jb) and (ib|ja)
//for opposite spin you just have to store for a>b because the integral is (ia|jb) - (ib|ja)
//now this class is made by just considering integrals that are smaller than threshhold
std::map<std::pair<short,short>, std::multimap<float, std::pair<short,short>, compAbs > > sameSpin;
std::map<std::pair<short,short>, std::multimap<float, std::pair<short,short>, compAbs > > oppositeSpin;
MatrixXd Singles;
double epsilon;
double zero ;
twoIntHeatBath(double epsilon_) :zero(0.0),epsilon(fabs(epsilon_)) {}
//the orbs contain all orbitals used to make the ij pair above
//typically these can be all orbitals of the problem or just the active space ones
//ab will typically contain all orbitals(norbs)
void constructClass(std::vector<int>& orbs, twoInt& I2, oneInt& I1, int norbs) {
for (int i=0; i<orbs.size(); i++)
for (int j=0;j<=i;j++) {
std::pair<short,short> IJ=make_pair(i,j);
//sameSpin[IJ]=std::map<double, std::pair<int,int> >();
//oppositeSpin[IJ]=std::map<double, std::pair<int,int> >();
for (int a=0; a<norbs; a++)
for (int b=0; b<norbs; b++) {
//opposite spin
if (fabs(I2(2*i, 2*a, 2*j, 2*b)) > epsilon)
oppositeSpin[IJ].insert(pair<float, std::pair<short,short> >(I2(2*i, 2*a, 2*j, 2*b), make_pair(a,b)));
//samespin
if (a>=b && fabs(I2(2*i,2*a,2*j,2*b) - I2(2*i,2*b,2*j,2*a)) > epsilon) {
sameSpin[IJ].insert(pair<float, std::pair<short,short> >( I2(2*i,2*a,2*j,2*b) - I2(2*i,2*b,2*j,2*a), make_pair(a,b)));
//sameSpin[IJ][fabs(I2(2*i,2*a,2*j,2*b) - I2(2*i,2*b,2*j,2*a))] = make_pair<int,int>(a,b);
}
}
} // ij
Singles = MatrixXd::Zero(2*norbs, 2*norbs);
for (int i=0; i<2*norbs; i++)
for (int a=0; a<2*norbs; a++) {
Singles(i,a) = std::abs(I1(i,a));
for (int j=0; j<2*norbs; j++) {
if (fabs(Singles(i,a)) < fabs(I2(i,a,j,j) - I2(i, j, j, a)))
Singles(i,a) = std::abs(I2(i,a,j,j) - I2(i, j, j, a));
}
}
} // end constructClass
};
class twoIntHeatBathSHM {
public:
float* sameSpinIntegrals;
float* oppositeSpinIntegrals;
size_t* startingIndicesSameSpin;
size_t* startingIndicesOppositeSpin;
short* sameSpinPairs;
short* oppositeSpinPairs;
double* singleExcitation;
MatrixXd Singles;
double epsilon;
twoIntHeatBathSHM(double epsilon_) : epsilon(fabs(epsilon_)) {}
void constructClass(int norbs, twoIntHeatBath& I2) ;
};
#ifdef Complex
void readSOCIntegrals(oneInt& I1soc, int norbs, string fileprefix, schedule& schd);
void readGTensorIntegrals(vector<oneInt>& I1soc, int norbs, string fileprefix);
#endif
int readNorbs(string fcidump);
void readIntegrals(
string fcidump,
twoInt& I2, oneInt& I1,
int& nelec, int& norbs, double& coreE,
std::vector<int>& irrep);
#endif