1 |
chuckv |
4 |
module goddard_module |
2 |
|
|
use definitions, ONLY : ndim,DP |
3 |
|
|
use simulation |
4 |
|
|
use parameter |
5 |
|
|
use second_deriv |
6 |
|
|
use force_utilities, ONLY : wrap,check,save_nlist |
7 |
|
|
#ifdef MPI |
8 |
|
|
use mpi_module |
9 |
|
|
#endif |
10 |
|
|
|
11 |
|
|
real( kind = DP ),allocatable, dimension(:) :: c |
12 |
|
|
real( kind = DP ),allocatable, dimension(:,:) :: BigD |
13 |
|
|
real( kind = DP ),allocatable, dimension(:,:) :: alpha |
14 |
|
|
real( kind = DP ),allocatable, dimension(:,:) :: m |
15 |
|
|
real( kind = DP ),allocatable, dimension(:,:) :: n |
16 |
|
|
real( kind = DP ),allocatable, dimension(:,:) :: rcutg |
17 |
|
|
real( kind = DP ),allocatable, dimension(:,:) :: vpair_rcut |
18 |
|
|
real( kind = DP ),allocatable, dimension(:,:) :: rho_rcut |
19 |
|
|
|
20 |
|
|
!! private force arrays |
21 |
|
|
|
22 |
|
|
|
23 |
|
|
|
24 |
|
|
|
25 |
|
|
|
26 |
|
|
|
27 |
|
|
private :: c, BigD, alpha, m, n, vpair_rcut |
28 |
|
|
public :: rcutg |
29 |
|
|
!private :: allocate_goddard_module, mass_weight |
30 |
|
|
|
31 |
|
|
public :: calc_goddard_dens,calc_goddard_forces,initialize_goddard |
32 |
|
|
!public :: deallocate_goddard_module |
33 |
|
|
|
34 |
|
|
contains |
35 |
|
|
|
36 |
|
|
|
37 |
|
|
subroutine allocate_goddard_module(n_size_atype) |
38 |
|
|
integer, intent(in) :: n_size_atype |
39 |
|
|
|
40 |
|
|
|
41 |
|
|
allocate(c(n_size_atype)) |
42 |
|
|
allocate(BigD(n_size_atype,n_size_atype)) |
43 |
|
|
allocate(alpha(n_size_atype,n_size_atype)) |
44 |
|
|
allocate(m(n_size_atype,n_size_atype)) |
45 |
|
|
allocate(n(n_size_atype,n_size_atype)) |
46 |
|
|
allocate(rcutg(n_size_atype,n_size_atype)) |
47 |
|
|
allocate(vpair_rcut(n_size_atype,n_size_atype)) |
48 |
|
|
allocate(rho_rcut(n_size_atype,n_size_atype)) |
49 |
|
|
|
50 |
|
|
end subroutine allocate_goddard_module |
51 |
|
|
|
52 |
|
|
subroutine deallocate_goddard_module |
53 |
|
|
deallocate(c) |
54 |
|
|
deallocate(BigD) |
55 |
|
|
deallocate(alpha) |
56 |
|
|
deallocate(m) |
57 |
|
|
deallocate(n) |
58 |
|
|
deallocate(rcutg) |
59 |
|
|
deallocate(vpair_rcut) |
60 |
|
|
deallocate(rho_rcut) |
61 |
|
|
|
62 |
|
|
end subroutine deallocate_goddard_module |
63 |
|
|
|
64 |
|
|
subroutine calc_goddard_dens(update_nlist) |
65 |
|
|
|
66 |
|
|
! include 'headers/sizes.h' |
67 |
|
|
|
68 |
|
|
|
69 |
|
|
real( kind = DP ) ptmp |
70 |
|
|
integer :: i, j, atype1, atype2, nlist, jbeg, jend, jnab |
71 |
|
|
real( kind = DP ) :: rxi, ryi, rzi, rxij, ryij, rzij, rijsq |
72 |
|
|
real( kind = DP ) :: aij, mij, rcij, r |
73 |
|
|
integer :: j_start |
74 |
|
|
integer :: tag_i,tag_j |
75 |
|
|
logical, intent(inout) :: update_nlist |
76 |
|
|
|
77 |
|
|
#ifndef MPI |
78 |
|
|
integer :: nrow |
79 |
|
|
integer :: ncol |
80 |
|
|
|
81 |
|
|
nrow = natoms - 1 |
82 |
|
|
ncol = natoms |
83 |
|
|
#endif |
84 |
|
|
|
85 |
|
|
|
86 |
|
|
#ifdef MPI |
87 |
|
|
rho_row = 0.0E0_DP |
88 |
|
|
rho_col = 0.0E0_DP |
89 |
|
|
j_start = 1 |
90 |
|
|
#else |
91 |
|
|
rho = 0.0E0_DP |
92 |
|
|
#endif |
93 |
|
|
|
94 |
|
|
if (update_nlist) then |
95 |
|
|
|
96 |
|
|
! save current configuration, contruct neighbor list, |
97 |
|
|
! and calculate forces |
98 |
|
|
call save_nlist() |
99 |
|
|
|
100 |
|
|
nlist = 0 |
101 |
|
|
|
102 |
|
|
do i = 1, nrow |
103 |
|
|
point(i) = nlist + 1 |
104 |
|
|
#ifdef MPI |
105 |
|
|
tag_i = tag_row(i) |
106 |
|
|
rxi = q_row(1,i) |
107 |
|
|
ryi = q_row(2,i) |
108 |
|
|
rzi = q_row(3,i) |
109 |
|
|
#else |
110 |
|
|
j_start = i + 1 |
111 |
|
|
rxi = q(1,i) |
112 |
|
|
ryi = q(2,i) |
113 |
|
|
rzi = q(3,i) |
114 |
|
|
#endif |
115 |
|
|
inner: do j = j_start, ncol |
116 |
|
|
#ifdef MPI |
117 |
|
|
tag_j = tag_col(j) |
118 |
|
|
if (newtons_thrd) then |
119 |
|
|
if (tag_i <= tag_j) then |
120 |
|
|
if (mod(tag_i + tag_j,2) == 0) cycle inner |
121 |
|
|
else |
122 |
|
|
if (mod(tag_i + tag_j,2) == 1) cycle inner |
123 |
|
|
endif |
124 |
|
|
|
125 |
|
|
endif |
126 |
|
|
|
127 |
|
|
|
128 |
|
|
|
129 |
|
|
rxij = wrap(rxi - q_col(1,j), 1) |
130 |
|
|
ryij = wrap(ryi - q_col(2,j), 2) |
131 |
|
|
rzij = wrap(rzi - q_col(3,j), 3) |
132 |
|
|
|
133 |
|
|
#else |
134 |
|
|
rxij = wrap(rxi - q(1,j),1) |
135 |
|
|
ryij = wrap(ryi - q(2,j),2) |
136 |
|
|
rzij = wrap(rzi - q(3,j),3) |
137 |
|
|
#endif |
138 |
|
|
rijsq = rxij*rxij + ryij*ryij + rzij*rzij |
139 |
|
|
|
140 |
|
|
#ifdef MPI |
141 |
|
|
if (rijsq <= rlstsq .AND. & |
142 |
|
|
tag_j /= tag_i) then |
143 |
|
|
#else |
144 |
|
|
|
145 |
|
|
if (rijsq < rlstsq) then |
146 |
|
|
#endif |
147 |
|
|
|
148 |
|
|
nlist = nlist + 1 |
149 |
|
|
list(nlist) = j |
150 |
|
|
|
151 |
|
|
! rcutsq is the largest of the pairwise cutoffs |
152 |
|
|
|
153 |
|
|
if (rijsq < rcutsq) then |
154 |
|
|
|
155 |
|
|
! go ahead and do the atomic interactions: |
156 |
|
|
|
157 |
|
|
#ifdef MPI |
158 |
|
|
atype1 = ident_row(i) |
159 |
|
|
atype2 = ident_col(j) |
160 |
|
|
#else |
161 |
|
|
atype1 = ident(i) |
162 |
|
|
atype2 = ident(j) |
163 |
|
|
#endif |
164 |
|
|
aij = alpha(atype1, atype2) |
165 |
|
|
mij = m(atype1, atype2) |
166 |
|
|
rcij = rcutg(atype1, atype2) |
167 |
|
|
r = dsqrt(rijsq) |
168 |
|
|
|
169 |
|
|
!rcij is the pairwise cutoff radius |
170 |
|
|
if (r <= rcij) then |
171 |
|
|
|
172 |
|
|
ptmp = (aij/r)**mij |
173 |
|
|
#ifdef MPI |
174 |
|
|
rho_row(i) = rho_row(i) + ptmp |
175 |
|
|
rho_col(i) = rho_col(i) + ptmp |
176 |
|
|
#else |
177 |
|
|
rho(i) = rho(i) + ptmp |
178 |
|
|
rho(j) = rho(j) + ptmp |
179 |
|
|
#endif |
180 |
|
|
|
181 |
|
|
endif |
182 |
|
|
endif |
183 |
|
|
endif |
184 |
|
|
enddo inner |
185 |
|
|
enddo |
186 |
|
|
#ifdef MP |
187 |
|
|
point(nrow + 1) = nlist + 1 |
188 |
|
|
#else |
189 |
|
|
point(natoms) = nlist + 1 |
190 |
|
|
#endif |
191 |
|
|
|
192 |
|
|
else |
193 |
|
|
! use the list to find the neighbors |
194 |
|
|
do i = 1, natoms -1 |
195 |
|
|
JBEG = POINT(I) |
196 |
|
|
JEND = POINT(I+1) - 1 |
197 |
|
|
! check thiat molecule i has neighbors |
198 |
|
|
if (jbeg <= jend) then |
199 |
|
|
#ifdef MPI |
200 |
|
|
rxi = q_row(1,i) |
201 |
|
|
ryi = q_row(2,i) |
202 |
|
|
rzi = q_row(3,i) |
203 |
|
|
#else |
204 |
|
|
rxi = q(1,i) |
205 |
|
|
ryi = q(2,i) |
206 |
|
|
rzi = q(3,i) |
207 |
|
|
#endif |
208 |
|
|
do jnab = jbeg, jend |
209 |
|
|
j = list(jnab) |
210 |
|
|
#ifdef MPI |
211 |
|
|
rxij = wrap(rxi - q_col(1,j), 1) |
212 |
|
|
ryij = wrap(ryi - q_col(2,j), 2) |
213 |
|
|
rzij = wrap(rzi - q_col(3,j), 3) |
214 |
|
|
#else |
215 |
|
|
rxij = wrap(rxi - q(1,j), 1) |
216 |
|
|
ryij = wrap(ryi - q(2,j), 2) |
217 |
|
|
rzij = wrap(rzi - q(3,j), 3) |
218 |
|
|
#endif |
219 |
|
|
rijsq = rxij*rxij + ryij*ryij + rzij*rzij |
220 |
|
|
|
221 |
|
|
if (rijsq .lt. rcutsq) then |
222 |
|
|
! go ahead and do the atomic interactions: |
223 |
|
|
#ifdef MPI |
224 |
|
|
atype1 = ident_row(i) |
225 |
|
|
atype2 = ident_col(j) |
226 |
|
|
#else |
227 |
|
|
atype1 = ident(i) |
228 |
|
|
atype2 = ident(j) |
229 |
|
|
#endif |
230 |
|
|
|
231 |
|
|
aij = alpha(atype1, atype2) |
232 |
|
|
mij = m(atype1, atype2) |
233 |
|
|
rcij = rcutg(atype1, atype2) |
234 |
|
|
r = dsqrt(rijsq) |
235 |
|
|
|
236 |
|
|
!rcij is the pairwise cutoff radius |
237 |
|
|
|
238 |
|
|
if (r.le.rcij) then |
239 |
|
|
|
240 |
|
|
ptmp = (aij/r)**mij |
241 |
|
|
|
242 |
|
|
#ifdef MPI |
243 |
|
|
rho_row(i) = rho_row(i) + ptmp |
244 |
|
|
rho_col(j) = rho_col(j) + ptmp |
245 |
|
|
#else |
246 |
|
|
rho(i) = rho(i) + ptmp |
247 |
|
|
rho(j) = rho(j) + ptmp |
248 |
|
|
#endif |
249 |
|
|
endif |
250 |
|
|
endif |
251 |
|
|
enddo |
252 |
|
|
endif |
253 |
|
|
enddo |
254 |
|
|
endif |
255 |
|
|
|
256 |
|
|
#ifdef MPI |
257 |
|
|
!! communicate densities |
258 |
|
|
call scatter(rho_row,rho,plan_row) |
259 |
|
|
if (newtons_thrd) then |
260 |
|
|
call scatter(rho_col,rho_tmp,plan_col) |
261 |
|
|
do i = 1, nlocal |
262 |
|
|
rho(i) = rho(i) + rho_tmp(i) |
263 |
|
|
end do |
264 |
|
|
endif |
265 |
|
|
#endif |
266 |
|
|
|
267 |
|
|
|
268 |
|
|
|
269 |
|
|
return |
270 |
|
|
end subroutine calc_goddard_dens |
271 |
|
|
|
272 |
|
|
subroutine calc_goddard_forces(nmflag,pot) |
273 |
|
|
|
274 |
|
|
! include 'headers/sizes.h' |
275 |
|
|
|
276 |
|
|
|
277 |
|
|
#ifdef MPI |
278 |
|
|
real( kind = DP ), dimension(nlocal) :: frho |
279 |
|
|
real( kind = DP ), dimension(nrow) :: frho_row |
280 |
|
|
real( kind = DP ), dimension(ncol) :: frho_col |
281 |
|
|
|
282 |
|
|
real( kind = DP ), dimension(nlocal) :: dfrhodrho |
283 |
|
|
real( kind = DP ), dimension(nrow) :: dfrhodrho_row |
284 |
|
|
real( kind = DP ), dimension(ncol) :: dfrhodrho_col |
285 |
|
|
|
286 |
|
|
real( kind = DP ), dimension(nlocal) :: d2frhodrhodrho |
287 |
|
|
real( kind = DP ), dimension(nrow) :: d2frhodrhodrho_row |
288 |
|
|
real( kind = DP ), dimension(ncol) :: d2frhodrhodrho_col |
289 |
|
|
real( kind = DP ), dimension(ndim,ncol) :: efr |
290 |
|
|
|
291 |
|
|
real( kind = DP ) :: pot_local, pot_phi_row, pot_Frho, pot_phi, pot_row |
292 |
|
|
|
293 |
|
|
#else |
294 |
|
|
real( kind = DP ), dimension(natoms) :: frho |
295 |
|
|
real( kind = DP ), dimension(natoms) :: dfrhodrho |
296 |
|
|
real( kind = DP ), dimension(natoms) :: d2frhodrhodrho |
297 |
|
|
real( kind = DP ), dimension(ndim,natoms) :: efr |
298 |
|
|
#endif |
299 |
|
|
|
300 |
|
|
|
301 |
|
|
real( kind = DP ), intent(out), optional :: pot |
302 |
|
|
real( kind = DP ) :: aij, mij, rcij, nij, dij, vcij, vptmp, dudr, ftmp |
303 |
|
|
real( kind = DP ) :: drhodr, dvpdr, drdx1, d2vpdrdr, d2rhodrdr, d2 |
304 |
|
|
real( kind = DP ) :: kt1, kt2, kt3, ktmp |
305 |
|
|
|
306 |
|
|
integer i, j, dim, atype1, atype2, idim, jdim, dim2, idim2, jdim2 |
307 |
|
|
integer jbeg, jend, jnab |
308 |
|
|
real( kind = DP ) :: rxij, ryij, rzij, rxi, ryi, rzi, rijsq, r |
309 |
|
|
|
310 |
|
|
|
311 |
|
|
integer :: nlist |
312 |
|
|
logical, intent(in) :: nmflag |
313 |
|
|
logical :: do_pot |
314 |
|
|
|
315 |
|
|
|
316 |
|
|
#ifndef MPI |
317 |
|
|
integer :: nrow |
318 |
|
|
integer :: ncol |
319 |
|
|
|
320 |
|
|
|
321 |
|
|
nrow = natoms - 1 |
322 |
|
|
ncol = natoms |
323 |
|
|
#endif |
324 |
|
|
|
325 |
|
|
do_pot = .false. |
326 |
|
|
if (present(pot)) do_pot = .true. |
327 |
|
|
|
328 |
|
|
|
329 |
|
|
#ifdef MPI |
330 |
|
|
f_row = 0.0E0_DP |
331 |
|
|
f_col = 0.0E0_DP |
332 |
|
|
|
333 |
|
|
pot_phi_row = 0.0E0_DP |
334 |
|
|
pot_phi = 0.0E0_DP |
335 |
|
|
pot_Frho = 0.0E0_DP |
336 |
|
|
pot_local = 0.0E0_DP |
337 |
|
|
pot_row = 0.0E0_DP |
338 |
|
|
|
339 |
|
|
e_row = 0.0E0_DP |
340 |
|
|
e_col = 0.0E0_DP |
341 |
|
|
e_tmp = 0.0E0_DP |
342 |
|
|
#else |
343 |
|
|
if (do_pot) pot = 0.0E0_DP |
344 |
|
|
f = 0.0E0_DP |
345 |
|
|
frho = 0.0E0_DP |
346 |
|
|
dfrhodrho = 0.0E0_DP |
347 |
|
|
d2frhodrhodrho = 0.0E0_DP |
348 |
|
|
efr = 0.0E0_DP |
349 |
|
|
#endif |
350 |
|
|
|
351 |
|
|
do i = 1, nlocal |
352 |
|
|
atype1 = ident(i) |
353 |
|
|
frho(i) = -c(atype1)*BigD(atype1,atype1)*dsqrt(rho(i)) |
354 |
|
|
dfrhodrho(i) = 0.5E0_DP*frho(i)/rho(i) |
355 |
|
|
d2frhodrhodrho(i) = -0.5E0_DP*dfrhodrho(i)/rho(i) |
356 |
|
|
#ifndef MPI |
357 |
|
|
if (do_pot) pot = pot + frho(i) |
358 |
|
|
#endif |
359 |
|
|
enddo |
360 |
|
|
|
361 |
|
|
#ifdef MPI |
362 |
|
|
!! communicate f(rho) and derivatives |
363 |
|
|
|
364 |
|
|
call gather(frho,frho_row,plan_row) |
365 |
|
|
call gather(dfrhodrho,dfrhodrho_row,plan_row) |
366 |
|
|
call gather(frho,frho_col,plan_col) |
367 |
|
|
call gather(dfrhodrho,dfrhodrho_col,plan_col) |
368 |
|
|
|
369 |
|
|
if (nmflag) then |
370 |
|
|
call gather(d2frhodrhodrho,d2frhodrhodrho_row,plan_row) |
371 |
|
|
call gather(d2frhodrhodrho,d2frhodrhodrho_col,plan_col) |
372 |
|
|
endif |
373 |
|
|
#endif |
374 |
|
|
|
375 |
|
|
|
376 |
|
|
|
377 |
|
|
do i = 1, nrow |
378 |
|
|
JBEG = POINT(i) |
379 |
|
|
JEND = POINT(i+1) - 1 |
380 |
|
|
! check thiat molecule i has neighbors |
381 |
|
|
if (jbeg <= jend) then |
382 |
|
|
#ifdef MPI |
383 |
|
|
atype1 = ident_row(i) |
384 |
|
|
rxi = q_row(1,i) |
385 |
|
|
ryi = q_row(2,i) |
386 |
|
|
rzi = q_row(3,i) |
387 |
|
|
#else |
388 |
|
|
atype1 = ident(i) |
389 |
|
|
rxi = q(1,i) |
390 |
|
|
ryi = q(2,i) |
391 |
|
|
rzi = q(3,i) |
392 |
|
|
#endif |
393 |
|
|
do jnab = jbeg, jend |
394 |
|
|
j = list(jnab) |
395 |
|
|
#ifdef MPI |
396 |
|
|
rxij = wrap(rxi - q_col(1,j), 1) |
397 |
|
|
ryij = wrap(ryi - q_col(2,j), 2) |
398 |
|
|
rzij = wrap(rzi - q_col(3,j), 3) |
399 |
|
|
#else |
400 |
|
|
rxij = wrap(rxi - q(1,j), 1) |
401 |
|
|
ryij = wrap(ryi - q(2,j), 2) |
402 |
|
|
rzij = wrap(rzi - q(3,j), 3) |
403 |
|
|
#endif |
404 |
|
|
rijsq = rxij*rxij + ryij*ryij + rzij*rzij |
405 |
|
|
|
406 |
|
|
if (rijsq .lt. rcutsq) then |
407 |
|
|
#ifdef MPI |
408 |
|
|
atype2 = ident_col(j) |
409 |
|
|
#else |
410 |
|
|
atype2 = ident(j) |
411 |
|
|
#endif |
412 |
|
|
r = dsqrt(rijsq) |
413 |
|
|
rcij = rcutg(atype1, atype2) |
414 |
|
|
|
415 |
|
|
if (r <= rcij) then |
416 |
|
|
|
417 |
|
|
! distance is within interaction region: |
418 |
|
|
efr(1,j) = -rxij |
419 |
|
|
efr(2,j) = -ryij |
420 |
|
|
efr(3,j) = -rzij |
421 |
|
|
|
422 |
|
|
dij = BigD(atype1,atype2) |
423 |
|
|
aij = alpha(atype1,atype2) |
424 |
|
|
nij = n(atype1,atype2) |
425 |
|
|
mij = m(atype1, atype2) |
426 |
|
|
vcij = vpair_rcut(atype1,atype2) |
427 |
|
|
|
428 |
|
|
vptmp = dij*((aij/r)**nij) |
429 |
|
|
#ifdef MPI |
430 |
|
|
e_row(i) = e_row(i) + (vptmp + vcij)*0.5 |
431 |
|
|
e_col(i) = e_col(i) + (vptmp + vcij)*0.5 |
432 |
|
|
#else |
433 |
|
|
if (do_pot) pot = pot + vptmp + vcij |
434 |
|
|
#endif |
435 |
|
|
|
436 |
|
|
dvpdr = -nij*vptmp/r |
437 |
|
|
d2vpdrdr = -dvpdr*(nij+1.0E0_DP)/r |
438 |
|
|
|
439 |
|
|
drhodr = -mij*((aij/r)**mij)/r |
440 |
|
|
d2rhodrdr = -drhodr*(mij+1.0E0_DP)/r |
441 |
|
|
|
442 |
|
|
#ifdef MPI |
443 |
|
|
dudr = drhodr*(dfrhodrho_row(i)+dfrhodrho_col(j)) & |
444 |
|
|
+ dvpdr |
445 |
|
|
|
446 |
|
|
d2 = d2vpdrdr + d2rhodrdr*(dfrhodrho_row(i)+dfrhodrho_col(j)) & |
447 |
|
|
+ drhodr*drhodr* & |
448 |
|
|
(d2frhodrhodrho_row(i)+d2frhodrhodrho_col(j)) |
449 |
|
|
#else |
450 |
|
|
dudr = drhodr*(dfrhodrho(i)+dfrhodrho(j)) & |
451 |
|
|
+ dvpdr |
452 |
|
|
|
453 |
|
|
d2 = d2vpdrdr + d2rhodrdr*(dfrhodrho(i)+dfrhodrho(j)) & |
454 |
|
|
+ drhodr*drhodr* & |
455 |
|
|
(d2frhodrhodrho(i)+d2frhodrhodrho(j)) |
456 |
|
|
#endif |
457 |
|
|
|
458 |
|
|
do dim = 1, 3 |
459 |
|
|
|
460 |
|
|
drdx1 = efr(dim,j) / r |
461 |
|
|
ftmp = dudr * drdx1 |
462 |
|
|
#ifdef MPI |
463 |
|
|
f_col(dim,j) = f_col(dim,j) - ftmp |
464 |
|
|
f_row(dim,i) = f_row(dim,i) + ftmp |
465 |
|
|
#else |
466 |
|
|
f(dim,j) = f(dim,j) - ftmp |
467 |
|
|
f(dim,i) = f(dim,i) + ftmp |
468 |
|
|
#endif |
469 |
|
|
|
470 |
|
|
|
471 |
|
|
if (nmflag) then |
472 |
|
|
idim = 3 * (i-1) + dim |
473 |
|
|
jdim = 3 * (j-1) + dim |
474 |
|
|
|
475 |
|
|
do dim2 = 1, 3 |
476 |
|
|
|
477 |
|
|
kt1 = d2 * efr(dim,j) * efr(dim2,j)/r/r |
478 |
|
|
kt2 = - dudr * efr(dim,j) * efr(dim2,j)/r/r/r |
479 |
|
|
|
480 |
|
|
if (dim.eq.dim2) then |
481 |
|
|
kt3 = dudr / r |
482 |
|
|
else |
483 |
|
|
kt3 = 0.0E0_DP |
484 |
|
|
endif |
485 |
|
|
|
486 |
|
|
! The factor of 2 below is to compensate for |
487 |
|
|
! overcounting. |
488 |
|
|
! Mass weighting is done separately... |
489 |
|
|
|
490 |
|
|
ktmp = (kt1+kt2+kt3)/2.0E0_DP |
491 |
|
|
idim2 = 3 * (i-1) + dim2 |
492 |
|
|
jdim2 = 3 * (j-1) + dim2 |
493 |
|
|
|
494 |
|
|
d(idim, idim2) = d(idim,idim2) + ktmp |
495 |
|
|
d(idim2, idim) = d(idim2,idim) + ktmp |
496 |
|
|
|
497 |
|
|
d(idim, jdim2) = d(idim,jdim2) - ktmp |
498 |
|
|
d(idim2, jdim) = d(idim2,jdim) - ktmp |
499 |
|
|
|
500 |
|
|
d(jdim, idim2) = d(jdim,idim2) - ktmp |
501 |
|
|
d(jdim2, idim) = d(jdim2,idim) - ktmp |
502 |
|
|
|
503 |
|
|
d(jdim, jdim2) = d(jdim,jdim2) + ktmp |
504 |
|
|
d(jdim2, jdim) = d(jdim2,jdim) + ktmp |
505 |
|
|
|
506 |
|
|
enddo |
507 |
|
|
endif |
508 |
|
|
enddo |
509 |
|
|
|
510 |
|
|
endif |
511 |
|
|
endif |
512 |
|
|
enddo |
513 |
|
|
endif |
514 |
|
|
enddo |
515 |
|
|
|
516 |
|
|
#ifdef MPI |
517 |
|
|
!!distribute forces |
518 |
|
|
call scatter(f_row,f,plan_row3) |
519 |
|
|
if (newtons_thrd) then |
520 |
|
|
call scatter(f_col,f_tmp,plan_col3) |
521 |
|
|
do i = 1,nlocal |
522 |
|
|
do dim = 1,3 |
523 |
|
|
f(dim,i) = f(dim,i) + f_tmp(dim,i) |
524 |
|
|
end do |
525 |
|
|
end do |
526 |
|
|
endif |
527 |
|
|
|
528 |
|
|
|
529 |
|
|
if (do_pot) then |
530 |
|
|
! scatter/gather pot_row into the members of my column |
531 |
|
|
call scatter(e_row,e_tmp,plan_row) |
532 |
|
|
|
533 |
|
|
! scatter/gather pot_local into all other procs |
534 |
|
|
! add resultant to get total pot |
535 |
|
|
do i = 1, nlocal |
536 |
|
|
pot_local = pot_local + frho(i) + e_tmp(i) |
537 |
|
|
enddo |
538 |
|
|
if (newtons_thrd) then |
539 |
|
|
e_tmp = 0.0E0_DP |
540 |
|
|
call scatter(e_col,e_tmp,plan_col) |
541 |
|
|
do i = 1, nlocal |
542 |
|
|
pot_local = pot_local + e_tmp(i) |
543 |
|
|
enddo |
544 |
|
|
endif |
545 |
|
|
endif |
546 |
|
|
#endif |
547 |
|
|
|
548 |
|
|
|
549 |
|
|
|
550 |
|
|
if (nmflag) then |
551 |
|
|
call mass_weight() |
552 |
|
|
endif |
553 |
|
|
|
554 |
|
|
return |
555 |
|
|
end subroutine calc_goddard_forces |
556 |
|
|
|
557 |
|
|
subroutine initialize_goddard() |
558 |
|
|
use model_module, ONLY: get_max_atype |
559 |
|
|
include 'headers/atom.h' |
560 |
|
|
|
561 |
|
|
|
562 |
|
|
integer n_gatypes |
563 |
|
|
parameter (n_gatypes = 8) |
564 |
|
|
integer gatype(n_gatypes), atype1, atype2, i, j |
565 |
|
|
integer :: n_size_atypes |
566 |
|
|
|
567 |
|
|
n_size_atypes = get_max_atype() |
568 |
|
|
|
569 |
|
|
call allocate_goddard_module(n_size_atypes) |
570 |
|
|
|
571 |
|
|
|
572 |
|
|
gatype(1) = Ni_atom |
573 |
|
|
gatype(2) = Cu_atom |
574 |
|
|
gatype(3) = Rh_atom |
575 |
|
|
gatype(4) = Pd_atom |
576 |
|
|
gatype(5) = Ag_atom |
577 |
|
|
gatype(6) = Ir_atom |
578 |
|
|
gatype(7) = Pt_atom |
579 |
|
|
gatype(8) = Au_atom |
580 |
|
|
|
581 |
|
|
! first set up the primaries: |
582 |
|
|
|
583 |
|
|
c(Ni_atom) = 84.745E0_DP |
584 |
|
|
c(Cu_atom) = 84.843E0_DP |
585 |
|
|
c(Rh_atom) = 305.499E0_DP |
586 |
|
|
c(Pd_atom) = 148.205E0_DP |
587 |
|
|
c(Ag_atom) = 96.524E0_DP |
588 |
|
|
c(Ir_atom) = 224.815E0_DP |
589 |
|
|
c(Pt_atom) = 71.336E0_DP |
590 |
|
|
c(Au_atom) = 53.581E0_DP |
591 |
|
|
|
592 |
|
|
m(Ni_atom,Ni_atom) = 5.0E0_DP |
593 |
|
|
m(Cu_atom,Cu_atom) = 5.0E0_DP |
594 |
|
|
m(Rh_atom,Rh_atom) = 5.0E0_DP |
595 |
|
|
m(Pd_atom,Pd_atom) = 6.0E0_DP |
596 |
|
|
m(Ag_atom,Ag_atom) = 6.0E0_DP |
597 |
|
|
m(Ir_atom,Ir_atom) = 6.0E0_DP |
598 |
|
|
m(Pt_atom,Pt_atom) = 7.0E0_DP |
599 |
|
|
m(Au_atom,Au_atom) = 8.0E0_DP |
600 |
|
|
|
601 |
|
|
n(Ni_atom,Ni_atom) = 10.0E0_DP |
602 |
|
|
n(Cu_atom,Cu_atom) = 10.0E0_DP |
603 |
|
|
n(Rh_atom,Rh_atom) = 13.0E0_DP |
604 |
|
|
n(Pd_atom,Pd_atom) = 12.0E0_DP |
605 |
|
|
n(Ag_atom,Ag_atom) = 11.0E0_DP |
606 |
|
|
n(Ir_atom,Ir_atom) = 13.0E0_DP |
607 |
|
|
n(Pt_atom,Pt_atom) = 11.0E0_DP |
608 |
|
|
n(Au_atom,Au_atom) = 11.0E0_DP |
609 |
|
|
|
610 |
|
|
alpha(Ni_atom,Ni_atom) = 3.5157E0_DP |
611 |
|
|
alpha(Cu_atom,Cu_atom) = 3.6030E0_DP |
612 |
|
|
alpha(Rh_atom,Rh_atom) = 3.7984E0_DP |
613 |
|
|
alpha(Pd_atom,Pd_atom) = 3.8813E0_DP |
614 |
|
|
alpha(Ag_atom,Ag_atom) = 4.0691E0_DP |
615 |
|
|
alpha(Ir_atom,Ir_atom) = 3.8344E0_DP |
616 |
|
|
alpha(Pt_atom,Pt_atom) = 3.9163E0_DP |
617 |
|
|
alpha(Au_atom,Au_atom) = 4.0651E0_DP |
618 |
|
|
|
619 |
|
|
BigD(Ni_atom,Ni_atom) = 7.3767E0_DP*0.02306054E0_DP |
620 |
|
|
BigD(Cu_atom,Cu_atom) = 5.7921E0_DP*0.02306054E0_DP |
621 |
|
|
BigD(Rh_atom,Rh_atom) = 2.4612E0_DP*0.02306054E0_DP |
622 |
|
|
BigD(Pd_atom,Pd_atom) = 3.2864E0_DP*0.02306054E0_DP |
623 |
|
|
BigD(Ag_atom,Ag_atom) = 3.9450E0_DP*0.02306054E0_DP |
624 |
|
|
BigD(Ir_atom,Ir_atom) = 3.7674E0_DP*0.02306054E0_DP |
625 |
|
|
BigD(Pt_atom,Pt_atom) = 9.7894E0_DP*0.02306054E0_DP |
626 |
|
|
BigD(Au_atom,Au_atom) = 7.8052E0_DP*0.02306054E0_DP |
627 |
|
|
|
628 |
|
|
! then the secondaries |
629 |
|
|
|
630 |
|
|
do i = 1, n_gatypes-1 |
631 |
|
|
atype1 = gatype(i) |
632 |
|
|
do j = i+1, n_gatypes |
633 |
|
|
atype2 = gatype(j) |
634 |
|
|
|
635 |
|
|
BigD(atype1, atype2) = dsqrt(BigD(atype1,atype1)* & |
636 |
|
|
BigD(atype2,atype2)) |
637 |
|
|
|
638 |
|
|
BigD(atype2, atype1) = BigD(atype1, atype2) |
639 |
|
|
|
640 |
|
|
m(atype1,atype2) = 0.5E0_DP*(m(atype1,atype1)+m(atype2,atype2)) |
641 |
|
|
m(atype2,atype1) = m(atype1,atype2) |
642 |
|
|
|
643 |
|
|
n(atype1,atype2) = 0.5E0_DP*(n(atype1,atype1)+n(atype2,atype2)) |
644 |
|
|
n(atype2,atype1) = n(atype1,atype2) |
645 |
|
|
|
646 |
|
|
alpha(atype1,atype2) = 0.5E0_DP*(alpha(atype1,atype1) + & |
647 |
|
|
alpha(atype2,atype2)) |
648 |
|
|
alpha(atype2,atype1) = alpha(atype1,atype2) |
649 |
|
|
|
650 |
|
|
enddo |
651 |
|
|
enddo |
652 |
|
|
|
653 |
|
|
do i = 1, n_gatypes |
654 |
|
|
atype1 = gatype(i) |
655 |
|
|
do j = 1, n_gatypes |
656 |
|
|
atype2 = gatype(j) |
657 |
|
|
|
658 |
|
|
rcutg(atype1,atype2) = 2.0E0_DP*alpha(atype1,atype2) |
659 |
|
|
|
660 |
|
|
vpair_rcut(atype1,atype2) = BigD(atype1,atype2)* & |
661 |
|
|
(alpha(atype1,atype2)/rcutg(atype1,atype2))**n(atype1,atype2) |
662 |
|
|
rho_rcut(atype1, atype2) = & |
663 |
|
|
(alpha(atype1,atype2)/rcutg(atype1,atype2))**m(atype1,atype2) |
664 |
|
|
enddo |
665 |
|
|
enddo |
666 |
|
|
|
667 |
|
|
return |
668 |
|
|
end subroutine initialize_goddard |
669 |
|
|
|
670 |
|
|
|
671 |
|
|
|
672 |
|
|
|
673 |
|
|
|
674 |
|
|
|
675 |
|
|
subroutine mass_weight() |
676 |
|
|
integer ia, ja, dim, dim2, idim, idim2 |
677 |
|
|
real( kind = DP ) :: mt, m1, m2, wt |
678 |
|
|
|
679 |
|
|
|
680 |
|
|
do ia = 1, natoms |
681 |
|
|
m1 = mass(ia) |
682 |
|
|
do ja = 1, natoms |
683 |
|
|
m2 = mass(ja) |
684 |
|
|
wt = 1.0E0_DP/dsqrt(m1*m2) |
685 |
|
|
do dim = 1, 3 |
686 |
|
|
idim = 3 * (ia-1) + dim |
687 |
|
|
do dim2 = 1, 3 |
688 |
|
|
idim2 = 3 * (ja-1) + dim2 |
689 |
|
|
d(idim,idim2) = d(idim,idim2)*wt |
690 |
|
|
enddo |
691 |
|
|
enddo |
692 |
|
|
enddo |
693 |
|
|
enddo |
694 |
|
|
|
695 |
|
|
end subroutine mass_weight |
696 |
|
|
|
697 |
|
|
end module goddard_module |