-
Notifications
You must be signed in to change notification settings - Fork 0
/
vgenpPatk.f
194 lines (192 loc) · 5.21 KB
/
vgenpPatk.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
c***********************************************************************
SUBROUTINE VGENP(ISTATE,RDIST,VDIST,dVdR,d2VdR2,IDAT)
c***********************************************************************
C Fit of the Ar-Ar potential
C K. Patkowski and K. Szalewicz, to be submitted (2010).
C
C Input:
C r - distance in Angstrom
C
C Output:
C value - interaction potential in cm-1
C deriv - analytic first derivative of potential in cm-1/bohr
c***********************************************************************
INCLUDE 'arrsizes.h'
INTEGER I,IDAT,ISTATE
REAL*8 RDIST(8),VDIST(8),BETADIST,dVdR(8),d2VdR2(8),ar2pot,ar2der
c
BETADIST= 0.d0
DO I= 1,8
VDIST(I)= ar2pot(RDIST(I))
d2VdR2(I)= 0.0
dVdR(I)= ar2der(RDIST(I))
ENDDO
cc write(6,*) RDIST,' ',VDIST,' ',dRdV
return
end
c***********************************************************************
c***********************************************************************
c*********************************************************************
real*8 function ar2pot(r)
implicit real*8 (a-h,o-z)
dimension cn(6)
data cn /64.288984d0,1514.86211d0,50240.0d0,1898195.d0,
$ 86445426.d0,4619452502.d0/
data ifirst /1/
data coef /127641878.945519894d0/
data coef2/-26138949.621478189d0/
data coef3/-115672346.174201056d0/
data coef4/2064381.526204719d0/
data coef5/-58371.409016267d0/
data alpha /1.553386357296d0/
data ddd /2.393847610341d0/
data conv /219474.63d0/
data rjoin /1.3d0/
data e0/324.d0/
data a0/-84134425.1811d0/
data b0/33397004.2848d0/
data c0/-4443000.0d0/
data alpha0/2.3577d0/
data beta0/-1.2756d0/
save
c
nc=6
if(ifirst.eq.1) then
do nn=1,nc
cn(nn) = cn(nn)*conv
enddo
e0=e0*conv
ifirst = 0
endif
if (r.ge.rjoin) then
rpt = r/0.529177209d0
term = (coef+coef2*rpt+coef3/rpt+coef4*rpt*rpt
$ +coef5*rpt*rpt*rpt)*dexp(-alpha*rpt)
sumc=0.0d0
do nn=1,nc
sumc=sumc+d(2*nn+4,ddd,rpt)*cn(nn)/rpt**(2*nn+4)
enddo
asy = sumc
term = term - asy
ar2pot = term
else
rpt = r/0.529177209d0
ar2pot=dexp(-alpha0*rpt-beta0*rpt*rpt)
ar2pot=ar2pot*(e0/rpt+a0+b0*rpt+c0*rpt*rpt)
end if
return
end
cc
real*8 function ar2der(r)
implicit real*8 (a-h,o-z)
dimension cn(6)
data cn /64.288984d0,1514.86211d0,50240.0d0,1898195.d0,
$ 86445426.d0,4619452502.d0/
data ifirst /1/
data coef /127641878.945519894d0/
data coef2/-26138949.621478189d0/
data coef3/-115672346.174201056d0/
data coef4/2064381.526204719d0/
data coef5/-58371.409016267d0/
data alpha /1.553386357296d0/
data ddd /2.393847610341d0/
data conv /219474.63d0/
data rjoin /1.3d0/
data e0/324.d0/
data a0/-84134425.1811d0/
data b0/33397004.2848d0/
data c0/-4443000.0d0/
data alpha0/2.3577d0/
data beta0/-1.2756d0/
save
c
nc=6
if(ifirst.eq.1) then
do nn=1,nc
cn(nn) = cn(nn)*conv
enddo
e0=e0*conv
ifirst = 0
endif
if (r.ge.rjoin) then
rpt = r/0.529177209d0
expfac = exp(-alpha*rpt)
dercoe = -alpha
der= coef * dercoe * expfac
der = der + coef2 * (1.0 + rpt*dercoe) * expfac
der = der + coef3 * (-1.0/rpt**2 + dercoe/rpt)*expfac
der = der + coef4 * (2.d0*rpt + rpt*rpt*dercoe) * expfac
der = der + coef5 * (3.d0*rpt*rpt + rpt*rpt*rpt*dercoe) * expfac
sumcd=0.0d0
do nn=1,nc
nnn=2*nn+4
sumcd = sumcd + cn(nn)*(-nnn/rpt**(nnn+1)
& + ddd*exp(-ddd*rpt)*d1(nnn,ddd,rpt)/rpt**nnn
& - exp(-ddd*rpt)*da(nnn,ddd,rpt) )
enddo
der = der - sumcd
ar2der = der
else
rpt = r/0.529177209d0
e1=dexp(-alpha0*rpt-beta0*rpt*rpt)
de1=-alpha0-2.d0*beta0*rpt
e2=(e0/rpt+a0+b0*rpt+c0*rpt*rpt)
de2=-e0/(rpt*rpt)+b0+2.0*c0*rpt
ar2der=e1*de2+e2*de1*e1
end if
return
end
cc
double precision function d(n,beta,r)
c
c calculate the damping factor
c
implicit real*8 (a-h,o-z)
br=beta*r
sum=1.0d0
term=1.0d0
ncn=n
do i=1,ncn
term=term*br/i
sum=sum+term
enddo
d=1.0d0 - dexp(-br)*sum
c write(6,*) n,beta,r,d
return
end
cc
double precision function d1(n,beta,r)
c
c calculate part of the derivative of the damping factor
c
implicit real*8 (a-h,o-z)
br=beta*r
sum=1.0d0
term=1.0d0
ncn=n
do i=1,ncn
term=term*br/i
sum=sum+term
enddo
d1= sum
c write(6,*) n,beta,r,d
return
end
cc
double precision function da(n,beta,r)
c
c calculate part of the derivative of the damping factor
c
implicit real*8 (a-h,o-z)
br=beta*r
sum=-n/r**(n+1)
term=1.0d0
ncn=n
do i=1,ncn-1
term=term*br/i
sum=sum+term*(i-n)/r**(n+1)
enddo
da = sum
c write(6,*) n,beta,r,d
return
end