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