(git:d18deda)
Loading...
Searching...
No Matches
pair_potential_types.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!> Teodoro Laino [Teo] 11.2005 : Reorganizing the structures to optimize
11!> memory management
12!> \author CJM
13! **************************************************************************************************
15
16 USE kinds, ONLY: default_path_length,&
18 dp
26#include "./base/base_uses.f90"
27
28 IMPLICIT NONE
29
30 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'pair_potential_types'
31
32 PRIVATE
33 ! when adding a new nonbonded potential please update also the list_pot
34 ! used for the linear scaling screening of potential calculation
35 INTEGER, PUBLIC, PARAMETER :: multi_type = -1, &
36 nn_type = 0, &
37 lj_type = 1, &
38 lj_charmm_type = 2, &
39 ft_type = 3, &
40 wl_type = 4, &
41 gw_type = 5, &
42 ip_type = 6, &
43 ea_type = 7, &
44 b4_type = 8, &
45 bm_type = 9, &
46 gp_type = 10, &
47 tersoff_type = 11, &
48 ftd_type = 12, &
49 siepmann_type = 13, &
50 gal_type = 14, &
51 quip_type = 15, &
52 nequip_type = 16, &
53 allegro_type = 17, &
54 gal21_type = 18, &
55 tab_type = 19, &
56 deepmd_type = 20
57
58 INTEGER, PUBLIC, PARAMETER, DIMENSION(21) :: list_pot = (/nn_type, &
59 lj_type, &
61 ft_type, &
62 wl_type, &
63 gw_type, &
64 ip_type, &
65 ea_type, &
66 b4_type, &
67 bm_type, &
68 gp_type, &
70 ftd_type, &
72 gal_type, &
73 quip_type, &
76 gal21_type, &
77 tab_type, &
79
80 ! Shell model
81 INTEGER, PUBLIC, PARAMETER :: nosh_nosh = 0, &
82 nosh_sh = 1, &
83 sh_sh = 2
84
85 INTEGER, PUBLIC, PARAMETER, DIMENSION(3) :: list_sh_type = (/nosh_nosh, nosh_sh, sh_sh/)
86
87 ! Single Spline generation info
88 REAL(kind=dp), PARAMETER, PUBLIC :: not_initialized = -huge(0.0_dp)
89 INTEGER, PARAMETER, DIMENSION(2), PUBLIC :: do_potential_single_allocation = (/lj_type, lj_charmm_type/)
90 INTEGER, PARAMETER, DIMENSION(2), PUBLIC :: no_potential_single_allocation = (/-huge(0), -huge(0)/)
91 INTEGER, DIMENSION(2), PUBLIC :: potential_single_allocation
92
94
99
100 PUBLIC :: pair_potential_pp_create, &
103
104 PUBLIC :: pair_potential_p_type, &
106
107 PUBLIC :: ft_pot_type, &
109 eam_pot_type, &
116 gal_pot_type, &
119
121 PUBLIC :: compare_pot
122
123! **************************************************************************************************
125 REAL(kind=dp), DIMENSION(2:15) :: a = 0.0_dp
126 REAL(kind=dp) :: rcore = 0.0_dp
127 REAL(kind=dp) :: m = 0.0_dp
128 REAL(kind=dp) :: b = 0.0_dp
129 END TYPE ipbv_pot_type
130
131! **************************************************************************************************
132 TYPE lj_pot_type
133 REAL(kind=dp) :: epsilon = 0.0_dp
134 REAL(kind=dp) :: sigma6 = 0.0_dp
135 REAL(kind=dp) :: sigma12 = 0.0_dp
136 END TYPE lj_pot_type
137
138! **************************************************************************************************
140 REAL(kind=dp) :: a = 0.0_dp
141 REAL(kind=dp) :: b = 0.0_dp
142 REAL(kind=dp) :: c = 0.0_dp
143 REAL(kind=dp) :: d = 0.0_dp
144 END TYPE ft_pot_type
145
146! **************************************************************************************************
147 TYPE ftd_pot_type
148 REAL(kind=dp) :: a = 0.0_dp
149 REAL(kind=dp) :: b = 0.0_dp
150 REAL(kind=dp) :: c = 0.0_dp
151 REAL(kind=dp) :: d = 0.0_dp
152 REAL(kind=dp), DIMENSION(2) :: bd = 0.0_dp
153 END TYPE ftd_pot_type
154
155! **************************************************************************************************
156 TYPE williams_pot_type
157 REAL(kind=dp) :: a = 0.0_dp
158 REAL(kind=dp) :: b = 0.0_dp
159 REAL(kind=dp) :: c = 0.0_dp
160 END TYPE williams_pot_type
161
162! **************************************************************************************************
163 TYPE goodwin_pot_type
164 REAL(kind=dp) :: vr0 = 0.0_dp
165 REAL(kind=dp) :: m = 0.0_dp, mc = 0.0_dp
166 REAL(kind=dp) :: d = 0.0_dp, dc = 0.0_dp
167 END TYPE goodwin_pot_type
168
169! **************************************************************************************************
171 CHARACTER(LEN=default_path_length) :: eam_file_name = ""
172 INTEGER :: npoints = 0
173 REAL(kind=dp) :: drar = 0.0_dp, drhoar = 0.0_dp, acutal = 0.0_dp
174 REAL(kind=dp), POINTER, DIMENSION(:) :: rho => null(), phi => null(), frho => null(), rhoval => null(), rval => null()
175 REAL(kind=dp), POINTER, DIMENSION(:) :: rhop => null(), phip => null(), frhop => null()
176 END TYPE eam_pot_type
177
178! **************************************************************************************************
180 CHARACTER(LEN=default_path_length) :: deepmd_file_name = 'NULL'
181 INTEGER :: atom_deepmd_type = 0
182 END TYPE deepmd_pot_type
183
184! **************************************************************************************************
186 CHARACTER(LEN=default_path_length) :: quip_file_name = ""
187 CHARACTER(LEN=1024) :: init_args = ""
188 CHARACTER(LEN=1024) :: calc_args = ""
189 END TYPE quip_pot_type
190
191! **************************************************************************************************
193 CHARACTER(LEN=default_path_length) :: nequip_file_name = 'NULL', nequip_version = 'NULL', &
194 unit_coords = 'NULL', unit_forces = 'NULL', &
195 unit_energy = 'NULL', unit_cell = 'NULL'
196 CHARACTER(LEN=100), DIMENSION(:), ALLOCATABLE :: type_names_torch
197 REAL(kind=dp) :: rcutsq = 0.0_dp, unit_coords_val = 1.0_dp, &
198 unit_forces_val = 1.0_dp, unit_energy_val = 1.0_dp, &
199 unit_cell_val = 1.0_dp
200 LOGICAL :: do_nequip_sp = .false.
201 END TYPE nequip_pot_type
202
203! **************************************************************************************************
205 CHARACTER(LEN=default_path_length) :: allegro_file_name = 'NULL', unit_cell = 'NULL', &
206 nequip_version = 'NULL', unit_coords = 'NULL', &
207 unit_forces = 'NULL', unit_energy = 'NULL'
208
209 CHARACTER(LEN=100), DIMENSION(:), ALLOCATABLE :: type_names_torch
210
211 REAL(kind=dp) :: rcutsq = 0.0_dp, unit_coords_val = 1.0_dp, &
212 unit_forces_val = 1.0_dp, unit_cell_val = 1.0_dp, &
213 unit_energy_val = 1.0_dp
214 LOGICAL :: do_allegro_sp = .false.
215 END TYPE allegro_pot_type
216
217! **************************************************************************************************
218 TYPE buck4ran_pot_type
219 REAL(kind=dp) :: a = 0.0_dp
220 REAL(kind=dp) :: b = 0.0_dp
221 REAL(kind=dp) :: c = 0.0_dp
222 REAL(kind=dp) :: r1 = 0.0_dp
223 REAL(kind=dp) :: r2 = 0.0_dp
224 REAL(kind=dp) :: r3 = 0.0_dp
225 INTEGER :: npoly1 = 0, npoly2 = 0
226 REAL(kind=dp), DIMENSION(0:10) :: poly1 = 0.0_dp
227 REAL(kind=dp), DIMENSION(0:10) :: poly2 = 0.0_dp
228 END TYPE buck4ran_pot_type
229
230! **************************************************************************************************
231 TYPE buckmorse_pot_type
232 REAL(kind=dp) :: f0 = 0.0_dp
233 REAL(kind=dp) :: a1 = 0.0_dp
234 REAL(kind=dp) :: a2 = 0.0_dp
235 REAL(kind=dp) :: b1 = 0.0_dp
236 REAL(kind=dp) :: b2 = 0.0_dp
237 REAL(kind=dp) :: c = 0.0_dp
238 REAL(kind=dp) :: d = 0.0_dp
239 REAL(kind=dp) :: r0 = 0.0_dp
240 REAL(kind=dp) :: beta = 0.0_dp
241 END TYPE buckmorse_pot_type
242
243! **************************************************************************************************
244 TYPE gp_pot_type
245 INTEGER :: myid = 0
246 CHARACTER(LEN=default_path_length) :: potential = ""
247 CHARACTER(LEN=default_string_length), &
248 POINTER, DIMENSION(:) :: parameters => null(), units => null()
249 CHARACTER(LEN=default_string_length) :: variables = ""
250 REAL(kind=dp), DIMENSION(:), POINTER :: values => null()
251 END TYPE gp_pot_type
252
253! **************************************************************************************************
255 ! Get this stuff from the PRB V38, N14 9902 (1988) by Tersoff
256 REAL(kind=dp) :: a = 0.0_dp
257 REAL(kind=dp) :: b = 0.0_dp
258 REAL(kind=dp) :: lambda1 = 0.0_dp
259 REAL(kind=dp) :: lambda2 = 0.0_dp
260 REAL(kind=dp) :: alpha = 0.0_dp
261 REAL(kind=dp) :: beta = 0.0_dp
262 REAL(kind=dp) :: n = 0.0_dp
263 REAL(kind=dp) :: c = 0.0_dp
264 REAL(kind=dp) :: d = 0.0_dp
265 REAL(kind=dp) :: h = 0.0_dp
266 REAL(kind=dp) :: lambda3 = 0.0_dp
267 REAL(kind=dp) :: bigr = 0.0_dp ! Used to be R = Rij + D
268 REAL(kind=dp) :: bigd = 0.0_dp ! Used to be D = Rij - D
269 REAL(kind=dp) :: rcutsq = 0.0_dp ! Always set to (bigR+bigD)^2
270 END TYPE tersoff_pot_type
271
272! **************************************************************************************************
274 REAL(kind=dp) :: b = 0.0_dp
275 REAL(kind=dp) :: d = 0.0_dp
276 REAL(kind=dp) :: e = 0.0_dp
277 REAL(kind=dp) :: f = 0.0_dp
278 REAL(kind=dp) :: beta = 0.0_dp
279 REAL(kind=dp) :: rcutsq = 0.0_dp
280 LOGICAL :: allow_oh_formation = .false.
281 LOGICAL :: allow_h3o_formation = .false.
282 LOGICAL :: allow_o_formation = .false.
283 END TYPE siepmann_pot_type
284
285! **************************************************************************************************
287 CHARACTER(LEN=2) :: met1 = ""
288 CHARACTER(LEN=2) :: met2 = ""
289 REAL(kind=dp) :: epsilon = 0.0_dp
290 REAL(kind=dp) :: bxy = 0.0_dp
291 REAL(kind=dp) :: bz = 0.0_dp
292 REAL(kind=dp) :: r1 = 0.0_dp
293 REAL(kind=dp) :: r2 = 0.0_dp
294 REAL(kind=dp) :: a1 = 0.0_dp
295 REAL(kind=dp) :: a2 = 0.0_dp
296 REAL(kind=dp) :: a3 = 0.0_dp
297 REAL(kind=dp) :: a4 = 0.0_dp
298 REAL(kind=dp) :: a = 0.0_dp
299 REAL(kind=dp) :: b = 0.0_dp
300 REAL(kind=dp) :: c = 0.0_dp
301 REAL(kind=dp), POINTER, DIMENSION(:) :: gcn => null()
302 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: n_vectors
303 REAL(kind=dp) :: rcutsq = 0.0_dp
304 LOGICAL :: express = .false.
305 END TYPE gal_pot_type
306
307! **************************************************************************************************
308
310 CHARACTER(LEN=2) :: met1 = ""
311 CHARACTER(LEN=2) :: met2 = ""
312 REAL(kind=dp) :: epsilon1 = 0.0_dp
313 REAL(kind=dp) :: epsilon2 = 0.0_dp
314 REAL(kind=dp) :: epsilon3 = 0.0_dp
315 REAL(kind=dp) :: bxy1 = 0.0_dp
316 REAL(kind=dp) :: bxy2 = 0.0_dp
317 REAL(kind=dp) :: bz1 = 0.0_dp
318 REAL(kind=dp) :: bz2 = 0.0_dp
319 REAL(kind=dp) :: r1 = 0.0_dp
320 REAL(kind=dp) :: r2 = 0.0_dp
321 REAL(kind=dp) :: a11 = 0.0_dp
322 REAL(kind=dp) :: a12 = 0.0_dp
323 REAL(kind=dp) :: a13 = 0.0_dp
324 REAL(kind=dp) :: a21 = 0.0_dp
325 REAL(kind=dp) :: a22 = 0.0_dp
326 REAL(kind=dp) :: a23 = 0.0_dp
327 REAL(kind=dp) :: a31 = 0.0_dp
328 REAL(kind=dp) :: a32 = 0.0_dp
329 REAL(kind=dp) :: a33 = 0.0_dp
330 REAL(kind=dp) :: a41 = 0.0_dp
331 REAL(kind=dp) :: a42 = 0.0_dp
332 REAL(kind=dp) :: a43 = 0.0_dp
333 REAL(kind=dp) :: ao1 = 0.0_dp
334 REAL(kind=dp) :: ao2 = 0.0_dp
335 REAL(kind=dp) :: bo1 = 0.0_dp
336 REAL(kind=dp) :: bo2 = 0.0_dp
337 REAL(kind=dp) :: c = 0.0_dp
338 REAL(kind=dp) :: ah1 = 0.0_dp
339 REAL(kind=dp) :: ah2 = 0.0_dp
340 REAL(kind=dp) :: bh1 = 0.0_dp
341 REAL(kind=dp) :: bh2 = 0.0_dp
342 REAL(kind=dp), POINTER, DIMENSION(:) :: gcn => null()
343 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: n_vectors
344 REAL(kind=dp) :: rcutsq = 0.0_dp
345 LOGICAL :: express = .false.
346 END TYPE gal21_pot_type
347
348! **************************************************************************************************
349
351 CHARACTER(LEN=default_path_length) :: tabpot_file_name = ""
352 INTEGER :: npoints = 0, index = 0
353 REAL(kind=dp) :: dr = 0.0_dp, rcut = 0.0_dp
354 REAL(kind=dp), POINTER, DIMENSION(:) :: r => null(), e => null(), f => null()
355 END TYPE tab_pot_type
356
357! **************************************************************************************************
358
359 TYPE pot_set_type
360 REAL(kind=dp) :: rmin = 0.0_dp, rmax = 0.0_dp
361 TYPE(ipbv_pot_type), POINTER :: ipbv => null()
362 TYPE(gp_pot_type), POINTER :: gp => null()
363 TYPE(lj_pot_type), POINTER :: lj => null()
364 TYPE(ft_pot_type), POINTER :: ft => null()
365 TYPE(williams_pot_type), POINTER :: willis => null()
366 TYPE(goodwin_pot_type), POINTER :: goodwin => null()
367 TYPE(eam_pot_type), POINTER :: eam => null()
368 TYPE(quip_pot_type), POINTER :: quip => null()
369 TYPE(nequip_pot_type), POINTER :: nequip => null()
370 TYPE(allegro_pot_type), POINTER :: allegro => null()
371 TYPE(deepmd_pot_type), POINTER :: deepmd => null()
372 TYPE(buck4ran_pot_type), POINTER :: buck4r => null()
373 TYPE(buckmorse_pot_type), POINTER :: buckmo => null()
374 TYPE(tersoff_pot_type), POINTER :: tersoff => null()
375 TYPE(siepmann_pot_type), POINTER :: siepmann => null()
376 TYPE(gal_pot_type), POINTER :: gal => null()
377 TYPE(gal21_pot_type), POINTER :: gal21 => null()
378 TYPE(ftd_pot_type), POINTER :: ftd => null()
379 TYPE(tab_pot_type), POINTER :: tab => null()
380 END TYPE pot_set_type
381
382! **************************************************************************************************
384 REAL(kind=dp) :: rcutsq = 0.0_dp
385 REAL(kind=dp) :: e_fac = 0.0_dp
386 REAL(kind=dp) :: e_fcc = 0.0_dp
387 REAL(kind=dp) :: e_fcs = 0.0_dp
388 REAL(kind=dp) :: e_fsc = 0.0_dp
389 REAL(kind=dp) :: z1 = 0.0_dp
390 REAL(kind=dp) :: z2 = 0.0_dp
391 REAL(kind=dp), DIMENSION(0:5) :: zbl_poly = 0.0_dp
392 REAL(kind=dp), DIMENSION(2) :: zbl_rcut = 0.0_dp
393 LOGICAL :: undef = .false., & ! non-bonding interaction not defined
394 no_mb = .false., & ! no many-body potential
395 no_pp = .false. ! no pair (=two-body) potential
396 INTEGER :: shell_type = 0
397 CHARACTER(LEN=default_string_length) :: at1 = ""
398 CHARACTER(LEN=default_string_length) :: at2 = ""
399 INTEGER, POINTER, DIMENSION(:) :: TYPE => null()
400 TYPE(pot_set_type), POINTER, DIMENSION(:) :: set => null()
401 TYPE(spline_data_p_type), POINTER, DIMENSION(:) :: pair_spline_data => null()
402 TYPE(spline_factor_type), POINTER :: spl_f => null()
404
405! **************************************************************************************************
406 TYPE pair_potential_type
407 TYPE(pair_potential_single_type), POINTER :: pot => null()
408 END TYPE pair_potential_type
409
410! **************************************************************************************************
412 TYPE(pair_potential_type), DIMENSION(:), POINTER :: pot => null()
413 END TYPE pair_potential_p_type
414
415! **************************************************************************************************
417 TYPE(pair_potential_type), DIMENSION(:, :), POINTER :: pot => null()
419
420CONTAINS
421
422! **************************************************************************************************
423!> \brief compare two different potentials
424!> \param pot1 ...
425!> \param pot2 ...
426!> \param compare ...
427!> \author Teodoro Laino [teo] 05.2006
428! **************************************************************************************************
429 SUBROUTINE compare_pot(pot1, pot2, compare)
430 TYPE(pair_potential_single_type), POINTER :: pot1, pot2
431 LOGICAL, INTENT(OUT) :: compare
432
433 INTEGER :: i
434 LOGICAL :: mycompare
435
436 compare = .false.
437 ! Preliminary checks
438
439 cpassert(ASSOCIATED(pot1%type))
440 cpassert(ASSOCIATED(pot2%type))
441 IF (SIZE(pot1%type) /= SIZE(pot2%type)) RETURN
442 IF (any(pot1%type /= pot2%type)) RETURN
443
444 ! Checking the real values of parameters
445 cpassert(ASSOCIATED(pot1%set))
446 cpassert(ASSOCIATED(pot2%set))
447 DO i = 1, SIZE(pot1%type)
448 mycompare = .false.
449 SELECT CASE (pot1%type(i))
450 CASE (lj_type, lj_charmm_type)
451 IF ((pot1%set(i)%lj%epsilon == pot2%set(i)%lj%epsilon) .AND. &
452 (pot1%set(i)%lj%sigma6 == pot2%set(i)%lj%sigma6) .AND. &
453 (pot1%set(i)%lj%sigma12 == pot2%set(i)%lj%sigma12)) mycompare = .true.
454 CASE (wl_type)
455 IF ((pot1%set(i)%willis%a == pot2%set(i)%willis%a) .AND. &
456 (pot1%set(i)%willis%b == pot2%set(i)%willis%b) .AND. &
457 (pot1%set(i)%willis%c == pot2%set(i)%willis%c)) mycompare = .true.
458 CASE (gw_type)
459 IF ((pot1%set(i)%goodwin%vr0 == pot2%set(i)%goodwin%vr0) .AND. &
460 (pot1%set(i)%goodwin%m == pot2%set(i)%goodwin%m) .AND. &
461 (pot1%set(i)%goodwin%mc == pot2%set(i)%goodwin%mc) .AND. &
462 (pot1%set(i)%goodwin%d == pot2%set(i)%goodwin%d) .AND. &
463 (pot1%set(i)%goodwin%dc == pot2%set(i)%goodwin%dc)) mycompare = .true.
464 CASE (ea_type)
465 ! Compare only if EAM have the same number of points
466 IF (pot1%set(i)%eam%npoints == pot2%set(i)%eam%npoints) THEN
467 IF ((pot1%set(i)%eam%drar == pot2%set(i)%eam%drar) .AND. &
468 (pot1%set(i)%eam%drhoar == pot2%set(i)%eam%drhoar) .AND. &
469 (pot1%set(i)%eam%acutal == pot2%set(i)%eam%acutal) .AND. &
470 (sum(abs(pot1%set(i)%eam%rho - pot2%set(i)%eam%rho)) == 0.0_dp) .AND. &
471 (sum(abs(pot1%set(i)%eam%phi - pot2%set(i)%eam%phi)) == 0.0_dp) .AND. &
472 (sum(abs(pot1%set(i)%eam%frho - pot2%set(i)%eam%frho)) == 0.0_dp) .AND. &
473 (sum(abs(pot1%set(i)%eam%rhoval - pot2%set(i)%eam%rhoval)) == 0.0_dp) .AND. &
474 (sum(abs(pot1%set(i)%eam%rval - pot2%set(i)%eam%rval)) == 0.0_dp) .AND. &
475 (sum(abs(pot1%set(i)%eam%rhop - pot2%set(i)%eam%rhop)) == 0.0_dp) .AND. &
476 (sum(abs(pot1%set(i)%eam%phip - pot2%set(i)%eam%phip)) == 0.0_dp) .AND. &
477 (sum(abs(pot1%set(i)%eam%frhop - pot2%set(i)%eam%frhop)) == 0.0_dp)) mycompare = .true.
478 END IF
479 CASE (deepmd_type)
480 IF ((pot1%set(i)%deepmd%deepmd_file_name == pot2%set(i)%deepmd%deepmd_file_name) .AND. &
481 (pot1%set(i)%deepmd%atom_deepmd_type == pot2%set(i)%deepmd%atom_deepmd_type)) mycompare = .true.
482 CASE (quip_type)
483 IF ((pot1%set(i)%quip%quip_file_name == pot2%set(i)%quip%quip_file_name) .AND. &
484 (pot1%set(i)%quip%init_args == pot2%set(i)%quip%init_args) .AND. &
485 (pot1%set(i)%quip%calc_args == pot2%set(i)%quip%calc_args)) mycompare = .true.
486 CASE (nequip_type)
487 IF ((pot1%set(i)%nequip%nequip_file_name == pot2%set(i)%nequip%nequip_file_name) .AND. &
488 (pot1%set(i)%nequip%unit_coords == pot2%set(i)%nequip%unit_coords) .AND. &
489 (pot1%set(i)%nequip%unit_forces == pot2%set(i)%nequip%unit_forces) .AND. &
490 (pot1%set(i)%nequip%unit_energy == pot2%set(i)%nequip%unit_energy) .AND. &
491 (pot1%set(i)%nequip%unit_cell == pot2%set(i)%nequip%unit_cell)) mycompare = .true.
492 CASE (allegro_type)
493 IF ((pot1%set(i)%allegro%allegro_file_name == pot2%set(i)%allegro%allegro_file_name) .AND. &
494 (pot1%set(i)%allegro%unit_coords == pot2%set(i)%allegro%unit_coords) .AND. &
495 (pot1%set(i)%allegro%unit_forces == pot2%set(i)%allegro%unit_forces) .AND. &
496 (pot1%set(i)%allegro%unit_energy == pot2%set(i)%allegro%unit_energy) .AND. &
497 (pot1%set(i)%allegro%unit_cell == pot2%set(i)%allegro%unit_cell)) mycompare = .true.
498 CASE (ft_type)
499 IF ((pot1%set(i)%ft%A == pot2%set(i)%ft%A) .AND. &
500 (pot1%set(i)%ft%B == pot2%set(i)%ft%B) .AND. &
501 (pot1%set(i)%ft%C == pot2%set(i)%ft%C) .AND. &
502 (pot1%set(i)%ft%D == pot2%set(i)%ft%D)) mycompare = .true.
503 CASE (ftd_type)
504 IF ((pot1%set(i)%ftd%A == pot2%set(i)%ftd%A) .AND. &
505 (pot1%set(i)%ftd%B == pot2%set(i)%ftd%B) .AND. &
506 (pot1%set(i)%ftd%C == pot2%set(i)%ftd%C) .AND. &
507 (pot1%set(i)%ftd%D == pot2%set(i)%ftd%D) .AND. &
508 (all(pot1%set(i)%ftd%BD(:) == pot2%set(i)%ftd%BD(:)))) mycompare = .true.
509 CASE (ip_type)
510 IF ((sum(abs(pot1%set(i)%ipbv%a - pot2%set(i)%ipbv%a)) == 0.0_dp) .AND. &
511 (pot1%set(i)%ipbv%rcore == pot2%set(i)%ipbv%rcore) .AND. &
512 (pot1%set(i)%ipbv%m == pot2%set(i)%ipbv%m) .AND. &
513 (pot1%set(i)%ipbv%b == pot2%set(i)%ipbv%b)) mycompare = .true.
514 CASE (tersoff_type)
515 IF ((pot1%set(i)%tersoff%A == pot2%set(i)%tersoff%A) .AND. &
516 (pot1%set(i)%tersoff%B == pot2%set(i)%tersoff%B) .AND. &
517 (pot1%set(i)%tersoff%lambda1 == pot2%set(i)%tersoff%lambda1) .AND. &
518 (pot1%set(i)%tersoff%lambda2 == pot2%set(i)%tersoff%lambda2) .AND. &
519 (pot1%set(i)%tersoff%alpha == pot2%set(i)%tersoff%alpha) .AND. &
520 (pot1%set(i)%tersoff%beta == pot2%set(i)%tersoff%beta) .AND. &
521 (pot1%set(i)%tersoff%n == pot2%set(i)%tersoff%n) .AND. &
522 (pot1%set(i)%tersoff%c == pot2%set(i)%tersoff%c) .AND. &
523 (pot1%set(i)%tersoff%d == pot2%set(i)%tersoff%d) .AND. &
524 (pot1%set(i)%tersoff%h == pot2%set(i)%tersoff%h) .AND. &
525 (pot1%set(i)%tersoff%lambda3 == pot2%set(i)%tersoff%lambda3) .AND. &
526 (pot1%set(i)%tersoff%rcutsq == pot2%set(i)%tersoff%rcutsq) .AND. &
527 (pot1%set(i)%tersoff%bigR == pot2%set(i)%tersoff%bigR) .AND. &
528 (pot1%set(i)%tersoff%bigD == pot2%set(i)%tersoff%bigD)) mycompare = .true.
529 CASE (siepmann_type)
530 IF ((pot1%set(i)%siepmann%B == pot2%set(i)%siepmann%B) .AND. &
531 (pot1%set(i)%siepmann%D == pot2%set(i)%siepmann%D) .AND. &
532 (pot1%set(i)%siepmann%E == pot2%set(i)%siepmann%E) .AND. &
533 (pot1%set(i)%siepmann%F == pot2%set(i)%siepmann%F) .AND. &
534 (pot1%set(i)%siepmann%beta == pot2%set(i)%siepmann%beta) .AND. &
535 (pot1%set(i)%siepmann%rcutsq == pot2%set(i)%siepmann%rcutsq) .AND. &
536 (pot1%set(i)%siepmann%allow_oh_formation .EQV. &
537 pot2%set(i)%siepmann%allow_oh_formation) .AND. &
538 (pot1%set(i)%siepmann%allow_o_formation .EQV. &
539 pot2%set(i)%siepmann%allow_o_formation) .AND. &
540 (pot1%set(i)%siepmann%allow_h3o_formation .EQV. &
541 pot2%set(i)%siepmann%allow_h3o_formation)) mycompare = .true.
542 CASE (gal_type)
543 IF ((pot1%set(i)%gal%epsilon == pot2%set(i)%gal%epsilon) .AND. &
544 (pot1%set(i)%gal%bxy == pot2%set(i)%gal%bxy) .AND. &
545 (pot1%set(i)%gal%bz == pot2%set(i)%gal%bz) .AND. &
546 (pot1%set(i)%gal%r1 == pot2%set(i)%gal%r1) .AND. &
547 (pot1%set(i)%gal%r2 == pot2%set(i)%gal%r2) .AND. &
548 (pot1%set(i)%gal%a1 == pot2%set(i)%gal%a1) .AND. &
549 (pot1%set(i)%gal%a2 == pot2%set(i)%gal%a2) .AND. &
550 (pot1%set(i)%gal%a3 == pot2%set(i)%gal%a3) .AND. &
551 (pot1%set(i)%gal%a4 == pot2%set(i)%gal%a4) .AND. &
552 (pot1%set(i)%gal%a == pot2%set(i)%gal%a) .AND. &
553 (pot1%set(i)%gal%b == pot2%set(i)%gal%b) .AND. &
554 (pot1%set(i)%gal%c == pot2%set(i)%gal%c) .AND. &
555 (pot1%set(i)%gal%express .EQV. &
556 pot2%set(i)%gal%express) .AND. &
557 (pot1%set(i)%gal%rcutsq == pot2%set(i)%gal%rcutsq)) mycompare = .true.
558 CASE (gal21_type)
559 IF ((pot1%set(i)%gal21%epsilon1 == pot2%set(i)%gal21%epsilon1) .AND. &
560 (pot1%set(i)%gal21%epsilon2 == pot2%set(i)%gal21%epsilon2) .AND. &
561 (pot1%set(i)%gal21%epsilon3 == pot2%set(i)%gal21%epsilon3) .AND. &
562 (pot1%set(i)%gal21%bxy1 == pot2%set(i)%gal21%bxy1) .AND. &
563 (pot1%set(i)%gal21%bxy2 == pot2%set(i)%gal21%bxy1) .AND. &
564 (pot1%set(i)%gal21%bz1 == pot2%set(i)%gal21%bz1) .AND. &
565 (pot1%set(i)%gal21%bz2 == pot2%set(i)%gal21%bz2) .AND. &
566 (pot1%set(i)%gal21%r1 == pot2%set(i)%gal21%r1) .AND. &
567 (pot1%set(i)%gal21%r2 == pot2%set(i)%gal21%r2) .AND. &
568 (pot1%set(i)%gal21%a11 == pot2%set(i)%gal21%a11) .AND. &
569 (pot1%set(i)%gal21%a12 == pot2%set(i)%gal21%a12) .AND. &
570 (pot1%set(i)%gal21%a13 == pot2%set(i)%gal21%a13) .AND. &
571 (pot1%set(i)%gal21%a21 == pot2%set(i)%gal21%a21) .AND. &
572 (pot1%set(i)%gal21%a22 == pot2%set(i)%gal21%a22) .AND. &
573 (pot1%set(i)%gal21%a23 == pot2%set(i)%gal21%a23) .AND. &
574 (pot1%set(i)%gal21%a31 == pot2%set(i)%gal21%a31) .AND. &
575 (pot1%set(i)%gal21%a32 == pot2%set(i)%gal21%a32) .AND. &
576 (pot1%set(i)%gal21%a33 == pot2%set(i)%gal21%a33) .AND. &
577 (pot1%set(i)%gal21%a41 == pot2%set(i)%gal21%a41) .AND. &
578 (pot1%set(i)%gal21%a42 == pot2%set(i)%gal21%a42) .AND. &
579 (pot1%set(i)%gal21%a43 == pot2%set(i)%gal21%a43) .AND. &
580 (pot1%set(i)%gal21%AO1 == pot2%set(i)%gal21%AO1) .AND. &
581 (pot1%set(i)%gal21%AO2 == pot2%set(i)%gal21%AO2) .AND. &
582 (pot1%set(i)%gal21%BO1 == pot2%set(i)%gal21%BO1) .AND. &
583 (pot1%set(i)%gal21%BO2 == pot2%set(i)%gal21%BO2) .AND. &
584 (pot1%set(i)%gal21%c == pot2%set(i)%gal21%c) .AND. &
585 (pot1%set(i)%gal21%AH1 == pot2%set(i)%gal21%AH1) .AND. &
586 (pot1%set(i)%gal21%AH2 == pot2%set(i)%gal21%AH2) .AND. &
587 (pot1%set(i)%gal21%BH1 == pot2%set(i)%gal21%BH1) .AND. &
588 (pot1%set(i)%gal21%BH2 == pot2%set(i)%gal21%BH2) .AND. &
589 (pot1%set(i)%gal21%express .EQV. &
590 pot2%set(i)%gal21%express) .AND. &
591 (pot1%set(i)%gal21%rcutsq == pot2%set(i)%gal21%rcutsq)) mycompare = .true.
592
593 END SELECT
594 mycompare = mycompare .AND. &
595 (pot1%set(i)%rmin == pot2%set(i)%rmin) .AND. (pot1%set(i)%rmax == pot2%set(i)%rmax)
596 IF ((mycompare) .AND. (i == 1)) compare = .true.
597 compare = compare .AND. mycompare
598 END DO
599
600 END SUBROUTINE compare_pot
601
602! **************************************************************************************************
603!> \brief Creates the potential parameter type
604!> \param potparm ...
605!> \param nset ...
606!> \author Teodoro Laino [teo] 11.2005
607! **************************************************************************************************
608 SUBROUTINE pair_potential_single_create(potparm, nset)
609 TYPE(pair_potential_single_type), POINTER :: potparm
610 INTEGER, INTENT(IN), OPTIONAL :: nset
611
612 INTEGER :: i, lnset
613
614 cpassert(.NOT. ASSOCIATED(potparm))
615 ALLOCATE (potparm)
616 lnset = 1
617 IF (PRESENT(nset)) lnset = nset
618 ! Standard allocation to size 1
619 ALLOCATE (potparm%type(lnset))
620 ALLOCATE (potparm%set(lnset))
621 NULLIFY (potparm%spl_f, &
622 potparm%pair_spline_data)
623 DO i = 1, lnset
624 potparm%set(i)%rmin = not_initialized
625 potparm%set(i)%rmax = not_initialized
626 NULLIFY (potparm%set(i)%ipbv, &
627 potparm%set(i)%lj, &
628 potparm%set(i)%gp, &
629 potparm%set(i)%ft, &
630 potparm%set(i)%willis, &
631 potparm%set(i)%goodwin, &
632 potparm%set(i)%eam, &
633 potparm%set(i)%quip, &
634 potparm%set(i)%nequip, &
635 potparm%set(i)%allegro, &
636 potparm%set(i)%deepmd, &
637 potparm%set(i)%buck4r, &
638 potparm%set(i)%buckmo, &
639 potparm%set(i)%tersoff, &
640 potparm%set(i)%siepmann, &
641 potparm%set(i)%gal, &
642 potparm%set(i)%gal21, &
643 potparm%set(i)%ftd, &
644 potparm%set(i)%tab)
645 END DO
646 CALL pair_potential_single_clean(potparm)
647 END SUBROUTINE pair_potential_single_create
648
649! **************************************************************************************************
650!> \brief Cleans the potential parameter type
651!> \param potparm ...
652!> \author unknown
653! **************************************************************************************************
654 SUBROUTINE pair_potential_single_clean(potparm)
655 TYPE(pair_potential_single_type), POINTER :: potparm
656
657 INTEGER :: i
658
659 potparm%type = nn_type
660 potparm%shell_type = nosh_nosh
661 potparm%undef = .true.
662 potparm%no_pp = .false.
663 potparm%no_mb = .false.
664 potparm%at1 = 'NULL'
665 potparm%at2 = 'NULL'
666 potparm%rcutsq = 0.0_dp
667 IF (ASSOCIATED(potparm%pair_spline_data)) &
668 CALL spline_data_p_release(potparm%pair_spline_data)
669 IF (ASSOCIATED(potparm%spl_f)) &
670 CALL spline_factor_release(potparm%spl_f)
671
672 DO i = 1, SIZE(potparm%type)
673 potparm%set(i)%rmin = not_initialized
674 potparm%set(i)%rmax = not_initialized
675 CALL pair_potential_lj_clean(potparm%set(i)%lj)
676 CALL pair_potential_williams_clean(potparm%set(i)%willis)
677 CALL pair_potential_goodwin_clean(potparm%set(i)%goodwin)
678 CALL pair_potential_eam_clean(potparm%set(i)%eam)
679 CALL pair_potential_quip_clean(potparm%set(i)%quip)
680 CALL pair_potential_nequip_clean(potparm%set(i)%nequip)
681 CALL pair_potential_allegro_clean(potparm%set(i)%allegro)
682 CALL pair_potential_deepmd_clean(potparm%set(i)%deepmd)
683 CALL pair_potential_buck4r_clean(potparm%set(i)%buck4r)
684 CALL pair_potential_buckmo_clean(potparm%set(i)%buckmo)
685 CALL pair_potential_bmhft_clean(potparm%set(i)%ft)
686 CALL pair_potential_bmhftd_clean(potparm%set(i)%ftd)
687 CALL pair_potential_ipbv_clean(potparm%set(i)%ipbv)
688 CALL pair_potential_gp_clean(potparm%set(i)%gp)
689 CALL pair_potential_tersoff_clean(potparm%set(i)%tersoff)
690 CALL pair_potential_siepmann_clean(potparm%set(i)%siepmann)
691 CALL pair_potential_gal_clean(potparm%set(i)%gal)
692 CALL pair_potential_gal21_clean(potparm%set(i)%gal21)
693 CALL pair_potential_tab_clean(potparm%set(i)%tab)
694 END DO
695 END SUBROUTINE pair_potential_single_clean
696
697! **************************************************************************************************
698!> \brief Copy two potential parameter type
699!> \param potparm_source ...
700!> \param potparm_dest ...
701!> \author Teodoro Laino [teo] 11.2005
702! **************************************************************************************************
703 SUBROUTINE pair_potential_single_copy(potparm_source, potparm_dest)
704 TYPE(pair_potential_single_type), POINTER :: potparm_source, potparm_dest
705
706 INTEGER :: i
707
708 cpassert(ASSOCIATED(potparm_source))
709 IF (.NOT. ASSOCIATED(potparm_dest)) THEN
710 CALL pair_potential_single_create(potparm_dest, SIZE(potparm_source%type))
711 ELSE
712 CALL pair_potential_single_clean(potparm_dest)
713 END IF
714 potparm_dest%type = potparm_source%type
715 potparm_dest%shell_type = potparm_source%shell_type
716 potparm_dest%undef = potparm_source%undef
717 potparm_dest%no_mb = potparm_source%no_mb
718 potparm_dest%no_pp = potparm_source%no_pp
719 potparm_dest%at1 = potparm_source%at1
720 potparm_dest%at2 = potparm_source%at2
721 potparm_dest%rcutsq = potparm_source%rcutsq
722 IF (ASSOCIATED(potparm_source%pair_spline_data)) THEN
723 CALL spline_data_p_copy(potparm_source%pair_spline_data, potparm_dest%pair_spline_data)
724 END IF
725
726 IF (ASSOCIATED(potparm_source%spl_f)) THEN
727 CALL spline_factor_copy(potparm_source%spl_f, potparm_dest%spl_f)
728 END IF
729
730 DO i = 1, SIZE(potparm_source%type)
731 potparm_dest%set(i)%rmin = potparm_source%set(i)%rmin
732 potparm_dest%set(i)%rmax = potparm_source%set(i)%rmax
733 CALL pair_potential_lj_copy(potparm_source%set(i)%lj, potparm_dest%set(i)%lj)
734 CALL pair_potential_williams_copy(potparm_source%set(i)%willis, potparm_dest%set(i)%willis)
735 CALL pair_potential_goodwin_copy(potparm_source%set(i)%goodwin, potparm_dest%set(i)%goodwin)
736 CALL pair_potential_eam_copy(potparm_source%set(i)%eam, potparm_dest%set(i)%eam)
737 CALL pair_potential_quip_copy(potparm_source%set(i)%quip, potparm_dest%set(i)%quip)
738 CALL pair_potential_nequip_copy(potparm_source%set(i)%nequip, potparm_dest%set(i)%nequip)
739 CALL pair_potential_allegro_copy(potparm_source%set(i)%allegro, potparm_dest%set(i)%allegro)
740 CALL pair_potential_deepmd_copy(potparm_source%set(i)%deepmd, potparm_dest%set(i)%deepmd)
741 CALL pair_potential_bmhft_copy(potparm_source%set(i)%ft, potparm_dest%set(i)%ft)
742 CALL pair_potential_bmhftd_copy(potparm_source%set(i)%ftd, potparm_dest%set(i)%ftd)
743 CALL pair_potential_ipbv_copy(potparm_source%set(i)%ipbv, potparm_dest%set(i)%ipbv)
744 CALL pair_potential_buck4r_copy(potparm_source%set(i)%buck4r, potparm_dest%set(i)%buck4r)
745 CALL pair_potential_buckmo_copy(potparm_source%set(i)%buckmo, potparm_dest%set(i)%buckmo)
746 CALL pair_potential_gp_copy(potparm_source%set(i)%gp, potparm_dest%set(i)%gp)
747 CALL pair_potential_tersoff_copy(potparm_source%set(i)%tersoff, potparm_dest%set(i)%tersoff)
748 CALL pair_potential_siepmann_copy(potparm_source%set(i)%siepmann, potparm_dest%set(i)%siepmann)
749 CALL pair_potential_gal_copy(potparm_source%set(i)%gal, potparm_dest%set(i)%gal)
750 CALL pair_potential_gal21_copy(potparm_source%set(i)%gal21, potparm_dest%set(i)%gal21)
751 CALL pair_potential_tab_copy(potparm_source%set(i)%tab, potparm_dest%set(i)%tab)
752 END DO
753 END SUBROUTINE pair_potential_single_copy
754
755! **************************************************************************************************
756!> \brief Add potential parameter type to an existing potential parameter type
757!> Used in case of multiple_potential definition
758!> \param potparm_source ...
759!> \param potparm_dest ...
760!> \author Teodoro Laino [teo] 11.2005
761! **************************************************************************************************
762 SUBROUTINE pair_potential_single_add(potparm_source, potparm_dest)
763 TYPE(pair_potential_single_type), POINTER :: potparm_source, potparm_dest
764
765 INTEGER :: i, j, size_dest, size_source
766 LOGICAL :: allocate_new, check
767 TYPE(pair_potential_single_type), POINTER :: potparm_tmp
768
769 cpassert(ASSOCIATED(potparm_source))
770 ! At this level we expect all splines types
771 ! be not allocated.. No sense add splines at this level.. in case fail!
772 check = (.NOT. ASSOCIATED(potparm_source%pair_spline_data)) .AND. &
773 (.NOT. ASSOCIATED(potparm_source%spl_f))
774 cpassert(check)
775 check = (.NOT. ASSOCIATED(potparm_dest%pair_spline_data)) .AND. &
776 (.NOT. ASSOCIATED(potparm_dest%spl_f))
777 cpassert(check)
778 ! Increase the size of the destination potparm (in case) and copy the new data
779 size_source = SIZE(potparm_source%type)
780 allocate_new = .NOT. ASSOCIATED(potparm_dest)
781 IF (.NOT. allocate_new) THEN
782 size_dest = SIZE(potparm_dest%type)
783 IF (size_dest == 1) THEN
784 check = (ASSOCIATED(potparm_dest%set(1)%lj)) .OR. &
785 (ASSOCIATED(potparm_dest%set(1)%willis)) .OR. &
786 (ASSOCIATED(potparm_dest%set(1)%goodwin)) .OR. &
787 (ASSOCIATED(potparm_dest%set(1)%eam)) .OR. &
788 (ASSOCIATED(potparm_dest%set(1)%quip)) .OR. &
789 (ASSOCIATED(potparm_dest%set(1)%nequip)) .OR. &
790 (ASSOCIATED(potparm_dest%set(1)%allegro)) .OR. &
791 (ASSOCIATED(potparm_dest%set(1)%deepmd)) .OR. &
792 (ASSOCIATED(potparm_dest%set(1)%ft)) .OR. &
793 (ASSOCIATED(potparm_dest%set(1)%ftd)) .OR. &
794 (ASSOCIATED(potparm_dest%set(1)%ipbv)) .OR. &
795 (ASSOCIATED(potparm_dest%set(1)%buck4r)) .OR. &
796 (ASSOCIATED(potparm_dest%set(1)%buckmo)) .OR. &
797 (ASSOCIATED(potparm_dest%set(1)%gp)) .OR. &
798 (ASSOCIATED(potparm_dest%set(1)%tersoff)) .OR. &
799 (ASSOCIATED(potparm_dest%set(1)%siepmann)) .OR. &
800 (ASSOCIATED(potparm_dest%set(1)%gal)) .OR. &
801 (ASSOCIATED(potparm_dest%set(1)%gal)) .OR. &
802 (ASSOCIATED(potparm_dest%set(1)%tab))
803 IF (.NOT. check) THEN
804 allocate_new = .true.
805 CALL pair_potential_single_release(potparm_dest)
806 END IF
807 END IF
808 END IF
809 IF (allocate_new) THEN
810 size_dest = 0
811 CALL pair_potential_single_create(potparm_dest, size_source)
812 potparm_dest%shell_type = potparm_source%shell_type
813 potparm_dest%undef = potparm_source%undef
814 potparm_dest%no_mb = potparm_source%no_mb
815 potparm_dest%no_pp = potparm_source%no_pp
816 potparm_dest%at1 = potparm_source%at1
817 potparm_dest%at2 = potparm_source%at2
818 potparm_dest%rcutsq = potparm_source%rcutsq
819 ELSE
820 size_dest = SIZE(potparm_dest%type)
821 NULLIFY (potparm_tmp)
822 CALL pair_potential_single_copy(potparm_dest, potparm_tmp)
823 CALL pair_potential_single_release(potparm_dest)
824 CALL pair_potential_single_create(potparm_dest, size_dest + size_source)
825 ! Copy back original informations..
826 potparm_dest%shell_type = potparm_tmp%shell_type
827 potparm_dest%undef = potparm_tmp%undef
828 potparm_dest%no_mb = potparm_tmp%no_mb
829 potparm_dest%no_pp = potparm_tmp%no_pp
830 potparm_dest%at1 = potparm_tmp%at1
831 potparm_dest%at2 = potparm_tmp%at2
832 potparm_dest%rcutsq = potparm_tmp%rcutsq
833 DO i = 1, size_dest
834 potparm_dest%type(i) = potparm_tmp%type(i)
835 potparm_dest%set(i)%rmin = potparm_tmp%set(i)%rmin
836 potparm_dest%set(i)%rmax = potparm_tmp%set(i)%rmax
837 CALL pair_potential_lj_copy(potparm_tmp%set(i)%lj, potparm_dest%set(i)%lj)
838 CALL pair_potential_williams_copy(potparm_tmp%set(i)%willis, potparm_dest%set(i)%willis)
839 CALL pair_potential_goodwin_copy(potparm_tmp%set(i)%goodwin, potparm_dest%set(i)%goodwin)
840 CALL pair_potential_eam_copy(potparm_tmp%set(i)%eam, potparm_dest%set(i)%eam)
841 CALL pair_potential_quip_copy(potparm_tmp%set(i)%quip, potparm_dest%set(i)%quip)
842 CALL pair_potential_nequip_copy(potparm_tmp%set(i)%nequip, potparm_dest%set(i)%nequip)
843 CALL pair_potential_allegro_copy(potparm_tmp%set(i)%allegro, potparm_dest%set(i)%allegro)
844 CALL pair_potential_deepmd_copy(potparm_tmp%set(i)%deepmd, potparm_dest%set(i)%deepmd)
845 CALL pair_potential_bmhft_copy(potparm_tmp%set(i)%ft, potparm_dest%set(i)%ft)
846 CALL pair_potential_bmhftd_copy(potparm_tmp%set(i)%ftd, potparm_dest%set(i)%ftd)
847 CALL pair_potential_ipbv_copy(potparm_tmp%set(i)%ipbv, potparm_dest%set(i)%ipbv)
848 CALL pair_potential_buck4r_copy(potparm_tmp%set(i)%buck4r, potparm_dest%set(i)%buck4r)
849 CALL pair_potential_buckmo_copy(potparm_tmp%set(i)%buckmo, potparm_dest%set(i)%buckmo)
850 CALL pair_potential_gp_copy(potparm_tmp%set(i)%gp, potparm_dest%set(i)%gp)
851 CALL pair_potential_tersoff_copy(potparm_tmp%set(i)%tersoff, potparm_dest%set(i)%tersoff)
852 CALL pair_potential_siepmann_copy(potparm_tmp%set(i)%siepmann, potparm_dest%set(i)%siepmann)
853 CALL pair_potential_gal_copy(potparm_tmp%set(i)%gal, potparm_dest%set(i)%gal)
854 CALL pair_potential_gal21_copy(potparm_tmp%set(i)%gal21, potparm_dest%set(i)%gal21)
855 CALL pair_potential_tab_copy(potparm_tmp%set(i)%tab, potparm_dest%set(i)%tab)
856 END DO
857 CALL pair_potential_single_release(potparm_tmp)
858 END IF
859 ! Further check with main option with source and dest (already filled with few informations)
860 check = (potparm_dest%shell_type == potparm_source%shell_type) .AND. &
861 (potparm_dest%undef .EQV. potparm_source%undef) .AND. &
862 (potparm_dest%no_mb .EQV. potparm_source%no_mb) .AND. &
863 (potparm_dest%no_pp .EQV. potparm_source%no_pp) .AND. &
864 (potparm_dest%at1 == potparm_source%at1) .AND. &
865 (potparm_dest%at2 == potparm_source%at2) .AND. &
866 (potparm_dest%rcutsq == potparm_source%rcutsq)
867 cpassert(check)
868 ! Now copy the new pair_potential type
869 DO i = size_dest + 1, size_dest + size_source
870 j = i - size_dest
871 potparm_dest%type(i) = potparm_source%type(j)
872 potparm_dest%set(i)%rmin = potparm_source%set(j)%rmin
873 potparm_dest%set(i)%rmax = potparm_source%set(j)%rmax
874 CALL pair_potential_lj_copy(potparm_source%set(j)%lj, potparm_dest%set(i)%lj)
875 CALL pair_potential_williams_copy(potparm_source%set(j)%willis, potparm_dest%set(i)%willis)
876 CALL pair_potential_goodwin_copy(potparm_source%set(j)%goodwin, potparm_dest%set(i)%goodwin)
877 CALL pair_potential_eam_copy(potparm_source%set(j)%eam, potparm_dest%set(i)%eam)
878 CALL pair_potential_quip_copy(potparm_source%set(j)%quip, potparm_dest%set(i)%quip)
879 CALL pair_potential_nequip_copy(potparm_source%set(j)%nequip, potparm_dest%set(i)%nequip)
880 CALL pair_potential_allegro_copy(potparm_source%set(j)%allegro, potparm_dest%set(i)%allegro)
881 CALL pair_potential_deepmd_copy(potparm_source%set(j)%deepmd, potparm_dest%set(i)%deepmd)
882 CALL pair_potential_bmhft_copy(potparm_source%set(j)%ft, potparm_dest%set(i)%ft)
883 CALL pair_potential_bmhftd_copy(potparm_source%set(j)%ftd, potparm_dest%set(i)%ftd)
884 CALL pair_potential_ipbv_copy(potparm_source%set(j)%ipbv, potparm_dest%set(i)%ipbv)
885 CALL pair_potential_buck4r_copy(potparm_source%set(j)%buck4r, potparm_dest%set(i)%buck4r)
886 CALL pair_potential_buckmo_copy(potparm_source%set(j)%buckmo, potparm_dest%set(i)%buckmo)
887 CALL pair_potential_gp_copy(potparm_source%set(j)%gp, potparm_dest%set(i)%gp)
888 CALL pair_potential_tersoff_copy(potparm_source%set(j)%tersoff, potparm_dest%set(i)%tersoff)
889 CALL pair_potential_siepmann_copy(potparm_source%set(j)%siepmann, potparm_dest%set(i)%siepmann)
890 CALL pair_potential_gal_copy(potparm_source%set(j)%gal, potparm_dest%set(i)%gal)
891 CALL pair_potential_gal21_copy(potparm_source%set(j)%gal21, potparm_dest%set(i)%gal21)
892 CALL pair_potential_tab_copy(potparm_source%set(j)%tab, potparm_dest%set(i)%tab)
893 END DO
894 END SUBROUTINE pair_potential_single_add
895
896! **************************************************************************************************
897!> \brief Release Data-structure that constains potential parameters of a single pair
898!> \param potparm ...
899!> \author Teodoro Laino [Teo] 11.2005
900! **************************************************************************************************
901 SUBROUTINE pair_potential_single_release(potparm)
902 TYPE(pair_potential_single_type), POINTER :: potparm
903
904 INTEGER :: i
905
906 cpassert(ASSOCIATED(potparm))
907 CALL spline_data_p_release(potparm%pair_spline_data)
908 CALL spline_factor_release(potparm%spl_f)
909 DO i = 1, SIZE(potparm%type)
910 CALL pair_potential_ipbv_release(potparm%set(i)%ipbv)
911 CALL pair_potential_lj_release(potparm%set(i)%lj)
912 CALL pair_potential_bmhft_release(potparm%set(i)%ft)
913 CALL pair_potential_bmhftd_release(potparm%set(i)%ftd)
914 CALL pair_potential_williams_release(potparm%set(i)%willis)
915 CALL pair_potential_goodwin_release(potparm%set(i)%goodwin)
916 CALL pair_potential_eam_release(potparm%set(i)%eam)
917 CALL pair_potential_quip_release(potparm%set(i)%quip)
918 CALL pair_potential_nequip_release(potparm%set(i)%nequip)
919 CALL pair_potential_allegro_release(potparm%set(i)%allegro)
920 CALL pair_potential_deepmd_release(potparm%set(i)%deepmd)
921 CALL pair_potential_buck4r_release(potparm%set(i)%buck4r)
922 CALL pair_potential_buckmo_release(potparm%set(i)%buckmo)
923 CALL pair_potential_gp_release(potparm%set(i)%gp)
924 CALL pair_potential_tersoff_release(potparm%set(i)%tersoff)
925 CALL pair_potential_siepmann_release(potparm%set(i)%siepmann)
926 CALL pair_potential_gal_release(potparm%set(i)%gal)
927 CALL pair_potential_gal21_release(potparm%set(i)%gal21)
928 CALL pair_potential_tab_release(potparm%set(i)%tab)
929 END DO
930 DEALLOCATE (potparm%type)
931 DEALLOCATE (potparm%set)
932 DEALLOCATE (potparm)
933 END SUBROUTINE pair_potential_single_release
934
935! **************************************************************************************************
936!> \brief Data-structure that constains potential parameters
937!> \param potparm ...
938!> \param nkinds ...
939!> \author unknown
940! **************************************************************************************************
941 SUBROUTINE pair_potential_pp_create(potparm, nkinds)
942 TYPE(pair_potential_pp_type), POINTER :: potparm
943 INTEGER, INTENT(IN) :: nkinds
944
945 INTEGER :: i, j
946
947 cpassert(.NOT. ASSOCIATED(potparm))
948 ALLOCATE (potparm)
949 ALLOCATE (potparm%pot(nkinds, nkinds))
950 DO i = 1, nkinds
951 DO j = 1, nkinds
952 NULLIFY (potparm%pot(i, j)%pot)
953 END DO
954 END DO
955 ! Use no-redundancy in the potential definition
956 DO i = 1, nkinds
957 DO j = i, nkinds
958 CALL pair_potential_single_create(potparm%pot(i, j)%pot)
959 potparm%pot(j, i)%pot => potparm%pot(i, j)%pot
960 END DO
961 END DO
962 END SUBROUTINE pair_potential_pp_create
963
964! **************************************************************************************************
965!> \brief Release Data-structure that constains potential parameters
966!> \param potparm ...
967!> \par History
968!> Teodoro Laino [Teo] 11.2005 : Reorganizing the structures to optimize
969!> memory management
970!> \author unknown
971! **************************************************************************************************
972 SUBROUTINE pair_potential_pp_release(potparm)
973 TYPE(pair_potential_pp_type), POINTER :: potparm
974
975 INTEGER :: i, j
976
977 IF (ASSOCIATED(potparm)) THEN
978 IF (ASSOCIATED(potparm%pot)) THEN
979 DO i = 1, SIZE(potparm%pot, 1)
980 DO j = i, SIZE(potparm%pot, 2)
981 CALL pair_potential_single_release(potparm%pot(i, j)%pot)
982 NULLIFY (potparm%pot(j, i)%pot)
983 END DO
984 END DO
985 DEALLOCATE (potparm%pot)
986 END IF
987 DEALLOCATE (potparm)
988 END IF
989 NULLIFY (potparm)
990 END SUBROUTINE pair_potential_pp_release
991
992! **************************************************************************************************
993!> \brief Data-structure that constains potential parameters
994!> \param potparm ...
995!> \param ndim ...
996!> \param ub ...
997!> \param lb ...
998!> \author unknown
999! **************************************************************************************************
1000 SUBROUTINE pair_potential_p_create(potparm, ndim, ub, lb)
1001 TYPE(pair_potential_p_type), POINTER :: potparm
1002 INTEGER, INTENT(IN), OPTIONAL :: ndim, ub, lb
1003
1004 INTEGER :: i, loc_lb, loc_ub
1005
1006 cpassert(.NOT. ASSOCIATED(potparm))
1007 ALLOCATE (potparm)
1008 IF (PRESENT(ndim)) THEN
1009 loc_lb = 1
1010 loc_ub = ndim
1011 ALLOCATE (potparm%pot(loc_lb:loc_ub))
1012 IF (PRESENT(lb) .OR. PRESENT(ub)) THEN
1013 cpabort("")
1014 END IF
1015 ELSE IF (PRESENT(lb) .AND. PRESENT(ub)) THEN
1016 loc_lb = lb
1017 loc_ub = ub
1018 ALLOCATE (potparm%pot(loc_lb:loc_ub))
1019 IF (PRESENT(ndim)) THEN
1020 cpabort("")
1021 END IF
1022 ELSE
1023 cpabort("")
1024 END IF
1025 DO i = loc_lb, loc_ub
1026 NULLIFY (potparm%pot(i)%pot)
1027 CALL pair_potential_single_create(potparm%pot(i)%pot)
1028 END DO
1029 END SUBROUTINE pair_potential_p_create
1030
1031! **************************************************************************************************
1032!> \brief Release Data-structure that constains potential parameters
1033!> \param potparm ...
1034!> \par History
1035!> Teodoro Laino [Teo] 11.2005 : Reorganizing the structures to optimize
1036!> memory management
1037!> \author unknown
1038! **************************************************************************************************
1039 SUBROUTINE pair_potential_p_release(potparm)
1040 TYPE(pair_potential_p_type), POINTER :: potparm
1041
1042 INTEGER :: i
1043
1044 IF (ASSOCIATED(potparm)) THEN
1045 IF (ASSOCIATED(potparm%pot)) THEN
1046 DO i = 1, SIZE(potparm%pot)
1047 CALL pair_potential_single_release(potparm%pot(i)%pot)
1048 END DO
1049 DEALLOCATE (potparm%pot)
1050 END IF
1051 DEALLOCATE (potparm)
1052 END IF
1053 NULLIFY (potparm)
1054 END SUBROUTINE pair_potential_p_release
1055
1056! **************************************************************************************************
1057!> \brief Copy structures between two pair_potential_p_type
1058!> \param source ...
1059!> \param dest ...
1060!> \param istart ...
1061!> \param iend ...
1062!> \author Teodoro Laino [Teo] 11.2005
1063! **************************************************************************************************
1064 SUBROUTINE pair_potential_p_copy(source, dest, istart, iend)
1065 TYPE(pair_potential_p_type), POINTER :: source, dest
1066 INTEGER, INTENT(IN), OPTIONAL :: istart, iend
1067
1068 INTEGER :: i, l_end, l_start
1069
1070 cpassert(ASSOCIATED(source))
1071 cpassert(ASSOCIATED(dest))
1072 l_start = lbound(source%pot, 1)
1073 l_end = ubound(source%pot, 1)
1074 IF (PRESENT(istart)) l_start = istart
1075 IF (PRESENT(iend)) l_end = iend
1076 DO i = l_start, l_end
1077 IF (.NOT. ASSOCIATED(source%pot(i)%pot)) &
1078 CALL pair_potential_single_create(source%pot(i)%pot)
1079 CALL pair_potential_single_copy(source%pot(i)%pot, dest%pot(i)%pot)
1080 END DO
1081 END SUBROUTINE pair_potential_p_copy
1082
1083! **************************************************************************************************
1084!> \brief Cleans the potential parameter type
1085!> \param p ...
1086!> \param lb1_new ...
1087!> \param ub1_new ...
1088!> \param lj ...
1089!> \param lj_charmm ...
1090!> \param williams ...
1091!> \param goodwin ...
1092!> \param eam ...
1093!> \param quip ...
1094!> \param nequip ...
1095!> \param allegro ...
1096!> \param bmhft ...
1097!> \param bmhftd ...
1098!> \param ipbv ...
1099!> \param buck4r ...
1100!> \param buckmo ...
1101!> \param gp ...
1102!> \param tersoff ...
1103!> \param siepmann ...
1104!> \param gal ...
1105!> \param gal21 ...
1106!> \param tab ...
1107!> \param deepmd ...
1108!> \author Teodoro Laino [Teo] 11.2005
1109! **************************************************************************************************
1110 SUBROUTINE pair_potential_reallocate(p, lb1_new, ub1_new, lj, lj_charmm, williams, goodwin, eam, &
1111 quip, nequip, allegro, bmhft, bmhftd, ipbv, buck4r, buckmo, &
1112 gp, tersoff, siepmann, gal, gal21, tab, deepmd)
1113 TYPE(pair_potential_p_type), POINTER :: p
1114 INTEGER, INTENT(IN) :: lb1_new, ub1_new
1115 LOGICAL, INTENT(IN), OPTIONAL :: lj, lj_charmm, williams, goodwin, eam, quip, nequip, &
1116 allegro, bmhft, bmhftd, ipbv, buck4r, buckmo, gp, tersoff, siepmann, gal, gal21, tab, &
1117 deepmd
1118
1119 INTEGER :: i, ipot, lb1_old, std_dim, ub1_old
1120 LOGICAL :: check, lallegro, lbmhft, lbmhftd, lbuck4r, lbuckmo, ldeepmd, leam, lgal, lgal21, &
1121 lgoodwin, lgp, lipbv, llj, llj_charmm, lnequip, lquip, lsiepmann, ltab, ltersoff, &
1122 lwilliams
1123 TYPE(pair_potential_p_type), POINTER :: work
1124
1125 NULLIFY (work)
1126 ipot = 0
1127 llj = .false.; IF (PRESENT(lj)) llj = lj
1128 llj_charmm = .false.; IF (PRESENT(lj_charmm)) llj_charmm = lj_charmm
1129 lwilliams = .false.; IF (PRESENT(williams)) lwilliams = williams
1130 lgoodwin = .false.; IF (PRESENT(goodwin)) lgoodwin = goodwin
1131 leam = .false.; IF (PRESENT(eam)) leam = eam
1132 lquip = .false.; IF (PRESENT(quip)) lquip = quip
1133 lnequip = .false.; IF (PRESENT(nequip)) lnequip = nequip
1134 lallegro = .false.; IF (PRESENT(allegro)) lallegro = allegro
1135 ldeepmd = .false.; IF (PRESENT(deepmd)) ldeepmd = deepmd
1136 lbmhft = .false.; IF (PRESENT(bmhft)) lbmhft = bmhft
1137 lbmhftd = .false.; IF (PRESENT(bmhftd)) lbmhftd = bmhftd
1138 lipbv = .false.; IF (PRESENT(ipbv)) lipbv = ipbv
1139 lbuck4r = .false.; IF (PRESENT(buck4r)) lbuck4r = buck4r
1140 lbuckmo = .false.; IF (PRESENT(buckmo)) lbuckmo = buckmo
1141 lgp = .false.; IF (PRESENT(gp)) lgp = gp
1142 ltersoff = .false.; IF (PRESENT(tersoff)) ltersoff = tersoff
1143 lsiepmann = .false.; IF (PRESENT(siepmann)) lsiepmann = siepmann
1144 lgal = .false.; IF (PRESENT(gal)) lgal = gal
1145 lgal21 = .false.; IF (PRESENT(gal21)) lgal21 = gal21
1146 ltab = .false.; IF (PRESENT(tab)) ltab = tab
1147
1148 IF (llj) THEN
1149 ipot = lj_type
1150 check = .NOT. (llj_charmm .OR. lwilliams .OR. lgoodwin .OR. leam .OR. lquip .OR. lnequip .OR. lallegro &
1151 .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lbuckmo .OR. lgp .OR. ltersoff &
1152 .OR. lsiepmann .OR. lgal .OR. lgal21 .OR. ltab .OR. ldeepmd)
1153 cpassert(check)
1154 END IF
1155 IF (llj_charmm) THEN
1156 ipot = lj_charmm_type
1157 check = .NOT. (llj .OR. lwilliams .OR. lgoodwin .OR. leam .OR. lquip .OR. lnequip .OR. lallegro &
1158 .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lbuckmo .OR. lgp .OR. ltersoff &
1159 .OR. lsiepmann .OR. lgal .OR. lgal21 .OR. ltab .OR. ldeepmd)
1160 cpassert(check)
1161 END IF
1162 IF (lwilliams) THEN
1163 ipot = wl_type
1164 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. leam .OR. lquip .OR. lnequip .OR. lallegro &
1165 .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lbuckmo .OR. lgp .OR. ltersoff &
1166 .OR. lsiepmann .OR. lgal .OR. lgal21 .OR. ltab .OR. ldeepmd)
1167 cpassert(check)
1168 END IF
1169 IF (lgoodwin) THEN
1170 ipot = gw_type
1171 check = .NOT. (llj .OR. llj_charmm .OR. lwilliams .OR. leam .OR. lquip .OR. lnequip .OR. lallegro &
1172 .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lbuckmo .OR. lgp .OR. ltersoff &
1173 .OR. lsiepmann .OR. lgal .OR. lgal21 .OR. ltab .OR. ldeepmd)
1174 cpassert(check)
1175 END IF
1176 IF (leam) THEN
1177 ipot = ea_type
1178 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. lquip .OR. lnequip .OR. lallegro &
1179 .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lbuckmo .OR. lgp .OR. ltersoff &
1180 .OR. lsiepmann .OR. lgal .OR. lgal21 .OR. ltab .OR. ldeepmd)
1181 cpassert(check)
1182 END IF
1183 IF (lquip) THEN
1184 ipot = quip_type
1185 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. leam .OR. lnequip .OR. lallegro &
1186 .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lbuckmo .OR. lgp .OR. ltersoff &
1187 .OR. lsiepmann .OR. lgal .OR. lgal21 .OR. ltab .OR. ldeepmd)
1188 cpassert(check)
1189 END IF
1190 IF (lnequip) THEN
1191 ipot = nequip_type
1192 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. leam .OR. lquip .OR. lallegro &
1193 .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lbuckmo .OR. lgp .OR. ltersoff &
1194 .OR. lsiepmann .OR. lgal .OR. lgal21 .OR. ltab .OR. ldeepmd)
1195 cpassert(check)
1196 END IF
1197 IF (lallegro) THEN
1198 ipot = allegro_type
1199 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. leam .OR. lquip .OR. lnequip &
1200 .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lbuckmo .OR. lgp .OR. ltersoff &
1201 .OR. lsiepmann .OR. lgal .OR. lgal21 .OR. ltab .OR. ldeepmd)
1202 cpassert(check)
1203 END IF
1204 IF (ldeepmd) THEN
1205 ipot = deepmd_type
1206 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. leam .OR. lquip .OR. lnequip &
1207 .OR. lallegro .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lbuckmo .OR. lgp &
1208 .OR. ltersoff .OR. lsiepmann .OR. lgal .OR. lgal21 .OR. ltab)
1209 cpassert(check)
1210 END IF
1211 IF (lbmhft) THEN
1212 ipot = ft_type
1213 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. leam .OR. lquip .OR. lnequip &
1214 .OR. lallegro .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lbuckmo .OR. lgp .OR. ltersoff &
1215 .OR. lsiepmann .OR. lgal .OR. lgal21 .OR. ltab .OR. ldeepmd)
1216 cpassert(check)
1217 END IF
1218 IF (lbmhftd) THEN
1219 ipot = ftd_type
1220 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. leam .OR. lquip .OR. lnequip &
1221 .OR. lallegro .OR. lbmhft .OR. lipbv .OR. lbuck4r .OR. lbuckmo .OR. lgp .OR. ltersoff &
1222 .OR. lsiepmann .OR. lgal .OR. lgal21 .OR. ltab .OR. ldeepmd)
1223 cpassert(check)
1224 END IF
1225 IF (lipbv) THEN
1226 ipot = ip_type
1227 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. leam .OR. lquip .OR. lnequip &
1228 .OR. lallegro .OR. lbmhft .OR. lbmhftd .OR. lbuck4r .OR. lbuckmo .OR. lgp .OR. ltersoff &
1229 .OR. lsiepmann .OR. lgal .OR. lgal21 .OR. ltab .OR. ldeepmd)
1230 cpassert(check)
1231 END IF
1232 IF (lbuck4r) THEN
1233 ipot = b4_type
1234 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. leam .OR. lquip .OR. lnequip &
1235 .OR. lallegro .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuckmo .OR. lgp .OR. ltersoff &
1236 .OR. lsiepmann .OR. lgal .OR. lgal21 .OR. ltab .OR. ldeepmd)
1237 cpassert(check)
1238 END IF
1239 IF (lbuckmo) THEN
1240 ipot = bm_type
1241 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. leam .OR. lquip .OR. lnequip &
1242 .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lgp .OR. ltersoff &
1243 .OR. lsiepmann .OR. lgal .OR. lgal21 .OR. ltab .OR. ldeepmd)
1244 cpassert(check)
1245 END IF
1246 IF (ltersoff) THEN
1247 ipot = tersoff_type
1248 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. leam .OR. lquip .OR. lnequip &
1249 .OR. lallegro .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lgp .OR. lbuckmo &
1250 .OR. lsiepmann .OR. lgal .OR. lgal21 .OR. ltab .OR. ldeepmd)
1251 cpassert(check)
1252 END IF
1253 IF (lsiepmann) THEN
1254 ipot = siepmann_type
1255 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. leam .OR. lquip .OR. lnequip &
1256 .OR. lallegro .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lgp .OR. lbuckmo &
1257 .OR. ltersoff .OR. lgal .OR. lgal21 .OR. ltab .OR. ldeepmd)
1258 cpassert(check)
1259 END IF
1260 IF (lgal) THEN
1261 ipot = gal_type
1262 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. leam .OR. lquip .OR. lnequip &
1263 .OR. lallegro .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lgp .OR. lbuckmo &
1264 .OR. ltersoff .OR. lsiepmann .OR. lgal21 .OR. ltab .OR. ldeepmd)
1265 cpassert(check)
1266 END IF
1267 IF (lgal21) THEN
1268 ipot = gal21_type
1269 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. leam .OR. lquip .OR. lnequip &
1270 .OR. lallegro .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lgp .OR. lbuckmo &
1271 .OR. ltersoff .OR. lsiepmann .OR. lgal .OR. ltab .OR. ldeepmd)
1272 cpassert(check)
1273 END IF
1274 IF (lgp) THEN
1275 ipot = gp_type
1276 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. leam .OR. lquip .OR. lnequip &
1277 .OR. lallegro .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lgal21 .OR. lbuckmo &
1278 .OR. ltersoff .OR. lsiepmann .OR. lgal .OR. ltab .OR. ldeepmd)
1279 cpassert(check)
1280 END IF
1281 IF (ltab) THEN
1282 ipot = tab_type
1283 check = .NOT. (llj .OR. llj_charmm .OR. lgoodwin .OR. lwilliams .OR. leam .OR. lquip .OR. lnequip &
1284 .OR. lallegro .OR. lbmhft .OR. lbmhftd .OR. lipbv .OR. lbuck4r .OR. lgp .OR. lgal21 &
1285 .OR. lbuckmo .OR. ltersoff .OR. lsiepmann .OR. lgal)
1286 cpassert(check)
1287 END IF
1288
1289 lb1_old = 0
1290 ub1_old = 0
1291 IF (ASSOCIATED(p)) THEN
1292 lb1_old = lbound(p%pot, 1)
1293 ub1_old = ubound(p%pot, 1)
1294 CALL pair_potential_p_create(work, lb=lb1_old, ub=ub1_old)
1295 CALL pair_potential_p_copy(p, work)
1297 END IF
1298
1299 CALL pair_potential_p_create(p, lb=lb1_new, ub=ub1_new)
1300 IF (ASSOCIATED(work)) THEN
1301 CALL pair_potential_p_copy(work, p, istart=lb1_old, iend=ub1_old)
1302 END IF
1303 std_dim = 1
1304 DO i = ub1_old + 1, ub1_new
1305 check = (SIZE(p%pot(i)%pot%type) == std_dim) .AND. (SIZE(p%pot(i)%pot%type) == std_dim)
1306 cpassert(check)
1307 p%pot(i)%pot%type = nn_type
1308 p%pot(i)%pot%shell_type = nosh_nosh
1309 p%pot(i)%pot%undef = .true.
1310 p%pot(i)%pot%no_mb = .false.
1311 p%pot(i)%pot%no_pp = .false.
1312 p%pot(i)%pot%at1 = 'NULL'
1313 p%pot(i)%pot%at2 = 'NULL'
1314 p%pot(i)%pot%set(std_dim)%rmin = not_initialized
1315 p%pot(i)%pot%set(std_dim)%rmax = not_initialized
1316 SELECT CASE (ipot)
1317 CASE (lj_type, lj_charmm_type)
1318 CALL pair_potential_lj_create(p%pot(i)%pot%set(std_dim)%lj)
1319 CASE (wl_type)
1320 CALL pair_potential_williams_create(p%pot(i)%pot%set(std_dim)%willis)
1321 CASE (gw_type)
1322 CALL pair_potential_goodwin_create(p%pot(i)%pot%set(std_dim)%goodwin)
1323 CASE (ea_type)
1324 CALL pair_potential_eam_create(p%pot(i)%pot%set(std_dim)%eam)
1325 CASE (quip_type)
1326 CALL pair_potential_quip_create(p%pot(i)%pot%set(std_dim)%quip)
1327 CASE (nequip_type)
1328 CALL pair_potential_nequip_create(p%pot(i)%pot%set(std_dim)%nequip)
1329 CASE (allegro_type)
1330 CALL pair_potential_allegro_create(p%pot(i)%pot%set(std_dim)%allegro)
1331 CASE (deepmd_type)
1332 CALL pair_potential_deepmd_create(p%pot(i)%pot%set(std_dim)%deepmd)
1333 CASE (ft_type)
1334 CALL pair_potential_bmhft_create(p%pot(i)%pot%set(std_dim)%ft)
1335 CASE (ftd_type)
1336 CALL pair_potential_bmhftd_create(p%pot(i)%pot%set(std_dim)%ftd)
1337 CASE (ip_type)
1338 CALL pair_potential_ipbv_create(p%pot(i)%pot%set(std_dim)%ipbv)
1339 CASE (b4_type)
1340 CALL pair_potential_buck4r_create(p%pot(i)%pot%set(std_dim)%buck4r)
1341 CASE (bm_type)
1342 CALL pair_potential_buckmo_create(p%pot(i)%pot%set(std_dim)%buckmo)
1343 CASE (gp_type)
1344 CALL pair_potential_gp_create(p%pot(i)%pot%set(std_dim)%gp)
1345 CASE (tersoff_type)
1346 CALL pair_potential_tersoff_create(p%pot(i)%pot%set(std_dim)%tersoff)
1347 CASE (siepmann_type)
1348 CALL pair_potential_siepmann_create(p%pot(i)%pot%set(std_dim)%siepmann)
1349 CASE (gal_type)
1350 CALL pair_potential_gal_create(p%pot(i)%pot%set(std_dim)%gal)
1351 CASE (gal21_type)
1352 CALL pair_potential_gal21_create(p%pot(i)%pot%set(std_dim)%gal21)
1353 CASE (tab_type)
1354 CALL pair_potential_tab_create(p%pot(i)%pot%set(std_dim)%tab)
1355 END SELECT
1356 NULLIFY (p%pot(i)%pot%spl_f)
1357 NULLIFY (p%pot(i)%pot%pair_spline_data)
1358 END DO
1359
1360 IF (ASSOCIATED(work)) CALL pair_potential_p_release(work)
1361 END SUBROUTINE pair_potential_reallocate
1362
1363! **************************************************************************************************
1364!> \brief Creates the generic potential type
1365!> \param gp ...
1366!> \author Teodoro Laino [teo] 11.2005
1367! **************************************************************************************************
1368 SUBROUTINE pair_potential_gp_create(gp)
1369 TYPE(gp_pot_type), POINTER :: gp
1370
1371 cpassert(.NOT. ASSOCIATED(gp))
1372 ALLOCATE (gp)
1373 NULLIFY (gp%parameters)
1374 NULLIFY (gp%values)
1375 CALL pair_potential_gp_clean(gp)
1376 END SUBROUTINE pair_potential_gp_create
1377
1378! **************************************************************************************************
1379!> \brief Copy two generic potential type
1380!> \param gp_source ...
1381!> \param gp_dest ...
1382!> \author Teodoro Laino [teo] 11.2005
1383! **************************************************************************************************
1384 SUBROUTINE pair_potential_gp_copy(gp_source, gp_dest)
1385 TYPE(gp_pot_type), POINTER :: gp_source, gp_dest
1386
1387 INTEGER :: idim
1388
1389 IF (.NOT. ASSOCIATED(gp_source)) RETURN
1390 IF (ASSOCIATED(gp_dest)) CALL pair_potential_gp_release(gp_dest)
1391 CALL pair_potential_gp_create(gp_dest)
1392 gp_dest%myid = gp_source%myid
1393 gp_dest%potential = gp_source%potential
1394 gp_dest%variables = gp_source%variables
1395 IF (ASSOCIATED(gp_source%parameters)) THEN
1396 idim = SIZE(gp_source%parameters)
1397 ALLOCATE (gp_dest%parameters(idim))
1398 gp_dest%parameters = gp_source%parameters
1399 END IF
1400 IF (ASSOCIATED(gp_source%values)) THEN
1401 idim = SIZE(gp_source%values)
1402 ALLOCATE (gp_dest%values(idim))
1403 gp_dest%values = gp_source%values
1404 END IF
1405 END SUBROUTINE pair_potential_gp_copy
1406
1407! **************************************************************************************************
1408!> \brief Cleans the generic potential type
1409!> \param gp ...
1410!> \author Teodoro Laino [teo] 11.2005
1411! **************************************************************************************************
1412 SUBROUTINE pair_potential_gp_clean(gp)
1413 TYPE(gp_pot_type), POINTER :: gp
1414
1415 IF (.NOT. ASSOCIATED(gp)) RETURN
1416 gp%myid = 0
1417 gp%potential = ""
1418 gp%variables = ""
1419 IF (ASSOCIATED(gp%values)) THEN
1420 DEALLOCATE (gp%values)
1421 END IF
1422 IF (ASSOCIATED(gp%parameters)) THEN
1423 DEALLOCATE (gp%parameters)
1424 END IF
1425 END SUBROUTINE pair_potential_gp_clean
1426
1427! **************************************************************************************************
1428!> \brief Destroys the generic potential type
1429!> \param gp ...
1430!> \author Teodoro Laino [teo] 11.2005
1431! **************************************************************************************************
1432 SUBROUTINE pair_potential_gp_release(gp)
1433 TYPE(gp_pot_type), POINTER :: gp
1434
1435 IF (ASSOCIATED(gp)) THEN
1436 IF (ASSOCIATED(gp%parameters)) THEN
1437 DEALLOCATE (gp%parameters)
1438 END IF
1439 IF (ASSOCIATED(gp%values)) THEN
1440 DEALLOCATE (gp%values)
1441 END IF
1442 DEALLOCATE (gp)
1443 END IF
1444 NULLIFY (gp)
1445 END SUBROUTINE pair_potential_gp_release
1446
1447! **************************************************************************************************
1448!> \brief Cleans the LJ potential type
1449!> \param lj ...
1450!> \author Teodoro Laino [teo] 11.2005
1451! **************************************************************************************************
1453 TYPE(lj_pot_type), POINTER :: lj
1454
1455 cpassert(.NOT. ASSOCIATED(lj))
1456 ALLOCATE (lj)
1457 CALL pair_potential_lj_clean(lj)
1458 END SUBROUTINE pair_potential_lj_create
1459
1460! **************************************************************************************************
1461!> \brief Copy two LJ potential type
1462!> \param lj_source ...
1463!> \param lj_dest ...
1464!> \author Teodoro Laino [teo] 11.2005
1465! **************************************************************************************************
1466 SUBROUTINE pair_potential_lj_copy(lj_source, lj_dest)
1467 TYPE(lj_pot_type), POINTER :: lj_source, lj_dest
1468
1469 IF (.NOT. ASSOCIATED(lj_source)) RETURN
1470 IF (ASSOCIATED(lj_dest)) CALL pair_potential_lj_release(lj_dest)
1471 CALL pair_potential_lj_create(lj_dest)
1472 lj_dest%epsilon = lj_source%epsilon
1473 lj_dest%sigma6 = lj_source%sigma6
1474 lj_dest%sigma12 = lj_source%sigma12
1475 END SUBROUTINE pair_potential_lj_copy
1476
1477! **************************************************************************************************
1478!> \brief Creates the LJ potential type
1479!> \param lj ...
1480!> \author Teodoro Laino [teo] 11.2005
1481! **************************************************************************************************
1482 SUBROUTINE pair_potential_lj_clean(lj)
1483 TYPE(lj_pot_type), POINTER :: lj
1484
1485 IF (.NOT. ASSOCIATED(lj)) RETURN
1486 lj%epsilon = 0.0_dp
1487 lj%sigma6 = 0.0_dp
1488 lj%sigma12 = 0.0_dp
1489 END SUBROUTINE pair_potential_lj_clean
1490
1491! **************************************************************************************************
1492!> \brief Destroys the LJ potential type
1493!> \param lj ...
1494!> \author Teodoro Laino [teo] 11.2005
1495! **************************************************************************************************
1496 SUBROUTINE pair_potential_lj_release(lj)
1497 TYPE(lj_pot_type), POINTER :: lj
1498
1499 IF (ASSOCIATED(lj)) THEN
1500 DEALLOCATE (lj)
1501 END IF
1502 NULLIFY (lj)
1503 END SUBROUTINE pair_potential_lj_release
1504
1505! **************************************************************************************************
1506!> \brief Creates the WILLIAMS potential type
1507!> \param willis ...
1508!> \author Teodoro Laino [teo] 11.2005
1509! **************************************************************************************************
1510 SUBROUTINE pair_potential_williams_create(willis)
1511 TYPE(williams_pot_type), POINTER :: willis
1512
1513 cpassert(.NOT. ASSOCIATED(willis))
1514 ALLOCATE (willis)
1515 CALL pair_potential_williams_clean(willis)
1516 END SUBROUTINE pair_potential_williams_create
1517
1518! **************************************************************************************************
1519!> \brief Copy two WILLIAMS potential type
1520!> \param willis_source ...
1521!> \param willis_dest ...
1522!> \author Teodoro Laino [teo] 11.2005
1523! **************************************************************************************************
1524 SUBROUTINE pair_potential_williams_copy(willis_source, willis_dest)
1525 TYPE(williams_pot_type), POINTER :: willis_source, willis_dest
1526
1527 IF (.NOT. ASSOCIATED(willis_source)) RETURN
1528 IF (ASSOCIATED(willis_dest)) CALL pair_potential_williams_release(willis_dest)
1529 CALL pair_potential_williams_create(willis_dest)
1530 willis_dest%a = willis_source%a
1531 willis_dest%b = willis_source%b
1532 willis_dest%c = willis_source%c
1533 END SUBROUTINE pair_potential_williams_copy
1534
1535! **************************************************************************************************
1536!> \brief Creates the WILLIAMS potential type
1537!> \param willis ...
1538!> \author Teodoro Laino [teo] 11.2005
1539! **************************************************************************************************
1540 SUBROUTINE pair_potential_williams_clean(willis)
1541 TYPE(williams_pot_type), POINTER :: willis
1542
1543 IF (.NOT. ASSOCIATED(willis)) RETURN
1544 willis%a = 0.0_dp
1545 willis%b = 0.0_dp
1546 willis%c = 0.0_dp
1547 END SUBROUTINE pair_potential_williams_clean
1548
1549! **************************************************************************************************
1550!> \brief Destroys the WILLIAMS potential type
1551!> \param willis ...
1552!> \author Teodoro Laino [teo] 11.2005
1553! **************************************************************************************************
1554 SUBROUTINE pair_potential_williams_release(willis)
1555 TYPE(williams_pot_type), POINTER :: willis
1556
1557 IF (ASSOCIATED(willis)) THEN
1558 DEALLOCATE (willis)
1559 END IF
1560 NULLIFY (willis)
1561 END SUBROUTINE pair_potential_williams_release
1562
1563! **************************************************************************************************
1564!> \brief Creates the GOODWIN potential type
1565!> \param goodwin ...
1566!> \author Teodoro Laino [teo] 11.2005
1567! **************************************************************************************************
1568 SUBROUTINE pair_potential_goodwin_create(goodwin)
1569 TYPE(goodwin_pot_type), POINTER :: goodwin
1570
1571 cpassert(.NOT. ASSOCIATED(goodwin))
1572 ALLOCATE (goodwin)
1573 CALL pair_potential_goodwin_clean(goodwin)
1574 END SUBROUTINE pair_potential_goodwin_create
1575
1576! **************************************************************************************************
1577!> \brief Copy two GOODWIN potential type
1578!> \param goodwin_source ...
1579!> \param goodwin_dest ...
1580!> \author Teodoro Laino [teo] 11.2005
1581! **************************************************************************************************
1582 SUBROUTINE pair_potential_goodwin_copy(goodwin_source, goodwin_dest)
1583 TYPE(goodwin_pot_type), POINTER :: goodwin_source, goodwin_dest
1584
1585 IF (.NOT. ASSOCIATED(goodwin_source)) RETURN
1586 IF (ASSOCIATED(goodwin_dest)) CALL pair_potential_goodwin_release(goodwin_dest)
1587 CALL pair_potential_goodwin_create(goodwin_dest)
1588 goodwin_dest%vr0 = goodwin_source%vr0
1589 goodwin_dest%d = goodwin_source%d
1590 goodwin_dest%dc = goodwin_source%dc
1591 goodwin_dest%m = goodwin_source%m
1592 goodwin_dest%mc = goodwin_source%mc
1593 END SUBROUTINE pair_potential_goodwin_copy
1594
1595! **************************************************************************************************
1596!> \brief Creates the GOODWIN potential type
1597!> \param goodwin ...
1598!> \author Teodoro Laino [teo] 11.2005
1599! **************************************************************************************************
1600 SUBROUTINE pair_potential_goodwin_clean(goodwin)
1601 TYPE(goodwin_pot_type), POINTER :: goodwin
1602
1603 IF (.NOT. ASSOCIATED(goodwin)) RETURN
1604 goodwin%vr0 = 0.0_dp
1605 goodwin%d = 0.0_dp
1606 goodwin%dc = 0.0_dp
1607 goodwin%m = 0.0_dp
1608 goodwin%mc = 0.0_dp
1609 END SUBROUTINE pair_potential_goodwin_clean
1610
1611! **************************************************************************************************
1612!> \brief Destroys the GOODWIN potential type
1613!> \param goodwin ...
1614!> \author Teodoro Laino [teo] 11.2005
1615! **************************************************************************************************
1616 SUBROUTINE pair_potential_goodwin_release(goodwin)
1617 TYPE(goodwin_pot_type), POINTER :: goodwin
1618
1619 IF (ASSOCIATED(goodwin)) THEN
1620 DEALLOCATE (goodwin)
1621 END IF
1622 NULLIFY (goodwin)
1623 END SUBROUTINE pair_potential_goodwin_release
1624
1625! **************************************************************************************************
1626!> \brief Creates the EAM potential type
1627!> \param eam ...
1628!> \author Teodoro Laino [teo] 11.2005
1629! **************************************************************************************************
1630 SUBROUTINE pair_potential_eam_create(eam)
1631 TYPE(eam_pot_type), POINTER :: eam
1632
1633 cpassert(.NOT. ASSOCIATED(eam))
1634 ALLOCATE (eam)
1635 NULLIFY (eam%rho, eam%phi, eam%frho, eam%rhoval, eam%rval, &
1636 eam%rhop, eam%phip, eam%frhop)
1637 CALL pair_potential_eam_clean(eam)
1638 END SUBROUTINE pair_potential_eam_create
1639
1640! **************************************************************************************************
1641!> \brief Copy two EAM potential type
1642!> \param eam_source ...
1643!> \param eam_dest ...
1644!> \author Teodoro Laino [teo] 11.2005
1645! **************************************************************************************************
1646 SUBROUTINE pair_potential_eam_copy(eam_source, eam_dest)
1647 TYPE(eam_pot_type), POINTER :: eam_source, eam_dest
1648
1649 IF (.NOT. ASSOCIATED(eam_source)) RETURN
1650 IF (ASSOCIATED(eam_dest)) CALL pair_potential_eam_release(eam_dest)
1651 CALL pair_potential_eam_create(eam_dest)
1652 eam_dest%eam_file_name = eam_source%eam_file_name
1653 eam_dest%drar = eam_source%drar
1654 eam_dest%drhoar = eam_source%drhoar
1655 eam_dest%acutal = eam_source%acutal
1656 eam_dest%npoints = eam_source%npoints
1657 ! Allocate arrays with the proper size
1658 CALL reallocate(eam_dest%rho, 1, eam_dest%npoints)
1659 CALL reallocate(eam_dest%rhop, 1, eam_dest%npoints)
1660 CALL reallocate(eam_dest%phi, 1, eam_dest%npoints)
1661 CALL reallocate(eam_dest%phip, 1, eam_dest%npoints)
1662 CALL reallocate(eam_dest%frho, 1, eam_dest%npoints)
1663 CALL reallocate(eam_dest%frhop, 1, eam_dest%npoints)
1664 CALL reallocate(eam_dest%rval, 1, eam_dest%npoints)
1665 CALL reallocate(eam_dest%rhoval, 1, eam_dest%npoints)
1666 eam_dest%rho = eam_source%rho
1667 eam_dest%phi = eam_source%phi
1668 eam_dest%frho = eam_source%frho
1669 eam_dest%rhoval = eam_source%rhoval
1670 eam_dest%rval = eam_source%rval
1671 eam_dest%rhop = eam_source%rhop
1672 eam_dest%phip = eam_source%phip
1673 eam_dest%frhop = eam_source%frhop
1674 END SUBROUTINE pair_potential_eam_copy
1675
1676! **************************************************************************************************
1677!> \brief Creates the EAM potential type
1678!> \param eam ...
1679!> \author Teodoro Laino [teo] 11.2005
1680! **************************************************************************************************
1681 SUBROUTINE pair_potential_eam_clean(eam)
1682 TYPE(eam_pot_type), POINTER :: eam
1683
1684 IF (.NOT. ASSOCIATED(eam)) RETURN
1685 eam%eam_file_name = 'NULL'
1686 eam%drar = 0.0_dp
1687 eam%drhoar = 0.0_dp
1688 eam%acutal = 0.0_dp
1689 eam%npoints = 0
1690 CALL reallocate(eam%rho, 1, eam%npoints)
1691 CALL reallocate(eam%rhop, 1, eam%npoints)
1692 CALL reallocate(eam%phi, 1, eam%npoints)
1693 CALL reallocate(eam%phip, 1, eam%npoints)
1694 CALL reallocate(eam%frho, 1, eam%npoints)
1695 CALL reallocate(eam%frhop, 1, eam%npoints)
1696 CALL reallocate(eam%rval, 1, eam%npoints)
1697 CALL reallocate(eam%rhoval, 1, eam%npoints)
1698 END SUBROUTINE pair_potential_eam_clean
1699
1700! **************************************************************************************************
1701!> \brief Destroys the EAM potential type
1702!> \param eam ...
1703!> \author Teodoro Laino [teo] 11.2005
1704! **************************************************************************************************
1705 SUBROUTINE pair_potential_eam_release(eam)
1706 TYPE(eam_pot_type), POINTER :: eam
1707
1708 IF (ASSOCIATED(eam)) THEN
1709 IF (ASSOCIATED(eam%rho)) THEN
1710 DEALLOCATE (eam%rho)
1711 END IF
1712 IF (ASSOCIATED(eam%rhop)) THEN
1713 DEALLOCATE (eam%rhop)
1714 END IF
1715 IF (ASSOCIATED(eam%phi)) THEN
1716 DEALLOCATE (eam%phi)
1717 END IF
1718 IF (ASSOCIATED(eam%phip)) THEN
1719 DEALLOCATE (eam%phip)
1720 END IF
1721 IF (ASSOCIATED(eam%frho)) THEN
1722 DEALLOCATE (eam%frho)
1723 END IF
1724 IF (ASSOCIATED(eam%frhop)) THEN
1725 DEALLOCATE (eam%frhop)
1726 END IF
1727 IF (ASSOCIATED(eam%rval)) THEN
1728 DEALLOCATE (eam%rval)
1729 END IF
1730 IF (ASSOCIATED(eam%rhoval)) THEN
1731 DEALLOCATE (eam%rhoval)
1732 END IF
1733 DEALLOCATE (eam)
1734 END IF
1735 END SUBROUTINE pair_potential_eam_release
1736
1737! **************************************************************************************************
1738!> \brief Creates the DEEPMD potential type
1739!> \param deepmd ...
1740!> \author Yongbin Zhuang 07.2019
1741! **************************************************************************************************
1742 SUBROUTINE pair_potential_deepmd_create(deepmd)
1743 TYPE(deepmd_pot_type), POINTER :: deepmd
1744
1745 cpassert(.NOT. ASSOCIATED(deepmd))
1746 ALLOCATE (deepmd)
1747 END SUBROUTINE pair_potential_deepmd_create
1748
1749! **************************************************************************************************
1750!> \brief Copy two DEEPMD potential type
1751!> \param deepmd_source ...
1752!> \param deepmd_dest ...
1753!> \author Yongbin Zhuang 07.2019
1754! **************************************************************************************************
1755 SUBROUTINE pair_potential_deepmd_copy(deepmd_source, deepmd_dest)
1756 TYPE(deepmd_pot_type), POINTER :: deepmd_source, deepmd_dest
1757
1758 IF (.NOT. ASSOCIATED(deepmd_source)) RETURN
1759 NULLIFY (deepmd_dest)
1760 IF (ASSOCIATED(deepmd_dest)) CALL pair_potential_deepmd_release(deepmd_dest)
1761 CALL pair_potential_deepmd_create(deepmd_dest)
1762 deepmd_dest = deepmd_source
1763 END SUBROUTINE pair_potential_deepmd_copy
1764
1765! **************************************************************************************************
1766!> \brief CLEAN the DEEPMD potential type
1767!> \param deepmd ...
1768!> \author Yongbin Zhuang 07.2019
1769! **************************************************************************************************
1770 SUBROUTINE pair_potential_deepmd_clean(deepmd)
1771 TYPE(deepmd_pot_type), POINTER :: deepmd
1772
1773 IF (.NOT. ASSOCIATED(deepmd)) RETURN
1774 deepmd = deepmd_pot_type()
1775 END SUBROUTINE pair_potential_deepmd_clean
1776
1777! **************************************************************************************************
1778!> \brief Destroys the DEEPMD potential type
1779!> \param deepmd ...
1780!> \author Yongbin Zhuang 07.2019
1781! **************************************************************************************************
1782 SUBROUTINE pair_potential_deepmd_release(deepmd)
1783 TYPE(deepmd_pot_type), POINTER :: deepmd
1784
1785 IF (ASSOCIATED(deepmd)) THEN
1786 DEALLOCATE (deepmd)
1787 END IF
1788 END SUBROUTINE pair_potential_deepmd_release
1789
1790! **************************************************************************************************
1791!> \brief Creates the QUIP potential type
1792!> \param quip ...
1793!> \author Teodoro Laino [teo] 11.2005
1794! **************************************************************************************************
1795 SUBROUTINE pair_potential_quip_create(quip)
1796 TYPE(quip_pot_type), POINTER :: quip
1797
1798 cpassert(.NOT. ASSOCIATED(quip))
1799 ALLOCATE (quip)
1800 quip%quip_file_name = ""
1801 quip%init_args = ""
1802 quip%calc_args = ""
1803 CALL pair_potential_quip_clean(quip)
1804 END SUBROUTINE pair_potential_quip_create
1805
1806! **************************************************************************************************
1807!> \brief Copy two QUIP potential type
1808!> \param quip_source ...
1809!> \param quip_dest ...
1810!> \author Teodoro Laino [teo] 11.2005
1811! **************************************************************************************************
1812 SUBROUTINE pair_potential_quip_copy(quip_source, quip_dest)
1813 TYPE(quip_pot_type), POINTER :: quip_source, quip_dest
1814
1815 IF (.NOT. ASSOCIATED(quip_source)) RETURN
1816 IF (ASSOCIATED(quip_dest)) CALL pair_potential_quip_release(quip_dest)
1817 CALL pair_potential_quip_create(quip_dest)
1818 quip_dest%quip_file_name = quip_source%quip_file_name
1819 quip_dest%init_args = quip_source%init_args
1820 quip_dest%calc_args = quip_source%calc_args
1821 END SUBROUTINE pair_potential_quip_copy
1822
1823! **************************************************************************************************
1824!> \brief Creates the QUIP potential type
1825!> \param quip ...
1826!> \author Teodoro Laino [teo] 11.2005
1827! **************************************************************************************************
1828 SUBROUTINE pair_potential_quip_clean(quip)
1829 TYPE(quip_pot_type), POINTER :: quip
1830
1831 IF (.NOT. ASSOCIATED(quip)) RETURN
1832 quip%quip_file_name = 'NULL'
1833 quip%init_args = ''
1834 quip%calc_args = ''
1835 END SUBROUTINE pair_potential_quip_clean
1836
1837! **************************************************************************************************
1838!> \brief Destroys the QUIP potential type
1839!> \param quip ...
1840!> \author Teodoro Laino [teo] 11.2005
1841! **************************************************************************************************
1842 SUBROUTINE pair_potential_quip_release(quip)
1843 TYPE(quip_pot_type), POINTER :: quip
1844
1845 IF (ASSOCIATED(quip)) THEN
1846 DEALLOCATE (quip)
1847 END IF
1848 END SUBROUTINE pair_potential_quip_release
1849
1850! **************************************************************************************************
1851!> \brief Creates the NEQUIP potential type
1852!> \param nequip ...
1853!> \author Gabriele Tocci 2023
1854! **************************************************************************************************
1855 SUBROUTINE pair_potential_nequip_create(nequip)
1856 TYPE(nequip_pot_type), POINTER :: nequip
1857
1858 cpassert(.NOT. ASSOCIATED(nequip))
1859 ALLOCATE (nequip)
1860 END SUBROUTINE pair_potential_nequip_create
1861
1862! **************************************************************************************************
1863!> \brief Copy two NEQUIP potential type
1864!> \param nequip_source ...
1865!> \param nequip_dest ...
1866!> \author Gabriele Tocci 2023
1867! **************************************************************************************************
1868 SUBROUTINE pair_potential_nequip_copy(nequip_source, nequip_dest)
1869 TYPE(nequip_pot_type), POINTER :: nequip_source, nequip_dest
1870
1871 IF (.NOT. ASSOCIATED(nequip_source)) RETURN
1872 IF (ASSOCIATED(nequip_dest)) CALL pair_potential_nequip_release(nequip_dest)
1873 CALL pair_potential_nequip_create(nequip_dest)
1874 nequip_dest = nequip_source
1875
1876 END SUBROUTINE pair_potential_nequip_copy
1877
1878! **************************************************************************************************
1879!> \brief Creates the NEQUIP potential type
1880!> \param nequip ...
1881!> \author Gabriele Tocci 2023
1882! **************************************************************************************************
1883 SUBROUTINE pair_potential_nequip_clean(nequip)
1884 TYPE(nequip_pot_type), POINTER :: nequip
1885
1886 IF (.NOT. ASSOCIATED(nequip)) RETURN
1887 nequip = nequip_pot_type()
1888
1889 END SUBROUTINE pair_potential_nequip_clean
1890
1891! **************************************************************************************************
1892!> \brief Destroys the NEQUIP potential type
1893!> \param nequip ...
1894!> \author Gabriele Tocci 2023
1895! **************************************************************************************************
1896 SUBROUTINE pair_potential_nequip_release(nequip)
1897 TYPE(nequip_pot_type), POINTER :: nequip
1898
1899 IF (ASSOCIATED(nequip)) THEN
1900 DEALLOCATE (nequip)
1901 END IF
1902 END SUBROUTINE pair_potential_nequip_release
1903
1904! **************************************************************************************************
1905!> \brief Creates the ALLEGRO potential type
1906!> \param allegro ...
1907!> \author Gabriele Tocci 2023
1908! **************************************************************************************************
1909 SUBROUTINE pair_potential_allegro_create(allegro)
1910 TYPE(allegro_pot_type), POINTER :: allegro
1911
1912 cpassert(.NOT. ASSOCIATED(allegro))
1913 ALLOCATE (allegro)
1914 END SUBROUTINE pair_potential_allegro_create
1915
1916! **************************************************************************************************
1917!> \brief Copy two ALLEGRO potential type
1918!> \param allegro_source ...
1919!> \param allegro_dest ...
1920!> \author Gabriele Tocci 2023
1921! **************************************************************************************************
1922 SUBROUTINE pair_potential_allegro_copy(allegro_source, allegro_dest)
1923 TYPE(allegro_pot_type), POINTER :: allegro_source, allegro_dest
1924
1925 IF (.NOT. ASSOCIATED(allegro_source)) RETURN
1926 IF (ASSOCIATED(allegro_dest)) CALL pair_potential_allegro_release(allegro_dest)
1927 CALL pair_potential_allegro_create(allegro_dest)
1928 allegro_dest = allegro_source
1929 END SUBROUTINE pair_potential_allegro_copy
1930
1931! **************************************************************************************************
1932!> \brief Creates the ALLEGRO potential type
1933!> \param allegro ...
1934!> \author Gabriele Tocci 2023
1935! **************************************************************************************************
1936 SUBROUTINE pair_potential_allegro_clean(allegro)
1937 TYPE(allegro_pot_type), POINTER :: allegro
1938
1939 IF (.NOT. ASSOCIATED(allegro)) RETURN
1940 allegro = allegro_pot_type()
1941
1942 END SUBROUTINE pair_potential_allegro_clean
1943
1944! **************************************************************************************************
1945!> \brief Destroys the ALLEGRO potential type
1946!> \param allegro ...
1947!> \author Gabriele Tocci 2023
1948! **************************************************************************************************
1949 SUBROUTINE pair_potential_allegro_release(allegro)
1950 TYPE(allegro_pot_type), POINTER :: allegro
1951
1952 IF (ASSOCIATED(allegro)) THEN
1953 DEALLOCATE (allegro)
1954 END IF
1955 END SUBROUTINE pair_potential_allegro_release
1956
1957! **************************************************************************************************
1958!> \brief Creates the BMHFT (TOSI-FUMI) potential type
1959!> \param ft ...
1960!> \author Teodoro Laino [teo] 11.2005
1961! **************************************************************************************************
1962 SUBROUTINE pair_potential_bmhft_create(ft)
1963 TYPE(ft_pot_type), POINTER :: ft
1964
1965 cpassert(.NOT. ASSOCIATED(ft))
1966 ALLOCATE (ft)
1967 CALL pair_potential_bmhft_clean(ft)
1968 END SUBROUTINE pair_potential_bmhft_create
1969
1970! **************************************************************************************************
1971!> \brief Copy two BMHFT (TOSI-FUMI) potential type
1972!> \param ft_source ...
1973!> \param ft_dest ...
1974!> \author Teodoro Laino [teo] 11.2005
1975! **************************************************************************************************
1976 SUBROUTINE pair_potential_bmhft_copy(ft_source, ft_dest)
1977 TYPE(ft_pot_type), POINTER :: ft_source, ft_dest
1978
1979 IF (.NOT. ASSOCIATED(ft_source)) RETURN
1980 IF (ASSOCIATED(ft_dest)) CALL pair_potential_bmhft_release(ft_dest)
1981 CALL pair_potential_bmhft_create(ft_dest)
1982 ft_dest%A = ft_source%A
1983 ft_dest%B = ft_source%B
1984 ft_dest%C = ft_source%C
1985 ft_dest%D = ft_source%D
1986 END SUBROUTINE pair_potential_bmhft_copy
1987
1988! **************************************************************************************************
1989!> \brief Creates the BMHFT (TOSI-FUMI) potential type
1990!> \param ft ...
1991!> \author Teodoro Laino [teo] 11.2005
1992! **************************************************************************************************
1993 SUBROUTINE pair_potential_bmhft_clean(ft)
1994 TYPE(ft_pot_type), POINTER :: ft
1995
1996 IF (.NOT. ASSOCIATED(ft)) RETURN
1997 ft%A = 0.0_dp
1998 ft%B = 0.0_dp
1999 ft%C = 0.0_dp
2000 ft%D = 0.0_dp
2001 END SUBROUTINE pair_potential_bmhft_clean
2002
2003! **************************************************************************************************
2004!> \brief Destroys the BMHFT potential type
2005!> \param ft ...
2006!> \author Teodoro Laino [teo] 11.2005
2007! **************************************************************************************************
2008 SUBROUTINE pair_potential_bmhft_release(ft)
2009 TYPE(ft_pot_type), POINTER :: ft
2010
2011 IF (ASSOCIATED(ft)) THEN
2012 DEALLOCATE (ft)
2013 END IF
2014 NULLIFY (ft)
2015 END SUBROUTINE pair_potential_bmhft_release
2016
2017! **************************************************************************************************
2018!> \brief Creates the BMHFTD (damped TOSI-FUMI) potential type
2019!> \param ftd ...
2020!> \author Mathieu Salanne 05.2010
2021! **************************************************************************************************
2022 SUBROUTINE pair_potential_bmhftd_create(ftd)
2023 TYPE(ftd_pot_type), POINTER :: ftd
2024
2025 cpassert(.NOT. ASSOCIATED(ftd))
2026 ALLOCATE (ftd)
2027 CALL pair_potential_bmhftd_clean(ftd)
2028 END SUBROUTINE pair_potential_bmhftd_create
2029
2030! **************************************************************************************************
2031!> \brief Copy two BMHFTD (Damped TOSI-FUMI) potential type
2032!> \param ftd_source ...
2033!> \param ftd_dest ...
2034!> \author Mathieu Salanne 05.2010
2035! **************************************************************************************************
2036 SUBROUTINE pair_potential_bmhftd_copy(ftd_source, ftd_dest)
2037 TYPE(ftd_pot_type), POINTER :: ftd_source, ftd_dest
2038
2039 IF (.NOT. ASSOCIATED(ftd_source)) RETURN
2040 IF (ASSOCIATED(ftd_dest)) CALL pair_potential_bmhftd_release(ftd_dest)
2041 CALL pair_potential_bmhftd_create(ftd_dest)
2042 ftd_dest%A = ftd_source%A
2043 ftd_dest%B = ftd_source%B
2044 ftd_dest%C = ftd_source%C
2045 ftd_dest%D = ftd_source%D
2046 ftd_dest%BD = ftd_source%BD
2047 END SUBROUTINE pair_potential_bmhftd_copy
2048
2049! **************************************************************************************************
2050!> \brief Cleans the BMHFTD (damped TOSI-FUMI) potential type
2051!> \param ftd ...
2052!> \author Mathieu Salanne
2053! **************************************************************************************************
2054 SUBROUTINE pair_potential_bmhftd_clean(ftd)
2055 TYPE(ftd_pot_type), POINTER :: ftd
2056
2057 IF (.NOT. ASSOCIATED(ftd)) RETURN
2058 ftd%A = 0.0_dp
2059 ftd%B = 0.0_dp
2060 ftd%C = 0.0_dp
2061 ftd%D = 0.0_dp
2062 ftd%BD = 0.0_dp
2063 END SUBROUTINE pair_potential_bmhftd_clean
2064
2065! **************************************************************************************************
2066!> \brief Destroys the BMHFTD potential type
2067!> \param ftd ...
2068!> \author Mathieu Salanne 05.2010
2069! **************************************************************************************************
2070 SUBROUTINE pair_potential_bmhftd_release(ftd)
2071 TYPE(ftd_pot_type), POINTER :: ftd
2072
2073 IF (ASSOCIATED(ftd)) THEN
2074 DEALLOCATE (ftd)
2075 END IF
2076 NULLIFY (ftd)
2077 END SUBROUTINE pair_potential_bmhftd_release
2078
2079! **************************************************************************************************
2080!> \brief Creates the IPBV potential type
2081!> \param ipbv ...
2082!> \author Teodoro Laino [teo] 11.2005
2083! **************************************************************************************************
2084 SUBROUTINE pair_potential_ipbv_create(ipbv)
2085 TYPE(ipbv_pot_type), POINTER :: ipbv
2086
2087 cpassert(.NOT. ASSOCIATED(ipbv))
2088 ALLOCATE (ipbv)
2089 CALL pair_potential_ipbv_clean(ipbv)
2090 END SUBROUTINE pair_potential_ipbv_create
2091
2092! **************************************************************************************************
2093!> \brief Copy two IPBV potential type
2094!> \param ipbv_source ...
2095!> \param ipbv_dest ...
2096!> \author Teodoro Laino [teo] 11.2005
2097! **************************************************************************************************
2098 SUBROUTINE pair_potential_ipbv_copy(ipbv_source, ipbv_dest)
2099 TYPE(ipbv_pot_type), POINTER :: ipbv_source, ipbv_dest
2100
2101 IF (.NOT. ASSOCIATED(ipbv_source)) RETURN
2102 IF (ASSOCIATED(ipbv_dest)) CALL pair_potential_ipbv_release(ipbv_dest)
2103 CALL pair_potential_ipbv_create(ipbv_dest)
2104 ipbv_dest%a = ipbv_source%a
2105 ipbv_dest%rcore = ipbv_source%rcore
2106 ipbv_dest%b = ipbv_source%b
2107 ipbv_dest%m = ipbv_source%m
2108 END SUBROUTINE pair_potential_ipbv_copy
2109
2110! **************************************************************************************************
2111!> \brief Creates the IPBV potential type
2112!> \param ipbv ...
2113!> \author Teodoro Laino [teo] 11.2005
2114! **************************************************************************************************
2115 SUBROUTINE pair_potential_ipbv_clean(ipbv)
2116 TYPE(ipbv_pot_type), POINTER :: ipbv
2117
2118 IF (.NOT. ASSOCIATED(ipbv)) RETURN
2119 ipbv%a = 0.0_dp
2120 ipbv%rcore = 0.0_dp
2121 ipbv%b = 0.0_dp
2122 ipbv%m = 0.0_dp
2123 END SUBROUTINE pair_potential_ipbv_clean
2124
2125! **************************************************************************************************
2126!> \brief Destroys the IPBV potential type
2127!> \param ipbv ...
2128!> \author Teodoro Laino [teo] 11.2005
2129! **************************************************************************************************
2130 SUBROUTINE pair_potential_ipbv_release(ipbv)
2131 TYPE(ipbv_pot_type), POINTER :: ipbv
2132
2133 IF (ASSOCIATED(ipbv)) THEN
2134 DEALLOCATE (ipbv)
2135 END IF
2136 NULLIFY (ipbv)
2137 END SUBROUTINE pair_potential_ipbv_release
2138
2139! **************************************************************************************************
2140!> \brief Creates the Buckingham 4 ranges potential type
2141!> \param buck4r ...
2142!> \author MI 10.2006
2143! **************************************************************************************************
2144 SUBROUTINE pair_potential_buck4r_create(buck4r)
2145 TYPE(buck4ran_pot_type), POINTER :: buck4r
2146
2147 cpassert(.NOT. ASSOCIATED(buck4r))
2148 ALLOCATE (buck4r)
2149 CALL pair_potential_buck4r_clean(buck4r)
2150 END SUBROUTINE pair_potential_buck4r_create
2151
2152! **************************************************************************************************
2153!> \brief Copy two Buckingham 4 ranges potential type
2154!> \param buck4r_source ...
2155!> \param buck4r_dest ...
2156!> \author MI 10.2006
2157! **************************************************************************************************
2158 SUBROUTINE pair_potential_buck4r_copy(buck4r_source, buck4r_dest)
2159 TYPE(buck4ran_pot_type), POINTER :: buck4r_source, buck4r_dest
2160
2161 IF (.NOT. ASSOCIATED(buck4r_source)) RETURN
2162 IF (ASSOCIATED(buck4r_dest)) CALL pair_potential_buck4r_release(buck4r_dest)
2163 CALL pair_potential_buck4r_create(buck4r_dest)
2164 buck4r_dest%a = buck4r_source%a
2165 buck4r_dest%b = buck4r_source%b
2166 buck4r_dest%c = buck4r_source%c
2167 buck4r_dest%r1 = buck4r_source%r1
2168 buck4r_dest%r2 = buck4r_source%r2
2169 buck4r_dest%r3 = buck4r_source%r3
2170 buck4r_dest%poly1 = buck4r_source%poly1
2171 buck4r_dest%poly2 = buck4r_source%poly2
2172 buck4r_dest%npoly1 = buck4r_source%npoly1
2173 buck4r_dest%npoly2 = buck4r_source%npoly2
2174 END SUBROUTINE pair_potential_buck4r_copy
2175
2176! **************************************************************************************************
2177!> \brief Creates the Buckingham 4 ranges potential type
2178!> \param buck4r ...
2179!> \author MI 10.2006
2180! **************************************************************************************************
2181 SUBROUTINE pair_potential_buck4r_clean(buck4r)
2182 TYPE(buck4ran_pot_type), POINTER :: buck4r
2183
2184 IF (.NOT. ASSOCIATED(buck4r)) RETURN
2185 buck4r%a = 0.0_dp
2186 buck4r%b = 0.0_dp
2187 buck4r%c = 0.0_dp
2188 buck4r%r1 = 0.0_dp
2189 buck4r%r2 = 0.0_dp
2190 buck4r%r3 = 0.0_dp
2191 buck4r%poly1 = 0.0_dp
2192 buck4r%npoly1 = 0
2193 buck4r%poly2 = 0.0_dp
2194 buck4r%npoly2 = 0
2195 END SUBROUTINE pair_potential_buck4r_clean
2196
2197! **************************************************************************************************
2198!> \brief Destroys the Buckingham 4 ranges potential type
2199!> \param buck4r ...
2200!> \author MI 10.2006
2201! **************************************************************************************************
2202 SUBROUTINE pair_potential_buck4r_release(buck4r)
2203 TYPE(buck4ran_pot_type), POINTER :: buck4r
2204
2205 IF (ASSOCIATED(buck4r)) THEN
2206 DEALLOCATE (buck4r)
2207 END IF
2208 NULLIFY (buck4r)
2209 END SUBROUTINE pair_potential_buck4r_release
2210
2211! **************************************************************************************************
2212!> \brief Creates the Buckingham plus Morse potential type
2213!> \param buckmo ...
2214!> \author MI 10.2006
2215! **************************************************************************************************
2216 SUBROUTINE pair_potential_buckmo_create(buckmo)
2217 TYPE(buckmorse_pot_type), POINTER :: buckmo
2218
2219 cpassert(.NOT. ASSOCIATED(buckmo))
2220 ALLOCATE (buckmo)
2221 CALL pair_potential_buckmo_clean(buckmo)
2222 END SUBROUTINE pair_potential_buckmo_create
2223
2224! **************************************************************************************************
2225!> \brief Copy two Buckingham plus Morse potential type
2226!> \param buckmo_source ...
2227!> \param buckmo_dest ...
2228!> \author MI 10.2006
2229! **************************************************************************************************
2230 SUBROUTINE pair_potential_buckmo_copy(buckmo_source, buckmo_dest)
2231 TYPE(buckmorse_pot_type), POINTER :: buckmo_source, buckmo_dest
2232
2233 IF (.NOT. ASSOCIATED(buckmo_source)) RETURN
2234 IF (ASSOCIATED(buckmo_dest)) CALL pair_potential_buckmo_release(buckmo_dest)
2235 CALL pair_potential_buckmo_create(buckmo_dest)
2236 buckmo_dest%f0 = buckmo_source%f0
2237 buckmo_dest%a1 = buckmo_source%a1
2238 buckmo_dest%a2 = buckmo_source%a2
2239 buckmo_dest%b1 = buckmo_source%b1
2240 buckmo_dest%b2 = buckmo_source%b2
2241 buckmo_dest%c = buckmo_source%c
2242 buckmo_dest%d = buckmo_source%d
2243 buckmo_dest%r0 = buckmo_source%r0
2244 buckmo_dest%beta = buckmo_source%beta
2245 END SUBROUTINE pair_potential_buckmo_copy
2246
2247! **************************************************************************************************
2248!> \brief Creates the Buckingham plus Morse potential type
2249!> \param buckmo ...
2250!> \author MI 10.2006
2251! **************************************************************************************************
2252 SUBROUTINE pair_potential_buckmo_clean(buckmo)
2253 TYPE(buckmorse_pot_type), POINTER :: buckmo
2254
2255 IF (.NOT. ASSOCIATED(buckmo)) RETURN
2256 buckmo%f0 = 0.0_dp
2257 buckmo%a1 = 0.0_dp
2258 buckmo%a2 = 0.0_dp
2259 buckmo%b1 = 0.0_dp
2260 buckmo%b2 = 0.0_dp
2261 buckmo%c = 0.0_dp
2262 buckmo%d = 0.0_dp
2263 buckmo%r0 = 0.0_dp
2264 buckmo%beta = 0.0_dp
2265 END SUBROUTINE pair_potential_buckmo_clean
2266
2267! **************************************************************************************************
2268!> \brief Destroys the Buckingham plus Morse potential type
2269!> \param buckmo ...
2270!> \author MI 10.2006
2271! **************************************************************************************************
2272 SUBROUTINE pair_potential_buckmo_release(buckmo)
2273 TYPE(buckmorse_pot_type), POINTER :: buckmo
2274
2275 IF (ASSOCIATED(buckmo)) THEN
2276 DEALLOCATE (buckmo)
2277 END IF
2278 NULLIFY (buckmo)
2279 END SUBROUTINE pair_potential_buckmo_release
2280
2281! **************************************************************************************************
2282!> \brief Creates the Tersoff potential type
2283!> (Tersoff, J. PRB 39(8), 5566, 1989)
2284!> \param tersoff ...
2285! **************************************************************************************************
2286 SUBROUTINE pair_potential_tersoff_create(tersoff)
2287 TYPE(tersoff_pot_type), POINTER :: tersoff
2288
2289 cpassert(.NOT. ASSOCIATED(tersoff))
2290 ALLOCATE (tersoff)
2291 CALL pair_potential_tersoff_clean(tersoff)
2292 END SUBROUTINE pair_potential_tersoff_create
2293
2294! **************************************************************************************************
2295!> \brief Copy two Tersoff potential type
2296!> (Tersoff, J. PRB 39(8), 5566, 1989)
2297!> \param tersoff_source ...
2298!> \param tersoff_dest ...
2299! **************************************************************************************************
2300 SUBROUTINE pair_potential_tersoff_copy(tersoff_source, tersoff_dest)
2301 TYPE(tersoff_pot_type), POINTER :: tersoff_source, tersoff_dest
2302
2303 IF (.NOT. ASSOCIATED(tersoff_source)) RETURN
2304 IF (ASSOCIATED(tersoff_dest)) CALL pair_potential_tersoff_release(tersoff_dest)
2305 CALL pair_potential_tersoff_create(tersoff_dest)
2306 tersoff_dest%A = tersoff_source%A
2307 tersoff_dest%B = tersoff_source%B
2308 tersoff_dest%lambda1 = tersoff_source%lambda1
2309 tersoff_dest%lambda2 = tersoff_source%lambda2
2310 tersoff_dest%alpha = tersoff_source%alpha
2311 tersoff_dest%beta = tersoff_source%beta
2312 tersoff_dest%n = tersoff_source%n
2313 tersoff_dest%c = tersoff_source%c
2314 tersoff_dest%d = tersoff_source%d
2315 tersoff_dest%h = tersoff_source%h
2316 tersoff_dest%lambda3 = tersoff_source%lambda3
2317 tersoff_dest%bigR = tersoff_source%bigR
2318 tersoff_dest%bigD = tersoff_source%bigD
2319 tersoff_dest%rcutsq = tersoff_source%rcutsq
2320 END SUBROUTINE pair_potential_tersoff_copy
2321
2322! **************************************************************************************************
2323!> \brief Creates the Tersoff potential type
2324!> (Tersoff, J. PRB 39(8), 5566, 1989)
2325!> \param tersoff ...
2326! **************************************************************************************************
2327 SUBROUTINE pair_potential_tersoff_clean(tersoff)
2328 TYPE(tersoff_pot_type), POINTER :: tersoff
2329
2330 IF (.NOT. ASSOCIATED(tersoff)) RETURN
2331 tersoff%A = 0.0_dp
2332 tersoff%B = 0.0_dp
2333 tersoff%lambda1 = 0.0_dp
2334 tersoff%lambda2 = 0.0_dp
2335 tersoff%alpha = 0.0_dp
2336 tersoff%beta = 0.0_dp
2337 tersoff%n = 0.0_dp
2338 tersoff%c = 0.0_dp
2339 tersoff%d = 0.0_dp
2340 tersoff%h = 0.0_dp
2341 tersoff%lambda3 = 0.0_dp
2342 tersoff%bigR = 0.0_dp
2343 tersoff%bigD = 0.0_dp
2344 tersoff%rcutsq = 0.0_dp
2345 END SUBROUTINE pair_potential_tersoff_clean
2346
2347! **************************************************************************************************
2348!> \brief Destroys the Tersoff
2349!> (Tersoff, J. PRB 39(8), 5566, 1989)
2350!> \param tersoff ...
2351! **************************************************************************************************
2352 SUBROUTINE pair_potential_tersoff_release(tersoff)
2353 TYPE(tersoff_pot_type), POINTER :: tersoff
2354
2355 IF (ASSOCIATED(tersoff)) THEN
2356 DEALLOCATE (tersoff)
2357 END IF
2358 NULLIFY (tersoff)
2359 END SUBROUTINE pair_potential_tersoff_release
2360
2361! **************************************************************************************************
2362!> \brief Creates the Siepmann-Sprik potential type
2363!> (Siepmann and Sprik, J. Chem. Phys. 102(1) 511, 1995)
2364!> \param siepmann ...
2365! **************************************************************************************************
2366 SUBROUTINE pair_potential_siepmann_create(siepmann)
2367 TYPE(siepmann_pot_type), POINTER :: siepmann
2368
2369 cpassert(.NOT. ASSOCIATED(siepmann))
2370 ALLOCATE (siepmann)
2371 CALL pair_potential_siepmann_clean(siepmann)
2372 END SUBROUTINE pair_potential_siepmann_create
2373! **************************************************************************************************
2374!> \brief Copy two Siepmann potential type
2375!> (Siepmann and Sprik, J. Chem. Phys. 102(1) 511, 1995)
2376!> \param siepmann_source ...
2377!> \param siepmann_dest ...
2378! **************************************************************************************************
2379 SUBROUTINE pair_potential_siepmann_copy(siepmann_source, siepmann_dest)
2380 TYPE(siepmann_pot_type), POINTER :: siepmann_source, siepmann_dest
2381
2382 IF (.NOT. ASSOCIATED(siepmann_source)) RETURN
2383 IF (ASSOCIATED(siepmann_dest)) CALL pair_potential_siepmann_release(siepmann_dest)
2384 CALL pair_potential_siepmann_create(siepmann_dest)
2385 siepmann_dest%B = siepmann_source%B
2386 siepmann_dest%D = siepmann_source%D
2387 siepmann_dest%E = siepmann_source%E
2388 siepmann_dest%F = siepmann_source%F
2389 siepmann_dest%beta = siepmann_source%beta
2390 siepmann_dest%rcutsq = siepmann_source%rcutsq
2391 siepmann_dest%allow_oh_formation = siepmann_source%allow_oh_formation
2392 siepmann_dest%allow_h3o_formation = siepmann_source%allow_h3o_formation
2393 siepmann_dest%allow_o_formation = siepmann_source%allow_o_formation
2394
2395 END SUBROUTINE pair_potential_siepmann_copy
2396
2397! **************************************************************************************************
2398!> \brief Creates the Siepmann-Sprik potential type
2399!> (Siepmann and Sprik, J. Chem. Phys. 102(1) 511, 1995)
2400!> \param siepmann ...
2401! **************************************************************************************************
2402 SUBROUTINE pair_potential_siepmann_clean(siepmann)
2403 TYPE(siepmann_pot_type), POINTER :: siepmann
2404
2405 IF (.NOT. ASSOCIATED(siepmann)) RETURN
2406 siepmann%B = 0.0_dp
2407 siepmann%D = 0.0_dp
2408 siepmann%E = 0.0_dp
2409 siepmann%F = 0.0_dp
2410 siepmann%beta = 0.0_dp
2411 siepmann%rcutsq = 0.0_dp
2412 siepmann%allow_oh_formation = .false.
2413 siepmann%allow_h3o_formation = .false.
2414 siepmann%allow_o_formation = .false.
2415
2416 END SUBROUTINE pair_potential_siepmann_clean
2417
2418! **************************************************************************************************
2419!> \brief Destroys the Siepmann-Sprik potential
2420!> (Siepmann and Sprik, J. Chem. Phys. 102(1) 511, 1995)
2421!> \param siepmann ...
2422! **************************************************************************************************
2423 SUBROUTINE pair_potential_siepmann_release(siepmann)
2424 TYPE(siepmann_pot_type), POINTER :: siepmann
2425
2426 IF (ASSOCIATED(siepmann)) THEN
2427 DEALLOCATE (siepmann)
2428 END IF
2429 NULLIFY (siepmann)
2430 END SUBROUTINE pair_potential_siepmann_release
2431
2432! **************************************************************************************************
2433!> \brief Creates the GAL19 potential type
2434!> (??)
2435!> \param gal ...
2436! **************************************************************************************************
2437 SUBROUTINE pair_potential_gal_create(gal)
2438 TYPE(gal_pot_type), POINTER :: gal
2439
2440 cpassert(.NOT. ASSOCIATED(gal))
2441 ALLOCATE (gal)
2442 CALL pair_potential_gal_clean(gal)
2443 END SUBROUTINE pair_potential_gal_create
2444
2445! **************************************************************************************************
2446!> \brief Copy two GAL potential type
2447!> (??)
2448!> \param gal_source ...
2449!> \param gal_dest ...
2450! **************************************************************************************************
2451 SUBROUTINE pair_potential_gal_copy(gal_source, gal_dest)
2452 TYPE(gal_pot_type), POINTER :: gal_source, gal_dest
2453
2454 IF (.NOT. ASSOCIATED(gal_source)) RETURN
2455 IF (ASSOCIATED(gal_dest)) CALL pair_potential_gal_release(gal_dest)
2456 CALL pair_potential_gal_create(gal_dest)
2457 gal_dest%met1 = gal_source%met1
2458 gal_dest%met2 = gal_source%met2
2459 gal_dest%epsilon = gal_source%epsilon
2460 gal_dest%bxy = gal_source%bxy
2461 gal_dest%bz = gal_source%bz
2462 gal_dest%r1 = gal_source%r1
2463 gal_dest%r2 = gal_source%r2
2464 gal_dest%a1 = gal_source%a1
2465 gal_dest%a2 = gal_source%a2
2466 gal_dest%a3 = gal_source%a3
2467 gal_dest%a4 = gal_source%a4
2468 gal_dest%a = gal_source%a
2469 gal_dest%b = gal_source%b
2470 gal_dest%c = gal_source%c
2471 ALLOCATE (gal_dest%gcn(SIZE(gal_source%gcn)))
2472 gal_dest%gcn = gal_source%gcn
2473 gal_dest%express = gal_source%express
2474 gal_dest%rcutsq = gal_source%rcutsq
2475
2476 END SUBROUTINE pair_potential_gal_copy
2477
2478! **************************************************************************************************
2479!> \brief Creates the GAL19 potential type
2480!> (??)
2481!> \param gal ...
2482! **************************************************************************************************
2483 SUBROUTINE pair_potential_gal_clean(gal)
2484 TYPE(gal_pot_type), POINTER :: gal
2485
2486 IF (.NOT. ASSOCIATED(gal)) RETURN
2487 gal%epsilon = 0.0_dp
2488 gal%bxy = 0.0_dp
2489 gal%bz = 0.0_dp
2490 gal%r1 = 0.0_dp
2491 gal%r2 = 0.0_dp
2492 gal%a1 = 0.0_dp
2493 gal%a2 = 0.0_dp
2494 gal%a3 = 0.0_dp
2495 gal%a4 = 0.0_dp
2496 gal%a = 0.0_dp
2497 gal%b = 0.0_dp
2498 gal%c = 0.0_dp
2499 gal%rcutsq = 0.0_dp
2500 gal%express = .false.
2501
2502 END SUBROUTINE pair_potential_gal_clean
2503
2504! **************************************************************************************************
2505!> \brief Destroys the GAL19 potential
2506!> (??)
2507!> \param gal ...
2508! **************************************************************************************************
2509 SUBROUTINE pair_potential_gal_release(gal)
2510 TYPE(gal_pot_type), POINTER :: gal
2511
2512 IF (ASSOCIATED(gal)) THEN
2513 DEALLOCATE (gal%gcn)
2514 DEALLOCATE (gal)
2515 END IF
2516 NULLIFY (gal)
2517 END SUBROUTINE pair_potential_gal_release
2518
2519! **************************************************************************************************
2520!> \brief Creates the GAL21 potential type
2521!> (??)
2522!> \param gal21 ...
2523! **************************************************************************************************
2524 SUBROUTINE pair_potential_gal21_create(gal21)
2525 TYPE(gal21_pot_type), POINTER :: gal21
2526
2527 cpassert(.NOT. ASSOCIATED(gal21))
2528 ALLOCATE (gal21)
2529 CALL pair_potential_gal21_clean(gal21)
2530 END SUBROUTINE pair_potential_gal21_create
2531
2532! **************************************************************************************************
2533!> \brief Copy two GAL21 potential type
2534!> (??)
2535!> \param gal21_source ...
2536!> \param gal21_dest ...
2537! **************************************************************************************************
2538 SUBROUTINE pair_potential_gal21_copy(gal21_source, gal21_dest)
2539 TYPE(gal21_pot_type), POINTER :: gal21_source, gal21_dest
2540
2541 IF (.NOT. ASSOCIATED(gal21_source)) RETURN
2542 IF (ASSOCIATED(gal21_dest)) CALL pair_potential_gal21_release(gal21_dest)
2543 CALL pair_potential_gal21_create(gal21_dest)
2544 gal21_dest%met1 = gal21_source%met1
2545 gal21_dest%met2 = gal21_source%met2
2546 gal21_dest%epsilon1 = gal21_source%epsilon1
2547 gal21_dest%epsilon2 = gal21_source%epsilon2
2548 gal21_dest%epsilon3 = gal21_source%epsilon3
2549 gal21_dest%bxy1 = gal21_source%bxy1
2550 gal21_dest%bxy2 = gal21_source%bxy2
2551 gal21_dest%bz1 = gal21_source%bz1
2552 gal21_dest%bz2 = gal21_source%bz2
2553 gal21_dest%r1 = gal21_source%r1
2554 gal21_dest%r2 = gal21_source%r2
2555 gal21_dest%a11 = gal21_source%a11
2556 gal21_dest%a12 = gal21_source%a12
2557 gal21_dest%a13 = gal21_source%a13
2558 gal21_dest%a21 = gal21_source%a21
2559 gal21_dest%a22 = gal21_source%a22
2560 gal21_dest%a23 = gal21_source%a23
2561 gal21_dest%a31 = gal21_source%a31
2562 gal21_dest%a32 = gal21_source%a32
2563 gal21_dest%a33 = gal21_source%a33
2564 gal21_dest%a41 = gal21_source%a41
2565 gal21_dest%a42 = gal21_source%a42
2566 gal21_dest%a43 = gal21_source%a43
2567 gal21_dest%AO1 = gal21_source%AO1
2568 gal21_dest%AO2 = gal21_source%AO2
2569 gal21_dest%BO1 = gal21_source%BO1
2570 gal21_dest%BO2 = gal21_source%BO2
2571 gal21_dest%c = gal21_source%c
2572 gal21_dest%AH1 = gal21_source%AH1
2573 gal21_dest%AH2 = gal21_source%AH2
2574 gal21_dest%BH1 = gal21_source%BH1
2575 gal21_dest%BH2 = gal21_source%BH2
2576 ALLOCATE (gal21_dest%gcn(SIZE(gal21_source%gcn)))
2577 gal21_dest%gcn = gal21_source%gcn
2578 gal21_dest%express = gal21_source%express
2579 gal21_dest%rcutsq = gal21_source%rcutsq
2580
2581 END SUBROUTINE pair_potential_gal21_copy
2582
2583! **************************************************************************************************
2584!> \brief Creates the GAL21 potential type
2585!> (??)
2586!> \param gal21 ...
2587! **************************************************************************************************
2588 SUBROUTINE pair_potential_gal21_clean(gal21)
2589 TYPE(gal21_pot_type), POINTER :: gal21
2590
2591 IF (.NOT. ASSOCIATED(gal21)) RETURN
2592 gal21%epsilon1 = 0.0_dp
2593 gal21%epsilon2 = 0.0_dp
2594 gal21%epsilon3 = 0.0_dp
2595 gal21%bxy1 = 0.0_dp
2596 gal21%bxy2 = 0.0_dp
2597 gal21%bz1 = 0.0_dp
2598 gal21%bz2 = 0.0_dp
2599 gal21%r1 = 0.0_dp
2600 gal21%r2 = 0.0_dp
2601 gal21%a11 = 0.0_dp
2602 gal21%a12 = 0.0_dp
2603 gal21%a13 = 0.0_dp
2604 gal21%a21 = 0.0_dp
2605 gal21%a22 = 0.0_dp
2606 gal21%a23 = 0.0_dp
2607 gal21%a31 = 0.0_dp
2608 gal21%a32 = 0.0_dp
2609 gal21%a33 = 0.0_dp
2610 gal21%a41 = 0.0_dp
2611 gal21%a42 = 0.0_dp
2612 gal21%a43 = 0.0_dp
2613 gal21%AO1 = 0.0_dp
2614 gal21%AO2 = 0.0_dp
2615 gal21%BO1 = 0.0_dp
2616 gal21%BO2 = 0.0_dp
2617 gal21%c = 0.0_dp
2618 gal21%AH1 = 0.0_dp
2619 gal21%AH2 = 0.0_dp
2620 gal21%BH1 = 0.0_dp
2621 gal21%BH2 = 0.0_dp
2622 gal21%rcutsq = 0.0_dp
2623 gal21%express = .false.
2624
2625 END SUBROUTINE pair_potential_gal21_clean
2626
2627! **************************************************************************************************
2628!> \brief Destroys the GAL21 potential
2629!> (??)
2630!> \param gal21 ...
2631! **************************************************************************************************
2632 SUBROUTINE pair_potential_gal21_release(gal21)
2633 TYPE(gal21_pot_type), POINTER :: gal21
2634
2635 IF (ASSOCIATED(gal21)) THEN
2636 DEALLOCATE (gal21%gcn)
2637 DEALLOCATE (gal21)
2638 END IF
2639 NULLIFY (gal21)
2640 END SUBROUTINE pair_potential_gal21_release
2641
2642! **************************************************************************************************
2643!> \brief Creates the TABPOT potential type
2644!> \param tab ...
2645!> \author Alex Mironenko, Da Teng 2019-2022
2646! **************************************************************************************************
2647 SUBROUTINE pair_potential_tab_create(tab)
2648 TYPE(tab_pot_type), POINTER :: tab
2649
2650 cpassert(.NOT. ASSOCIATED(tab))
2651 ALLOCATE (tab)
2652 NULLIFY (tab%r, tab%e, tab%f)
2653 CALL pair_potential_tab_clean(tab)
2654 END SUBROUTINE pair_potential_tab_create
2655
2656! **************************************************************************************************
2657!> \brief Copy two TABPOT potential type
2658!> \param tab_source ...
2659!> \param tab_dest ...
2660! **************************************************************************************************
2661 SUBROUTINE pair_potential_tab_copy(tab_source, tab_dest)
2662 TYPE(tab_pot_type), POINTER :: tab_source, tab_dest
2663
2664 IF (.NOT. ASSOCIATED(tab_source)) RETURN
2665 IF (ASSOCIATED(tab_dest)) CALL pair_potential_tab_release(tab_dest)
2666 CALL pair_potential_tab_create(tab_dest)
2667 tab_dest%tabpot_file_name = tab_source%tabpot_file_name
2668 tab_dest%dr = tab_source%dr
2669 tab_dest%rcut = tab_source%rcut
2670 tab_dest%npoints = tab_source%npoints
2671 tab_dest%index = tab_source%index
2672 ! Allocate arrays with the proper size
2673 CALL reallocate(tab_dest%r, 1, tab_dest%npoints)
2674 CALL reallocate(tab_dest%e, 1, tab_dest%npoints)
2675 CALL reallocate(tab_dest%f, 1, tab_dest%npoints)
2676 tab_dest%r = tab_source%r
2677 tab_dest%e = tab_source%e
2678 tab_dest%f = tab_source%f
2679 END SUBROUTINE pair_potential_tab_copy
2680
2681! **************************************************************************************************
2682!> \brief Creates the TABPOT potential type
2683!> \param tab ...
2684! **************************************************************************************************
2685 SUBROUTINE pair_potential_tab_clean(tab)
2686 TYPE(tab_pot_type), POINTER :: tab
2687
2688 IF (.NOT. ASSOCIATED(tab)) RETURN
2689 tab%tabpot_file_name = 'NULL'
2690 tab%dr = 0.0_dp
2691 tab%rcut = 0.0_dp
2692 tab%npoints = 0
2693 tab%index = 0
2694 CALL reallocate(tab%r, 1, tab%npoints)
2695 CALL reallocate(tab%e, 1, tab%npoints)
2696 CALL reallocate(tab%f, 1, tab%npoints)
2697
2698 END SUBROUTINE pair_potential_tab_clean
2699
2700! **************************************************************************************************
2701!> \brief Destroys the TABPOT potential type
2702!> \param tab ...
2703! **************************************************************************************************
2704 SUBROUTINE pair_potential_tab_release(tab)
2705 TYPE(tab_pot_type), POINTER :: tab
2706
2707 IF (ASSOCIATED(tab)) THEN
2708 IF (ASSOCIATED(tab%r)) THEN
2709 DEALLOCATE (tab%r)
2710 END IF
2711 IF (ASSOCIATED(tab%e)) THEN
2712 DEALLOCATE (tab%e)
2713 END IF
2714 IF (ASSOCIATED(tab%f)) THEN
2715 DEALLOCATE (tab%f)
2716 END IF
2717 DEALLOCATE (tab)
2718 END IF
2719 END SUBROUTINE pair_potential_tab_release
2720
2721END MODULE pair_potential_types
2722
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_string_length
Definition kinds.F:57
integer, parameter, public default_path_length
Definition kinds.F:58
Utility routines for the memory handling.
integer, parameter, public sh_sh
integer, parameter, public nosh_nosh
integer, dimension(21), parameter, public list_pot
integer, parameter, public lj_charmm_type
integer, parameter, public allegro_type
integer, parameter, public bm_type
integer, parameter, public gal_type
subroutine, public pair_potential_pp_release(potparm)
Release Data-structure that constains potential parameters.
integer, parameter, public nequip_type
integer, parameter, public wl_type
integer, parameter, public ft_type
integer, parameter, public tab_type
integer, parameter, public ftd_type
integer, parameter, public ip_type
subroutine, public pair_potential_p_release(potparm)
Release Data-structure that constains potential parameters.
integer, parameter, public lj_type
integer, parameter, public deepmd_type
subroutine, public pair_potential_single_copy(potparm_source, potparm_dest)
Copy two potential parameter type.
integer, parameter, public nn_type
integer, parameter, public multi_type
integer, parameter, public quip_type
integer, parameter, public gp_type
subroutine, public pair_potential_single_add(potparm_source, potparm_dest)
Add potential parameter type to an existing potential parameter type Used in case of multiple_potenti...
integer, parameter, public siepmann_type
integer, parameter, public nosh_sh
subroutine, public pair_potential_single_clean(potparm)
Cleans the potential parameter type.
subroutine, public pair_potential_lj_create(lj)
Cleans the LJ potential type.
integer, dimension(2), parameter, public do_potential_single_allocation
subroutine, public compare_pot(pot1, pot2, compare)
compare two different potentials
integer, parameter, public gw_type
subroutine, public pair_potential_reallocate(p, lb1_new, ub1_new, lj, lj_charmm, williams, goodwin, eam, quip, nequip, allegro, bmhft, bmhftd, ipbv, buck4r, buckmo, gp, tersoff, siepmann, gal, gal21, tab, deepmd)
Cleans the potential parameter type.
real(kind=dp), parameter, public not_initialized
subroutine, public pair_potential_pp_create(potparm, nkinds)
Data-structure that constains potential parameters.
integer, dimension(3), parameter, public list_sh_type
integer, dimension(2), parameter, public no_potential_single_allocation
integer, parameter, public b4_type
integer, parameter, public gal21_type
integer, dimension(2), public potential_single_allocation
integer, parameter, public ea_type
integer, parameter, public tersoff_type
routines for handling splines_types
subroutine, public spline_data_p_copy(spl_p_source, spl_p_dest)
Copy Data-structure of spline_data_p_type.
subroutine, public spline_factor_release(spline_factor)
releases spline_factor
subroutine, public spline_data_p_release(spl_p)
releases spline_data_p
subroutine, public spline_factor_copy(spline_factor_source, spline_factor_dest)
releases spline_factor