ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/nano_mpi/src/force_glue_module.F90
Revision: 4
Committed: Mon Jun 10 17:18:36 2002 UTC (22 years ago) by chuckv
File size: 28604 byte(s)
Log Message:
Import Root

File Contents

# Content
1 module glue_module
2 use definitions, ONLY : ndim,DP
3 use parameter
4 use simulation
5 use second_deriv
6 use force_utilities, ONLY : check,wrap,save_nlist
7 use status, ONLY: error,info,warning
8 #ifdef MPI
9 use mpi_module
10 #endif
11
12
13
14 integer, parameter :: n_glatypes = 2
15
16 integer , allocatable, dimension(:) :: glatype
17 real( kind = DP ), allocatable, dimension(:,:) :: dpars ! density parameters
18 real( kind = DP ), allocatable, dimension(:,:) :: gpars ! glue parameters
19 real( kind = DP ), allocatable, dimension(:,:) :: ppars ! pair parameters
20
21 !! arrays for force calculations
22
23
24
25 private :: n_glatypes,glatype,dpars,gpars,ppars
26 private :: allocate_glue_module
27 private :: get_ppars, get_gpars
28 public :: get_dpars
29 private :: uu,v2
30 public :: deallocate_glue_module,calc_glue_dens,calc_glue_forces
31 public :: initialize_glue
32 private :: mass_weight
33
34
35 contains
36
37
38 subroutine allocate_glue_module(n_size)
39 integer, intent(in) :: n_size
40 allocate(glatype(n_size))
41 allocate(dpars(n_size,13))
42 allocate(gpars(n_size,13))
43 allocate(ppars(n_size,20))
44 end subroutine allocate_glue_module
45
46 subroutine deallocate_glue_module
47 deallocate(glatype)
48 deallocate(dpars)
49 deallocate(gpars)
50 deallocate(ppars)
51 end subroutine deallocate_glue_module
52
53 subroutine calc_glue_dens(update_nlist)
54
55 ! include 'headers/sizes.h'
56
57 real( kind = DP ) :: ptmp
58 integer :: i, j, atype1, atype2, nlist, jbeg, jend, jnab
59 integer :: j_start
60 integer :: tag_i,tag_j
61 real( kind = DP ) :: rxi, ryi, rzi, rxij, ryij, rzij, rijsq
62 real( kind = DP ), dimension(13) :: g_dpars
63 real( kind = DP ) :: r, drho, d2rho
64 logical, intent(inout) :: update_nlist
65
66
67 #ifndef MPI
68 integer :: nrow
69 integer :: ncol
70
71 nrow = natoms - 1
72 ncol = natoms
73 #endif
74
75 #ifdef MPI
76 rho_row = 0.0E0_DP
77 rho_col = 0.0E0_DP
78 j_start = 1
79
80 call gather(q,q_row,plan_row3)
81 call gather(q,q_col,plan_col3)
82 #endif
83
84 rho = 0.0E0_DP
85
86
87 if (update_nlist) then
88
89 ! save current configuration, contruct neighbor list,
90 ! and calculate forces
91 call save_nlist()
92
93 nlist = 0
94
95 do i = 1, nrow
96 point(i) = nlist + 1
97
98 #ifdef MPI
99 tag_i = tag_row(i)
100 atype1 = ident_row(i)
101 rxi = q_row(1,i)
102 ryi = q_row(2,i)
103 rzi = q_row(3,i)
104 #else
105 j_start = i + 1
106
107 rxi = q(1,i)
108 ryi = q(2,i)
109 rzi = q(3,i)
110 #endif
111 inner: do j = j_start, ncol
112
113
114
115
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 rxij = wrap(rxi - q_col(1,j), 1)
127 ryij = wrap(ryi - q_col(2,j), 2)
128 rzij = wrap(rzi - q_col(3,j), 3)
129 #else
130 rxij = wrap(rxi - q(1,j), 1)
131 ryij = wrap(ryi - q(2,j), 2)
132 rzij = wrap(rzi - q(3,j), 3)
133 #endif
134 rijsq = rxij*rxij + ryij*ryij + rzij*rzij
135
136 #ifdef MPI
137 if (rijsq <= rlstsq .AND. &
138 tag_j /= tag_i) then
139 #else
140 if (rijsq < rlstsq) then
141 #endif
142 nlist = nlist + 1
143 list(nlist) = j
144
145
146 if (rijsq < rcutsq) then
147
148 r = dsqrt(rijsq)
149 #ifdef MPI
150 atype1 = ident_row(i)
151 #else
152 atype1 = ident(i)
153 #endif
154 call get_dpars(atype1, g_dpars)
155 call rh(r, ptmp, drho, d2rho, g_dpars)
156
157 ! density at site j depends on type of atom at site i
158 #ifdef MPI
159 rho_col(j) = rho_col(j) + ptmp
160 #else
161 rho(j) = rho(j) + ptmp
162 #endif
163
164 #ifdef MPI
165 atype2 = ident_col(j)
166 #else
167 atype2 = ident(j)
168 #endif
169 call get_dpars(atype2, g_dpars)
170 call rh(r, ptmp, drho, d2rho, g_dpars)
171
172 ! density at site i depends on type of atom at site j
173 #ifdef MPI
174 rho_row(i) = rho_row(i) + ptmp
175 #else
176 rho(i) = rho(i) + ptmp
177 #endif
178 endif
179 endif
180 enddo inner
181 enddo
182
183 #ifdef MPI
184 point(nrow + 1) = nlist + 1
185 #else
186 point(natoms) = nlist + 1
187 #endif
188
189 else
190 ! use the list to find the neighbors
191 do i = 1, nrow
192 JBEG = POINT(I)
193 JEND = POINT(I+1) - 1
194
195 ! check thiat molecule i has neighbors
196 if (jbeg <= jend) then
197
198 #ifdef MPI
199 rxi = q_row(1,i)
200 ryi = q_row(2,i)
201 rzi = q_row(3,i)
202 #else
203 rxi = q(1,i)
204 ryi = q(2,i)
205 rzi = q(3,i)
206 #endif
207
208 do jnab = jbeg, jend
209 j = list(jnab)
210
211 #ifdef MPI
212 rxij = wrap(rxi - q_col(1,j), 1)
213 ryij = wrap(ryi - q_col(2,j), 2)
214 rzij = wrap(rzi - q_col(3,j), 3)
215 #else
216 rxij = wrap(rxi - q(1,j), 1)
217 ryij = wrap(ryi - q(2,j), 2)
218 rzij = wrap(rzi - q(3,j), 3)
219 #endif
220 rijsq = rxij*rxij + ryij*ryij + rzij*rzij
221
222 if (rijsq < rcutsq) then
223
224 r = dsqrt(rijsq)
225 #ifdef MPI
226 atype1 = ident_row(i)
227 atype2 = ident_col(j)
228 #else
229 atype1 = ident(i)
230 atype2 = ident(j)
231 #endif
232
233 !! get density for each atom site
234 call get_dpars(atype1, g_dpars)
235 call rh(r, ptmp, drho, d2rho, g_dpars)
236 #ifdef MPI
237 rho_col(j) = rho_col(j) + ptmp
238 #else
239 ! density at site j depends on type of atom at site i
240 rho(j) = rho(j) + ptmp
241 #endif
242
243 call get_dpars(atype2, g_dpars)
244 call rh(r, ptmp, drho, d2rho, g_dpars)
245 #ifdef MPI
246 rho_row(i) = rho_row(i) + ptmp
247 #else
248 ! density at site i depends on type of atom at site j
249 rho(i) = rho(i) + ptmp
250 #endif
251 endif
252 enddo
253 endif
254 enddo
255 endif
256
257 #ifdef MPI
258 !! communicate densities
259 call scatter(rho_row,rho,plan_row)
260 if (newtons_thrd) then
261 call scatter(rho_col,rho_tmp,plan_col)
262 do i = 1, nlocal
263 rho(i) = rho(i) + rho_tmp(i)
264 end do
265 endif
266 #endif
267
268 return
269 end subroutine calc_glue_dens
270
271 subroutine calc_glue_forces(nmflag,pot)
272
273
274 #ifdef MPI
275 real( kind = DP ), dimension(nlocal) :: frho
276 real( kind = DP ), dimension(nrow) :: frho_row
277 real( kind = DP ), dimension(ncol) :: frho_col
278
279 real( kind = DP ), dimension(nlocal) :: dfrhodrho
280 real( kind = DP ), dimension(nrow) :: dfrhodrho_row
281 real( kind = DP ), dimension(ncol) :: dfrhodrho_col
282
283 real( kind = DP ), dimension(nlocal) :: d2frhodrhodrho
284 real( kind = DP ), dimension(nrow) :: d2frhodrhodrho_row
285 real( kind = DP ), dimension(ncol) :: d2frhodrhodrho_col
286 real( kind = DP ), dimension(3,ncol) :: efr
287
288 real( kind = DP ) :: pot_local, pot_phi_row, pot_Frho, pot_phi, pot_row
289
290 #else
291 real( kind = DP ), dimension(natoms) :: frho
292 real( kind = DP ), dimension(natoms) :: dfrhodrho
293 real( kind = DP ), dimension(natoms) :: d2frhodrhodrho
294 real( kind = DP ), dimension(3,natoms) :: efr
295 #endif
296
297 real( kind = DP ),intent(out), optional :: pot
298 real( kind = DP ) :: vptmp, dudr, ftmp
299
300 real( kind = DP ) :: g_gpars(13), g_dpars(13), g_ppars(20)
301 real( kind = DP ) :: u, u1, u2, phab, rci, rcj
302 real( kind = DP ) :: rha, drha, d2rha, pha, dpha, d2pha
303 real( kind = DP ) :: rhb, drhb, d2rhb, phb, dphb, d2phb
304
305
306
307
308 real( kind = DP ) :: drhoidr, drhojdr, d2rhoidrdr, d2rhojdrdr
309 real( kind = DP ) :: dvpdr, drdx1, d2vpdrdr, d2
310 real( kind = DP ) :: kt1, kt2, kt3, ktmp
311
312 integer :: i, j, dim, atype1, atype2, idim, jdim, dim2, idim2, jdim2
313 integer :: jbeg, jend, jnab
314 real( kind = DP ) :: rxij, ryij, rzij, rxi, ryi, rzi, rijsq, r
315
316
317 integer :: nlist
318
319 logical, intent(in) :: nmflag
320 logical :: do_pot
321
322 #ifndef MPI
323 integer :: nrow
324 integer :: ncol
325
326 nrow = natoms - 1
327 ncol = natoms
328 #endif
329
330
331 do_pot = .false.
332 if (present(pot)) do_pot = .true.
333
334 #ifndef MPI
335 if (do_pot) pot = 0.0E0_DP
336 f = 0.0E0_DP
337 e = 0.0E0_DP
338 #else
339 f_row = 0.0E0_DP
340 f_col = 0.0E0_DP
341
342 pot_phi_row = 0.0E0_DP
343 pot_phi = 0.0E0_DP
344 pot_Frho = 0.0E0_DP
345 pot_local = 0.0E0_DP
346 pot_row = 0.0E0_DP
347
348 e_row = 0.0E0_DP
349 e_col = 0.0E0_DP
350 e_tmp = 0.0E0_DP
351 #endif
352
353 ! get functional for electron density
354 ! MPI we calculate this locally then
355 do i = 1, nlocal
356 atype1 = ident(i)
357
358 call get_gpars(atype1, g_gpars)
359 call uu(rho(i), u, u1, u2, g_gpars)
360
361 frho(i) = u
362 dfrhodrho(i) = u1
363 d2frhodrhodrho(i) = u2
364 #ifndef MPI
365 if (do_pot) pot = pot + u
366 #endif
367 enddo
368
369 #ifdef MPI
370 !! communicate f(rho) and derivatives
371
372 call gather(frho,frho_row,plan_row)
373 call gather(dfrhodrho,dfrhodrho_row,plan_row)
374 call gather(frho,frho_col,plan_col)
375 call gather(dfrhodrho,dfrhodrho_col,plan_col)
376
377 if (nmflag) then
378 call gather(d2frhodrhodrho,d2frhodrhodrho_row,plan_row)
379 call gather(d2frhodrhodrho,d2frhodrhodrho_col,plan_col)
380 endif
381 #endif
382
383 do i = 1, nrow
384 JBEG = POINT(i)
385 JEND = POINT(i+1) - 1
386 ! check thiat molecule i has neighbors
387 if (jbeg .le. jend) then
388 #ifdef MPI
389 atype1 = ident_row(i)
390 rxi = q_row(1,i)
391 ryi = q_row(2,i)
392 rzi = q_row(3,i)
393 #else
394 atype1 = ident(i)
395 rxi = q(1,i)
396 ryi = q(2,i)
397 rzi = q(3,i)
398 #endif
399 do jnab = jbeg, jend
400 j = list(jnab)
401
402 #ifdef MPI
403 rxij = wrap(rxi - q_col(1,j), 1)
404 ryij = wrap(ryi - q_col(2,j), 2)
405 rzij = wrap(rzi - q_col(3,j), 3)
406 #else
407 rxij = wrap(rxi - q(1,j), 1)
408 ryij = wrap(ryi - q(2,j), 2)
409 rzij = wrap(rzi - q(3,j), 3)
410 #endif
411 rijsq = rxij*rxij + ryij*ryij + rzij*rzij
412
413 if (rijsq .lt. rcutsq) then
414
415
416 r = dsqrt(rijsq)
417
418 efr(1,j) = -rxij
419 efr(2,j) = -ryij
420 efr(3,j) = -rzij
421
422 #ifdef MPI
423 atype1 = ident_row(i)
424 #else
425 atype1 = ident(i)
426 #endif
427
428 call get_dpars(atype1, g_dpars)
429 rci = g_dpars(3)
430 call rh(r, rha, drha, d2rha, g_dpars)
431 call get_ppars(atype1, g_ppars)
432 call v2(r, pha, dpha, d2pha, g_ppars)
433
434
435 #ifdef MPI
436 atype2 = ident_col(j)
437 #else
438 atype2 = ident(j)
439 #endif
440
441 call get_dpars(atype2, g_dpars)
442 rcj = g_dpars(3)
443 call rh(r, rhb, drhb, d2rhb, g_dpars)
444 call get_ppars(atype2, g_ppars)
445 call v2(r, phb, dphb, d2phb, g_ppars)
446
447
448 phab = 0.0E0_DP
449 dvpdr = 0.0E0_DP
450 d2vpdrdr = 0.0E0_DP
451
452 if (r.lt.rci) then
453 phab = phab + 0.5E0_DP*(rhb/rha)*pha
454 dvpdr = dvpdr + 0.5E0_DP*((rhb/rha)*dpha + &
455 pha*((drhb/rha) - (rhb*drha/rha/rha)))
456 d2vpdrdr = d2vpdrdr + 0.5E0_DP*((rhb/rha)*d2pha + &
457 2.0E0_DP*dpha*((drhb/rha) - (rhb*drha/rha/rha)) + &
458 pha*((d2rhb/rha) - 2.0E0_DP*(drhb*drha/rha/rha) + &
459 (2.0E0_DP*rhb*drha*drha/rha/rha/rha) - (rhb*d2rha/rha/rha)))
460 endif
461
462 if (r.lt.rcj) then
463 phab = phab + 0.5E0_DP*(rha/rhb)*phb
464 dvpdr = dvpdr + 0.5E0_DP*((rha/rhb)*dphb + &
465 phb*((drha/rhb) - (rha*drhb/rhb/rhb)))
466 d2vpdrdr = d2vpdrdr + 0.5E0_DP*((rha/rhb)*d2phb + &
467 2.0E0_DP*dphb*((drha/rhb) - (rha*drhb/rhb/rhb)) + &
468 phb*((d2rha/rhb) - 2.0E0_DP*(drha*drhb/rhb/rhb) + &
469 (2.0E0_DP*rha*drhb*drhb/rhb/rhb/rhb) - (rha*d2rhb/rhb/rhb)))
470 endif
471
472 !! add to the total potential energy
473 #ifdef MPI
474 e_row(i) = e_row(i) + phab*0.5
475 e_col(i) = e_col(i) + phab*0.5
476 #else
477 if (do_pot) pot = pot + phab
478 #endif
479
480 drhoidr = drha
481 drhojdr = drhb
482
483 d2rhoidrdr = d2rha
484 d2rhojdrdr = d2rhb
485 #ifdef MPI
486 dudr = drhojdr*dfrhodrho_row(i)+drhoidr*dfrhodrho_col(j) &
487 + dvpdr
488
489 if (nmflag) then
490 d2 = d2vpdrdr + &
491 d2rhoidrdr*dfrhodrho_col(j) + &
492 d2rhojdrdr*dfrhodrho_row(i) + &
493 drhoidr*drhoidr*d2frhodrhodrho_col(j) + &
494 drhojdr*drhojdr*d2frhodrhodrho_row(i)
495 endif
496 #else
497
498 dudr = drhojdr*dfrhodrho(i)+drhoidr*dfrhodrho(j) &
499 + dvpdr
500
501 if (nmflag) then
502 d2 = d2vpdrdr + &
503 d2rhoidrdr*dfrhodrho(j) + &
504 d2rhojdrdr*dfrhodrho(i) + &
505 drhoidr*drhoidr*d2frhodrhodrho(j) + &
506 drhojdr*drhojdr*d2frhodrhodrho(i)
507 endif
508 #endif
509
510 do dim = 1, 3
511
512 drdx1 = efr(dim,j) / r
513 ftmp = dudr * drdx1
514 #ifdef MPI
515 f_col(dim,j) = f_col(dim,j) - ftmp
516 f_row(dim,i) = f_row(dim,i) + ftmp
517 #else
518 f(dim,j) = f(dim,j) - ftmp
519 f(dim,i) = f(dim,i) + ftmp
520 #endif
521 if (nmflag) then
522 idim = 3 * (i-1) + dim
523 jdim = 3 * (j-1) + dim
524
525 do dim2 = 1, 3
526
527 kt1 = d2 * efr(dim,j) * efr(dim2,j)/r/r
528 kt2 = - dudr * efr(dim,j) * efr(dim2,j)/r/r/r
529
530 if (dim.eq.dim2) then
531 kt3 = dudr / r
532 else
533 kt3 = 0.0E0_DP
534 endif
535
536 ! The factor of 2 below is to compensate for
537 ! overcounting.
538 ! Mass weighting is done separately...
539
540 ktmp = (kt1+kt2+kt3)/2.0E0_DP
541 idim2 = 3 * (i-1) + dim2
542 jdim2 = 3 * (j-1) + dim2
543
544 d(idim, idim2) = d(idim,idim2) + ktmp
545 d(idim2, idim) = d(idim2,idim) + ktmp
546
547 d(idim, jdim2) = d(idim,jdim2) - ktmp
548 d(idim2, jdim) = d(idim2,jdim) - ktmp
549
550 d(jdim, idim2) = d(jdim,idim2) - ktmp
551 d(jdim2, idim) = d(jdim2,idim) - ktmp
552
553 d(jdim, jdim2) = d(jdim,jdim2) + ktmp
554 d(jdim2, jdim) = d(jdim2,jdim) + ktmp
555
556 enddo
557 endif
558 enddo
559
560 endif
561 enddo
562 endif
563 enddo
564
565
566
567 #ifdef MPI
568 !!distribute forces
569 call scatter(f_row,f,plan_row3)
570 if (newtons_thrd) then
571 call scatter(f_col,f_tmp,plan_col3)
572 do i = 1,nlocal
573 do dim = 1,3
574 f(dim,i) = f(dim,i) + f_tmp(dim,i)
575 end do
576 end do
577 endif
578
579
580 if (do_pot) then
581 ! scatter/gather pot_row into the members of my column
582 call scatter(e_row,e_tmp,plan_row)
583
584 ! scatter/gather pot_local into all other procs
585 ! add resultant to get total pot
586 do i = 1, nlocal
587 pot_local = pot_local + frho(i) + e_tmp(i)
588 enddo
589 if (newtons_thrd) then
590 e_tmp = 0.0E0_DP
591 call scatter(e_col,e_tmp,plan_col)
592 do i = 1, nlocal
593 pot_local = pot_local + e_tmp(i)
594 enddo
595 endif
596
597
598 call mpi_reduce(pot_local,pot,1,mpi_double_precision, &
599 mpi_sum,0,mpi_comm_world,mpi_err)
600 endif
601 #endif
602
603 if (nmflag) then
604 call mass_weight()
605 endif
606
607
608 return
609 end subroutine calc_glue_forces
610
611 subroutine initialize_glue()
612 use model_module
613 include 'headers/atom.h'
614 ! Order of the dpars array:
615 ! 1 , 2 , 3 , 4 , 5 , 6 , 7 , 8 , 9 , 10 , 11 , 12 , 13
616 ! RRD,RRB,RRC,RHOD,RHOA,R1I,R2I,R3I,R1II,R2II,R3II,R2III,R3III
617 !
618 ! Order of the gpars array:
619 ! 1 , 2 , 3 , 4 , 5 , 6 , 7 , 8 , 9 , 10 , 11 , 12 , 13
620 ! DB, UB,DSW,B0I,B1I,B2I,B3I,B4I,B2II,B3II,B4II,B2III,B3III
621 ! Order of the ppars array:
622 !
623 ! 1,2, 3, 4 , 5 , 6 , 7 , 8 , 9 , 10, 11 , 12 , 13 , 14 , 15 , 16 , 17 ,
624 ! D,A,RC,PHI1,PHI2,A0I,A1I,A2I,A3I,A4I,A0II,A1II,A2II,A3II,A4II,A5II,A6II,
625 !
626 ! 18 , 19 , 20
627 ! A3III,A4III,A5III
628 !
629 ! units are Angstroms for the distances, kcal/mol for the potentials
630 ! the r coefficients are unchanged since they just give the density.
631
632
633
634
635 call allocate_glue_module(n_glatypes)
636
637 glatype(1) = Au_atom
638 glatype(2) = Pb_atom
639
640 dpars(1,1) = 0.2878207442141723D+01
641 dpars(1,2) = 0.3500000000000000D+01
642 dpars(1,3) = 0.3900000000000000D+01
643 dpars(1,4) = 0.1000000000000000D+01
644 dpars(1,5) = 0.0000000000000000D+00
645 dpars(1,6) = -0.6800000000000000D+00
646 dpars(1,7) = 0.7500000000000000D+00
647 dpars(1,8) = -0.1333333333333333D+01
648 dpars(1,9) = -0.6800000000000000D+00
649 dpars(1,10) = 0.7500000000000000D+00
650 dpars(1,11) = -0.1527241171296038D+01
651 dpars(1,12) = 0.5578188675490974D+01
652 dpars(1,13) = 0.6132971688727435D+01
653
654 dpars(2,1) = 0.3471540742235355D+01
655 dpars(2,2) = 0.4909500000000000D+01
656 dpars(2,3) = 0.5503000000000000D+01
657 dpars(2,4) = 0.8500000000000000D+00
658 dpars(2,5) = 0.3000000000000000D+00
659 dpars(2,6) = -0.3800000000000000D+00
660 dpars(2,7) = 0.2450000000000000D+00
661 dpars(2,8) = -0.5166666666666667D+00
662 dpars(2,9) = -0.3800000000000000D+00
663 dpars(2,10) = 0.2450000000000000D+00
664 dpars(2,11) = -0.1715828759323021D+00
665 dpars(2,12) = 0.1308624193974091D+01
666 dpars(2,13) = 0.7699032959082642D+00
667
668 gpars(1,1) = 0.1200000000000000D+02
669 gpars(1,2) = -0.3300000000000000D+01*23.06054E0_DP
670 gpars(1,3) = 0.9358157767784574D+01
671 gpars(1,4) = -0.2793388616771698D+01*23.06054E0_DP
672 gpars(1,5) = -0.3419999999999999D+00*23.06054E0_DP
673 gpars(1,6) = 0.3902327808424106D-01*23.06054E0_DP
674 gpars(1,7) = 0.7558829951858879D-02*23.06054E0_DP
675 gpars(1,8) = 0.3090472511796849D-03*23.06054E0_DP
676 gpars(1,9) = 0.8618226772941980D-01*23.06054E0_DP
677 gpars(1,10) = 0.4341701445034724D-02*23.06054E0_DP
678 gpars(1,11) = -0.3044398779375916D-03*23.06054E0_DP
679 gpars(1,12) = 0.8618226772941980D-01*23.06054E0_DP
680 gpars(1,13) = 0.4325981467602070D-02*23.06054E0_DP
681
682 gpars(2,1) = 0.1200000000000000D+02
683 gpars(2,2) = -0.1850000000000000D+01*23.06054E0_DP
684 gpars(2,3) = 0.9081433382788202D+01
685 gpars(2,4) = -0.1536837858573856D+01*23.06054E0_DP
686 gpars(2,5) = -0.1850000000000000D+00*23.06054E0_DP
687 gpars(2,6) = -0.1515954156009047D-01*23.06054E0_DP
688 gpars(2,7) = -0.1478056600250295D-02*23.06054E0_DP
689 gpars(2,8) = 0.0000000000000000D+00*23.06054E0_DP
690 gpars(2,9) = 0.1526631307568407D-01*23.06054E0_DP
691 gpars(2,10) = -0.1820707358322264D-01*23.06054E0_DP
692 gpars(2,11) = -0.3714503267605675D-02*23.06054E0_DP
693 gpars(2,12) = 0.1526631307568407D-01*23.06054E0_DP
694 gpars(2,13) = 0.3239697893498678D-01*23.06054E0_DP
695
696 ppars(1,1) = 0.2878207442141723D+01
697 ppars(1,2) = 0.4070400000000000D+01
698 ppars(1,3) = 0.3700000000000000D+01
699 ppars(1,4) = -0.8000000000000000D-01
700 ppars(1,5) = 0.0000000000000000D+00
701 ppars(1,6) = -0.8000000000000000D-01*23.06054E0_DP
702 ppars(1,7) = 0.0000000000000000D+00*23.06054E0_DP
703 ppars(1,8) = 0.7619231375231362D+00*23.06054E0_DP
704 ppars(1,9) = -0.8333333333333333D+00*23.06054E0_DP
705 ppars(1,10) = -0.1211483464993159D+00*23.06054E0_DP
706 ppars(1,11) = -0.8000000000000000D-01*23.06054E0_DP
707 ppars(1,12) = 0.0000000000000000D+00*23.06054E0_DP
708 ppars(1,13) = 0.7619231375231362D+00*23.06054E0_DP
709 ppars(1,14) = -0.8333333333333333D+00*23.06054E0_DP
710 ppars(1,15) = -0.1096009851140349D+01*23.06054E0_DP
711 ppars(1,16) = 0.2158417178555998D+01*23.06054E0_DP
712 ppars(1,17) = -0.9128915709636862D+00*23.06054E0_DP
713 ppars(1,18) = 0.0000000000000000D+00*23.06054E0_DP
714 ppars(1,19) = 0.0000000000000000D+00*23.06054E0_DP
715 ppars(1,20) = 0.0000000000000000D+00*23.06054E0_DP
716
717 ppars(2,1) = 0.3471540742235355D+01
718 ppars(2,2) = 0.4909500000000000D+01
719 ppars(2,3) = 0.4230000000000000D+01
720 ppars(2,4) = -0.3000000000000000D-01
721 ppars(2,5) = 0.0000000000000000D+00
722 ppars(2,6) = -0.3000000000000000D-01*23.06054E0_DP
723 ppars(2,7) = 0.0000000000000000D+00*23.06054E0_DP
724 ppars(2,8) = 0.1102661976296813D+00*23.06054E0_DP
725 ppars(2,9) = -0.8166666666666667D+00*23.06054E0_DP
726 ppars(2,10) = 0.4072201316247714D-01*23.06054E0_DP
727 ppars(2,11) = -0.3000000000000000D-01*23.06054E0_DP
728 ppars(2,12) = 0.0000000000000000D+00*23.06054E0_DP
729 ppars(2,13) = 0.1102661976296813D+00*23.06054E0_DP
730 ppars(2,14) = -0.8166666666666667D+00*23.06054E0_DP
731 ppars(2,15) = 0.3439976422630956D+01*23.06054E0_DP
732 ppars(2,16) = -0.5105760527431719D+01*23.06054E0_DP
733 ppars(2,17) = 0.2448028237231130D+01*23.06054E0_DP
734 ppars(2,18) = 0.0000000000000000D+00*23.06054E0_DP
735 ppars(2,19) = 0.0000000000000000D+00*23.06054E0_DP
736 ppars(2,20) = 0.0000000000000000D+00*23.06054E0_DP
737
738 return
739 end subroutine initialize_glue
740
741 subroutine get_dpars(atype, atmp)
742
743
744
745
746 real( kind = DP ), dimension(13) :: atmp(13)
747 integer :: atype, i, j
748
749 do i = 1, n_glatypes
750 if (atype.eq.glatype(i)) then
751 do j = 1, 13
752 atmp(j) = dpars(i,j)
753 enddo
754 return
755 endif
756 enddo
757 call error('GET_DPARS','Unknown atom type for force field')
758 end subroutine get_dpars
759
760 subroutine get_gpars(atype, atmp)
761
762
763 double precision atmp(13)
764 integer atype, i, j
765
766 do i = 1, n_glatypes
767 if (atype.eq.glatype(i)) then
768 do j = 1, 13
769 atmp(j) = gpars(i,j)
770 enddo
771 return
772 endif
773 enddo
774
775
776 call error('GET_GPARS','Unknown atom type for force field')
777
778 end subroutine get_gpars
779
780 subroutine get_ppars(atype, atmp)
781
782
783
784 double precision atmp(20)
785 integer atype, i, j
786
787 do i = 1, n_glatypes
788 if (atype.eq.glatype(i)) then
789 do j = 1, 20
790 atmp(j) = ppars(i,j)
791 enddo
792 return
793 endif
794 enddo
795
796 call error('GET_GPARS','Unknown atom type for force field')
797
798 end subroutine get_ppars
799
800 subroutine rh(r, rho, drho, d2rho, g_dpars)
801
802 ! returns the density function rho(r) and its first two derivatives
803 ! at distance r.
804 !
805 ! Order of the dpars array:
806 ! 1 , 2 , 3 , 4 , 5 , 6 , 7 , 8 , 9 , 10 , 11 , 12 , 13
807 ! RRD,RRB,RRC,RHOD,RHOA,R1I,R2I,R3I,R1II,R2II,R3II,R2III,R3III
808
809
810 double precision r, rho, drho, d2rho,g_dpars(13)
811 double precision x
812
813 if (r.ge.g_dpars(3)) then
814 ! after cut radius it is all zero :
815 rho = 0.E0_DP
816 drho = 0.E0_DP
817 d2rho = 0.E0_DP
818 else if (r.ge.g_dpars(2)) then
819 ! region iii (release) with a spline :
820 x = r - g_dpars(3)
821 rho = (x**2) * (g_dpars(12) + x*g_dpars(13))
822 drho = x * (2.E0_DP*g_dpars(12) + x*3.E0_DP*g_dpars(13))
823 d2rho = 2.E0_DP*g_dpars(12) + x*6.E0_DP*g_dpars(13)
824 else if (r.ge.g_dpars(1)) then
825 ! region ii (sustain) with a spline :
826 x = r - g_dpars(1)
827 rho = g_dpars(4) + x*(g_dpars(9) + x*(g_dpars(10) + x*g_dpars(11)))
828 drho = g_dpars(9) + x*(2.E0_DP*g_dpars(10) + x*3.E0_DP*g_dpars(11))
829 d2rho = 2.E0_DP*g_dpars(10) + x*6.E0_DP*g_dpars(11)
830 else
831 ! region i (decay) with a spline :
832 x = r - g_dpars(1)
833 rho = g_dpars(4) + x*(g_dpars(6) + x*(g_dpars(7) + x*g_dpars(8)))
834 drho = g_dpars(6) + x*(2.E0_DP*g_dpars(7) + x*3.E0_DP*g_dpars(8))
835 d2rho = 2.E0_DP*g_dpars(7) + x*6.E0_DP*g_dpars(8)
836 endif
837
838 return
839 end subroutine rh
840
841 subroutine uu(dens, u, u1, u2, g_gpars)
842
843 ! returns the function u(n) and its two first derivatives at n=dens
844 !
845 ! Order of the gpars array:
846 ! 1 , 2 , 3 , 4 , 5 , 6 , 7 , 8 , 9 , 10 , 11 , 12 , 13
847 ! DB, UB,DSW,B0I,B1I,B2I,B3I,B4I,B2II,B3II,B4II,B2III,B3III
848
849 double precision dens, u, u1, u2, g_gpars(13)
850 double precision x
851
852 if (dens.gt.g_gpars(1)) then
853 ! region iii
854 x = dens - g_gpars(1)
855 u = g_gpars(2) + (x**2) * (g_gpars(12) + x*g_gpars(13))
856 u1 = x*(2.E0_DP*g_gpars(12) + x*3.E0_DP*g_gpars(13))
857 u2 = 2.E0_DP*g_gpars(12) + x*6.E0_DP*g_gpars(13)
858 else if (dens.gt.g_gpars(3)) then
859 ! region ii
860 x = dens - g_gpars(1)
861 u = g_gpars(2) + (x**2) * (g_gpars(9) + &
862 x*(g_gpars(10) + x*g_gpars(11)))
863 u1 = x*(2.E0_DP*g_gpars(9) + x*(3.E0_DP*g_gpars(10) + x*4.E0_DP*g_gpars(11)))
864 u2 = 2.E0_DP*g_gpars(9) + x*(6.E0_DP*g_gpars(10) + x*12.E0_DP*g_gpars(11))
865 else
866 ! region i
867 x = dens - g_gpars(3)
868 u = g_gpars(4) + x*(g_gpars(5) + x*(g_gpars(6) + &
869 x*(g_gpars(7) + x*g_gpars(8))))
870 u1 = g_gpars(5) + x*(2.E0_DP*g_gpars(6) + x*(3.E0_DP*g_gpars(7) &
871 + x*4.E0_DP*g_gpars(8)))
872 u2 = 2.E0_DP*g_gpars(6) + x*(6.E0_DP*g_gpars(7) + x*12.E0_DP*g_gpars(8))
873 endif
874 return
875 end subroutine uu
876
877 subroutine v2(r, phi, dphi, d2phi, g_ppars)
878
879 ! returns the potential and its first two derivatives at distance r
880 !
881 ! Order of the ppars array:
882 !
883 ! 1,2, 3, 4 , 5 , 6 , 7 , 8 , 9 , 10, 11 , 12 , 13 , 14 , 15 , 16 , 17 ,
884 ! D,A,RC,PHI1,PHI2,A0I,A1I,A2I,A3I,A4I,A0II,A1II,A2II,A3II,A4II,A5II,A6II,
885 !
886 ! 18 , 19 , 20
887 ! A3III,A4III,A5III
888
889 double precision r, phi, dphi, d2phi, g_ppars(20)
890 double precision x
891
892 if (r.ge.g_ppars(3)) then
893 phi = 0.E0_DP
894 dphi = 0.E0_DP
895 d2phi = 0.E0_DP
896 else if (r.ge.g_ppars(2)) then
897 ! region iii, after second neighbours.
898 ! this works only if rc.gt.a, i.e. when the potential is
899 ! a true second neighbours potential ; otherwise control
900 ! passes to region ii.
901 x = r - g_ppars(3)
902 phi = (x**3) * (g_ppars(20)*x**2 + g_ppars(19)*x + g_ppars(18))
903 dphi = (x**2) * (5.E0_DP*g_ppars(20)*x**2 + 4.E0_DP*g_ppars(19)*x + &
904 3.E0_DP*g_ppars(18))
905 d2phi = x * (20.E0_DP*g_ppars(20)*x**2 + 12.E0_DP*g_ppars(19)*x &
906 + 6.E0_DP*g_ppars(18))
907 else if (r.ge.g_ppars(1)) then
908 ! region ii, between first and second neighbours.
909 x = r - g_ppars(1)
910 phi = g_ppars(11) + x*(g_ppars(12) + x*(g_ppars(13) + &
911 x*(g_ppars(14) + x*(g_ppars(15) + x*(g_ppars(16) &
912 + x*g_ppars(17))))))
913 dphi = g_ppars(12) + x*(2.E0_DP*g_ppars(13) + x*(3.E0_DP*g_ppars(14) + &
914 x*(4.E0_DP*g_ppars(15) + x*(5.E0_DP*g_ppars(16) + x*6.E0_DP*g_ppars(17)))))
915 d2phi = 2.E0_DP*g_ppars(13) + x*(6.E0_DP*g_ppars(14) + x*(12.E0_DP*g_ppars(15) + &
916 x*(20.E0_DP*g_ppars(16) + x*30.E0_DP*g_ppars(17))))
917 else
918 ! region i, before first neighbours.
919 x = r - g_ppars(1)
920 phi = g_ppars(6) + x*(g_ppars(7) + x*(g_ppars(8) + &
921 x*(g_ppars(9) + x*g_ppars(10))))
922 dphi = g_ppars(7) + x*(2.E0_DP*g_ppars(8) + x*(3.E0_DP*g_ppars(9) &
923 + x*4.E0_DP*g_ppars(10)))
924 d2phi = 2.E0_DP*(g_ppars(8) + x*(3.E0_DP*g_ppars(9) + &
925 x*6.E0_DP*g_ppars(10)))
926 endif
927 return
928 end subroutine v2
929
930 subroutine mass_weight()
931 integer ia, ja, dim, dim2, idim, idim2
932 real( kind = DP ) :: mt, m1, m2, wt
933
934
935 do ia = 1, natoms
936 m1 = mass(ia)
937 do ja = 1, natoms
938 m2 = mass(ja)
939 wt = 1.0E0_DP/dsqrt(m1*m2)
940 do dim = 1, 3
941 idim = 3 * (ia-1) + dim
942 do dim2 = 1, 3
943 idim2 = 3 * (ja-1) + dim2
944 d(idim,idim2) = d(idim,idim2)*wt
945 enddo
946 enddo
947 enddo
948 enddo
949
950 end subroutine mass_weight
951
952
953
954
955
956
957 end module glue_module