-
Notifications
You must be signed in to change notification settings - Fork 1
/
ci_blocks_el.f90
147 lines (107 loc) · 3.9 KB
/
ci_blocks_el.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
module ci_blocks_el
use globals
! use misc
use mat_element_el
implicit none
contains
! the following functions needed to construct the CI matrix
!
! main blocks:
! function main_2h1p (csf_type, a, nrow, ncol) result(mat)
!
! coupling blocks:
! function coupling_2h1p (a, type_1, type_2, nrow, ncol) result(mat)
!
! =================
! MAIN BLOCK
! =================
! -----------------
! 2h1p
! -----------------
function main_2h1p (csf_type, a, nrow, ncol) result(mat)
integer :: a
integer :: ov1, ov2, ov3, ov4, row, col
integer :: i, j
integer :: nrow, ncol
integer :: csf_type
double precision, dimension(nrow, ncol) :: mat
if (csf_type == 1) then !singlet I
do i = 1, n_config
ov1 = config_2h(i, 1)
ov2 = config_2h(i, 2)
do j = 1, n_config
ov3 = config_2h(j, 1)
ov4 = config_2h(j, 2)
mat(i, j) = mat_el_2h1p_sI (ov1, ov2, a, ov3, ov4, a)
end do
end do
else if (csf_type == 2) then ! singlet II
do i = h_f_min, h_f_max-1
row = i - h_f_min + 1
mat(row,row) = mat_el_2h1p_sII (i, a, i, a)
do j = i + 1, h_f_max
col = j - h_f_min + 1
mat(row,col) = mat_el_2h1p_sII (i, a, j, a)
mat(col,row) = mat(row,col)
end do
end do
mat(n_ov,n_ov) = mat_el_2h1p_sII (h_f_max, a, h_f_max, a)
else if (csf_type == 3) then ! triplet
do i = 1, n_config-1
ov1 = config_2h(i, 1)
ov2 = config_2h(i, 2)
mat(i, i) = mat_el_2h1p_t (ov1, ov2, a, ov1, ov2, a)
do j = i+1, n_config
ov3 = config_2h(j, 1)
ov4 = config_2h(j, 2)
mat(i, j) = mat_el_2h1p_t (ov1, ov2, a, ov3, ov4, a)
mat(j, i) = mat(i, j)
end do
end do
ov1 = config_2h(n_config,1)
ov2 = config_2h(n_config,2)
mat(n_config, n_config) = mat_el_2h1p_t (ov1, ov2, a, ov1, ov2, a)
end if
end function main_2h1p
! ====================
! COUPLING BLOCK
! ====================
function coupling_2h1p (a, type_1, type_2, nrow, ncol) result(mat)
integer :: a
integer :: ov1, ov2, ov3, ov4
integer :: i, j, row
integer :: nrow, ncol
integer :: type_1, type_2
double precision, dimension(nrow, ncol) :: mat
if ((type_1 == 1) .and. (type_2 == 2)) then
do i = 1, n_config
ov1 = config_2h(i, 1)
ov2 = config_2h(i, 2)
do j = 1, n_ov
ov3 = j + h_f_min - 1
mat(i, j) = mat_el_2h1p_sI_sII (ov1, ov2, a, ov3, a)
end do
end do
else if ((type_1 == 1) .and. (type_2 == 3)) then
do i = 1, n_config
ov1 = config_2h(i, 1)
ov2 = config_2h(i, 2)
mat(i, i) = mat_el_2h1p_sI_t (ov1, ov2, a, ov1, ov2, a)
do j = 1, n_config
ov3 = config_2h(j, 1)
ov4 = config_2h(j, 2)
mat(i, j) = mat_el_2h1p_sI_t (ov1, ov2, a, ov3, ov4, a)
end do
end do
else if ((type_1 == 2) .and. (type_2 == 3)) then
do i = 1, n_ov
ov1 = i + h_f_min - 1
do j = 1, n_config
ov2 = config_2h(j, 1)
ov3 = config_2h(j, 2)
mat(i, j) = mat_el_2h1p_sII_t (ov1, a, ov2, ov3, a)
end do
end do
end if
end function coupling_2h1p
end module ci_blocks_el