-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathrescale.irp.f
320 lines (281 loc) · 8.46 KB
/
rescale.irp.f
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
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
BEGIN_PROVIDER [ double precision, kappa ]
implicit none
BEGIN_DOC
! Constant in rescaling
END_DOC
kappa = 0.6d0
END_PROVIDER
BEGIN_PROVIDER [ double precision, kappa_inv ]
implicit none
BEGIN_DOC
! inverse of kappa
END_DOC
kappa_inv = 1.0d0 / kappa
END_PROVIDER
BEGIN_PROVIDER [ double precision, rescale_ee, (nelec_8, nelec) ]
implicit none
BEGIN_DOC
! R = (1 - exp(-kappa r))/kappa for electron-electron for $J_{ee}$
END_DOC
integer :: i, j
double precision :: x
do j = 1, nelec
do i = 1, j-1
x = (1.0d0 - dexp(-kappa * elec_dist(i, j))) * kappa_inv
rescale_ee(i, j) = x
rescale_ee(j, i) = x
enddo
rescale_ee(j, j) = 0.d0
enddo
END_PROVIDER
BEGIN_PROVIDER [ double precision, rescale_ee_deriv_e, (4, nelec, nelec) ]
implicit none
BEGIN_DOC
! R = (1 - exp(-kappa r))/kappa derived wrt x
! Dimensions 1-3 : dx, dy, dz
! Dimension 4 : d2x + d2y + d2z
END_DOC
integer :: i, j, ii
double precision :: f
do j = 1, nelec
do i = 1, nelec
f = 1.d0 - kappa*rescale_ee(i,j) ! == dexp(-kappa * elec_dist(i, j))
do ii = 1, 4
rescale_ee_deriv_e(ii, i, j) = elec_dist_deriv_e(ii, i, j)
end do
rescale_ee_deriv_e(4, i, j) = rescale_ee_deriv_e(4, i, j) + &
(-kappa * rescale_ee_deriv_e(1, i, j) * rescale_ee_deriv_e(1, i, j)) + &
(-kappa * rescale_ee_deriv_e(2, i, j) * rescale_ee_deriv_e(2, i, j)) + &
(-kappa * rescale_ee_deriv_e(3, i, j) * rescale_ee_deriv_e(3, i, j))
do ii = 1, 4
rescale_ee_deriv_e(ii, i, j) = rescale_ee_deriv_e(ii, i, j) &
* f
enddo
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [ double precision, rescale_en, (nelec_8, nnuc) ]
implicit none
BEGIN_DOC
! R = (1 - exp(-kappa r))/kappa for electron-nucleus for $J_{en}$
END_DOC
integer :: i, a
do a = 1, nnuc
do i = 1, nelec
rescale_en(i, a) = (1.0d0 - dexp(-kappa * elnuc_dist(i, a))) * kappa_inv
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [ double precision, rescale_en_deriv_e, (4, nelec, nnuc) ]
implicit none
BEGIN_DOC
! R = (1 - exp(-kappa r))/kappa derived wrt x
! Dimensions 1-3 : dx, dy, dz
! Dimension 4 : d2x + d2y + d2z
END_DOC
integer :: i, ii, a
double precision :: f
do a = 1, nnuc
do i = 1, nelec
f = 1.d0 - kappa*rescale_en(i,a) ! == dexp(-kappa * elnuc_dist(i, a))
do ii = 1, 4
rescale_en_deriv_e(ii, i, a) = elnuc_dist_deriv_e(ii, i, a)
end do
rescale_en_deriv_e(4, i, a) = rescale_en_deriv_e(4, i, a) + &
(-kappa * rescale_en_deriv_e(1, i, a) * rescale_en_deriv_e(1, i, a)) + &
(-kappa * rescale_en_deriv_e(2, i, a) * rescale_en_deriv_e(2, i, a)) + &
(-kappa * rescale_en_deriv_e(3, i, a) * rescale_en_deriv_e(3, i, a))
do ii = 1, 4
rescale_en_deriv_e(ii, i, a) = rescale_en_deriv_e(ii, i, a) &
* f
enddo
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [double precision, rescale_een_e, (nelec_8, nelec, 0:ncord)]
implicit none
BEGIN_DOC
! R = exp(-kappa r) for electron-electron for $J_{een}$
END_DOC
integer :: i, j, k, l
double precision :: x
double precision, parameter :: f = dexp(1.d0)
rescale_een_e(:, :, 0) = 1.d0
do l = 1, ncord
k=0
do j = 1, nelec
do i = 1, j-1
k = k+1
x = rescale_een_e_ij(k,l)
rescale_een_e(i, j, l) = x
rescale_een_e(j, i, l) = x
enddo
enddo
enddo
do l = 0, ncord
do j = 1, nelec
rescale_een_e(j, j, l) = 0.d0
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [double precision, rescale_een_e_ij, (nelec*(nelec-1)/2, 0:ncord)]
implicit none
BEGIN_DOC
! R = exp(-kappa r) for electron-electron for $J_{een}$
END_DOC
integer :: i, j, l,k
double precision :: x
double precision, parameter :: f = dexp(1.d0)
rescale_een_e_ij(:, 0) = 1.d0
k=0
do j = 1, nelec
do i = 1, j-1
k = k+1
rescale_een_e_ij(k, 1) = dexp(-kappa * elec_dist(i, j))
enddo
enddo
do l = 2, ncord
do k=1,(nelec*nelec-nelec)/2
rescale_een_e_ij(k, l) = rescale_een_e_ij(k, l-1) * rescale_een_e_ij(k, 1)
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [double precision, rescale_een_n, (nelec_8, nnuc, 0:ncord)]
implicit none
BEGIN_DOC
! R = exp(-kappa r) for electron-electron for $J_{een}$
END_DOC
integer :: i, a, l
double precision :: x
double precision, parameter :: f = dexp(1.d0)
rescale_een_n(:,:,0) = 1.d0
do a = 1, nnuc
do i = 1, nelec
rescale_een_n(i, a, 1) = dexp(-kappa * elnuc_dist(i, a))
enddo
enddo
do l = 2, ncord
do a = 1, nnuc
do i = 1, nelec
rescale_een_n(i, a, l) = rescale_een_n(i, a, l-1) * rescale_een_n(i, a, 1)
enddo
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [double precision, rescale_een_n_deriv_e, (nelec_8, 4, nnuc, 0:ncord)]
implicit none
BEGIN_DOC
! Derivative of the scaled distance J_{een} wrt R_{ia}
END_DOC
integer :: i, ii, j, l, a
double precision :: kappa_l
do l = 0, ncord
kappa_l = - dble(l) * kappa
do a = 1, nnuc
do i = 1, nelec
do ii = 1, 4
! r'(x) \lor r''(x)
rescale_een_n_deriv_e(i, ii, a, l) = &
kappa_l * elnuc_dist_deriv_e(ii, i, a)
!print *, "pp", ii, i, a, elnuc_dist_deriv_e(ii, i, a)
enddo
! \left(r''(x)+r'(x)^2\right)
rescale_een_n_deriv_e(i, 4, a, l) = rescale_een_n_deriv_e(i, 4, a, l) + &
rescale_een_n_deriv_e(i, 1, a, l) * rescale_een_n_deriv_e(i, 1, a, l) + &
rescale_een_n_deriv_e(i, 2, a, l) * rescale_een_n_deriv_e(i, 2, a, l) + &
rescale_een_n_deriv_e(i, 3, a, l) * rescale_een_n_deriv_e(i, 3, a, l)
! \times e^{r(x)}
do ii = 1, 4
rescale_een_n_deriv_e(i, ii, a, l) = &
rescale_een_n_deriv_e(i, ii, a, l) * rescale_een_n(i, a, l)
enddo
enddo
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [double precision, elnuc_dist_deriv_e, (4, nelec, nnuc)]
BEGIN_DOC
! Derivative of R_{ia} wrt x
! Dimensions 1-3 : dx, dy, dz
! Dimension 4 : d2x + d2y + d2z
END_DOC
implicit none
integer :: i, ii, a
double precision :: ria_inv
do a = 1, nnuc
do i = 1, nelec
ria_inv = 1.0d0 / elnuc_dist(i, a)
do ii = 1, 3
! \frac{x-x0}{\sqrt{c+(x-x0)^2}}
elnuc_dist_deriv_e(ii, i, a) = (elec_coord(i, ii) - nuc_coord(a, ii)) * ria_inv
end do
elnuc_dist_deriv_e(4, i, a) = 2.0d0 * ria_inv
end do
end do
END_PROVIDER
BEGIN_PROVIDER [double precision, rescale_een_e_deriv_e, (nelec_8, 4, nelec, 0:ncord)]
BEGIN_DOC
! Derivative of the scaled distance J_{een} wrt R_{ia}
END_DOC
implicit none
integer :: i, ii, j, l
double precision :: kappa_l
!TODO: Check if rescale_een_e_deriv_e(:,:,0) = 0.d0
do l = 0, ncord
kappa_l = - dble(l) * kappa
do j = 1, nelec
do i = 1, nelec
! r'(x) \lor r''(x)
do ii = 1, 4
rescale_een_e_deriv_e(i, ii, j, l) = &
kappa_l * elec_dist_deriv_e(ii, i, j)
enddo
! \left(r''(x)+r'(x)^2\right)
rescale_een_e_deriv_e(i, 4, j, l) = rescale_een_e_deriv_e(i, 4, j, l) + &
rescale_een_e_deriv_e(i, 1, j, l) * rescale_een_e_deriv_e(i, 1, j, l) + &
rescale_een_e_deriv_e(i, 2, j, l) * rescale_een_e_deriv_e(i, 2, j, l) + &
rescale_een_e_deriv_e(i, 3, j, l) * rescale_een_e_deriv_e(i, 3, j, l)
! \times e^{r(x)}
do ii = 1, 4
rescale_een_e_deriv_e(i, ii, j, l) = &
rescale_een_e_deriv_e(i, ii, j, l) * rescale_een_e(i, j, l)
enddo
enddo
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [double precision, rescale_een_e_deriv_e_t, (nelec_8, 4, nelec, 0:ncord)]
implicit none
BEGIN_DOC
! Transposed rescale_een_e_deriv_e
END_DOC
integer :: i,j,k,l
do l=0,ncord
do j=1,nelec
do i=1,nelec
rescale_een_e_deriv_e_t(j,1:4,i,l) = rescale_een_e_deriv_e(i,1:4,j,l)
enddo
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [double precision, elec_dist_deriv_e, (4, nelec, nelec)]
BEGIN_DOC
! Derivative of R_{ij} wrt x
! Dimensions 1-3 : dx, dy, dz
! Dimension 4 : d2x + d2y + d2z = 2/rij
END_DOC
implicit none
integer :: i, ii, j
double precision :: rij_inv
do j = 1, nelec
do i = 1, nelec
rij_inv = 1.0d0 / elec_dist(i, j)
do ii = 1, 3
! \frac{x-x0}{\sqrt{c+(x-x0)^2}}
elec_dist_deriv_e(ii, i, j) = (elec_coord(i, ii) - elec_coord(j, ii)) * rij_inv
end do
elec_dist_deriv_e(4, i, j) = 2.0d0 * rij_inv
end do
elec_dist_deriv_e(:, j, j) = 0.0d0
end do
END_PROVIDER