-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathbuild_ci_mat.f90
204 lines (156 loc) · 7.02 KB
/
build_ci_mat.f90
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
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
module build_ci_mat
use globals
use ci_blocks_el
use ci_blocks_pol
use mat_element_2h
implicit none
contains
subroutine build_ci_mat_init ()
integer :: i, j
integer :: h1, h2, p1, p2
CI_init_matrix(:,:) = 0.0d0
if (prop_type == "pol") then
do i = 1, ci_mat_in_size
h1 = config_1h1p(i,1)
p1 = config_1h1p(i,2)
do j = 1, ci_mat_in_size
h2 = config_1h1p(j,1)
p2 = config_1h1p(j,2)
CI_init_matrix(i, j) = mat_el_1h1p (h1, p1, h2, p2)
end do
end do
end if
end subroutine build_ci_mat_init
subroutine build_cimatrix(a)
integer, intent(in) :: a
integer :: offset_2h1p_t
integer :: offset_2h2p_a, offset_2h2p_b, offset_2h2p_1
CI_matrix(:,:) = 0.0d0
if (prop_type == "el") then
!--------------------------------------------------------
! Construction of the CI_matrix for a given virtual
! orbital k: only the upper triangular part of the
! matrix is created
!
! 1. main block singlet I
! 2. main block singlet II
! 3. main block triplet III
!--------------------------------------------------------
offset_2h1p_t = n_config + n_ov
if (gam == "singlet") then
CI_matrix(1:n_config, 1:n_config) = main_2h1p(1, a, n_config, n_config)
CI_matrix((n_config+1):ci_mat_size, (n_config+1):ci_mat_size) = &
& main_2h1p(2, a, n_ov, n_ov)
else if (gam == "triplet") then
CI_matrix = main_2h1p(3, a, n_config, n_config)
else
CI_matrix(1:n_config, 1:n_config) = main_2h1p(1, a, n_config, n_config)
CI_matrix((n_config+1):offset_2h1p_t, (n_config+1):offset_2h1p_t) =&
& main_2h1p(2, a, n_ov, n_ov)
CI_matrix((offset_2h1p_t+1):ci_mat_size, (offset_2h1p_t+1):ci_mat_size) =&
& main_2h1p(3, a, n_config, n_config)
end if
!--------------------------------------------------------
! 4. coupling block S I - S II
! 5. coupling block S I - T III
! 6. coupling block S II - T III
!--------------------------------------------------------
if (gam == "singlet") then
CI_matrix(1:n_config, (n_config+1):offset_2h1p_t) = &
& coupling_2h1p(a, 1, 2, n_config, n_ov)
else if (gam == "total") then
CI_matrix(1:n_config, (n_config+1):offset_2h1p_t) = &
& coupling_2h1p(a, 1, 2, n_config, n_ov)
CI_matrix(1:n_config, (offset_2h1p_t+1):ci_mat_size) = &
& coupling_2h1p(a, 1, 3, n_config, n_config)
CI_matrix((n_config+1):offset_2h1p_t, (offset_2h1p_t+1):ci_mat_size) = &
& coupling_2h1p(a, 2, 3, n_ov, n_config)
end if
else if (prop_type == "pol") then
!--------------------------------------------------------
! Construction of the CI_matrix for a given virtual
! orbital k: only the upper triangular part of the
! matrix is created
!
! 1. main block 1h1p
! 2. main block 2h2p | ii ab>
! 3. main block 2h2p | ij ab, A>
! 4. main block 2h2p | ij ab, B>
!--------------------------------------------------------
offset_2h2p_1 = n_ov + 1
offset_2h2p_a = 2*n_ov + 1
offset_2h2p_b = 2*n_ov + n_config + 1
CI_matrix(1:n_ov, 1:n_ov) = main_1h1p(a)
CI_matrix(offset_2h2p_1 : offset_2h2p_a-1, offset_2h2p_1 : offset_2h2p_a-1) =&
& main_2h2p (p_in, a, 1, n_ov, n_ov)
CI_matrix(offset_2h2p_a : offset_2h2p_b-1, offset_2h2p_a : offset_2h2p_b-1) =&
& main_2h2p (p_in, a, 2, n_config, n_config)
CI_matrix(offset_2h2p_b : ci_mat_size, offset_2h2p_b : ci_mat_size) =&
& main_2h2p (p_in, a, 3, n_config, n_config)
!--------------------------------------------------------
! 5. coupling block 1h1p 2h2p |ii ab>
! 6. coupling block 1h1p 2h2p |ij ab, A>
! 7. coupling block 1h1p 2h2p |ij ab, B>
! 8. coupling block 2h2p |ii ab> / |ij ab, A>
! 9. coupling block 2h2p |ii ab> / |ij ab, B>
! 10. coupling block 2h2p |ij ab, A> / |ij ab, B>
!--------------------------------------------------------
CI_matrix(1:n_ov, offset_2h2p_1 : offset_2h2p_a-1) = &
& coupling_1h1p_2h2p(p_in, a, 1, n_ov, n_ov)
CI_matrix(1:n_ov, offset_2h2p_a : offset_2h2p_b-1) = &
& coupling_1h1p_2h2p(p_in, a, 2, n_ov, n_config)
CI_matrix(1:n_ov, offset_2h2p_b : ci_mat_size) = &
& coupling_1h1p_2h2p(p_in, a, 3, n_ov, n_config)
CI_matrix(offset_2h2p_1 : offset_2h2p_a-1, offset_2h2p_a : offset_2h2p_b-1) = &
& coupling_2h2p_2h2p (p_in, a, 1, 2, n_config)
CI_matrix(offset_2h2p_1 : offset_2h2p_a-1, offset_2h2p_b : ci_mat_size) = &
& coupling_2h2p_2h2p (p_in, a, 1, 3, n_config)
CI_matrix(offset_2h2p_a : offset_2h2p_b-1, offset_2h2p_b : ci_mat_size) = &
& coupling_2h2p_2h2p (p_in, a, 2, 3, n_config)
end if
end subroutine build_cimatrix
subroutine build_ci2hmatrix (type, mat_2h, m_size)
integer, intent(in) :: type
integer, intent(in) :: m_size
integer :: i, j, h1, h2, h3, h4
double precision, dimension (m_size, m_size) :: mat_2h
mat_2h(:,:) = 0.d0
if (type == 0) then ! singlet
do i = 1, n_config
h1 = config_2h(i, 1)
h2 = config_2h(i, 2)
do j = 1, n_config
h3 = config_2h(j, 1)
h4 = config_2h(j, 2)
mat_2h(i, j) = mat_el_2h_s (1, h1, h2, h3, h4)
end do
end do
do i = 1, n_ov
h1 = h_f_min - 1 + i
do j = 1, n_ov
h2 = h_f_min - 1 + j
mat_2h(n_config + i, n_config + j) = mat_el_2h_s (2, h1, h1, h2, h2)
end do
end do
do i = 1, n_config
h1 = config_2h(i, 1)
h2 = config_2h(i, 2)
do j = 1, n_ov
h3 = h_f_min - 1 + j
mat_2h(i, n_config + j) = mat_el_2h_s (3, h3, h3, h1, h2)
mat_2h(n_config + j, i) = mat_el_2h_s (3, h3, h3, h1, h2)
end do
end do
else if (type == 1) then ! triplet
do i = 1, n_config
h1 = config_2h(i, 1)
h2 = config_2h(i, 2)
do j = 1, n_config
h3 = config_2h(j, 1)
h4 = config_2h(j, 2)
mat_2h(i, j) = mat_el_2h_t (h1, h2, h3, h4)
end do
end do
end if
end subroutine build_ci2hmatrix
end module build_ci_mat