(git:495eafe)
Loading...
Searching...
No Matches
manybody_potential.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \par History
10!> Efficient tersoff implementation and general "lifting" of manybody_potential module
11!> 12.2007 [tlaino] - Splitting manybody module : In this module we should only
12!> keep the main routines for computing energy and forces of
13!> manybody potentials. Each potential should have his own module!
14!> \author CJM, I-Feng W. Kuo, Teodoro Laino
15! **************************************************************************************************
17
19 USE cell_types, ONLY: cell_type
28 USE kinds, ONLY: dp
55 USE pair_potential_types, ONLY: &
60 USE util, ONLY: sort
61#include "./base/base_uses.f90"
62
63 IMPLICIT NONE
64
65 PRIVATE
66 PUBLIC :: energy_manybody
68 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'manybody_potential'
69
70CONTAINS
71
72! **************************************************************************************************
73!> \brief computes the embedding contribution to the energy
74!> \param fist_nonbond_env ...
75!> \param atomic_kind_set ...
76!> \param local_particles ...
77!> \param particle_set ...
78!> \param cell ...
79!> \param pot_manybody ...
80!> \param para_env ...
81!> \param mm_section ...
82!> \param use_virial ...
83!> \par History
84!> tlaino [2007] - New algorithm for tersoff potential
85!> \author CJM, I-Feng W. Kuo, Teodoro Laino
86! **************************************************************************************************
87 SUBROUTINE energy_manybody(fist_nonbond_env, atomic_kind_set, local_particles, &
88 particle_set, cell, pot_manybody, para_env, mm_section, use_virial)
89
90 TYPE(fist_nonbond_env_type), POINTER :: fist_nonbond_env
91 TYPE(atomic_kind_type), POINTER :: atomic_kind_set(:)
92 TYPE(distribution_1d_type), POINTER :: local_particles
93 TYPE(particle_type), POINTER :: particle_set(:)
94 TYPE(cell_type), POINTER :: cell
95 REAL(dp), INTENT(INOUT) :: pot_manybody
96 TYPE(mp_para_env_type), OPTIONAL, POINTER :: para_env
97 TYPE(section_vals_type), POINTER :: mm_section
98 LOGICAL, INTENT(IN) :: use_virial
99
100 CHARACTER(LEN=*), PARAMETER :: routinen = 'energy_manybody'
101
102 INTEGER :: atom_a, atom_b, handle, i, iend, ifirst, igrp, ikind, ilast, ilist, indexa, &
103 ipair, iparticle, iparticle_local, istart, iunique, jkind, junique, mpair, nkinds, &
104 nloc_size, npairs, nparticle, nparticle_local, nr_h3o, nr_o, nr_oh, nunique
105 INTEGER, DIMENSION(:), POINTER :: glob_loc_list_a, work_list
106 INTEGER, DIMENSION(:, :), POINTER :: glob_loc_list, list, sort_list
107 LOGICAL :: any_ace, any_allegro, any_deepmd, &
108 any_gal, any_gal21, any_nequip, &
109 any_siepmann, any_tersoff
110 REAL(kind=dp) :: drij, embed, pot_ace, pot_allegro, &
111 pot_deepmd, pot_loc, pot_nequip, qr, &
112 rab2_max, rij(3)
113 REAL(kind=dp), DIMENSION(3) :: cell_v, cvi
114 REAL(kind=dp), DIMENSION(:, :), POINTER :: glob_cell_v
115 REAL(kind=dp), POINTER :: fembed(:)
116 TYPE(eam_pot_type), POINTER :: eam
117 TYPE(eam_type), DIMENSION(:), POINTER :: eam_data
118 TYPE(fist_neighbor_type), POINTER :: nonbonded
119 TYPE(gal21_pot_type), POINTER :: gal21
120 TYPE(gal_pot_type), POINTER :: gal
121 TYPE(neighbor_kind_pairs_type), POINTER :: neighbor_kind_pair
122 TYPE(pair_potential_pp_type), POINTER :: potparm
123 TYPE(pair_potential_single_type), POINTER :: pot
124 TYPE(pos_type), DIMENSION(:), POINTER :: r_last_update_pbc
125 TYPE(siepmann_pot_type), POINTER :: siepmann
126 TYPE(tersoff_pot_type), POINTER :: tersoff
127
128 NULLIFY (eam, siepmann, tersoff, gal, gal21)
129 any_tersoff = .false.
130 any_siepmann = .false.
131 any_gal = .false.
132 any_gal21 = .false.
133 any_allegro = .false.
134 any_nequip = .false.
135 any_ace = .false.
136 any_deepmd = .false.
137 CALL timeset(routinen, handle)
138 CALL fist_nonbond_env_get(fist_nonbond_env, r_last_update_pbc=r_last_update_pbc, &
139 potparm=potparm, eam_data=eam_data)
140 ! EAM requires a single loop
141 DO ikind = 1, SIZE(atomic_kind_set)
142 pot => potparm%pot(ikind, ikind)%pot
143 DO i = 1, SIZE(pot%type)
144 IF (pot%type(i) /= ea_type) cycle
145 eam => pot%set(i)%eam
146 nparticle = SIZE(particle_set)
147 ALLOCATE (fembed(nparticle))
148 fembed(:) = 0._dp
149 cpassert(ASSOCIATED(eam_data))
150 ! computation of embedding function and energy
151 nparticle_local = local_particles%n_el(ikind)
152 DO iparticle_local = 1, nparticle_local
153 iparticle = local_particles%list(ikind)%array(iparticle_local)
154 indexa = int(eam_data(iparticle)%rho/eam%drhoar) + 1
155 IF (indexa > eam%npoints - 1) indexa = eam%npoints - 1
156 qr = eam_data(iparticle)%rho - eam%rhoval(indexa)
157
158 embed = eam%frho(indexa) + qr*eam%frhop(indexa)
159 fembed(iparticle) = eam%frhop(indexa) + qr*(eam%frhop(indexa + 1) - eam%frhop(indexa))/eam%drhoar
160
161 pot_manybody = pot_manybody + embed
162 END DO
163 ! communicate data
164 CALL para_env%sum(fembed)
165 DO iparticle = 1, nparticle
166 IF (particle_set(iparticle)%atomic_kind%kind_number == ikind) THEN
167 eam_data(iparticle)%f_embed = fembed(iparticle)
168 END IF
169 END DO
170
171 DEALLOCATE (fembed)
172 END DO
173 END DO
174 ! Other manybody potential
175 DO ikind = 1, SIZE(atomic_kind_set)
176 DO jkind = ikind, SIZE(atomic_kind_set)
177 pot => potparm%pot(ikind, jkind)%pot
178 any_tersoff = any_tersoff .OR. any(pot%type == tersoff_type)
179 any_nequip = any_nequip .OR. any(pot%type == nequip_type)
180 any_ace = any_ace .OR. any(pot%type == ace_type)
181 any_allegro = any_allegro .OR. any(pot%type == allegro_type)
182 any_deepmd = any_deepmd .OR. any(pot%type == deepmd_type)
183 any_siepmann = any_siepmann .OR. any(pot%type == siepmann_type)
184 any_gal = any_gal .OR. any(pot%type == gal_type)
185 any_gal21 = any_gal21 .OR. any(pot%type == gal21_type)
186 END DO
187 END DO
188 CALL fist_nonbond_env_get(fist_nonbond_env, nonbonded=nonbonded, natom_types=nkinds)
189 ! NEQUIP
190 IF (any_nequip) THEN
191 NULLIFY (glob_loc_list, glob_cell_v, glob_loc_list_a)
192 CALL nequip_energy_store_force_virial(nonbonded, particle_set, local_particles, cell, &
193 atomic_kind_set, potparm, r_last_update_pbc, &
194 pot_nequip, fist_nonbond_env, &
195 para_env, use_virial, nequip_type)
196 pot_manybody = pot_manybody + pot_nequip
197 END IF
198 ! ALLEGRO
199 IF (any_allegro) THEN
200 NULLIFY (glob_loc_list, glob_cell_v, glob_loc_list_a)
201 CALL nequip_energy_store_force_virial(nonbonded, particle_set, local_particles, cell, &
202 atomic_kind_set, potparm, r_last_update_pbc, &
203 pot_allegro, fist_nonbond_env, &
204 para_env, use_virial, allegro_type)
205 pot_manybody = pot_manybody + pot_allegro
206 END IF
207 ! ACE
208 IF (any_ace) THEN
209 CALL ace_energy_store_force_virial(particle_set, cell, atomic_kind_set, potparm, &
210 fist_nonbond_env, pot_ace)
211 pot_manybody = pot_manybody + pot_ace
212 END IF
213 ! DEEPMD
214 IF (any_deepmd) THEN
215 CALL deepmd_energy_store_force_virial(particle_set, cell, atomic_kind_set, potparm, &
216 fist_nonbond_env, pot_deepmd, para_env)
217 pot_manybody = pot_manybody + pot_deepmd
218 END IF
219
220 ! TERSOFF
221 IF (any_tersoff) THEN
222 NULLIFY (glob_loc_list, glob_cell_v, glob_loc_list_a)
223 CALL setup_tersoff_arrays(nonbonded, potparm, glob_loc_list, glob_cell_v, glob_loc_list_a, cell)
224 DO ilist = 1, nonbonded%nlists
225 neighbor_kind_pair => nonbonded%neighbor_kind_pairs(ilist)
226 npairs = neighbor_kind_pair%npairs
227 IF (npairs == 0) cycle
228 kind_group_loop: DO igrp = 1, neighbor_kind_pair%ngrp_kind
229 istart = neighbor_kind_pair%grp_kind_start(igrp)
230 iend = neighbor_kind_pair%grp_kind_end(igrp)
231 ikind = neighbor_kind_pair%ij_kind(1, igrp)
232 jkind = neighbor_kind_pair%ij_kind(2, igrp)
233 list => neighbor_kind_pair%list
234 cvi = neighbor_kind_pair%cell_vector
235 pot => potparm%pot(ikind, jkind)%pot
236 DO i = 1, SIZE(pot%type)
237 IF (pot%type(i) /= tersoff_type) cycle
238 rab2_max = pot%set(i)%tersoff%rcutsq
239 cell_v = matmul(cell%hmat, cvi)
240 pot => potparm%pot(ikind, jkind)%pot
241 tersoff => pot%set(i)%tersoff
242 npairs = iend - istart + 1
243 IF (npairs /= 0) THEN
244 ALLOCATE (sort_list(2, npairs), work_list(npairs))
245 sort_list = list(:, istart:iend)
246 ! Sort the list of neighbors, this increases the efficiency for single
247 ! potential contributions
248 CALL sort(sort_list(1, :), npairs, work_list)
249 DO ipair = 1, npairs
250 work_list(ipair) = sort_list(2, work_list(ipair))
251 END DO
252 sort_list(2, :) = work_list
253 ! find number of unique elements of array index 1
254 nunique = 1
255 DO ipair = 1, npairs - 1
256 IF (sort_list(1, ipair + 1) /= sort_list(1, ipair)) nunique = nunique + 1
257 END DO
258 ipair = 1
259 junique = sort_list(1, ipair)
260 ifirst = 1
261 DO iunique = 1, nunique
262 atom_a = junique
263 IF (glob_loc_list_a(ifirst) > atom_a) cycle
264 DO mpair = ifirst, SIZE(glob_loc_list_a)
265 IF (glob_loc_list_a(mpair) == atom_a) EXIT
266 END DO
267 ifirst = mpair
268 DO mpair = ifirst, SIZE(glob_loc_list_a)
269 IF (glob_loc_list_a(mpair) /= atom_a) EXIT
270 END DO
271 ilast = mpair - 1
272 nloc_size = 0
273 IF (ifirst /= 0) nloc_size = ilast - ifirst + 1
274 DO WHILE (ipair <= npairs)
275 IF (sort_list(1, ipair) /= junique) EXIT
276 atom_b = sort_list(2, ipair)
277 ! Energy terms
278 pot_loc = 0.0_dp
279 rij(:) = r_last_update_pbc(atom_b)%r(:) - r_last_update_pbc(atom_a)%r(:) + cell_v
280 drij = dot_product(rij, rij)
281 ipair = ipair + 1
282 IF (drij > rab2_max) cycle
283 drij = sqrt(drij)
284 CALL tersoff_energy(pot_loc, tersoff, r_last_update_pbc, atom_a, atom_b, nloc_size, &
285 glob_loc_list(:, ifirst:ilast), glob_cell_v(:, ifirst:ilast), cell_v, drij)
286 pot_manybody = pot_manybody + 0.5_dp*pot_loc
287 END DO
288 ifirst = ilast + 1
289 IF (ipair <= npairs) junique = sort_list(1, ipair)
290 END DO
291 DEALLOCATE (sort_list, work_list)
292 END IF
293 END DO
294 END DO kind_group_loop
295 END DO
296 CALL destroy_tersoff_arrays(glob_loc_list, glob_cell_v, glob_loc_list_a)
297 END IF
298
299 !SIEPMANN POTENTIAL
300 IF (any_siepmann) THEN
301 NULLIFY (glob_loc_list, glob_cell_v, glob_loc_list_a)
302 nr_oh = 0
303 nr_h3o = 0
304 nr_o = 0
305 CALL setup_siepmann_arrays(nonbonded, potparm, glob_loc_list, glob_cell_v, glob_loc_list_a, cell)
306 DO ilist = 1, nonbonded%nlists
307 neighbor_kind_pair => nonbonded%neighbor_kind_pairs(ilist)
308 npairs = neighbor_kind_pair%npairs
309 IF (npairs == 0) cycle
310 kind_group_loop_2: DO igrp = 1, neighbor_kind_pair%ngrp_kind
311 istart = neighbor_kind_pair%grp_kind_start(igrp)
312 iend = neighbor_kind_pair%grp_kind_end(igrp)
313 ikind = neighbor_kind_pair%ij_kind(1, igrp)
314 jkind = neighbor_kind_pair%ij_kind(2, igrp)
315 list => neighbor_kind_pair%list
316 cvi = neighbor_kind_pair%cell_vector
317 pot => potparm%pot(ikind, jkind)%pot
318 DO i = 1, SIZE(pot%type)
319 IF (pot%type(i) /= siepmann_type) cycle
320 rab2_max = pot%set(i)%siepmann%rcutsq
321 cell_v = matmul(cell%hmat, cvi)
322 pot => potparm%pot(ikind, jkind)%pot
323 siepmann => pot%set(i)%siepmann
324 npairs = iend - istart + 1
325 IF (npairs /= 0) THEN
326 ALLOCATE (sort_list(2, npairs), work_list(npairs))
327 sort_list = list(:, istart:iend)
328 ! Sort the list of neighbors, this increases the efficiency for single
329 ! potential contributions
330 CALL sort(sort_list(1, :), npairs, work_list)
331 DO ipair = 1, npairs
332 work_list(ipair) = sort_list(2, work_list(ipair))
333 END DO
334 sort_list(2, :) = work_list
335 ! find number of unique elements of array index 1
336 nunique = 1
337 DO ipair = 1, npairs - 1
338 IF (sort_list(1, ipair + 1) /= sort_list(1, ipair)) nunique = nunique + 1
339 END DO
340 ipair = 1
341 junique = sort_list(1, ipair)
342 ifirst = 1
343 DO iunique = 1, nunique
344 atom_a = junique
345 IF (glob_loc_list_a(ifirst) > atom_a) cycle
346 DO mpair = ifirst, SIZE(glob_loc_list_a)
347 IF (glob_loc_list_a(mpair) == atom_a) EXIT
348 END DO
349 ifirst = mpair
350 DO mpair = ifirst, SIZE(glob_loc_list_a)
351 IF (glob_loc_list_a(mpair) /= atom_a) EXIT
352 END DO
353 ilast = mpair - 1
354 nloc_size = 0
355 IF (ifirst /= 0) nloc_size = ilast - ifirst + 1
356 DO WHILE (ipair <= npairs)
357 IF (sort_list(1, ipair) /= junique) EXIT
358 atom_b = sort_list(2, ipair)
359 ! Energy terms
360 pot_loc = 0.0_dp
361 rij(:) = r_last_update_pbc(atom_b)%r(:) - r_last_update_pbc(atom_a)%r(:) + cell_v
362 drij = dot_product(rij, rij)
363 ipair = ipair + 1
364 IF (drij > rab2_max) cycle
365 drij = sqrt(drij)
366 CALL siepmann_energy(pot_loc, siepmann, r_last_update_pbc, atom_a, atom_b, nloc_size, &
367 glob_loc_list(:, ifirst:ilast), cell_v, cell, drij, &
368 particle_set, nr_oh, nr_h3o, nr_o)
369 pot_manybody = pot_manybody + pot_loc
370 END DO
371 ifirst = ilast + 1
372 IF (ipair <= npairs) junique = sort_list(1, ipair)
373 END DO
374 DEALLOCATE (sort_list, work_list)
375 END IF
376 END DO
377 END DO kind_group_loop_2
378 END DO
379 CALL destroy_siepmann_arrays(glob_loc_list, glob_cell_v, glob_loc_list_a)
380 CALL print_nr_ions_siepmann(nr_oh, mm_section, para_env, print_oh=.true., &
381 print_h3o=.false., print_o=.false.)
382 CALL print_nr_ions_siepmann(nr_h3o, mm_section, para_env, print_oh=.false., &
383 print_h3o=.true., print_o=.false.)
384 CALL print_nr_ions_siepmann(nr_o, mm_section, para_env, print_oh=.false., &
385 print_h3o=.false., print_o=.true.)
386 END IF
387
388 !GAL19 POTENTIAL
389 IF (any_gal) THEN
390 NULLIFY (glob_loc_list, glob_cell_v, glob_loc_list_a)
391 CALL setup_gal_arrays(nonbonded, potparm, glob_loc_list, glob_cell_v, glob_loc_list_a, cell)
392 DO ilist = 1, nonbonded%nlists
393 neighbor_kind_pair => nonbonded%neighbor_kind_pairs(ilist)
394 npairs = neighbor_kind_pair%npairs
395 IF (npairs == 0) cycle
396 kind_group_loop_3: DO igrp = 1, neighbor_kind_pair%ngrp_kind
397 istart = neighbor_kind_pair%grp_kind_start(igrp)
398 iend = neighbor_kind_pair%grp_kind_end(igrp)
399 ikind = neighbor_kind_pair%ij_kind(1, igrp)
400 jkind = neighbor_kind_pair%ij_kind(2, igrp)
401 list => neighbor_kind_pair%list
402 cvi = neighbor_kind_pair%cell_vector
403 pot => potparm%pot(ikind, jkind)%pot
404 DO i = 1, SIZE(pot%type)
405 IF (pot%type(i) /= gal_type) cycle
406 rab2_max = pot%set(i)%gal%rcutsq
407 cell_v = matmul(cell%hmat, cvi)
408 pot => potparm%pot(ikind, jkind)%pot
409 gal => pot%set(i)%gal
410 npairs = iend - istart + 1
411 IF (npairs /= 0) THEN
412 ALLOCATE (sort_list(2, npairs), work_list(npairs))
413 sort_list = list(:, istart:iend)
414 ! Sort the list of neighbors, this increases the efficiency for single
415 ! potential contributions
416 CALL sort(sort_list(1, :), npairs, work_list)
417 DO ipair = 1, npairs
418 work_list(ipair) = sort_list(2, work_list(ipair))
419 END DO
420 sort_list(2, :) = work_list
421 ! find number of unique elements of array index 1
422 nunique = 1
423 DO ipair = 1, npairs - 1
424 IF (sort_list(1, ipair + 1) /= sort_list(1, ipair)) nunique = nunique + 1
425 END DO
426 ipair = 1
427 junique = sort_list(1, ipair)
428 ifirst = 1
429 DO iunique = 1, nunique
430 atom_a = junique
431 IF (glob_loc_list_a(ifirst) > atom_a) cycle
432 DO mpair = ifirst, SIZE(glob_loc_list_a)
433 IF (glob_loc_list_a(mpair) == atom_a) EXIT
434 END DO
435 ifirst = mpair
436 DO mpair = ifirst, SIZE(glob_loc_list_a)
437 IF (glob_loc_list_a(mpair) /= atom_a) EXIT
438 END DO
439 ilast = mpair - 1
440 nloc_size = 0
441 IF (ifirst /= 0) nloc_size = ilast - ifirst + 1
442 DO WHILE (ipair <= npairs)
443 IF (sort_list(1, ipair) /= junique) EXIT
444 atom_b = sort_list(2, ipair)
445 ! Energy terms
446 pot_loc = 0.0_dp
447 rij(:) = r_last_update_pbc(atom_b)%r(:) - r_last_update_pbc(atom_a)%r(:) + cell_v
448 drij = dot_product(rij, rij)
449 ipair = ipair + 1
450 IF (drij > rab2_max) cycle
451 drij = sqrt(drij)
452 CALL gal_energy(pot_loc, gal, r_last_update_pbc, atom_a, atom_b, &
453 cell, particle_set, mm_section)
454
455 pot_manybody = pot_manybody + pot_loc
456 END DO
457 ifirst = ilast + 1
458 IF (ipair <= npairs) junique = sort_list(1, ipair)
459 END DO
460 DEALLOCATE (sort_list, work_list)
461 END IF
462 END DO
463 END DO kind_group_loop_3
464 END DO
465 CALL destroy_gal_arrays(glob_loc_list, glob_cell_v, glob_loc_list_a)
466 END IF
467
468 !GAL21 POTENTIAL
469 IF (any_gal21) THEN
470 NULLIFY (glob_loc_list, glob_cell_v, glob_loc_list_a)
471 CALL setup_gal21_arrays(nonbonded, potparm, glob_loc_list, glob_cell_v, glob_loc_list_a, cell)
472 DO ilist = 1, nonbonded%nlists
473 neighbor_kind_pair => nonbonded%neighbor_kind_pairs(ilist)
474 npairs = neighbor_kind_pair%npairs
475 IF (npairs == 0) cycle
476 kind_group_loop_5: DO igrp = 1, neighbor_kind_pair%ngrp_kind
477 istart = neighbor_kind_pair%grp_kind_start(igrp)
478 iend = neighbor_kind_pair%grp_kind_end(igrp)
479 ikind = neighbor_kind_pair%ij_kind(1, igrp)
480 jkind = neighbor_kind_pair%ij_kind(2, igrp)
481 list => neighbor_kind_pair%list
482 cvi = neighbor_kind_pair%cell_vector
483 pot => potparm%pot(ikind, jkind)%pot
484 DO i = 1, SIZE(pot%type)
485 IF (pot%type(i) /= gal21_type) cycle
486 rab2_max = pot%set(i)%gal21%rcutsq
487 cell_v = matmul(cell%hmat, cvi)
488 pot => potparm%pot(ikind, jkind)%pot
489 gal21 => pot%set(i)%gal21
490 npairs = iend - istart + 1
491 IF (npairs /= 0) THEN
492 ALLOCATE (sort_list(2, npairs), work_list(npairs))
493 sort_list = list(:, istart:iend)
494 ! Sort the list of neighbors, this increases the efficiency for single
495 ! potential contributions
496 CALL sort(sort_list(1, :), npairs, work_list)
497 DO ipair = 1, npairs
498 work_list(ipair) = sort_list(2, work_list(ipair))
499 END DO
500 sort_list(2, :) = work_list
501 ! find number of unique elements of array index 1
502 nunique = 1
503 DO ipair = 1, npairs - 1
504 IF (sort_list(1, ipair + 1) /= sort_list(1, ipair)) nunique = nunique + 1
505 END DO
506 ipair = 1
507 junique = sort_list(1, ipair)
508 ifirst = 1
509 DO iunique = 1, nunique
510 atom_a = junique
511 IF (glob_loc_list_a(ifirst) > atom_a) cycle
512 DO mpair = ifirst, SIZE(glob_loc_list_a)
513 IF (glob_loc_list_a(mpair) == atom_a) EXIT
514 END DO
515 ifirst = mpair
516 DO mpair = ifirst, SIZE(glob_loc_list_a)
517 IF (glob_loc_list_a(mpair) /= atom_a) EXIT
518 END DO
519 ilast = mpair - 1
520 nloc_size = 0
521 IF (ifirst /= 0) nloc_size = ilast - ifirst + 1
522 DO WHILE (ipair <= npairs)
523 IF (sort_list(1, ipair) /= junique) EXIT
524 atom_b = sort_list(2, ipair)
525 ! Energy terms
526 pot_loc = 0.0_dp
527 rij(:) = r_last_update_pbc(atom_b)%r(:) - r_last_update_pbc(atom_a)%r(:) + cell_v
528 drij = dot_product(rij, rij)
529 ipair = ipair + 1
530 IF (drij > rab2_max) cycle
531 drij = sqrt(drij)
532 CALL gal21_energy(pot_loc, gal21, r_last_update_pbc, atom_a, atom_b, &
533 cell, particle_set, mm_section)
534
535 pot_manybody = pot_manybody + pot_loc
536 END DO
537 ifirst = ilast + 1
538 IF (ipair <= npairs) junique = sort_list(1, ipair)
539 END DO
540 DEALLOCATE (sort_list, work_list)
541 END IF
542 END DO
543 END DO kind_group_loop_5
544 END DO
545 CALL destroy_gal21_arrays(glob_loc_list, glob_cell_v, glob_loc_list_a)
546 END IF
547
548 CALL timestop(handle)
549 END SUBROUTINE energy_manybody
550
551! **************************************************************************************************
552!> \brief ...
553!> \param fist_nonbond_env ...
554!> \param particle_set ...
555!> \param cell ...
556!> \param f_nonbond ...
557!> \param pv_nonbond ...
558!> \param use_virial ...
559!> \par History
560!> Fast implementation of the tersoff potential - [tlaino] 2007
561!> \author I-Feng W. Kuo, Teodoro Laino
562! **************************************************************************************************
563 SUBROUTINE force_nonbond_manybody(fist_nonbond_env, particle_set, cell, &
564 f_nonbond, pv_nonbond, use_virial)
565
566 TYPE(fist_nonbond_env_type), POINTER :: fist_nonbond_env
567 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
568 TYPE(cell_type), POINTER :: cell
569 REAL(kind=dp), DIMENSION(:, :), INTENT(INOUT) :: f_nonbond, pv_nonbond
570 LOGICAL, INTENT(IN) :: use_virial
571
572 CHARACTER(LEN=*), PARAMETER :: routinen = 'force_nonbond_manybody'
573
574 INTEGER :: atom_a, atom_b, handle, i, i_a, i_b, iend, ifirst, igrp, ikind, ilast, ilist, &
575 ipair, istart, iunique, jkind, junique, kind_a, kind_b, mpair, nkinds, nloc_size, npairs, &
576 nunique
577 INTEGER, ALLOCATABLE, DIMENSION(:, :) :: eam_kinds_index
578 INTEGER, DIMENSION(:), POINTER :: glob_loc_list_a, work_list
579 INTEGER, DIMENSION(:, :), POINTER :: glob_loc_list, list, sort_list
580 LOGICAL :: any_ace, any_allegro, any_deepmd, &
581 any_gal, any_gal21, any_nequip, &
582 any_siepmann, any_tersoff
583 REAL(kind=dp) :: f_eam, fac, fr(3), ptens11, ptens12, ptens13, ptens21, ptens22, ptens23, &
584 ptens31, ptens32, ptens33, rab(3), rab2, rab2_max, rtmp(3)
585 REAL(kind=dp), DIMENSION(3) :: cell_v, cvi
586 REAL(kind=dp), DIMENSION(:, :), POINTER :: glob_cell_v
587 TYPE(eam_pot_type), POINTER :: eam_a, eam_b
588 TYPE(eam_type), DIMENSION(:), POINTER :: eam_data
589 TYPE(fist_neighbor_type), POINTER :: nonbonded
590 TYPE(gal21_pot_type), POINTER :: gal21
591 TYPE(gal_pot_type), POINTER :: gal
592 TYPE(neighbor_kind_pairs_type), POINTER :: neighbor_kind_pair
593 TYPE(pair_potential_pp_type), POINTER :: potparm
594 TYPE(pair_potential_single_type), POINTER :: pot
595 TYPE(pos_type), DIMENSION(:), POINTER :: r_last_update_pbc
596 TYPE(siepmann_pot_type), POINTER :: siepmann
597 TYPE(tersoff_pot_type), POINTER :: tersoff
598
599 any_tersoff = .false.
600 any_allegro = .false.
601 any_nequip = .false.
602 any_siepmann = .false.
603 any_ace = .false.
604 any_deepmd = .false.
605 any_gal = .false.
606 any_gal21 = .false.
607 CALL timeset(routinen, handle)
608 NULLIFY (eam_a, eam_b, tersoff, siepmann, gal, gal21)
609
610 CALL fist_nonbond_env_get(fist_nonbond_env, nonbonded=nonbonded, potparm=potparm, &
611 natom_types=nkinds, eam_data=eam_data, r_last_update_pbc=r_last_update_pbc)
612
613 ! Initializing the potential energy, pressure tensor and force
614 IF (use_virial) THEN
615 ptens11 = 0.0_dp; ptens12 = 0.0_dp; ptens13 = 0.0_dp
616 ptens21 = 0.0_dp; ptens22 = 0.0_dp; ptens23 = 0.0_dp
617 ptens31 = 0.0_dp; ptens32 = 0.0_dp; ptens33 = 0.0_dp
618 END IF
619
620 nkinds = SIZE(potparm%pot, 1)
621 ALLOCATE (eam_kinds_index(nkinds, nkinds))
622 eam_kinds_index = -1
623 DO ikind = 1, nkinds
624 DO jkind = ikind, nkinds
625 DO i = 1, SIZE(potparm%pot(ikind, jkind)%pot%type)
626 IF (potparm%pot(ikind, jkind)%pot%type(i) == ea_type) THEN
627 ! At the moment we allow only 1 EAM per each kinds pair..
628 cpassert(eam_kinds_index(ikind, jkind) == -1)
629 cpassert(eam_kinds_index(jkind, ikind) == -1)
630 eam_kinds_index(ikind, jkind) = i
631 eam_kinds_index(jkind, ikind) = i
632 END IF
633 END DO
634 END DO
635 END DO
636 DO ikind = 1, nkinds
637 DO jkind = ikind, nkinds
638 any_ace = any_ace .OR. any(potparm%pot(ikind, jkind)%pot%type == ace_type)
639 END DO
640 END DO
641 ! ACE
642 IF (any_ace) &
643 CALL ace_add_force_virial(fist_nonbond_env, f_nonbond, pv_nonbond, use_virial)
644
645 DO ikind = 1, nkinds
646 DO jkind = ikind, nkinds
647 any_deepmd = any_deepmd .OR. any(potparm%pot(ikind, jkind)%pot%type == deepmd_type)
648 END DO
649 END DO
650 ! DEEPMD
651 IF (any_deepmd) &
652 CALL deepmd_add_force_virial(fist_nonbond_env, f_nonbond, pv_nonbond, use_virial)
653
654 ! NEQUIP
655 DO ikind = 1, nkinds
656 DO jkind = ikind, nkinds
657 any_nequip = any_nequip .OR. any(potparm%pot(ikind, jkind)%pot%type == nequip_type)
658 END DO
659 END DO
660 IF (any_nequip) THEN
661 CALL nequip_add_force_virial(fist_nonbond_env, f_nonbond, pv_nonbond, use_virial)
662 END IF
663
664 ! ALLEGRO
665 DO ikind = 1, nkinds
666 DO jkind = ikind, nkinds
667 any_allegro = any_allegro .OR. any(potparm%pot(ikind, jkind)%pot%type == allegro_type)
668 END DO
669 END DO
670 IF (any_allegro) THEN
671 CALL nequip_add_force_virial(fist_nonbond_env, f_nonbond, pv_nonbond, use_virial)
672 END IF
673
674 ! starting the force loop
675 DO ilist = 1, nonbonded%nlists
676 neighbor_kind_pair => nonbonded%neighbor_kind_pairs(ilist)
677 npairs = neighbor_kind_pair%npairs
678 IF (npairs == 0) cycle
679 kind_group_loop1: DO igrp = 1, neighbor_kind_pair%ngrp_kind
680 istart = neighbor_kind_pair%grp_kind_start(igrp)
681 iend = neighbor_kind_pair%grp_kind_end(igrp)
682 ikind = neighbor_kind_pair%ij_kind(1, igrp)
683 jkind = neighbor_kind_pair%ij_kind(2, igrp)
684 list => neighbor_kind_pair%list
685 cvi = neighbor_kind_pair%cell_vector
686 pot => potparm%pot(ikind, jkind)%pot
687 IF (pot%no_mb) cycle kind_group_loop1
688 rab2_max = pot%rcutsq
689 cell_v = matmul(cell%hmat, cvi)
690 any_tersoff = any_tersoff .OR. any(pot%type == tersoff_type)
691 any_siepmann = any_siepmann .OR. any(pot%type == siepmann_type)
692 any_ace = any_ace .OR. any(pot%type == ace_type)
693 any_deepmd = any_deepmd .OR. any(pot%type == deepmd_type)
694 any_gal = any_gal .OR. any(pot%type == gal_type)
695 any_gal21 = any_gal21 .OR. any(pot%type == gal21_type)
696 any_nequip = any_nequip .OR. any(pot%type == nequip_type)
697 any_allegro = any_allegro .OR. any(pot%type == allegro_type)
698 i = eam_kinds_index(ikind, jkind)
699 IF (i == -1) cycle kind_group_loop1
700 ! EAM
701 cpassert(ASSOCIATED(eam_data))
702 DO ipair = istart, iend
703 atom_a = list(1, ipair)
704 atom_b = list(2, ipair)
705 fac = 1.0_dp
706 IF (atom_a == atom_b) fac = 0.5_dp
707 kind_a = particle_set(atom_a)%atomic_kind%kind_number
708 kind_b = particle_set(atom_b)%atomic_kind%kind_number
709 i_a = eam_kinds_index(kind_a, kind_a)
710 i_b = eam_kinds_index(kind_b, kind_b)
711 eam_a => potparm%pot(kind_a, kind_a)%pot%set(i_a)%eam
712 eam_b => potparm%pot(kind_b, kind_b)%pot%set(i_b)%eam
713
714 !set this outside the potential type in case need multiple potentials
715 !Do everything necessary for EAM here
716 rab = r_last_update_pbc(atom_b)%r - r_last_update_pbc(atom_a)%r
717 rab = rab + cell_v
718 rab2 = rab(1)*rab(1) + rab(2)*rab(2) + rab(3)*rab(3)
719 IF (rab2 <= rab2_max) THEN
720 CALL get_force_eam(rab2, eam_a, eam_b, eam_data, atom_a, atom_b, f_eam)
721 f_eam = f_eam*fac
722
723 fr(1) = -f_eam*rab(1)
724 fr(2) = -f_eam*rab(2)
725 fr(3) = -f_eam*rab(3)
726 f_nonbond(1, atom_a) = f_nonbond(1, atom_a) - fr(1)
727 f_nonbond(2, atom_a) = f_nonbond(2, atom_a) - fr(2)
728 f_nonbond(3, atom_a) = f_nonbond(3, atom_a) - fr(3)
729
730 f_nonbond(1, atom_b) = f_nonbond(1, atom_b) + fr(1)
731 f_nonbond(2, atom_b) = f_nonbond(2, atom_b) + fr(2)
732 f_nonbond(3, atom_b) = f_nonbond(3, atom_b) + fr(3)
733 IF (use_virial) THEN
734 ptens11 = ptens11 + rab(1)*fr(1)
735 ptens21 = ptens21 + rab(2)*fr(1)
736 ptens31 = ptens31 + rab(3)*fr(1)
737 ptens12 = ptens12 + rab(1)*fr(2)
738 ptens22 = ptens22 + rab(2)*fr(2)
739 ptens32 = ptens32 + rab(3)*fr(2)
740 ptens13 = ptens13 + rab(1)*fr(3)
741 ptens23 = ptens23 + rab(2)*fr(3)
742 ptens33 = ptens33 + rab(3)*fr(3)
743 END IF
744 END IF
745 END DO
746 END DO kind_group_loop1
747 END DO
748 DEALLOCATE (eam_kinds_index)
749
750 ! Special way of handling the tersoff potential..
751 IF (any_tersoff) THEN
752 NULLIFY (glob_loc_list, glob_cell_v, glob_loc_list_a)
753 CALL setup_tersoff_arrays(nonbonded, potparm, glob_loc_list, glob_cell_v, glob_loc_list_a, cell)
754 DO ilist = 1, nonbonded%nlists
755 neighbor_kind_pair => nonbonded%neighbor_kind_pairs(ilist)
756 npairs = neighbor_kind_pair%npairs
757 IF (npairs == 0) cycle
758 kind_group_loop2: DO igrp = 1, neighbor_kind_pair%ngrp_kind
759 istart = neighbor_kind_pair%grp_kind_start(igrp)
760 iend = neighbor_kind_pair%grp_kind_end(igrp)
761 ikind = neighbor_kind_pair%ij_kind(1, igrp)
762 jkind = neighbor_kind_pair%ij_kind(2, igrp)
763 list => neighbor_kind_pair%list
764 cvi = neighbor_kind_pair%cell_vector
765 pot => potparm%pot(ikind, jkind)%pot
766
767 IF (pot%no_mb) cycle kind_group_loop2
768 rab2_max = pot%rcutsq
769 cell_v = matmul(cell%hmat, cvi)
770 DO i = 1, SIZE(pot%type)
771 ! TERSOFF
772 IF (pot%type(i) == tersoff_type) THEN
773 npairs = iend - istart + 1
774 tersoff => pot%set(i)%tersoff
775 ALLOCATE (sort_list(2, npairs), work_list(npairs))
776 sort_list = list(:, istart:iend)
777 ! Sort the list of neighbors, this increases the efficiency for single
778 ! potential contributions
779 CALL sort(sort_list(1, :), npairs, work_list)
780 DO ipair = 1, npairs
781 work_list(ipair) = sort_list(2, work_list(ipair))
782 END DO
783 sort_list(2, :) = work_list
784 ! find number of unique elements of array index 1
785 nunique = 1
786 DO ipair = 1, npairs - 1
787 IF (sort_list(1, ipair + 1) /= sort_list(1, ipair)) nunique = nunique + 1
788 END DO
789 ipair = 1
790 junique = sort_list(1, ipair)
791 ifirst = 1
792 DO iunique = 1, nunique
793 atom_a = junique
794 IF (glob_loc_list_a(ifirst) > atom_a) cycle
795 DO mpair = ifirst, SIZE(glob_loc_list_a)
796 IF (glob_loc_list_a(mpair) == atom_a) EXIT
797 END DO
798 ifirst = mpair
799 DO mpair = ifirst, SIZE(glob_loc_list_a)
800 IF (glob_loc_list_a(mpair) /= atom_a) EXIT
801 END DO
802 ilast = mpair - 1
803 nloc_size = 0
804 IF (ifirst /= 0) nloc_size = ilast - ifirst + 1
805 DO WHILE (ipair <= npairs)
806 IF (sort_list(1, ipair) /= junique) EXIT
807 atom_b = sort_list(2, ipair)
808 ! Derivative terms
809 rtmp = r_last_update_pbc(atom_b)%r(:) - r_last_update_pbc(atom_a)%r(:) + cell_v
810 ipair = ipair + 1
811 IF (dot_product(rtmp, rtmp) <= tersoff%rcutsq) THEN
812 CALL tersoff_forces(tersoff, r_last_update_pbc, cell_v, &
813 nloc_size, glob_loc_list(:, ifirst:ilast), glob_cell_v(:, ifirst:ilast), &
814 atom_a, atom_b, f_nonbond, pv_nonbond, use_virial, tersoff%rcutsq)
815 END IF
816 END DO
817 ifirst = ilast + 1
818 IF (ipair <= npairs) junique = sort_list(1, ipair)
819 END DO
820 DEALLOCATE (sort_list, work_list)
821 END IF
822 END DO
823 END DO kind_group_loop2
824 END DO
825 CALL destroy_tersoff_arrays(glob_loc_list, glob_cell_v, glob_loc_list_a)
826 END IF
827 ! Special way of handling the siepmann potential..
828 IF (any_siepmann) THEN
829 NULLIFY (glob_loc_list, glob_cell_v, glob_loc_list_a)
830 CALL setup_siepmann_arrays(nonbonded, potparm, glob_loc_list, glob_cell_v, glob_loc_list_a, cell)
831 DO ilist = 1, nonbonded%nlists
832 neighbor_kind_pair => nonbonded%neighbor_kind_pairs(ilist)
833 npairs = neighbor_kind_pair%npairs
834 IF (npairs == 0) cycle
835 kind_group_loop3: DO igrp = 1, neighbor_kind_pair%ngrp_kind
836 istart = neighbor_kind_pair%grp_kind_start(igrp)
837 iend = neighbor_kind_pair%grp_kind_end(igrp)
838 ikind = neighbor_kind_pair%ij_kind(1, igrp)
839 jkind = neighbor_kind_pair%ij_kind(2, igrp)
840 list => neighbor_kind_pair%list
841 cvi = neighbor_kind_pair%cell_vector
842 pot => potparm%pot(ikind, jkind)%pot
843
844 IF (pot%no_mb) cycle kind_group_loop3
845 rab2_max = pot%rcutsq
846 cell_v = matmul(cell%hmat, cvi)
847 DO i = 1, SIZE(pot%type)
848 ! SIEPMANN
849 IF (pot%type(i) == siepmann_type) THEN
850 npairs = iend - istart + 1
851 siepmann => pot%set(i)%siepmann
852 ALLOCATE (sort_list(2, npairs), work_list(npairs))
853 sort_list = list(:, istart:iend)
854 ! Sort the list of neighbors, this increases the efficiency for single
855 ! potential contributions
856 CALL sort(sort_list(1, :), npairs, work_list)
857 DO ipair = 1, npairs
858 work_list(ipair) = sort_list(2, work_list(ipair))
859 END DO
860 sort_list(2, :) = work_list
861 ! find number of unique elements of array index 1
862 nunique = 1
863 DO ipair = 1, npairs - 1
864 IF (sort_list(1, ipair + 1) /= sort_list(1, ipair)) nunique = nunique + 1
865 END DO
866 ipair = 1
867 junique = sort_list(1, ipair)
868 ifirst = 1
869 DO iunique = 1, nunique
870 atom_a = junique
871 IF (glob_loc_list_a(ifirst) > atom_a) cycle
872 DO mpair = ifirst, SIZE(glob_loc_list_a)
873 IF (glob_loc_list_a(mpair) == atom_a) EXIT
874 END DO
875 ifirst = mpair
876 DO mpair = ifirst, SIZE(glob_loc_list_a)
877 IF (glob_loc_list_a(mpair) /= atom_a) EXIT
878 END DO
879 ilast = mpair - 1
880 nloc_size = 0
881 IF (ifirst /= 0) nloc_size = ilast - ifirst + 1
882 DO WHILE (ipair <= npairs)
883 IF (sort_list(1, ipair) /= junique) EXIT
884 atom_b = sort_list(2, ipair)
885 ! Derivative terms
886 rtmp = r_last_update_pbc(atom_b)%r(:) - r_last_update_pbc(atom_a)%r(:) + cell_v
887 ipair = ipair + 1
888 IF (dot_product(rtmp, rtmp) <= siepmann%rcutsq) THEN
889 CALL siepmann_forces_v2(siepmann, r_last_update_pbc, cell_v, cell, &
890 atom_a, atom_b, f_nonbond, use_virial, siepmann%rcutsq, &
891 particle_set)
892 CALL siepmann_forces_v3(siepmann, r_last_update_pbc, cell_v, &
893 nloc_size, glob_loc_list(:, ifirst:ilast), &
894 atom_a, atom_b, f_nonbond, use_virial, siepmann%rcutsq, &
895 cell, particle_set)
896 END IF
897 END DO
898 ifirst = ilast + 1
899 IF (ipair <= npairs) junique = sort_list(1, ipair)
900 END DO
901 DEALLOCATE (sort_list, work_list)
902 END IF
903 END DO
904 END DO kind_group_loop3
905 END DO
906 CALL destroy_siepmann_arrays(glob_loc_list, glob_cell_v, glob_loc_list_a)
907 END IF
908
909 ! GAL19 potential..
910 IF (any_gal) THEN
911 NULLIFY (glob_loc_list, glob_cell_v, glob_loc_list_a)
912 CALL setup_gal_arrays(nonbonded, potparm, glob_loc_list, glob_cell_v, glob_loc_list_a, cell)
913 DO ilist = 1, nonbonded%nlists
914 neighbor_kind_pair => nonbonded%neighbor_kind_pairs(ilist)
915 npairs = neighbor_kind_pair%npairs
916 IF (npairs == 0) cycle
917 kind_group_loop4: DO igrp = 1, neighbor_kind_pair%ngrp_kind
918 istart = neighbor_kind_pair%grp_kind_start(igrp)
919 iend = neighbor_kind_pair%grp_kind_end(igrp)
920 ikind = neighbor_kind_pair%ij_kind(1, igrp)
921 jkind = neighbor_kind_pair%ij_kind(2, igrp)
922 list => neighbor_kind_pair%list
923 cvi = neighbor_kind_pair%cell_vector
924 pot => potparm%pot(ikind, jkind)%pot
925
926 IF (pot%no_mb) cycle kind_group_loop4
927 rab2_max = pot%rcutsq
928 cell_v = matmul(cell%hmat, cvi)
929 DO i = 1, SIZE(pot%type)
930 ! GAL19
931 IF (pot%type(i) == gal_type) THEN
932 npairs = iend - istart + 1
933 gal => pot%set(i)%gal
934 ALLOCATE (sort_list(2, npairs), work_list(npairs))
935 sort_list = list(:, istart:iend)
936 ! Sort the list of neighbors, this increases the efficiency for single
937 ! potential contributions
938 CALL sort(sort_list(1, :), npairs, work_list)
939 DO ipair = 1, npairs
940 work_list(ipair) = sort_list(2, work_list(ipair))
941 END DO
942 sort_list(2, :) = work_list
943 ! find number of unique elements of array index 1
944 nunique = 1
945 DO ipair = 1, npairs - 1
946 IF (sort_list(1, ipair + 1) /= sort_list(1, ipair)) nunique = nunique + 1
947 END DO
948 ipair = 1
949 junique = sort_list(1, ipair)
950 ifirst = 1
951 DO iunique = 1, nunique
952 atom_a = junique
953 IF (glob_loc_list_a(ifirst) > atom_a) cycle
954 DO mpair = ifirst, SIZE(glob_loc_list_a)
955 IF (glob_loc_list_a(mpair) == atom_a) EXIT
956 END DO
957 ifirst = mpair
958 DO mpair = ifirst, SIZE(glob_loc_list_a)
959 IF (glob_loc_list_a(mpair) /= atom_a) EXIT
960 END DO
961 ilast = mpair - 1
962 nloc_size = 0
963 IF (ifirst /= 0) nloc_size = ilast - ifirst + 1
964 DO WHILE (ipair <= npairs)
965 IF (sort_list(1, ipair) /= junique) EXIT
966 atom_b = sort_list(2, ipair)
967 ! Derivative terms
968 rtmp = r_last_update_pbc(atom_b)%r(:) - r_last_update_pbc(atom_a)%r(:) + cell_v
969 ipair = ipair + 1
970 IF (dot_product(rtmp, rtmp) <= gal%rcutsq) THEN
971 CALL gal_forces(gal, r_last_update_pbc, &
972 atom_a, atom_b, f_nonbond, use_virial, &
973 cell, particle_set)
974 END IF
975 END DO
976 ifirst = ilast + 1
977 IF (ipair <= npairs) junique = sort_list(1, ipair)
978 END DO
979 DEALLOCATE (sort_list, work_list)
980 END IF
981 END DO
982 END DO kind_group_loop4
983 END DO
984 CALL destroy_gal_arrays(glob_loc_list, glob_cell_v, glob_loc_list_a)
985 END IF
986
987 ! GAL21 potential..
988 IF (any_gal21) THEN
989 NULLIFY (glob_loc_list, glob_cell_v, glob_loc_list_a)
990 CALL setup_gal21_arrays(nonbonded, potparm, glob_loc_list, glob_cell_v, glob_loc_list_a, cell)
991 DO ilist = 1, nonbonded%nlists
992 neighbor_kind_pair => nonbonded%neighbor_kind_pairs(ilist)
993 npairs = neighbor_kind_pair%npairs
994 IF (npairs == 0) cycle
995 kind_group_loop6: DO igrp = 1, neighbor_kind_pair%ngrp_kind
996 istart = neighbor_kind_pair%grp_kind_start(igrp)
997 iend = neighbor_kind_pair%grp_kind_end(igrp)
998 ikind = neighbor_kind_pair%ij_kind(1, igrp)
999 jkind = neighbor_kind_pair%ij_kind(2, igrp)
1000 list => neighbor_kind_pair%list
1001 cvi = neighbor_kind_pair%cell_vector
1002 pot => potparm%pot(ikind, jkind)%pot
1003
1004 IF (pot%no_mb) cycle kind_group_loop6
1005 rab2_max = pot%rcutsq
1006 cell_v = matmul(cell%hmat, cvi)
1007 DO i = 1, SIZE(pot%type)
1008 ! GAL21
1009 IF (pot%type(i) == gal21_type) THEN
1010 npairs = iend - istart + 1
1011 gal21 => pot%set(i)%gal21
1012 ALLOCATE (sort_list(2, npairs), work_list(npairs))
1013 sort_list = list(:, istart:iend)
1014 ! Sort the list of neighbors, this increases the efficiency for single
1015 ! potential contributions
1016 CALL sort(sort_list(1, :), npairs, work_list)
1017 DO ipair = 1, npairs
1018 work_list(ipair) = sort_list(2, work_list(ipair))
1019 END DO
1020 sort_list(2, :) = work_list
1021 ! find number of unique elements of array index 1
1022 nunique = 1
1023 DO ipair = 1, npairs - 1
1024 IF (sort_list(1, ipair + 1) /= sort_list(1, ipair)) nunique = nunique + 1
1025 END DO
1026 ipair = 1
1027 junique = sort_list(1, ipair)
1028 ifirst = 1
1029 DO iunique = 1, nunique
1030 atom_a = junique
1031 IF (glob_loc_list_a(ifirst) > atom_a) cycle
1032 DO mpair = ifirst, SIZE(glob_loc_list_a)
1033 IF (glob_loc_list_a(mpair) == atom_a) EXIT
1034 END DO
1035 ifirst = mpair
1036 DO mpair = ifirst, SIZE(glob_loc_list_a)
1037 IF (glob_loc_list_a(mpair) /= atom_a) EXIT
1038 END DO
1039 ilast = mpair - 1
1040 nloc_size = 0
1041 IF (ifirst /= 0) nloc_size = ilast - ifirst + 1
1042 DO WHILE (ipair <= npairs)
1043 IF (sort_list(1, ipair) /= junique) EXIT
1044 atom_b = sort_list(2, ipair)
1045 ! Derivative terms
1046 rtmp = r_last_update_pbc(atom_b)%r(:) - r_last_update_pbc(atom_a)%r(:) + cell_v
1047 ipair = ipair + 1
1048 IF (dot_product(rtmp, rtmp) <= gal21%rcutsq) THEN
1049 CALL gal21_forces(gal21, r_last_update_pbc, &
1050 atom_a, atom_b, f_nonbond, pv_nonbond, use_virial, &
1051 cell, particle_set)
1052 END IF
1053 END DO
1054 ifirst = ilast + 1
1055 IF (ipair <= npairs) junique = sort_list(1, ipair)
1056 END DO
1057 DEALLOCATE (sort_list, work_list)
1058 END IF
1059 END DO
1060 END DO kind_group_loop6
1061 END DO
1062 CALL destroy_gal_arrays(glob_loc_list, glob_cell_v, glob_loc_list_a)
1063 END IF
1064
1065 IF (use_virial) THEN
1066 pv_nonbond(1, 1) = pv_nonbond(1, 1) + ptens11
1067 pv_nonbond(1, 2) = pv_nonbond(1, 2) + ptens12
1068 pv_nonbond(1, 3) = pv_nonbond(1, 3) + ptens13
1069 pv_nonbond(2, 1) = pv_nonbond(2, 1) + ptens21
1070 pv_nonbond(2, 2) = pv_nonbond(2, 2) + ptens22
1071 pv_nonbond(2, 3) = pv_nonbond(2, 3) + ptens23
1072 pv_nonbond(3, 1) = pv_nonbond(3, 1) + ptens31
1073 pv_nonbond(3, 2) = pv_nonbond(3, 2) + ptens32
1074 pv_nonbond(3, 3) = pv_nonbond(3, 3) + ptens33
1075 END IF
1076 CALL timestop(handle)
1077 END SUBROUTINE force_nonbond_manybody
1078
1079END MODULE manybody_potential
1080
static GRID_HOST_DEVICE double fac(const int i)
Factorial function, e.g. fac(5) = 5! = 120.
Definition grid_common.h:56
Define the atomic kind types and their sub types.
Handles all functions related to the CELL.
Definition cell_types.F:15
stores a lists of integer that are local to a processor. The idea is that these integers represent ob...
Define the neighbor list data types and the corresponding functionality.
subroutine, public fist_nonbond_env_get(fist_nonbond_env, potparm14, potparm, nonbonded, rlist_cut, rlist_lowsq, aup, lup, ei_scale14, vdw_scale14, shift_cutoff, do_electrostatics, r_last_update, r_last_update_pbc, rshell_last_update_pbc, rcore_last_update_pbc, cell_last_update, num_update, last_update, counter, natom_types, long_range_correction, ij_kind_full_fac, eam_data, nequip_data, deepmd_data, ace_data, charges)
sets a fist_nonbond_env
objects that represent the structure of input sections and the data contained in an input section
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
An array-based list which grows on demand. When the internal array is full, a new array of twice the ...
Definition list.F:24
subroutine, public ace_energy_store_force_virial(particle_set, cell, atomic_kind_set, potparm, fist_nonbond_env, pot_ace)
... >
subroutine, public ace_add_force_virial(fist_nonbond_env, force, pv_nonbond, use_virial)
...
subroutine, public deepmd_energy_store_force_virial(particle_set, cell, atomic_kind_set, potparm, fist_nonbond_env, pot_deepmd, para_env)
...
subroutine, public deepmd_add_force_virial(fist_nonbond_env, force, pv_nonbond, use_virial)
...
subroutine, public get_force_eam(rab2, eam_a, eam_b, eam_data, atom_a, atom_b, f_eam)
...
Implementation of the GAL21 potential.
subroutine, public destroy_gal21_arrays(glob_loc_list, glob_cell_v, glob_loc_list_a)
...
subroutine, public gal21_forces(gal21, r_last_update_pbc, iparticle, jparticle, f_nonbond, pv_nonbond, use_virial, cell, particle_set)
forces generated by the GAL2119 potential
subroutine, public setup_gal21_arrays(nonbonded, potparm, glob_loc_list, glob_cell_v, glob_loc_list_a, cell)
...
subroutine, public gal21_energy(pot_loc, gal21, r_last_update_pbc, iparticle, jparticle, cell, particle_set, mm_section)
Main part of the energy evaluation of GAL2119.
Implementation of the GAL19 potential.
subroutine, public destroy_gal_arrays(glob_loc_list, glob_cell_v, glob_loc_list_a)
...
subroutine, public setup_gal_arrays(nonbonded, potparm, glob_loc_list, glob_cell_v, glob_loc_list_a, cell)
...
subroutine, public gal_energy(pot_loc, gal, r_last_update_pbc, iparticle, jparticle, cell, particle_set, mm_section)
Main part of the energy evaluation of GAL19.
subroutine, public gal_forces(gal, r_last_update_pbc, iparticle, jparticle, f_nonbond, use_virial, cell, particle_set)
forces generated by the GAL19 potential
subroutine, public nequip_add_force_virial(fist_nonbond_env, f_nonbond, pv_nonbond, use_virial)
...
subroutine, public nequip_energy_store_force_virial(nonbonded, particle_set, local_particles, cell, atomic_kind_set, potparm, r_last_update_pbc, pot_total, fist_nonbond_env, para_env, use_virial, target_pot_type)
...
subroutine, public force_nonbond_manybody(fist_nonbond_env, particle_set, cell, f_nonbond, pv_nonbond, use_virial)
...
subroutine, public energy_manybody(fist_nonbond_env, atomic_kind_set, local_particles, particle_set, cell, pot_manybody, para_env, mm_section, use_virial)
computes the embedding contribution to the energy
implementation of dipole and three-body part of Siepmann-Sprik potential dipole term: 3rd term in Eq....
subroutine, public siepmann_energy(pot_loc, siepmann, r_last_update_pbc, atom_a, atom_b, nloc_size, full_loc_list, cell_v, cell, drij, particle_set, nr_oh, nr_h3o, nr_o)
energy of two-body dipole term and three-body term
subroutine, public setup_siepmann_arrays(nonbonded, potparm, glob_loc_list, glob_cell_v, glob_loc_list_a, cell)
...
subroutine, public siepmann_forces_v2(siepmann, r_last_update_pbc, cell_v, cell, iparticle, jparticle, f_nonbond, use_virial, rcutsq, particle_set)
forces generated by the dipole term
subroutine, public siepmann_forces_v3(siepmann, r_last_update_pbc, cell_v, n_loc_size, full_loc_list, iparticle, jparticle, f_nonbond, use_virial, rcutsq, cell, particle_set)
forces generated by the three-body term
subroutine, public print_nr_ions_siepmann(nr_ions, mm_section, para_env, print_oh, print_h3o, print_o)
prints the number of OH- ions or H3O+ ions near surface
subroutine, public destroy_siepmann_arrays(glob_loc_list, glob_cell_v, glob_loc_list_a)
...
subroutine, public setup_tersoff_arrays(nonbonded, potparm, glob_loc_list, glob_cell_v, glob_loc_list_a, cell)
...
subroutine, public tersoff_forces(tersoff, r_last_update_pbc, cell_v, n_loc_size, full_loc_list, loc_cell_v, iparticle, jparticle, f_nonbond, pv_nonbond, use_virial, rcutsq)
...
subroutine, public destroy_tersoff_arrays(glob_loc_list, glob_cell_v, glob_loc_list_a)
...
subroutine, public tersoff_energy(pot_loc, tersoff, r_last_update_pbc, atom_a, atom_b, nloc_size, full_loc_list, loc_cell_v, cell_v, drij)
...
Interface to the message passing library MPI.
integer, parameter, public allegro_type
integer, parameter, public gal_type
integer, parameter, public nequip_type
integer, parameter, public deepmd_type
integer, parameter, public siepmann_type
integer, parameter, public ace_type
integer, parameter, public gal21_type
integer, parameter, public ea_type
integer, parameter, public tersoff_type
Define the data structure for the particle information.
All kind of helpful little routines.
Definition util.F:14
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:60
structure to store local (to a processor) ordered lists of integers.
stores all the informations relevant to an mpi environment