root / src / egrad_chamfre.f90 @ 7
Historique | Voir | Annoter | Télécharger (25,84 ko)
1 |
!*************############### Attention ###############************* |
---|---|
2 |
!(a) Two parameters aa and br(3,3) must be provided in main program |
3 |
!(b) At first, z=0 for the outmost surface |
4 |
!(c) The units of all parameters are not atomic units, but angstrom, eV etc. |
5 |
!(d) Before using this routine, you must confirm that the energy is CONSISTENT with its derivative. |
6 |
! This check of consistency is VERY IMPORTANT for energy conservation in MD. |
7 |
!!! |
8 |
!! xa, ya,za present the positions of atoms. v is the energy of the system. dvdr present the force along the x,y,z directions. |
9 |
!! |
10 |
subroutine egrad_chamfre(Nat,V,Geom,Grad) |
11 |
! WAS xa,ya,za,v,dvdr) |
12 |
|
13 |
|
14 |
use Path_module, only : Lat_a,Lat_b,Lat_c,order,orderinv,renum |
15 |
use Io_module |
16 |
|
17 |
|
18 |
IMPLICIT NONE |
19 |
|
20 |
INTEGER(KINT), INTENT(IN) :: Nat |
21 |
REAL(KREAL), INTENT(OUT) :: V,grad(Nat*3) |
22 |
REAL(KREAL), INTENT(IN) :: Geom(Nat,3) |
23 |
|
24 |
! Bohr --> Angstr |
25 |
real(KREAL), parameter :: BOHR = 0.52917726D+00 |
26 |
|
27 |
!! PFL 2010 Nov. 28 -> |
28 |
! for now I keep this old F77 scheme. This has to be replaced by ALLOCATABLE |
29 |
! things |
30 |
! <- PFL 2010 Nov 28. |
31 |
REAL(KREAL) :: xa(100),ya(100),za(100),dvdr(100,3) |
32 |
REAL(KREAL) :: rshort,rlong,zshort,zlong |
33 |
REAL(KREAL), SAVE :: rnc,rnn,rhh,rch,rnh,ann(100),anc(100),anh(100),ach(100),ahh(100),a(100) |
34 |
REAL(KREAL) :: dv(100,10,3),df(100,10,3) |
35 |
REAL(KREAL) :: px(9000),py(9000),pz(9000),tmp,temp,rtemp,yf |
36 |
INTEGER(KINT) :: ka(9000),kx(9000),ky(9000),itx,ity,nf(100),mf(100,500) |
37 |
REAL(KREAL) :: xi(3),xj(3),xk(3),rij,rik,rjk,rji,rki,er,eb |
38 |
REAL(KREAL) :: fcut,tf,bodij,bodji,bod,cost,gtheta,dgdt,xt,t1,t2,t3 |
39 |
REAL(KREAL) :: zi,fcut2 |
40 |
INTEGER(KINT) :: np,nq,nh,nc |
41 |
INTEGER(KINT) :: i, j, nn, m, n, mid, n1, im, jm, n2,n3,n4, iat |
42 |
LOGICAL, SAVE :: First=.TRUE. |
43 |
LOGICAL :: debug, FExist |
44 |
|
45 |
common/zdis/zshort,zlong |
46 |
common/rdistanc/rshort,rlong |
47 |
|
48 |
|
49 |
INTERFACE |
50 |
function valid(string) result (isValid) |
51 |
CHARACTER(*), intent(in) :: string |
52 |
logical :: isValid |
53 |
END function VALID |
54 |
|
55 |
END INTERFACE |
56 |
|
57 |
!********************************************** |
58 |
|
59 |
|
60 |
debug=valid('egrad_chamfre') |
61 |
|
62 |
if (debug) Call Header('Entering Egrad_chamfre') |
63 |
|
64 |
v=0.0;dvdr=0.0 |
65 |
!******prepare the parameters of REBO potential |
66 |
!rshort=3.5 |
67 |
!rlong=3.9 |
68 |
zshort=3.5 |
69 |
zlong=4.5 |
70 |
|
71 |
rnn=3.2 !! to adjust |
72 |
rnc=3.5 !! to adjust |
73 |
rnh=3.0 !! to adjust |
74 |
rch=3.2 !! to adjust |
75 |
rhh=3.2 !! to adjust |
76 |
|
77 |
!rph=rlong |
78 |
|
79 |
if (First) THEN |
80 |
First=.FALSE. |
81 |
INQUIRE(FILE="parameter.dat",EXIST=FExist) |
82 |
if (.NOT.Fexist) THEN |
83 |
WRITE(*,*) "!!!!!!!!!!!! ERROR !!!!!!!!!!!!!!!!" |
84 |
WRITE(*,*) "!! egrad_chamfre !!!!" |
85 |
WRITE(*,*) "The file parameter.dat does not exists." |
86 |
WRITE(*,*) "!! egrad_chamfre !!!!" |
87 |
WRITE(*,*) "!!!!!!!!!!!! ERROR !!!!!!!!!!!!!!!!" |
88 |
STOP |
89 |
END IF |
90 |
open(1,file="parameter.dat") |
91 |
do i=1,5 !Ni-Ni |
92 |
read(1,*)ann(i) |
93 |
enddo |
94 |
|
95 |
do i=1,5 !Ni-C |
96 |
read(1,*)anc(i) |
97 |
enddo |
98 |
|
99 |
do i=1,5 !Ni-H |
100 |
read(1,*)anh(i) |
101 |
enddo |
102 |
|
103 |
do i=1,5 !H-H |
104 |
read(1,*)ahh(i) |
105 |
enddo |
106 |
do i=1,5 !C-H |
107 |
read(1,*)ach(i) |
108 |
enddo |
109 |
!********LONG INTERACTION ********** |
110 |
!do i=1,2 |
111 |
!read(1,*)a(i) |
112 |
!enddo |
113 |
close(1) |
114 |
END IF |
115 |
!************************************ |
116 |
np=45 |
117 |
nh=4 |
118 |
nc=1 |
119 |
nq=np+nh+nc |
120 |
!************************************ |
121 |
|
122 |
DO I=1,Nat |
123 |
Iat=Order(I) |
124 |
xa(I)=Geom(Iat,1) |
125 |
ya(I)=Geom(Iat,2) |
126 |
za(I)=Geom(Iat,3) |
127 |
END DO |
128 |
|
129 |
px=0.0;py=0.0;pz=0.0 |
130 |
nn=0 |
131 |
do m=-1,1 |
132 |
do n=-1,1 |
133 |
do i=1,nq |
134 |
if(m==0.and.n==0.and.i==1)mid=nn |
135 |
nn=nn+1 |
136 |
px(nn)=xa(i)+m*lat_a(1)+n*lat_b(1) |
137 |
py(nn)=ya(i)+m*lat_a(2)+n*lat_b(2) |
138 |
pz(nn)=za(i) |
139 |
ka(nn)=i;kx(nn)=m;ky(nn)=n |
140 |
enddo |
141 |
enddo |
142 |
enddo |
143 |
!***********************xyz file required for MS input |
144 |
open(123,file="111.xyz",ACCESS='APPEND') |
145 |
write(123,*)nn |
146 |
write(123,*) |
147 |
do i=1,nn |
148 |
if(ka(i)<(np+1))then |
149 |
write(123,10)'Ni',px(i),py(i),pz(i) |
150 |
elseif (ka(i)<np+5) then |
151 |
write(123,10)'H',px(i),py(i),pz(i) |
152 |
else |
153 |
write(123,10)'C',px(i),py(i),pz(i) |
154 |
10 format(a6,3f15.8) |
155 |
endif |
156 |
enddo |
157 |
close(123) |
158 |
!****************************************** |
159 |
do i=1,nq |
160 |
tmp=px(i);px(i)=px(i+mid);px(i+mid)=tmp |
161 |
tmp=py(i);py(i)=py(i+mid);py(i+mid)=tmp |
162 |
tmp=pz(i);pz(i)=pz(i+mid);pz(i+mid)=tmp |
163 |
itx=kx(i);kx(i)=kx(i+mid);kx(i+mid)=itx |
164 |
ity=ky(i);ky(i)=ky(i+mid);ky(i+mid)=ity |
165 |
enddo |
166 |
mid=0 |
167 |
!****************************************** |
168 |
nf=0;mf=0 |
169 |
do i=1,nq |
170 |
nf(i)=0 |
171 |
do j=1,nn |
172 |
temp=sqrt((px(i)-px(j))**2+(py(i)-py(j))**2+(pz(i)-pz(j))**2) |
173 |
! |
174 |
!rtemp=rlong |
175 |
if (ka(i+mid)<np+1.and.ka(j+mid)<np+1)then ! Ni-Ni |
176 |
|
177 |
rtemp=rnn |
178 |
endif |
179 |
if (ka(i+mid)<np+1.and.ka(j+mid)>np.and.ka(j+mid)<np+5)then ! Ni-H |
180 |
|
181 |
rtemp=rnh |
182 |
endif |
183 |
if (ka(i+mid)<np+1.and.ka(j+mid)>np+4) then ! Ni-C |
184 |
|
185 |
rtemp=rnc |
186 |
endif |
187 |
if (ka(i+mid)>np.and.ka(i+mid)<np+5.and.ka(j+mid)>np+4) then ! H-C |
188 |
|
189 |
rtemp=rch |
190 |
endif |
191 |
if (ka(i+mid)>np.and.ka(i+mid)<np+5.and.ka(j+mid)>np.and.ka(j+mid)<np+5)then ! H-H |
192 |
|
193 |
rtemp=rhh |
194 |
endif |
195 |
if (ka(i+mid)>np.and.ka(i+mid)<np+5.and.ka(j+mid)<np+1) then ! H-Ni |
196 |
|
197 |
rtemp=rnh |
198 |
endif |
199 |
if (ka(i+mid)>np+4.and.ka(j+mid)<np+1) then ! C-Ni |
200 |
|
201 |
rtemp=rnc |
202 |
endif |
203 |
if (ka(i+mid)>np+4.and.ka(j+mid)>np.and.ka(j+mid)<np+5) then ! C-H |
204 |
|
205 |
rtemp=rch |
206 |
endif |
207 |
|
208 |
if(temp>0.1.and.temp<rtemp)then |
209 |
nf(i)=nf(i)+1 |
210 |
mf(i,nf(i))=j |
211 |
endif |
212 |
! |
213 |
enddo |
214 |
!write(*,*)nf(i) |
215 |
enddo |
216 |
!************************************** |
217 |
!stop |
218 |
|
219 |
!*************************************Starting to calculate the sum |
220 |
yf=0.0;dvdr=0.0 |
221 |
do i=1,nq |
222 |
er=0.0;eb=0.0 |
223 |
dv=0.0;df=0.0 |
224 |
do j=1,nf(i) |
225 |
|
226 |
xi(1)=px(i); xi(2)=py(i); xi(3)=pz(i) |
227 |
xj(1)=px(mf(i,j));xj(2)=py(mf(i,j));xj(3)=pz(mf(i,j)) |
228 |
rij=sqrt((xj(1)-xi(1))**2.0+(xj(2)-xi(2))**2.0+(xj(3)-xi(3))**2.0) |
229 |
n1=ka(i)/(np+1);n2=ka(mf(i,j))/(np+1) !!! Ni 46 H 50 C |
230 |
n3=ka(i)/(np+5);n4=ka(mf(i,j))/(np+5) |
231 |
|
232 |
!**********************************Ni-Ni interaction |
233 |
! ********** |
234 |
if(n1==0.and.n3==0.and.n2==0.and.n4==0)then !Ni-Ni |
235 |
call fctNiNi(rnn,rij,fcut) |
236 |
call dfcutNiNi(rnn,rij,tf) |
237 |
eb=eb+ann(5)*exp(-2.0*ann(1)*(rij/ann(3)-1.0))*fcut |
238 |
er=er+ann(4)*exp(-1.0*ann(2)*(rij/ann(3)-1.0))*fcut |
239 |
|
240 |
!Calculation of FORCE |
241 |
do im=1,3 |
242 |
dv(i,1,im)=dv(i,1,im)+ann(4)*exp(-1.0*ann(2)*(rij/ann(3)-1.0))*fcut & |
243 |
*(-ann(2)/ann(3))*(xi(im)-xj(im))/rij |
244 |
dv(i,2,im)=dv(i,2,im)+ann(5)*exp(-2.0*ann(1)*(rij/ann(3)-1.0))*fcut |
245 |
dv(i,3,im)=dv(i,3,im)+ann(5)*exp(-2.0*ann(1)*(rij/ann(3)-1.0))*fcut & |
246 |
*(-2.0*ann(1)/ann(3))*(xi(im)-xj(im))/rij |
247 |
dv(ka(mf(i+mid,j)),1,im)=dv(ka(mf(i+mid,j)),1,im)+ann(4)*exp(-1.0*ann(2) & |
248 |
*(rij/ann(3)-1.0))*fcut*(-ann(2)/ann(3))*(xj(im)-xi(im))/rij |
249 |
dv(ka(mf(i+mid,j)),2,im)=dv(ka(mf(i+mid,j)),2,im)+ann(5)*exp(-2.0*ann(1) & |
250 |
*(rij/ann(3)-1.0))*fcut |
251 |
dv(ka(mf(i+mid,j)),3,im)=dv(ka(mf(i+mid,j)),3,im)+ann(5)*exp(-2.0*ann(1) & |
252 |
*(rij/ann(3)-1.0))*fcut*(-2.0*ann(1)/ann(3))*(xj(im)-xi(im))/rij |
253 |
|
254 |
df(i,1,im)=df(i,1,im)+ann(4)*exp(-1.0*ann(2)*(rij/ann(3)-1.0))*tf & |
255 |
*(xi(im)-xj(im))/rij |
256 |
df(i,2,im)=df(i,2,im)+ann(5)*exp(-2.0*ann(1)*(rij/ann(3)-1.0))*fcut |
257 |
df(i,3,im)=df(i,3,im)+ann(5)*exp(-2.0*ann(1)*(rij/ann(3)-1.0))*tf & |
258 |
*(xi(im)-xj(im))/rij |
259 |
df(ka(mf(i+mid,j)),1,im)=df(ka(mf(i+mid,j)),1,im)+ann(4)*exp(-1.0*ann(2) & |
260 |
*(rij/ann(3)-1.0))*tf*(xj(im)-xi(im))/rij |
261 |
df(ka(mf(i+mid,j)),2,im)=df(ka(mf(i+mid,j)),2,im)+ann(5)*exp(-2.0*ann(1) & |
262 |
*(rij/ann(3)-1.0))*fcut |
263 |
df(ka(mf(i+mid,j)),3,im)=df(ka(mf(i+mid,j)),3,im)+ann(5)*exp(-2.0*ann(1) & |
264 |
*(rij/ann(3)-1.0))*tf*(xj(im)-xi(im))/rij |
265 |
enddo |
266 |
! |
267 |
endif |
268 |
|
269 |
!**********************************Ni-H interaction |
270 |
! ********** |
271 |
if(n1==0.and.n3==0.and.n2==1.and.n4==0)then !Ni-H |
272 |
call fctNiH(rnh,rij,fcut) |
273 |
call dfcutNiH(rnh,rij,tf) |
274 |
eb=eb+anh(5)*exp(-2.0*anh(1)*(rij/anh(3)-1.0))*fcut |
275 |
er=er+anh(4)*exp(-1.0*anh(2)*(rij/anh(3)-1.0))*fcut |
276 |
|
277 |
!Calculation of FORCE |
278 |
do im=1,3 |
279 |
dv(i,1,im)=dv(i,1,im)+anh(4)*exp(-1.0*anh(2)*(rij/anh(3)-1.0))*fcut & |
280 |
*(-anh(2)/anh(3))*(xi(im)-xj(im))/rij |
281 |
dv(i,2,im)=dv(i,2,im)+anh(5)*exp(-2.0*anh(1)*(rij/anh(3)-1.0))*fcut |
282 |
dv(i,3,im)=dv(i,3,im)+anh(5)*exp(-2.0*anh(1)*(rij/anh(3)-1.0))*fcut & |
283 |
*(-2.0*anh(1)/anh(3))*(xi(im)-xj(im))/rij |
284 |
dv(ka(mf(i+mid,j)),1,im)=dv(ka(mf(i+mid,j)),1,im)+anh(4)* & |
285 |
exp(-1.0*anh(2)*(rij/anh(3)-1.0))*fcut*(-anh(2)/anh(3)) & |
286 |
*(xj(im)-xi(im))/rij |
287 |
dv(ka(mf(i+mid,j)),2,im)=dv(ka(mf(i+mid,j)),2,im)+anh(5)* & |
288 |
exp(-2.0*anh(1)*(rij/anh(3)-1.0))*fcut |
289 |
dv(ka(mf(i+mid,j)),3,im)=dv(ka(mf(i+mid,j)),3,im)+anh(5)* & |
290 |
exp(-2.0*anh(1)*(rij/anh(3)-1.0))*fcut*(-2.0*anh(1)/anh(3)) & |
291 |
*(xj(im)-xi(im))/rij |
292 |
|
293 |
df(i,1,im)=df(i,1,im)+anh(4)*exp(-1.0*anh(2)*(rij/anh(3)-1.0))*tf & |
294 |
*(xi(im)-xj(im))/rij |
295 |
df(i,2,im)=df(i,2,im)+anh(5)*exp(-2.0*anh(1)*(rij/anh(3)-1.0))*fcut |
296 |
df(i,3,im)=df(i,3,im)+anh(5)*exp(-2.0*anh(1)*(rij/anh(3)-1.0))*tf & |
297 |
*(xi(im)-xj(im))/rij |
298 |
df(ka(mf(i+mid,j)),1,im)=df(ka(mf(i+mid,j)),1,im)+anh(4)*exp(-1.0*anh(2) & |
299 |
*(rij/anh(3)-1.0))*tf*(xj(im)-xi(im))/rij |
300 |
df(ka(mf(i+mid,j)),2,im)=df(ka(mf(i+mid,j)),2,im)+anh(5)*exp(-2.0*anh(1) & |
301 |
*(rij/anh(3)-1.0))*fcut |
302 |
df(ka(mf(i+mid,j)),3,im)=df(ka(mf(i+mid,j)),3,im)+anh(5)*exp(-2.0*anh(1) & |
303 |
*(rij/anh(3)-1.0))*tf*(xj(im)-xi(im))/rij |
304 |
enddo |
305 |
! |
306 |
! |
307 |
endif |
308 |
|
309 |
!**********************************Ni-C interaction |
310 |
! ********** |
311 |
if(n1==0.and.n3==0.and.n2==1.and.n4==1)then !Ni-C |
312 |
call fctNiC(rnc,rij,fcut) |
313 |
call dfcutNiC(rnc,rij,tf) |
314 |
eb=eb+anc(5)*exp(-2.0*anc(1)*(rij/anc(3)-1.0))*fcut |
315 |
er=er+anc(4)*exp(-1.0*anc(2)*(rij/anc(3)-1.0))*fcut |
316 |
|
317 |
!Calculation of FORCE |
318 |
do im=1,3 |
319 |
dv(i,1,im)=dv(i,1,im)+anc(4)*exp(-1.0*anc(2)*(rij/anc(3)-1.0))*fcut & |
320 |
*(-anc(2)/anc(3))*(xi(im)-xj(im))/rij |
321 |
dv(i,2,im)=dv(i,2,im)+anc(5)*exp(-2.0*anc(1)*(rij/anc(3)-1.0))*fcut |
322 |
dv(i,3,im)=dv(i,3,im)+anc(5)*exp(-2.0*anc(1)*(rij/anc(3)-1.0))*fcut & |
323 |
*(-2.0*anc(1)/anc(3))*(xi(im)-xj(im))/rij |
324 |
dv(ka(mf(i+mid,j)),1,im)=dv(ka(mf(i+mid,j)),1,im)+anc(4)*exp(-1.0* & |
325 |
anc(2)*(rij/anc(3)-1.0))*fcut*(-anc(2)/anc(3))*(xj(im)-xi(im))/rij |
326 |
dv(ka(mf(i+mid,j)),2,im)=dv(ka(mf(i+mid,j)),2,im)+anc(5)*exp(-2.0* & |
327 |
anc(1)*(rij/anc(3)-1.0))*fcut |
328 |
dv(ka(mf(i+mid,j)),3,im)=dv(ka(mf(i+mid,j)),3,im)+anc(5)*exp(-2.0* & |
329 |
anc(1)*(rij/anc(3)-1.0))*fcut*(-2.0*anc(1)/anc(3))*(xj(im)-xi(im))/rij |
330 |
|
331 |
df(i,1,im)=df(i,1,im)+anc(4)*exp(-1.0*anc(2)*(rij/anc(3)-1.0))*tf & |
332 |
*(xi(im)-xj(im))/rij |
333 |
df(i,2,im)=df(i,2,im)+anc(5)*exp(-2.0*anc(1)*(rij/anc(3)-1.0))*fcut |
334 |
df(i,3,im)=df(i,3,im)+anc(5)*exp(-2.0*anc(1)*(rij/anc(3)-1.0))*tf & |
335 |
*(xi(im)-xj(im))/rij |
336 |
df(ka(mf(i+mid,j)),1,im)=df(ka(mf(i+mid,j)),1,im)+anc(4)*exp(-1.0* & |
337 |
anc(2)*(rij/anc(3)-1.0))*tf*(xj(im)-xi(im))/rij |
338 |
df(ka(mf(i+mid,j)),2,im)=df(ka(mf(i+mid,j)),2,im)+anc(5)*exp(-2.0* & |
339 |
anc(1)*(rij/anc(3)-1.0))*fcut |
340 |
df(ka(mf(i+mid,j)),3,im)=df(ka(mf(i+mid,j)),3,im)+anc(5)*exp(-2.0* & |
341 |
anc(1)*(rij/anc(3)-1.0))*tf*(xj(im)-xi(im))/rij |
342 |
enddo |
343 |
! |
344 |
endif |
345 |
!!**********************************H-Ni interaction |
346 |
! ********** |
347 |
if(n1==1.and.n3==0.and.n2==0.and.n4==0)then ! H-Ni |
348 |
call fctNiH(rnh,rij,fcut) |
349 |
call dfcutNiH(rnh,rij,tf) |
350 |
eb=eb+anh(5)*exp(-2.0*anh(1)*(rij/anh(3)-1.0))*fcut |
351 |
er=er+anh(4)*exp(-1.0*anh(2)*(rij/anh(3)-1.0))*fcut |
352 |
|
353 |
!Calculation of FORCE |
354 |
do im=1,3 |
355 |
dv(i,1,im)=dv(i,1,im)+anh(4)*exp(-1.0*anh(2)*(rij/anh(3)-1.0))*fcut & |
356 |
*(-anh(2)/anh(3))*(xi(im)-xj(im))/rij |
357 |
dv(i,2,im)=dv(i,2,im)+anh(5)*exp(-2.0*anh(1)*(rij/anh(3)-1.0))*fcut |
358 |
dv(i,3,im)=dv(i,3,im)+anh(5)*exp(-2.0*anh(1)*(rij/anh(3)-1.0))*fcut & |
359 |
*(-2.0*anh(1)/anh(3))*(xi(im)-xj(im))/rij |
360 |
dv(ka(mf(i+mid,j)),1,im)=dv(ka(mf(i+mid,j)),1,im)+anh(4)*exp(-1.0*anh(2) & |
361 |
*(rij/anh(3)-1.0))*fcut*(-anh(2)/anh(3))*(xj(im)-xi(im))/rij |
362 |
dv(ka(mf(i+mid,j)),2,im)=dv(ka(mf(i+mid,j)),2,im)+anh(5)*exp(-2.0* & |
363 |
anh(1)*(rij/anh(3)-1.0))*fcut |
364 |
dv(ka(mf(i+mid,j)),3,im)=dv(ka(mf(i+mid,j)),3,im)+anh(5)*exp(-2.0* & |
365 |
anh(1)*(rij/anh(3)-1.0))*fcut*(-2.0*anh(1)/anh(3))*(xj(im)-xi(im))/rij |
366 |
|
367 |
df(i,1,im)=df(i,1,im)+anh(4)*exp(-1.0*anh(2)*(rij/anh(3)-1.0))*tf & |
368 |
*(xi(im)-xj(im))/rij |
369 |
df(i,2,im)=df(i,2,im)+anh(5)*exp(-2.0*anh(1)*(rij/anh(3)-1.0))*fcut |
370 |
df(i,3,im)=df(i,3,im)+anh(5)*exp(-2.0*anh(1)*(rij/anh(3)-1.0))*tf & |
371 |
*(xi(im)-xj(im))/rij |
372 |
df(ka(mf(i+mid,j)),1,im)=df(ka(mf(i+mid,j)),1,im)+anh(4)*exp(-1.0* & |
373 |
anh(2)*(rij/anh(3)-1.0))*tf*(xj(im)-xi(im))/rij |
374 |
df(ka(mf(i+mid,j)),2,im)=df(ka(mf(i+mid,j)),2,im)+anh(5)*exp(-2.0* & |
375 |
anh(1)*(rij/anh(3)-1.0))*fcut |
376 |
df(ka(mf(i+mid,j)),3,im)=df(ka(mf(i+mid,j)),3,im)+anh(5)*exp(-2.0* & |
377 |
anh(1)*(rij/anh(3)-1.0))*tf*(xj(im)-xi(im))/rij |
378 |
enddo |
379 |
! |
380 |
endif |
381 |
!!**********************************H-H interaction |
382 |
! ********** |
383 |
if(n1==1.and.n3==0.and.n2==1.and.n4==0)then !H-H |
384 |
call fctHH(rhh,rij,fcut) |
385 |
call dfcutHH(rhh,rij,tf) |
386 |
eb=eb+ahh(5)*exp(-2.0*ahh(1)*(rij/ahh(3)-1.0))*fcut |
387 |
er=er+ahh(4)*exp(-1.0*ahh(2)*(rij/ahh(3)-1.0))*fcut |
388 |
|
389 |
!Calculation of FORCE |
390 |
do im=1,3 |
391 |
dv(i,1,im)=dv(i,1,im)+ahh(4)*exp(-1.0*ahh(2)*(rij/ahh(3)-1.0))*fcut & |
392 |
*(-ahh(2)/ahh(3))*(xi(im)-xj(im))/rij |
393 |
dv(i,2,im)=dv(i,2,im)+ahh(5)*exp(-2.0*ahh(1)*(rij/ahh(3)-1.0))*fcut |
394 |
dv(i,3,im)=dv(i,3,im)+ahh(5)*exp(-2.0*ahh(1)*(rij/ahh(3)-1.0))*fcut & |
395 |
*(-2.0*ahh(1)/ahh(3))*(xi(im)-xj(im))/rij |
396 |
dv(ka(mf(i+mid,j)),1,im)=dv(ka(mf(i+mid,j)),1,im)+ahh(4)*exp(-1.0* & |
397 |
ahh(2)*(rij/ahh(3)-1.0))*fcut*(-ahh(2)/ahh(3))*(xj(im)-xi(im))/rij |
398 |
dv(ka(mf(i+mid,j)),2,im)=dv(ka(mf(i+mid,j)),2,im)+ahh(5)*exp(-2.0* & |
399 |
ahh(1)*(rij/ahh(3)-1.0))*fcut |
400 |
dv(ka(mf(i+mid,j)),3,im)=dv(ka(mf(i+mid,j)),3,im)+ahh(5)*exp(-2.0* & |
401 |
ahh(1)*(rij/ahh(3)-1.0))*fcut*(-2.0*ahh(1)/ahh(3))*(xj(im)-xi(im))/rij |
402 |
|
403 |
df(i,1,im)=df(i,1,im)+ahh(4)*exp(-1.0*ahh(2)*(rij/ahh(3)-1.0))*tf & |
404 |
*(xi(im)-xj(im))/rij |
405 |
df(i,2,im)=df(i,2,im)+ahh(5)*exp(-2.0*ahh(1)*(rij/ahh(3)-1.0))*fcut |
406 |
df(i,3,im)=df(i,3,im)+ahh(5)*exp(-2.0*ahh(1)*(rij/ahh(3)-1.0))*tf & |
407 |
*(xi(im)-xj(im))/rij |
408 |
df(ka(mf(i+mid,j)),1,im)=df(ka(mf(i+mid,j)),1,im)+ahh(4)*exp(-1.0* & |
409 |
ahh(2)*(rij/ahh(3)-1.0))*tf*(xj(im)-xi(im))/rij |
410 |
df(ka(mf(i+mid,j)),2,im)=df(ka(mf(i+mid,j)),2,im)+ahh(5)*exp(-2.0* & |
411 |
ahh(1)*(rij/ahh(3)-1.0))*fcut |
412 |
df(ka(mf(i+mid,j)),3,im)=df(ka(mf(i+mid,j)),3,im)+ahh(5)*exp(-2.0* & |
413 |
ahh(1)*(rij/ahh(3)-1.0))*tf*(xj(im)-xi(im))/rij |
414 |
enddo |
415 |
! |
416 |
endif |
417 |
|
418 |
|
419 |
!!**********************************H-C interaction |
420 |
! ********** |
421 |
if(n1==1.and.n3==0.and.n2==1.and.n4==1)then !H-C |
422 |
call fctCH(rch,rij,fcut) |
423 |
call dfcutCH(rch,rij,tf) |
424 |
eb=eb+ach(5)*exp(-2.0*ach(1)*(rij/ach(3)-1.0))*fcut |
425 |
er=er+ach(4)*exp(-1.0*ach(2)*(rij/ach(3)-1.0))*fcut |
426 |
|
427 |
!Calculation of FORCE |
428 |
do im=1,3 |
429 |
dv(i,1,im)=dv(i,1,im)+ach(4)*exp(-1.0*ach(2)*(rij/ach(3)-1.0))*fcut & |
430 |
*(-ach(2)/ach(3))*(xi(im)-xj(im))/rij |
431 |
dv(i,2,im)=dv(i,2,im)+ach(5)*exp(-2.0*ach(1)*(rij/ach(3)-1.0))*fcut |
432 |
dv(i,3,im)=dv(i,3,im)+ach(5)*exp(-2.0*ach(1)*(rij/ach(3)-1.0))*fcut & |
433 |
*(-2.0*ach(1)/ach(3))*(xi(im)-xj(im))/rij |
434 |
dv(ka(mf(i+mid,j)),1,im)=dv(ka(mf(i+mid,j)),1,im)+ach(4)*exp(-1.0* & |
435 |
ach(2)*(rij/ach(3)-1.0))*fcut*(-ach(2)/ach(3))*(xj(im)-xi(im))/rij |
436 |
dv(ka(mf(i+mid,j)),2,im)=dv(ka(mf(i+mid,j)),2,im)+ach(5)*exp(-2.0* & |
437 |
ach(1)*(rij/ach(3)-1.0))*fcut |
438 |
dv(ka(mf(i+mid,j)),3,im)=dv(ka(mf(i+mid,j)),3,im)+ach(5)*exp(-2.0* & |
439 |
ach(1)*(rij/ach(3)-1.0))*fcut*(-2.0*ach(1)/ach(3))*(xj(im)-xi(im))/rij |
440 |
|
441 |
df(i,1,im)=df(i,1,im)+ach(4)*exp(-1.0*ach(2)*(rij/ach(3)-1.0))*tf & |
442 |
*(xi(im)-xj(im))/rij |
443 |
df(i,2,im)=df(i,2,im)+ach(5)*exp(-2.0*ach(1)*(rij/ach(3)-1.0))*fcut |
444 |
df(i,3,im)=df(i,3,im)+ach(5)*exp(-2.0*ach(1)*(rij/ach(3)-1.0))*tf & |
445 |
*(xi(im)-xj(im))/rij |
446 |
df(ka(mf(i+mid,j)),1,im)=df(ka(mf(i+mid,j)),1,im)+ach(4)*exp(-1.0* & |
447 |
ach(2)*(rij/ach(3)-1.0))*tf*(xj(im)-xi(im))/rij |
448 |
df(ka(mf(i+mid,j)),2,im)=df(ka(mf(i+mid,j)),2,im)+ach(5)*exp(-2.0* & |
449 |
ach(1)*(rij/ach(3)-1.0))*fcut |
450 |
df(ka(mf(i+mid,j)),3,im)=df(ka(mf(i+mid,j)),3,im)+ach(5)*exp(-2.0* & |
451 |
ach(1)*(rij/ach(3)-1.0))*tf*(xj(im)-xi(im))/rij |
452 |
enddo |
453 |
! |
454 |
endif |
455 |
!!**********************************C-Ni interaction |
456 |
! ********** |
457 |
if(n1==1.and.n3==1.and.n2==0.and.n4==0)then !C-Ni |
458 |
call fctNiC(rnc,rij,fcut) |
459 |
call dfcutNiC(rnc,rij,tf) |
460 |
eb=eb+anc(5)*exp(-2.0*anc(1)*(rij/anc(3)-1.0))*fcut |
461 |
er=er+anc(4)*exp(-1.0*anc(2)*(rij/anc(3)-1.0))*fcut |
462 |
|
463 |
!Calculation of FORCE |
464 |
do im=1,3 |
465 |
dv(i,1,im)=dv(i,1,im)+anc(4)*exp(-1.0*anc(2)*(rij/anc(3)-1.0))*fcut & |
466 |
*(-anc(2)/anc(3))*(xi(im)-xj(im))/rij |
467 |
dv(i,2,im)=dv(i,2,im)+anc(5)*exp(-2.0*anc(1)*(rij/anc(3)-1.0))*fcut |
468 |
dv(i,3,im)=dv(i,3,im)+anc(5)*exp(-2.0*anc(1)*(rij/anc(3)-1.0))*fcut & |
469 |
*(-2.0*anc(1)/anc(3))*(xi(im)-xj(im))/rij |
470 |
dv(ka(mf(i+mid,j)),1,im)=dv(ka(mf(i+mid,j)),1,im)+anc(4)*exp(-1.0* & |
471 |
anc(2)*(rij/anc(3)-1.0))*fcut*(-anc(2)/anc(3))*(xj(im)-xi(im))/rij |
472 |
dv(ka(mf(i+mid,j)),2,im)=dv(ka(mf(i+mid,j)),2,im)+anc(5)*exp(-2.0* & |
473 |
anc(1)*(rij/anc(3)-1.0))*fcut |
474 |
dv(ka(mf(i+mid,j)),3,im)=dv(ka(mf(i+mid,j)),3,im)+anc(5)*exp(-2.0* & |
475 |
anc(1)*(rij/anc(3)-1.0))*fcut*(-2.0*anc(1)/anc(3))*(xj(im)-xi(im))/rij |
476 |
|
477 |
df(i,1,im)=df(i,1,im)+anc(4)*exp(-1.0*anc(2)*(rij/anc(3)-1.0))*tf & |
478 |
*(xi(im)-xj(im))/rij |
479 |
df(i,2,im)=df(i,2,im)+anc(5)*exp(-2.0*anc(1)*(rij/anc(3)-1.0))*fcut |
480 |
df(i,3,im)=df(i,3,im)+anc(5)*exp(-2.0*anc(1)*(rij/anc(3)-1.0))*tf & |
481 |
*(xi(im)-xj(im))/rij |
482 |
df(ka(mf(i+mid,j)),1,im)=df(ka(mf(i+mid,j)),1,im)+anc(4)*exp(-1.0* & |
483 |
anc(2)*(rij/anc(3)-1.0))*tf*(xj(im)-xi(im))/rij |
484 |
df(ka(mf(i+mid,j)),2,im)=df(ka(mf(i+mid,j)),2,im)+anc(5)*exp(-2.0* & |
485 |
anc(1)*(rij/anc(3)-1.0))*fcut |
486 |
df(ka(mf(i+mid,j)),3,im)=df(ka(mf(i+mid,j)),3,im)+anc(5)*exp(-2.0* & |
487 |
anc(1)*(rij/anc(3)-1.0))*tf*(xj(im)-xi(im))/rij |
488 |
enddo |
489 |
! |
490 |
|
491 |
endif |
492 |
!!**********************************C-H interaction |
493 |
! ********** |
494 |
if(n1==1.and.n3==1.and.n2==1.and.n4==0)then !C-H |
495 |
call fctCH(rch,rij,fcut) |
496 |
call dfcutCH(rch,rij,tf) |
497 |
eb=eb+ach(5)*exp(-2.0*ach(1)*(rij/ach(3)-1.0))*fcut |
498 |
er=er+ach(4)*exp(-1.0*ach(2)*(rij/ach(3)-1.0))*fcut |
499 |
|
500 |
|
501 |
!Calculation of FORCE |
502 |
do im=1,3 |
503 |
dv(i,1,im)=dv(i,1,im)+ach(4)*exp(-1.0*ach(2)*(rij/ach(3)-1.0))*fcut & |
504 |
*(-ach(2)/ach(3))*(xi(im)-xj(im))/rij |
505 |
dv(i,2,im)=dv(i,2,im)+ach(5)*exp(-2.0*ach(1)*(rij/ach(3)-1.0))*fcut |
506 |
dv(i,3,im)=dv(i,3,im)+ach(5)*exp(-2.0*ach(1)*(rij/ach(3)-1.0))*fcut & |
507 |
*(-2.0*ach(1)/ach(3))*(xi(im)-xj(im))/rij |
508 |
dv(ka(mf(i+mid,j)),1,im)=dv(ka(mf(i+mid,j)),1,im)+ach(4)*exp(-1.0* & |
509 |
ach(2)*(rij/ach(3)-1.0))*fcut*(-ach(2)/ach(3))*(xj(im)-xi(im))/rij |
510 |
dv(ka(mf(i+mid,j)),2,im)=dv(ka(mf(i+mid,j)),2,im)+ach(5)*exp(-2.0* & |
511 |
ach(1)*(rij/ach(3)-1.0))*fcut |
512 |
dv(ka(mf(i+mid,j)),3,im)=dv(ka(mf(i+mid,j)),3,im)+ach(5)*exp(-2.0* & |
513 |
ach(1)*(rij/ach(3)-1.0))*fcut*(-2.0*ach(1)/ach(3))*(xj(im)-xi(im))/rij |
514 |
|
515 |
df(i,1,im)=df(i,1,im)+ach(4)*exp(-1.0*ach(2)*(rij/ach(3)-1.0))*tf & |
516 |
*(xi(im)-xj(im))/rij |
517 |
df(i,2,im)=df(i,2,im)+ach(5)*exp(-2.0*ach(1)*(rij/ach(3)-1.0))*fcut |
518 |
df(i,3,im)=df(i,3,im)+ach(5)*exp(-2.0*ach(1)*(rij/ach(3)-1.0))*tf & |
519 |
*(xi(im)-xj(im))/rij |
520 |
df(ka(mf(i+mid,j)),1,im)=df(ka(mf(i+mid,j)),1,im)+ach(4)*exp(-1.0* & |
521 |
ach(2)*(rij/ach(3)-1.0))*tf*(xj(im)-xi(im))/rij |
522 |
df(ka(mf(i+mid,j)),2,im)=df(ka(mf(i+mid,j)),2,im)+ach(5)*exp(-2.0* & |
523 |
ach(1)*(rij/ach(3)-1.0))*fcut |
524 |
df(ka(mf(i+mid,j)),3,im)=df(ka(mf(i+mid,j)),3,im)+ach(5)*exp(-2.0* & |
525 |
ach(1)*(rij/ach(3)-1.0))*tf*(xj(im)-xi(im))/rij |
526 |
enddo |
527 |
! |
528 |
|
529 |
endif |
530 |
!!***************************ONE end**************** |
531 |
|
532 |
enddo |
533 |
!!******** CALCULATION ENERGY |
534 |
! |
535 |
yf=yf+er-sqrt(eb) |
536 |
!!********* |
537 |
do im=1,3 |
538 |
dvdr(i,im)=dvdr(i,im)+dv(i,1,im)-0.5/(dv(i,2,im))**0.5*dv(i,3,im) |
539 |
dvdr(i,im)=dvdr(i,im)+df(i,1,im)-0.5/(df(i,2,im))**0.5*df(i,3,im) |
540 |
do j=1,nf(i+mid) |
541 |
jm=ka(mf(i+mid,j)) |
542 |
dvdr(jm,im)=dvdr(jm,im)+dv(jm,1,im)-0.5/(dv(i,2,im))**0.5*dv(jm,3,im) |
543 |
dvdr(jm,im)=dvdr(jm,im)+df(jm,1,im)-0.5/(df(i,2,im))**0.5*df(jm,3,im) |
544 |
enddo |
545 |
enddo |
546 |
! |
547 |
enddo |
548 |
!!*****************************TWO ENDS************************** |
549 |
!***************************************Long-ranged interaction |
550 |
!do i=nq-4,nq |
551 |
! zi=pz(i+mid) !-0.3*aa*br(3,3) |
552 |
! call fct2(zshort,zlong,zi,fcut2) |
553 |
! call dfcut2(zshort,zlong,zi,tf) |
554 |
! yf=yf+fcut2*(a(6)-a(7)/zi**2.0) |
555 |
|
556 |
! dvdr(i,3)=dvdr(i,3)+fcut2*(a(7)*2.0/zi**3.0)+tf*(a(6)-a(7)/zi**2.0) |
557 |
|
558 |
!enddo |
559 |
!******************************************************************* |
560 |
v=yf |
561 |
do im=1,3 |
562 |
do jm=1,nq |
563 |
dvdr(jm,im)=dvdr(jm,im) |
564 |
enddo |
565 |
enddo |
566 |
|
567 |
!PFL: We convert dvdr into Grad |
568 |
do jm=1,nq |
569 |
Iat=Jm |
570 |
if (renum) Iat=order(jm) |
571 |
Grad(3*Iat-2:3*Iat)=dvdr(jm,1:3) |
572 |
enddo |
573 |
|
574 |
if (debug) WRITE(*,*) "EGRAD_CHAMFRE, E=",V |
575 |
|
576 |
if (debug) call Header('Egrad_chamfre OVER') |
577 |
|
578 |
!********************************************ENDENDENDENDENDENDENDEND |
579 |
return |
580 |
end subroutine egrad_chamfre |
581 |
|
582 |
!****************************************** |
583 |
!**********PARAMETERS********************** |
584 |
subroutine fctNiNi(r2,r,fcut) |
585 |
|
586 |
IMPLICIT NONE |
587 |
integer, parameter :: KINT = kind(1) |
588 |
integer, parameter :: KREAL = kind(1.0d0) |
589 |
|
590 |
REAL(KREAL) :: r2,r,fcut,r1,pi |
591 |
pi=4.0*atan(1.0) |
592 |
r1=r2-0.5 |
593 |
|
594 |
if(r<=r1)then |
595 |
fcut=1.0 |
596 |
endif |
597 |
if(r>r1.and.r<=r2)then |
598 |
fcut=(1.0+cos(pi*(r-r1)/(r2-r1)))/2.0 |
599 |
endif |
600 |
if(r>r2)then |
601 |
fcut=0.0 |
602 |
endif |
603 |
return |
604 |
end subroutine fctNiNi |
605 |
!!************ |
606 |
subroutine dfcutNiNi(r2,r,df) |
607 |
IMPLICIT NONE |
608 |
integer, parameter :: KINT = kind(1) |
609 |
integer, parameter :: KREAL = kind(1.0d0) |
610 |
REAL(KREAL) :: r, r2,df,pi,r1 |
611 |
pi=4.0*atan(1.0) |
612 |
r1=r2-0.5 |
613 |
|
614 |
if(r<=r1)then |
615 |
df=0.0 |
616 |
endif |
617 |
if(r>r1.and.r<=r2)then |
618 |
df=-pi/2.0/(r2-r1)*sin(pi*(r-r1)/(r2-r1)) |
619 |
endif |
620 |
if(r>r2)then |
621 |
df=0.0 |
622 |
endif |
623 |
return |
624 |
end subroutine dfcutNiNi |
625 |
|
626 |
|
627 |
!!************ |
628 |
subroutine fctNiH(r2,r,fcut) |
629 |
IMPLICIT NONE |
630 |
integer, parameter :: KINT = kind(1) |
631 |
integer, parameter :: KREAL = kind(1.0d0) |
632 |
|
633 |
REAL(KREAL) ::r2,r,fcut |
634 |
REAL(KREAL) ::rshort,rlong,r1,pi |
635 |
|
636 |
rshort=2.7 |
637 |
rlong=3.0 |
638 |
pi=4.0*atan(1.0) |
639 |
r1=rshort |
640 |
r2=rlong |
641 |
if(r<=r1)then |
642 |
fcut=1.0 |
643 |
endif |
644 |
if(r>r1.and.r<=r2)then |
645 |
fcut=(1.0+cos(pi*(r-r1)/(r2-r1)))/2.0 |
646 |
endif |
647 |
if(r>r2)then |
648 |
fcut=0.0 |
649 |
endif |
650 |
return |
651 |
end subroutine fctNiH |
652 |
!!********** |
653 |
|
654 |
subroutine dfcutNiH(r2,r,df) |
655 |
|
656 |
|
657 |
IMPLICIT NONE |
658 |
integer, parameter :: KINT = kind(1) |
659 |
integer, parameter :: KREAL = kind(1.0d0) |
660 |
|
661 |
|
662 |
REAL(KREAL) ::r2,r,df |
663 |
REAL(KREAL) ::rshort,rlong,r1,pi |
664 |
|
665 |
rshort=2.7 |
666 |
rlong=3.0 |
667 |
pi=4.0*atan(1.0) |
668 |
r1=rshort |
669 |
r2=rlong |
670 |
if(r<=r1)then |
671 |
df=0.0 |
672 |
endif |
673 |
if(r>r1.and.r<=r2)then |
674 |
df=-pi/2.0/(r2-r1)*sin(pi*(r-r1)/(r2-r1)) |
675 |
endif |
676 |
if(r>r2)then |
677 |
df=0.0 |
678 |
endif |
679 |
return |
680 |
end subroutine dfcutNiH |
681 |
|
682 |
!!************ |
683 |
subroutine fctNiC(r2,r,fcut) |
684 |
|
685 |
IMPLICIT NONE |
686 |
integer, parameter :: KINT = kind(1) |
687 |
integer, parameter :: KREAL = kind(1.0d0) |
688 |
|
689 |
|
690 |
REAL(KREAL) ::r2,r,fcut |
691 |
REAL(KREAL) ::rshort,rlong,r1,pi |
692 |
|
693 |
rshort=3.2 |
694 |
rlong=3.5 |
695 |
pi=4.0*atan(1.0) |
696 |
r1=rshort |
697 |
r2=rlong |
698 |
if(r<=r1)then |
699 |
fcut=1.0 |
700 |
endif |
701 |
if(r>r1.and.r<=r2)then |
702 |
fcut=(1.0+cos(pi*(r-r1)/(r2-r1)))/2.0 |
703 |
endif |
704 |
if(r>r2)then |
705 |
fcut=0.0 |
706 |
endif |
707 |
return |
708 |
end subroutine fctNiC |
709 |
|
710 |
!!*********** |
711 |
subroutine dfcutNiC(r2,r,df) |
712 |
|
713 |
IMPLICIT NONE |
714 |
integer, parameter :: KINT = kind(1) |
715 |
integer, parameter :: KREAL = kind(1.0d0) |
716 |
|
717 |
|
718 |
REAL(KREAL) ::r2,r,df |
719 |
REAL(KREAL) ::rshort,rlong,r1,pi |
720 |
|
721 |
rshort=3.2 |
722 |
rlong=3.5 |
723 |
pi=4.0*atan(1.0) |
724 |
r1=rshort |
725 |
r2=rlong |
726 |
if(r<=r1)then |
727 |
df=0.0 |
728 |
endif |
729 |
if(r>r1.and.r<=r2)then |
730 |
df=-pi/2.0/(r2-r1)*sin(pi*(r-r1)/(r2-r1)) |
731 |
endif |
732 |
if(r>r2)then |
733 |
df=0.0 |
734 |
endif |
735 |
return |
736 |
end subroutine dfcutNiC |
737 |
|
738 |
!!************ |
739 |
subroutine fctCH(r2,r,fcut) |
740 |
|
741 |
|
742 |
IMPLICIT NONE |
743 |
integer, parameter :: KINT = kind(1) |
744 |
integer, parameter :: KREAL = kind(1.0d0) |
745 |
|
746 |
|
747 |
REAL(KREAL) ::r2,r,fcut |
748 |
REAL(KREAL) ::rshort,rlong,r1,pi |
749 |
|
750 |
rshort=2.1 |
751 |
rlong=3.2 |
752 |
pi=4.0*atan(1.0) |
753 |
r1=rshort |
754 |
r2=rlong |
755 |
if(r<=r1)then |
756 |
fcut=1.0 |
757 |
endif |
758 |
if(r>r1.and.r<=r2)then |
759 |
fcut=(1.0+cos(pi*(r-r1)/(r2-r1)))/2.0 |
760 |
endif |
761 |
if(r>r2)then |
762 |
fcut=0.0 |
763 |
endif |
764 |
return |
765 |
end subroutine fctCH |
766 |
!!***************** |
767 |
|
768 |
subroutine dfcutCH(r2,r,df) |
769 |
|
770 |
|
771 |
IMPLICIT NONE |
772 |
integer, parameter :: KINT = kind(1) |
773 |
integer, parameter :: KREAL = kind(1.0d0) |
774 |
|
775 |
|
776 |
REAL(KREAL) ::r2,r,df |
777 |
REAL(KREAL) ::rshort,rlong,r1,pi |
778 |
|
779 |
rshort=2.1 |
780 |
rlong=3.2 |
781 |
pi=4.0*atan(1.0) |
782 |
r1=rshort |
783 |
r2=rlong |
784 |
if(r<=r1)then |
785 |
df=0.0 |
786 |
endif |
787 |
if(r>r1.and.r<=r2)then |
788 |
df=-pi/2.0/(r2-r1)*sin(pi*(r-r1)/(r2-r1)) |
789 |
endif |
790 |
if(r>r2)then |
791 |
df=0.0 |
792 |
endif |
793 |
return |
794 |
end subroutine dfcutCH |
795 |
!!************ |
796 |
subroutine fctHH(r2,r,fcut) |
797 |
|
798 |
IMPLICIT NONE |
799 |
integer, parameter :: KINT = kind(1) |
800 |
integer, parameter :: KREAL = kind(1.0d0) |
801 |
|
802 |
|
803 |
REAL(KREAL) ::r2,r,fcut,r1,pi |
804 |
pi=4.0*atan(1.0) |
805 |
r1=r2-1.1 |
806 |
if(r<=r1)then |
807 |
fcut=1.0 |
808 |
endif |
809 |
if(r>r1.and.r<=r2)then |
810 |
fcut=(1.0+cos(pi*(r-r1)/(r2-r1)))/2.0 |
811 |
endif |
812 |
if(r>r2)then |
813 |
fcut=0.0 |
814 |
endif |
815 |
return |
816 |
end subroutine fctHH |
817 |
!!************ |
818 |
subroutine dfcutHH(r2,r,df) |
819 |
|
820 |
|
821 |
IMPLICIT NONE |
822 |
integer, parameter :: KINT = kind(1) |
823 |
integer, parameter :: KREAL = kind(1.0d0) |
824 |
|
825 |
|
826 |
REAL(KREAL) ::r2,r,df,r1,pi |
827 |
pi=4.0*atan(1.0) |
828 |
r1=r2-1.1 |
829 |
if(r<=r1)then |
830 |
df=0.0 |
831 |
endif |
832 |
if(r>r1.and.r<=r2)then |
833 |
df=-pi/2.0/(r2-r1)*sin(pi*(r-r1)/(r2-r1)) |
834 |
endif |
835 |
if(r>r2)then |
836 |
df=0.0 |
837 |
endif |
838 |
return |
839 |
end subroutine dfcutHH |
840 |
|
841 |
!!*********************LONG INTERACTION |
842 |
|
843 |
subroutine dfcut2(zshort,zlong,r,df) |
844 |
|
845 |
|
846 |
IMPLICIT NONE |
847 |
integer, parameter :: KINT = kind(1) |
848 |
integer, parameter :: KREAL = kind(1.0d0) |
849 |
|
850 |
|
851 |
REAL(KREAL) ::zshort,zlong,r,df,pi, r1, r2 |
852 |
|
853 |
pi=4.0*atan(1.0) |
854 |
|
855 |
r1=zshort |
856 |
r2=zlong |
857 |
|
858 |
if(r<=r1)then |
859 |
df=0.0 |
860 |
endif |
861 |
if(r>r1.and.r<=r2)then |
862 |
df=pi/2.0/(r2-r1)*sin(pi*(r-r1)/(r2-r1)) |
863 |
endif |
864 |
if(r>r2)then |
865 |
df=0.0 |
866 |
endif |
867 |
|
868 |
return |
869 |
end subroutine dfcut2 |
870 |
|
871 |
subroutine fct2(zshort,zlong,r,fcut) |
872 |
|
873 |
IMPLICIT NONE |
874 |
integer, parameter :: KINT = kind(1) |
875 |
integer, parameter :: KREAL = kind(1.0d0) |
876 |
|
877 |
|
878 |
REAL(KREAL) ::zshort,zlong,r,fcut,pi,r1,r2 |
879 |
|
880 |
pi=4.0*atan(1.0) |
881 |
|
882 |
r1=zshort |
883 |
r2=zlong |
884 |
|
885 |
if(r<=r1)then |
886 |
fcut=0.0 |
887 |
endif |
888 |
if(r>r1.and.r<=r2)then |
889 |
fcut=1.0-(1.0+cos(pi*(r-r1)/(r2-r1)))/2.0 |
890 |
endif |
891 |
if(r>r2)then |
892 |
fcut=1.0 |
893 |
endif |
894 |
|
895 |
return |
896 |
end subroutine fct2 |
897 |
|