(git:ec11232)
Loading...
Searching...
No Matches
external_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!> \brief Definition of the atomic potential types.
10!> \par History
11!> GT, 22.09.2002: added elp_potential_types
12!> \author Matthias Krack (04.07.2000)
13! **************************************************************************************************
15
16 USE ao_util, ONLY: exp_radius
17 USE bibliography, ONLY: goedecker1996,&
19 krack2000,&
20 krack2005,&
21 cite_reference
35 USE input_val_types, ONLY: val_get,&
37 USE kinds, ONLY: default_path_length,&
39 dp
40 USE mathconstants, ONLY: dfac,&
41 fac,&
42 pi,&
43 rootpi
44 USE mathlib, ONLY: symmetrize_matrix
47 USE orbital_pointers, ONLY: co,&
48 coset,&
50 nco,&
51 ncoset,&
52 nso
54 USE periodic_table, ONLY: ptable
55 USE string_utilities, ONLY: remove_word,&
57#include "../base/base_uses.f90"
58
59 IMPLICIT NONE
60
61 PRIVATE
62
63 ! Global parameters
64
65 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'external_potential_types'
66
67 ! Define the all-electron potential type
68 ! Literature: M. Krack and M. Parrinello,
69 ! Phys. Chem. Chem. Phys. 2, 2105 (2000)
71 !MK PRIVATE
72 CHARACTER(LEN=default_string_length) :: name = ""
73 CHARACTER(LEN=default_string_length), &
74 DIMENSION(2) :: description = ["All-electron potential ", &
75 "Krack, Parrinello, PCCP 2, 2105 (2000)"]
76 REAL(kind=dp) :: alpha_core_charge = 0.0_dp, &
77 ccore_charge = 0.0_dp, &
78 core_charge_radius = 0.0_dp, &
79 zeff = 0.0_dp, zeff_correction = 0.0_dp
80 INTEGER :: z = 0
81 INTEGER, DIMENSION(:), POINTER :: elec_conf => null()
82 END TYPE all_potential_type
83
84 ! Define the effective charge & inducible dipole potential type (for Fist)
86 PRIVATE
87 CHARACTER(LEN=default_string_length) :: name = ""
88 CHARACTER(LEN=default_string_length), &
89 DIMENSION(1) :: description = "Effective charge and inducible dipole potential"
90 REAL(kind=dp) :: apol = 0.0_dp, cpol = 0.0_dp, mm_radius = 0.0_dp, qeff = 0.0_dp, &
91 qmmm_corr_radius = 0.0_dp, qmmm_radius = 0.0_dp
92
93 END TYPE fist_potential_type
94
95 ! Local potential type
96 ! V(r) = SUM_i exp(0.5*(r/rci)**2) * ( C1i + C2i (r/rci)**2 + C3i (r/rci)**4 ...)
97 ! alpha = 0.5/rci**2
99 !PRIVATE
100 CHARACTER(LEN=default_string_length) :: name = ""
101 CHARACTER(LEN=default_string_length), &
102 DIMENSION(4) :: description = "Local short-range pseudopotential"
103 INTEGER :: ngau = 0, npol = 0
104 REAL(kind=dp) :: radius = 0.0_dp
105 REAL(kind=dp), DIMENSION(:), POINTER :: alpha => null()
106 REAL(kind=dp), DIMENSION(:, :), POINTER :: cval => null()
107 END TYPE local_potential_type
108
109 ! Define the GTH potential type
110 ! Literature: - S. Goedecker, M. Teter and J. Hutter,
111 ! Phys. Rev. B 54, 1703 (1996)
112 ! - C. Hartwigsen, S. Goedecker and J. Hutter,
113 ! Phys. Rev. B 58, 3641 (1998)
114 ! - M. Krack,
115 ! Theor. Chem. Acc. 114, 145 (2005)
117 CHARACTER(LEN=default_string_length) :: name = ""
118 CHARACTER(LEN=default_string_length) :: aliases = ""
119 CHARACTER(LEN=default_string_length), &
120 DIMENSION(4) :: description = ["Goedecker-Teter-Hutter pseudopotential", &
121 "Goedecker et al., PRB 54, 1703 (1996) ", &
122 "Hartwigsen et al., PRB 58, 3641 (1998)", &
123 "Krack, TCA 114, 145 (2005) "]
124 REAL(kind=dp) :: alpha_core_charge = 0.0_dp, &
125 alpha_ppl = 0.0_dp, &
126 ccore_charge = 0.0_dp, &
127 cerf_ppl = 0.0_dp, &
128 zeff = 0.0_dp, &
129 core_charge_radius = 0.0_dp, &
130 ppl_radius = 0.0_dp, &
131 ppnl_radius = 0.0_dp, &
132 zeff_correction = 0.0_dp
133 INTEGER :: lppnl = 0, &
134 lprj_ppnl_max = 0, &
135 nexp_ppl = 0, &
136 nppnl = 0, &
137 nprj_ppnl_max = 0, z = 0
138 REAL(kind=dp), DIMENSION(:), POINTER :: alpha_ppnl => null(), &
139 cexp_ppl => null()
140 INTEGER, DIMENSION(:), POINTER :: elec_conf => null()
141 ! Non-local projectors
142 INTEGER, DIMENSION(:), POINTER :: nprj_ppnl => null()
143 REAL(kind=dp), DIMENSION(:, :), POINTER :: cprj => null(), &
144 cprj_ppnl => null(), &
145 vprj_ppnl => null(), &
146 wprj_ppnl => null()
147 REAL(kind=dp), DIMENSION(:, :, :), POINTER :: hprj_ppnl => null(), &
148 kprj_ppnl => null()
149 ! Type extensions
150 ! Spin-orbit coupling (SOC) parameters
151 LOGICAL :: soc = .false.
152 ! NLCC
153 LOGICAL :: nlcc = .false.
154 INTEGER :: nexp_nlcc = 0
155 REAL(kind=dp), DIMENSION(:), POINTER :: alpha_nlcc => null()
156 INTEGER, DIMENSION(:), POINTER :: nct_nlcc => null()
157 REAL(kind=dp), DIMENSION(:, :), POINTER :: cval_nlcc => null()
158 ! LSD potential
159 LOGICAL :: lsdpot = .false.
160 INTEGER :: nexp_lsd = 0
161 REAL(kind=dp), DIMENSION(:), POINTER :: alpha_lsd => null()
162 INTEGER, DIMENSION(:), POINTER :: nct_lsd => null()
163 REAL(kind=dp), DIMENSION(:, :), POINTER :: cval_lsd => null()
164 ! Extended local potential
165 LOGICAL :: lpotextended = .false.
166 INTEGER :: nexp_lpot = 0
167 REAL(kind=dp), DIMENSION(:), POINTER :: alpha_lpot => null()
168 INTEGER, DIMENSION(:), POINTER :: nct_lpot => null()
169 REAL(kind=dp), DIMENSION(:, :), POINTER :: cval_lpot => null()
170 ! monovalent pseudopotential
171 LOGICAL :: monovalent = .false.
172 END TYPE gth_potential_type
173
175 CHARACTER(LEN=default_string_length) :: name = ""
176 CHARACTER(LEN=default_string_length) :: aliases = ""
177 CHARACTER(LEN=default_string_length), &
178 DIMENSION(4) :: description = ["Separable Gaussian pseudopotential ", &
179 "M. Pelissier, N. Komiha, J.P. Daudey, JCC, 9, 298 (1988)", &
180 "create from ", &
181 " "]
182 ! CHARGE
183 INTEGER :: z = 0
184 REAL(kind=dp) :: zeff = 0.0_dp, &
185 zeff_correction = 0.0_dp
186 REAL(kind=dp) :: alpha_core_charge = 0.0_dp, &
187 ccore_charge = 0.0_dp, &
188 core_charge_radius = 0.0_dp
189 REAL(kind=dp) :: ppl_radius = 0.0_dp, ppnl_radius = 0.0_dp
190 INTEGER, DIMENSION(:), POINTER :: elec_conf => null()
191 ! LOCAL
192 LOGICAL :: ecp_local = .false.
193 INTEGER :: n_local = 0
194 REAL(kind=dp), DIMENSION(:), POINTER :: a_local => null()
195 REAL(kind=dp), DIMENSION(:), POINTER :: c_local => null()
196 ! ECP local
197 INTEGER :: nloc = 0 ! # terms
198 INTEGER, DIMENSION(1:10) :: nrloc = 0 ! r**(n-2)
199 REAL(dp), DIMENSION(1:10) :: aloc = 0.0_dp ! coefficient
200 REAL(dp), DIMENSION(1:10) :: bloc = 0.0_dp ! exponent
201 ! ECP semi-local
202 LOGICAL :: ecp_semi_local = .false.
203 INTEGER :: sl_lmax = 0
204 INTEGER, DIMENSION(0:10) :: npot = 0 ! # terms
205 INTEGER, DIMENSION(1:15, 0:10) :: nrpot = 0 ! r**(n-2)
206 REAL(dp), DIMENSION(1:15, 0:10) :: apot = 0.0_dp ! coefficient
207 REAL(dp), DIMENSION(1:15, 0:10) :: bpot = 0.0_dp ! exponent
208 ! NON-LOCAL
209 INTEGER :: n_nonlocal = 0
210 INTEGER :: nppnl = 0
211 INTEGER :: lmax = -1
212 LOGICAL, DIMENSION(0:5) :: is_nonlocal = .false.
213 REAL(kind=dp), DIMENSION(:), POINTER :: a_nonlocal => null()
214 REAL(kind=dp), DIMENSION(:, :), POINTER :: h_nonlocal => null()
215 REAL(kind=dp), DIMENSION(:, :, :), POINTER :: c_nonlocal => null()
216 REAL(kind=dp), DIMENSION(:, :), POINTER :: cprj_ppnl => null()
217 REAL(kind=dp), DIMENSION(:), POINTER :: vprj_ppnl => null()
218 ! NLCC
219 LOGICAL :: has_nlcc = .false.
220 INTEGER :: n_nlcc = 0
221 REAL(kind=dp), DIMENSION(:), POINTER :: a_nlcc => null()
222 REAL(kind=dp), DIMENSION(:), POINTER :: c_nlcc => null()
223 END TYPE sgp_potential_type
224
225 TYPE all_potential_p_type
226 TYPE(all_potential_type), POINTER :: all_potential => null()
227 END TYPE all_potential_p_type
228
230 TYPE(gth_potential_type), POINTER :: gth_potential => null()
231 END TYPE gth_potential_p_type
232
233 TYPE local_potential_p_type
234 TYPE(local_potential_type), POINTER :: local_potential => null()
235 END TYPE local_potential_p_type
236
238 TYPE(sgp_potential_type), POINTER :: sgp_potential => null()
239 END TYPE sgp_potential_p_type
240
241 ! Public subroutines
242 PUBLIC :: allocate_potential, &
251
252 ! Public data types
253
254 PUBLIC :: all_potential_type, &
259 PUBLIC :: gth_potential_p_type, &
261
263 MODULE PROCEDURE allocate_all_potential, &
264 allocate_fist_potential, &
265 allocate_local_potential, &
266 allocate_gth_potential, &
267 allocate_sgp_potential
268 END INTERFACE
269
271 MODULE PROCEDURE deallocate_all_potential, &
272 deallocate_fist_potential, &
273 deallocate_local_potential, &
274 deallocate_sgp_potential, &
275 deallocate_gth_potential
276 END INTERFACE
277
279 MODULE PROCEDURE get_all_potential, &
280 get_fist_potential, &
281 get_local_potential, &
282 get_gth_potential, &
283 get_sgp_potential
284 END INTERFACE
285
287 MODULE PROCEDURE init_all_potential, &
288 init_gth_potential, &
289 init_sgp_potential
290 END INTERFACE
291
293 MODULE PROCEDURE read_all_potential, &
294 read_local_potential, &
295 read_gth_potential
296 END INTERFACE
297
299 MODULE PROCEDURE set_all_potential, &
300 set_fist_potential, &
301 set_local_potential, &
302 set_gth_potential, &
303 set_sgp_potential
304 END INTERFACE
305
307 MODULE PROCEDURE write_all_potential, &
308 write_local_potential, &
309 write_gth_potential, &
310 write_sgp_potential
311 END INTERFACE
312
314 MODULE PROCEDURE copy_all_potential, &
315 copy_gth_potential, &
316 copy_sgp_potential
317 END INTERFACE
318
319CONTAINS
320
321! **************************************************************************************************
322!> \brief Allocate an atomic all-electron potential data set.
323!> \param potential ...
324!> \date 25.07.2000,
325!> \author MK
326!> \version 1.0
327! **************************************************************************************************
328 SUBROUTINE allocate_all_potential(potential)
329 TYPE(all_potential_type), INTENT(INOUT), POINTER :: potential
330
331 IF (ASSOCIATED(potential)) CALL deallocate_potential(potential)
332
333 ALLOCATE (potential)
334
335 END SUBROUTINE allocate_all_potential
336
337! **************************************************************************************************
338!> \brief Allocate an effective charge and inducible dipole potential data set.
339!> \param potential ...
340!> \date 05.03.2010
341!> \author Toon.Verstraelen@gmail.com
342! **************************************************************************************************
343 SUBROUTINE allocate_fist_potential(potential)
344 TYPE(fist_potential_type), INTENT(INOUT), POINTER :: potential
345
346 IF (ASSOCIATED(potential)) CALL deallocate_potential(potential)
347
348 ALLOCATE (potential)
349
350 END SUBROUTINE allocate_fist_potential
351
352! **************************************************************************************************
353!> \brief Allocate an atomic local potential data set.
354!> \param potential ...
355!> \date 24.01.2014
356!> \author JGH
357!> \version 1.0
358! **************************************************************************************************
359 SUBROUTINE allocate_local_potential(potential)
360 TYPE(local_potential_type), INTENT(INOUT), POINTER :: potential
361
362 IF (ASSOCIATED(potential)) CALL deallocate_potential(potential)
363
364 ALLOCATE (potential)
365
366 END SUBROUTINE allocate_local_potential
367
368! **************************************************************************************************
369!> \brief Allocate an atomic GTH potential data set.
370!> \param potential ...
371!> \date 25.07.2000
372!> \author MK
373!> \version 1.0
374! **************************************************************************************************
375 SUBROUTINE allocate_gth_potential(potential)
376 TYPE(gth_potential_type), INTENT(INOUT), POINTER :: potential
377
378 IF (ASSOCIATED(potential)) CALL deallocate_potential(potential)
379
380 ALLOCATE (potential)
381
382 END SUBROUTINE allocate_gth_potential
383
384! **************************************************************************************************
385!> \brief Allocate an atomic SGP potential data set.
386!> \param potential ...
387!> \version 1.0
388! **************************************************************************************************
389 SUBROUTINE allocate_sgp_potential(potential)
390 TYPE(sgp_potential_type), INTENT(INOUT), POINTER :: potential
391
392 IF (ASSOCIATED(potential)) CALL deallocate_potential(potential)
393
394 ALLOCATE (potential)
395
396 END SUBROUTINE allocate_sgp_potential
397! **************************************************************************************************
398!> \brief Deallocate an atomic all-electron potential data set.
399!> \param potential ...
400!> \date 03.11.2000
401!> \author MK
402!> \version 1.0
403! **************************************************************************************************
404 SUBROUTINE deallocate_all_potential(potential)
405 TYPE(all_potential_type), POINTER :: potential
406
407 IF (.NOT. ASSOCIATED(potential)) THEN
408 cpabort("The pointer potential is not associated.")
409 END IF
410
411 DEALLOCATE (potential%elec_conf)
412 DEALLOCATE (potential)
413
414 END SUBROUTINE deallocate_all_potential
415
416! **************************************************************************************************
417!> \brief Deallocate an effective charge and inducible dipole potential data set.
418!> \param potential ...
419!> \date 05.03.2010
420!> \author Toon.Verstraelen@gmail.com
421! **************************************************************************************************
422 SUBROUTINE deallocate_fist_potential(potential)
423 TYPE(fist_potential_type), POINTER :: potential
424
425 IF (.NOT. ASSOCIATED(potential)) THEN
426 cpabort("The pointer potential is not associated.")
427 END IF
428
429 ! Nothing exciting here yet.
430 DEALLOCATE (potential)
431
432 END SUBROUTINE deallocate_fist_potential
433
434! **************************************************************************************************
435!> \brief Deallocate an atomic local potential data set.
436!> \param potential ...
437!> \date 24.01.2014
438!> \author JGH
439!> \version 1.0
440! **************************************************************************************************
441 SUBROUTINE deallocate_local_potential(potential)
442 TYPE(local_potential_type), POINTER :: potential
443
444 IF (.NOT. ASSOCIATED(potential)) THEN
445 cpabort("The pointer potential is not associated.")
446 END IF
447
448 IF (ASSOCIATED(potential%alpha)) THEN
449 DEALLOCATE (potential%alpha)
450 END IF
451 IF (ASSOCIATED(potential%cval)) THEN
452 DEALLOCATE (potential%cval)
453 END IF
454
455 DEALLOCATE (potential)
456
457 END SUBROUTINE deallocate_local_potential
458
459! **************************************************************************************************
460!> \brief Deallocate an atomic GTH potential data set.
461!> \param potential ...
462!> \date 03.11.2000
463!> \author MK
464!> \version 1.0
465! **************************************************************************************************
466 SUBROUTINE deallocate_gth_potential(potential)
467 TYPE(gth_potential_type), POINTER :: potential
468
469 IF (.NOT. ASSOCIATED(potential)) THEN
470 cpabort("The pointer potential is not associated.")
471 END IF
472
473 DEALLOCATE (potential%elec_conf)
474 ! Deallocate the parameters of the local part
475
476 IF (ASSOCIATED(potential%cexp_ppl)) THEN
477 DEALLOCATE (potential%cexp_ppl)
478 END IF
479
480 ! Deallocate the parameters of the non-local part
481 IF (ASSOCIATED(potential%alpha_ppnl)) THEN
482 DEALLOCATE (potential%alpha_ppnl)
483 DEALLOCATE (potential%cprj)
484 DEALLOCATE (potential%cprj_ppnl)
485 DEALLOCATE (potential%hprj_ppnl)
486 DEALLOCATE (potential%kprj_ppnl)
487 DEALLOCATE (potential%nprj_ppnl)
488 DEALLOCATE (potential%vprj_ppnl)
489 DEALLOCATE (potential%wprj_ppnl)
490 END IF
491
492 IF (ASSOCIATED(potential%alpha_lpot)) THEN
493 DEALLOCATE (potential%alpha_lpot)
494 DEALLOCATE (potential%nct_lpot)
495 DEALLOCATE (potential%cval_lpot)
496 END IF
497
498 IF (ASSOCIATED(potential%alpha_lsd)) THEN
499 DEALLOCATE (potential%alpha_lsd)
500 DEALLOCATE (potential%nct_lsd)
501 DEALLOCATE (potential%cval_lsd)
502 END IF
503
504 IF (ASSOCIATED(potential%alpha_nlcc)) THEN
505 DEALLOCATE (potential%alpha_nlcc)
506 DEALLOCATE (potential%nct_nlcc)
507 DEALLOCATE (potential%cval_nlcc)
508 END IF
509
510 DEALLOCATE (potential)
511
512 END SUBROUTINE deallocate_gth_potential
513
514! **************************************************************************************************
515!> \brief Deallocate an atomic SGP potential data set.
516!> \param potential ...
517! **************************************************************************************************
518 SUBROUTINE deallocate_sgp_potential(potential)
519 TYPE(sgp_potential_type), POINTER :: potential
520
521 IF (.NOT. ASSOCIATED(potential)) THEN
522 cpabort("The pointer potential is not associated.")
523 END IF
524
525 IF (ASSOCIATED(potential%elec_conf)) THEN
526 DEALLOCATE (potential%elec_conf)
527 END IF
528 IF (ASSOCIATED(potential%a_local)) THEN
529 DEALLOCATE (potential%a_local)
530 END IF
531 IF (ASSOCIATED(potential%c_local)) THEN
532 DEALLOCATE (potential%c_local)
533 END IF
534
535 IF (ASSOCIATED(potential%a_nonlocal)) THEN
536 DEALLOCATE (potential%a_nonlocal)
537 END IF
538 IF (ASSOCIATED(potential%h_nonlocal)) THEN
539 DEALLOCATE (potential%h_nonlocal)
540 END IF
541 IF (ASSOCIATED(potential%c_nonlocal)) THEN
542 DEALLOCATE (potential%c_nonlocal)
543 END IF
544 IF (ASSOCIATED(potential%cprj_ppnl)) THEN
545 DEALLOCATE (potential%cprj_ppnl)
546 END IF
547 IF (ASSOCIATED(potential%vprj_ppnl)) THEN
548 DEALLOCATE (potential%vprj_ppnl)
549 END IF
550
551 IF (ASSOCIATED(potential%a_nlcc)) THEN
552 DEALLOCATE (potential%a_nlcc)
553 END IF
554 IF (ASSOCIATED(potential%c_nlcc)) THEN
555 DEALLOCATE (potential%c_nlcc)
556 END IF
557
558 DEALLOCATE (potential)
559
560 END SUBROUTINE deallocate_sgp_potential
561
562! **************************************************************************************************
563!> \brief Get attributes of an all-electron potential data set.
564!> \param potential ...
565!> \param name ...
566!> \param alpha_core_charge ...
567!> \param ccore_charge ...
568!> \param core_charge_radius ...
569!> \param z ...
570!> \param zeff ...
571!> \param zeff_correction ...
572!> \param elec_conf ...
573!> \date 11.01.2002
574!> \author MK
575!> \version 1.0
576! **************************************************************************************************
577 SUBROUTINE get_all_potential(potential, name, alpha_core_charge, &
578 ccore_charge, core_charge_radius, z, zeff, &
579 zeff_correction, elec_conf)
580 TYPE(all_potential_type), INTENT(IN) :: potential
581 CHARACTER(LEN=default_string_length), &
582 INTENT(OUT), OPTIONAL :: name
583 REAL(KIND=dp), INTENT(OUT), OPTIONAL :: alpha_core_charge, ccore_charge, &
584 core_charge_radius
585 INTEGER, INTENT(OUT), OPTIONAL :: z
586 REAL(KIND=dp), INTENT(OUT), OPTIONAL :: zeff, zeff_correction
587 INTEGER, DIMENSION(:), OPTIONAL, POINTER :: elec_conf
588
589 IF (PRESENT(name)) name = potential%name
590 IF (PRESENT(alpha_core_charge)) &
591 alpha_core_charge = potential%alpha_core_charge
592 IF (PRESENT(ccore_charge)) ccore_charge = potential%ccore_charge
593 IF (PRESENT(core_charge_radius)) &
594 core_charge_radius = potential%core_charge_radius
595 IF (PRESENT(z)) z = potential%z
596 IF (PRESENT(zeff)) zeff = potential%zeff
597 IF (PRESENT(zeff_correction)) zeff_correction = potential%zeff_correction
598 IF (PRESENT(elec_conf)) elec_conf => potential%elec_conf
599
600 END SUBROUTINE get_all_potential
601
602! **************************************************************************************************
603!> \brief Get attributes of an effective point charge and inducible dipole
604!> potential.
605!> \param potential ...
606!> \param name ...
607!> \param apol ...
608!> \param cpol ...
609!> \param mm_radius ...
610!> \param qeff ...
611!> \param qmmm_corr_radius ...
612!> \param qmmm_radius ...
613!> \date 05.03-2010
614!> \author Toon.Verstraelen@UGent.be
615! **************************************************************************************************
616 ELEMENTAL SUBROUTINE get_fist_potential(potential, name, apol, cpol, mm_radius, qeff, &
617 qmmm_corr_radius, qmmm_radius)
618 TYPE(fist_potential_type), INTENT(IN) :: potential
619 CHARACTER(LEN=default_string_length), &
620 INTENT(OUT), OPTIONAL :: name
621 REAL(kind=dp), INTENT(OUT), OPTIONAL :: apol, cpol, mm_radius, qeff, &
622 qmmm_corr_radius, qmmm_radius
623
624 IF (PRESENT(name)) name = potential%name
625 IF (PRESENT(apol)) apol = potential%apol
626 IF (PRESENT(cpol)) cpol = potential%cpol
627 IF (PRESENT(mm_radius)) mm_radius = potential%mm_radius
628 IF (PRESENT(qeff)) qeff = potential%qeff
629 IF (PRESENT(qmmm_corr_radius)) qmmm_corr_radius = potential%qmmm_corr_radius
630 IF (PRESENT(qmmm_radius)) qmmm_radius = potential%qmmm_radius
631
632 END SUBROUTINE get_fist_potential
633
634! **************************************************************************************************
635!> \brief Get attributes of an atomic local potential data set.
636!> \param potential ...
637!> \param name ...
638!> \param ngau ...
639!> \param npol ...
640!> \param alpha ...
641!> \param cval ...
642!> \param radius ...
643!> \date 24.01.2014
644!> \author JGH
645!> \version 1.0
646! **************************************************************************************************
647 SUBROUTINE get_local_potential(potential, name, ngau, npol, alpha, cval, radius)
648 TYPE(local_potential_type), INTENT(IN) :: potential
649 CHARACTER(LEN=default_string_length), &
650 INTENT(OUT), OPTIONAL :: name
651 INTEGER, INTENT(OUT), OPTIONAL :: ngau, npol
652 REAL(kind=dp), DIMENSION(:), OPTIONAL, POINTER :: alpha
653 REAL(kind=dp), DIMENSION(:, :), OPTIONAL, POINTER :: cval
654 REAL(kind=dp), INTENT(OUT), OPTIONAL :: radius
655
656 IF (PRESENT(name)) name = potential%name
657 IF (PRESENT(ngau)) ngau = potential%ngau
658 IF (PRESENT(npol)) npol = potential%npol
659 IF (PRESENT(alpha)) alpha => potential%alpha
660 IF (PRESENT(cval)) cval => potential%cval
661 IF (PRESENT(radius)) radius = potential%radius
662
663 END SUBROUTINE get_local_potential
664
665! **************************************************************************************************
666!> \brief Get attributes of a GTH potential data set.
667!> \param potential ...
668!> \param name ...
669!> \param aliases ...
670!> \param alpha_core_charge ...
671!> \param alpha_ppl ...
672!> \param ccore_charge ...
673!> \param cerf_ppl ...
674!> \param core_charge_radius ...
675!> \param ppl_radius ...
676!> \param ppnl_radius ...
677!> \param lppnl ...
678!> \param lprj_ppnl_max ...
679!> \param nexp_ppl ...
680!> \param nppnl ...
681!> \param nprj_ppnl_max ...
682!> \param z ...
683!> \param zeff ...
684!> \param zeff_correction ...
685!> \param ppl_present ...
686!> \param ppnl_present ...
687!> \param soc_present ...
688!> \param alpha_ppnl ...
689!> \param cexp_ppl ...
690!> \param elec_conf ...
691!> \param nprj_ppnl ...
692!> \param cprj ...
693!> \param cprj_ppnl ...
694!> \param vprj_ppnl ...
695!> \param wprj_ppnl ...
696!> \param hprj_ppnl ...
697!> \param kprj_ppnl ...
698!> \param lpot_present ...
699!> \param nexp_lpot ...
700!> \param alpha_lpot ...
701!> \param nct_lpot ...
702!> \param cval_lpot ...
703!> \param lsd_present ...
704!> \param nexp_lsd ...
705!> \param alpha_lsd ...
706!> \param nct_lsd ...
707!> \param cval_lsd ...
708!> \param nlcc_present ...
709!> \param nexp_nlcc ...
710!> \param alpha_nlcc ...
711!> \param nct_nlcc ...
712!> \param cval_nlcc ...
713!> \param monovalent ...
714!> \date 11.01.2002
715!> \author MK
716!> \version 1.0
717! **************************************************************************************************
718 SUBROUTINE get_gth_potential(potential, name, aliases, alpha_core_charge, &
719 alpha_ppl, ccore_charge, cerf_ppl, &
720 core_charge_radius, ppl_radius, ppnl_radius, &
721 lppnl, lprj_ppnl_max, nexp_ppl, nppnl, &
722 nprj_ppnl_max, z, zeff, zeff_correction, &
723 ppl_present, ppnl_present, soc_present, &
724 alpha_ppnl, cexp_ppl, elec_conf, nprj_ppnl, cprj, &
725 cprj_ppnl, vprj_ppnl, wprj_ppnl, hprj_ppnl, kprj_ppnl, &
726 lpot_present, nexp_lpot, alpha_lpot, nct_lpot, cval_lpot, &
727 lsd_present, nexp_lsd, alpha_lsd, nct_lsd, cval_lsd, &
728 nlcc_present, nexp_nlcc, alpha_nlcc, nct_nlcc, cval_nlcc, &
729 monovalent)
730
731 TYPE(gth_potential_type), INTENT(IN) :: potential
732 CHARACTER(LEN=default_string_length), &
733 INTENT(OUT), OPTIONAL :: name, aliases
734 REAL(kind=dp), INTENT(OUT), OPTIONAL :: alpha_core_charge, alpha_ppl, &
735 ccore_charge, cerf_ppl, &
736 core_charge_radius, ppl_radius, &
737 ppnl_radius
738 INTEGER, INTENT(OUT), OPTIONAL :: lppnl, lprj_ppnl_max, nexp_ppl, nppnl, &
739 nprj_ppnl_max, z
740 REAL(kind=dp), INTENT(OUT), OPTIONAL :: zeff, zeff_correction
741 LOGICAL, INTENT(OUT), OPTIONAL :: ppl_present, ppnl_present, soc_present
742 REAL(kind=dp), DIMENSION(:), OPTIONAL, POINTER :: alpha_ppnl, cexp_ppl
743 INTEGER, DIMENSION(:), OPTIONAL, POINTER :: elec_conf, nprj_ppnl
744 REAL(kind=dp), DIMENSION(:, :), OPTIONAL, POINTER :: cprj, cprj_ppnl, vprj_ppnl, wprj_ppnl
745 REAL(kind=dp), DIMENSION(:, :, :), OPTIONAL, &
746 POINTER :: hprj_ppnl, kprj_ppnl
747 LOGICAL, INTENT(OUT), OPTIONAL :: lpot_present
748 INTEGER, INTENT(OUT), OPTIONAL :: nexp_lpot
749 REAL(kind=dp), DIMENSION(:), OPTIONAL, POINTER :: alpha_lpot
750 INTEGER, DIMENSION(:), OPTIONAL, POINTER :: nct_lpot
751 REAL(kind=dp), DIMENSION(:, :), OPTIONAL, POINTER :: cval_lpot
752 LOGICAL, INTENT(OUT), OPTIONAL :: lsd_present
753 INTEGER, INTENT(OUT), OPTIONAL :: nexp_lsd
754 REAL(kind=dp), DIMENSION(:), OPTIONAL, POINTER :: alpha_lsd
755 INTEGER, DIMENSION(:), OPTIONAL, POINTER :: nct_lsd
756 REAL(kind=dp), DIMENSION(:, :), OPTIONAL, POINTER :: cval_lsd
757 LOGICAL, INTENT(OUT), OPTIONAL :: nlcc_present
758 INTEGER, INTENT(OUT), OPTIONAL :: nexp_nlcc
759 REAL(kind=dp), DIMENSION(:), OPTIONAL, POINTER :: alpha_nlcc
760 INTEGER, DIMENSION(:), OPTIONAL, POINTER :: nct_nlcc
761 REAL(kind=dp), DIMENSION(:, :), OPTIONAL, POINTER :: cval_nlcc
762 LOGICAL, INTENT(OUT), OPTIONAL :: monovalent
763
764 IF (PRESENT(name)) name = potential%name
765 IF (PRESENT(aliases)) aliases = potential%aliases
766 IF (PRESENT(alpha_core_charge)) &
767 alpha_core_charge = potential%alpha_core_charge
768 IF (PRESENT(alpha_ppl)) alpha_ppl = potential%alpha_ppl
769 IF (PRESENT(ccore_charge)) ccore_charge = potential%ccore_charge
770 IF (PRESENT(cerf_ppl)) cerf_ppl = potential%cerf_ppl
771 IF (PRESENT(core_charge_radius)) &
772 core_charge_radius = potential%core_charge_radius
773 IF (PRESENT(ppl_radius)) ppl_radius = potential%ppl_radius
774 IF (PRESENT(ppnl_radius)) ppnl_radius = potential%ppnl_radius
775 IF (PRESENT(soc_present)) soc_present = potential%soc
776 IF (PRESENT(lppnl)) lppnl = potential%lppnl
777 IF (PRESENT(lprj_ppnl_max)) lprj_ppnl_max = potential%lprj_ppnl_max
778 IF (PRESENT(nexp_ppl)) nexp_ppl = potential%nexp_ppl
779 IF (PRESENT(nppnl)) nppnl = potential%nppnl
780 IF (PRESENT(nprj_ppnl_max)) nprj_ppnl_max = potential%nprj_ppnl_max
781 IF (PRESENT(z)) z = potential%z
782 IF (PRESENT(zeff)) zeff = potential%zeff
783 IF (PRESENT(zeff_correction)) zeff_correction = potential%zeff_correction
784 IF (PRESENT(ppl_present)) ppl_present = (potential%nexp_ppl > 0)
785 IF (PRESENT(ppnl_present)) ppnl_present = (potential%nppnl > 0)
786 IF (PRESENT(alpha_ppnl)) alpha_ppnl => potential%alpha_ppnl
787 IF (PRESENT(cexp_ppl)) cexp_ppl => potential%cexp_ppl
788 IF (PRESENT(elec_conf)) elec_conf => potential%elec_conf
789 IF (PRESENT(nprj_ppnl)) nprj_ppnl => potential%nprj_ppnl
790 IF (PRESENT(cprj)) cprj => potential%cprj
791 IF (PRESENT(cprj_ppnl)) cprj_ppnl => potential%cprj_ppnl
792 IF (PRESENT(hprj_ppnl)) hprj_ppnl => potential%hprj_ppnl
793 IF (PRESENT(kprj_ppnl)) kprj_ppnl => potential%kprj_ppnl
794 IF (PRESENT(vprj_ppnl)) vprj_ppnl => potential%vprj_ppnl
795 IF (PRESENT(wprj_ppnl)) wprj_ppnl => potential%wprj_ppnl
796
797 IF (PRESENT(lpot_present)) lpot_present = potential%lpotextended
798 IF (PRESENT(nexp_lpot)) nexp_lpot = potential%nexp_lpot
799 IF (PRESENT(alpha_lpot)) alpha_lpot => potential%alpha_lpot
800 IF (PRESENT(nct_lpot)) nct_lpot => potential%nct_lpot
801 IF (PRESENT(cval_lpot)) cval_lpot => potential%cval_lpot
802
803 IF (PRESENT(lsd_present)) lsd_present = potential%lsdpot
804 IF (PRESENT(nexp_lsd)) nexp_lsd = potential%nexp_lsd
805 IF (PRESENT(alpha_lsd)) alpha_lsd => potential%alpha_lsd
806 IF (PRESENT(nct_lsd)) nct_lsd => potential%nct_lsd
807 IF (PRESENT(cval_lsd)) cval_lsd => potential%cval_lsd
808
809 IF (PRESENT(nlcc_present)) nlcc_present = potential%nlcc
810 IF (PRESENT(nexp_nlcc)) nexp_nlcc = potential%nexp_nlcc
811 IF (PRESENT(alpha_nlcc)) alpha_nlcc => potential%alpha_nlcc
812 IF (PRESENT(nct_nlcc)) nct_nlcc => potential%nct_nlcc
813 IF (PRESENT(cval_nlcc)) cval_nlcc => potential%cval_nlcc
814
815 IF (PRESENT(monovalent)) monovalent = potential%monovalent
816
817 END SUBROUTINE get_gth_potential
818
819! **************************************************************************************************
820!> \brief ...
821!> \param potential ...
822!> \param name ...
823!> \param description ...
824!> \param aliases ...
825!> \param elec_conf ...
826!> \param z ...
827!> \param zeff ...
828!> \param zeff_correction ...
829!> \param alpha_core_charge ...
830!> \param ccore_charge ...
831!> \param core_charge_radius ...
832!> \param ppl_radius ...
833!> \param ppnl_radius ...
834!> \param ppl_present ...
835!> \param ppnl_present ...
836!> \param ppsl_present ...
837!> \param ecp_local ...
838!> \param n_local ...
839!> \param a_local ...
840!> \param c_local ...
841!> \param nloc ...
842!> \param nrloc ...
843!> \param aloc ...
844!> \param bloc ...
845!> \param ecp_semi_local ...
846!> \param sl_lmax ...
847!> \param npot ...
848!> \param nrpot ...
849!> \param apot ...
850!> \param bpot ...
851!> \param n_nonlocal ...
852!> \param nppnl ...
853!> \param lmax ...
854!> \param is_nonlocal ...
855!> \param a_nonlocal ...
856!> \param h_nonlocal ...
857!> \param c_nonlocal ...
858!> \param cprj_ppnl ...
859!> \param vprj_ppnl ...
860!> \param has_nlcc ...
861!> \param n_nlcc ...
862!> \param a_nlcc ...
863!> \param c_nlcc ...
864! **************************************************************************************************
865 SUBROUTINE get_sgp_potential(potential, name, description, aliases, elec_conf, &
866 z, zeff, zeff_correction, alpha_core_charge, &
867 ccore_charge, core_charge_radius, &
868 ppl_radius, ppnl_radius, ppl_present, ppnl_present, ppsl_present, &
869 ecp_local, n_local, a_local, c_local, &
870 nloc, nrloc, aloc, bloc, &
871 ecp_semi_local, sl_lmax, npot, nrpot, apot, bpot, &
872 n_nonlocal, nppnl, lmax, is_nonlocal, a_nonlocal, h_nonlocal, c_nonlocal, &
873 cprj_ppnl, vprj_ppnl, has_nlcc, n_nlcc, a_nlcc, c_nlcc)
874
875 TYPE(sgp_potential_type), INTENT(IN) :: potential
876 CHARACTER(LEN=default_string_length), &
877 INTENT(OUT), OPTIONAL :: name
878 CHARACTER(LEN=default_string_length), &
879 DIMENSION(4), INTENT(OUT), OPTIONAL :: description
880 CHARACTER(LEN=default_string_length), &
881 INTENT(OUT), OPTIONAL :: aliases
882 INTEGER, DIMENSION(:), OPTIONAL, POINTER :: elec_conf
883 INTEGER, INTENT(OUT), OPTIONAL :: z
884 REAL(KIND=dp), INTENT(OUT), OPTIONAL :: zeff, zeff_correction, &
885 alpha_core_charge, ccore_charge, &
886 core_charge_radius, ppl_radius, &
887 ppnl_radius
888 LOGICAL, INTENT(OUT), OPTIONAL :: ppl_present, ppnl_present, ppsl_present, &
889 ecp_local
890 INTEGER, INTENT(OUT), OPTIONAL :: n_local
891 REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: a_local, c_local
892 INTEGER, INTENT(OUT), OPTIONAL :: nloc
893 INTEGER, DIMENSION(1:10), INTENT(OUT), OPTIONAL :: nrloc
894 REAL(dp), DIMENSION(1:10), INTENT(OUT), OPTIONAL :: aloc, bloc
895 LOGICAL, INTENT(OUT), OPTIONAL :: ecp_semi_local
896 INTEGER, INTENT(OUT), OPTIONAL :: sl_lmax
897 INTEGER, DIMENSION(0:10), OPTIONAL :: npot
898 INTEGER, DIMENSION(1:15, 0:10), OPTIONAL :: nrpot
899 REAL(dp), DIMENSION(1:15, 0:10), OPTIONAL :: apot, bpot
900 INTEGER, INTENT(OUT), OPTIONAL :: n_nonlocal, nppnl, lmax
901 LOGICAL, DIMENSION(0:5), OPTIONAL :: is_nonlocal
902 REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: a_nonlocal
903 REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER :: h_nonlocal
904 REAL(KIND=dp), DIMENSION(:, :, :), OPTIONAL, &
905 POINTER :: c_nonlocal
906 REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER :: cprj_ppnl
907 REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: vprj_ppnl
908 LOGICAL, INTENT(OUT), OPTIONAL :: has_nlcc
909 INTEGER, INTENT(OUT), OPTIONAL :: n_nlcc
910 REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: a_nlcc, c_nlcc
911
912 IF (PRESENT(name)) name = potential%name
913 IF (PRESENT(aliases)) aliases = potential%aliases
914 IF (PRESENT(description)) description = potential%description
915
916 IF (PRESENT(elec_conf)) elec_conf => potential%elec_conf
917
918 IF (PRESENT(z)) z = potential%z
919 IF (PRESENT(zeff)) zeff = potential%zeff
920 IF (PRESENT(zeff_correction)) zeff_correction = potential%zeff_correction
921 IF (PRESENT(alpha_core_charge)) alpha_core_charge = potential%alpha_core_charge
922 IF (PRESENT(ccore_charge)) ccore_charge = potential%ccore_charge
923 IF (PRESENT(core_charge_radius)) core_charge_radius = potential%core_charge_radius
924
925 IF (PRESENT(ppl_radius)) ppl_radius = potential%ppl_radius
926 IF (PRESENT(ppnl_radius)) ppnl_radius = potential%ppnl_radius
927 IF (PRESENT(ppl_present)) THEN
928 ppl_present = (potential%nloc > 0 .OR. potential%n_local > 0)
929 END IF
930 IF (PRESENT(ppnl_present)) THEN
931 ppnl_present = any(potential%is_nonlocal)
932 END IF
933 IF (PRESENT(ppsl_present)) THEN
934 ppsl_present = potential%ecp_semi_local
935 END IF
936
937 IF (PRESENT(ecp_local)) ecp_local = potential%ecp_local
938 IF (PRESENT(n_local)) n_local = potential%n_local
939 IF (PRESENT(a_local)) a_local => potential%a_local
940 IF (PRESENT(c_local)) c_local => potential%c_local
941
942 IF (PRESENT(nloc)) nloc = potential%nloc
943 IF (PRESENT(nrloc)) nrloc = potential%nrloc
944 IF (PRESENT(aloc)) aloc = potential%aloc
945 IF (PRESENT(bloc)) bloc = potential%bloc
946
947 IF (PRESENT(ecp_semi_local)) ecp_semi_local = potential%ecp_semi_local
948 IF (PRESENT(sl_lmax)) sl_lmax = potential%sl_lmax
949 IF (PRESENT(npot)) npot = potential%npot
950 IF (PRESENT(nrpot)) nrpot = potential%nrpot
951 IF (PRESENT(apot)) apot = potential%apot
952 IF (PRESENT(bpot)) bpot = potential%bpot
953
954 IF (PRESENT(n_nonlocal)) n_nonlocal = potential%n_nonlocal
955 IF (PRESENT(nppnl)) nppnl = potential%nppnl
956 IF (PRESENT(lmax)) lmax = potential%lmax
957 IF (PRESENT(is_nonlocal)) is_nonlocal(:) = potential%is_nonlocal(:)
958 IF (PRESENT(a_nonlocal)) a_nonlocal => potential%a_nonlocal
959 IF (PRESENT(c_nonlocal)) c_nonlocal => potential%c_nonlocal
960 IF (PRESENT(h_nonlocal)) h_nonlocal => potential%h_nonlocal
961 IF (PRESENT(cprj_ppnl)) cprj_ppnl => potential%cprj_ppnl
962 IF (PRESENT(vprj_ppnl)) vprj_ppnl => potential%vprj_ppnl
963
964 IF (PRESENT(has_nlcc)) has_nlcc = potential%has_nlcc
965 IF (PRESENT(n_nlcc)) n_nlcc = potential%n_nlcc
966 IF (PRESENT(a_nlcc)) a_nlcc => potential%a_nlcc
967 IF (PRESENT(c_nlcc)) c_nlcc => potential%c_nlcc
968
969 END SUBROUTINE get_sgp_potential
970
971! **************************************************************************************************
972!> \brief Initialise the coefficients of the projectors of the non-local
973!> part of the GTH pseudopotential and the transformation matrices
974!> for Cartesian overlap integrals between the orbital basis
975!> functions and the projector functions.
976!> \param potential ...
977!> \date 16.10.2000
978!> \author MK
979!> \version 1.0
980! **************************************************************************************************
981 ELEMENTAL SUBROUTINE init_cprj_ppnl(potential)
982
983 TYPE(gth_potential_type), INTENT(INOUT) :: potential
984
985 INTEGER :: cpx, cpy, cpz, cx, cy, cz, ico, iprj, &
986 iprj_ppnl, l, lp, lprj_ppnl, nprj, px, &
987 py, pz
988 REAL(kind=dp) :: alpha_ppnl, cp
989
990 nprj = 0
991
992 DO l = 0, potential%lppnl
993 alpha_ppnl = potential%alpha_ppnl(l)
994 DO iprj_ppnl = 1, potential%nprj_ppnl(l)
995 lp = iprj_ppnl - 1
996 lprj_ppnl = l + 2*lp
997 cp = sqrt(2.0_dp**(2.0_dp*real(lprj_ppnl, dp) + 3.5_dp)* &
998 alpha_ppnl**(real(lprj_ppnl, dp) + 1.5_dp)/ &
999 (rootpi*dfac(2*lprj_ppnl + 1)))
1000 potential%cprj_ppnl(iprj_ppnl, l) = cp
1001 DO cx = 0, l
1002 DO cy = 0, l - cx
1003 cz = l - cx - cy
1004 iprj = nprj + co(cx, cy, cz)
1005 DO px = 0, lp
1006 DO py = 0, lp - px
1007 pz = lp - px - py
1008 cpx = cx + 2*px
1009 cpy = cy + 2*py
1010 cpz = cz + 2*pz
1011 ico = coset(cpx, cpy, cpz)
1012 potential%cprj(ico, iprj) = cp*fac(lp)/(fac(px)*fac(py)*fac(pz))
1013 END DO
1014 END DO
1015 END DO
1016 END DO
1017 nprj = nprj + nco(l)
1018 END DO
1019 END DO
1020
1021 END SUBROUTINE init_cprj_ppnl
1022
1023! **************************************************************************************************
1024!> \brief Initialise a GTH potential data set structure.
1025!> \param potential ...
1026!> \date 27.10.2000
1027!> \author MK
1028!> \version 1.0
1029! **************************************************************************************************
1030 SUBROUTINE init_gth_potential(potential)
1031
1032 TYPE(gth_potential_type), INTENT(IN), POINTER :: potential
1033
1034 IF (.NOT. ASSOCIATED(potential)) RETURN
1035
1036 IF (potential%nppnl > 0) THEN
1037
1038 ! Initialise the projector coefficients of the non-local part of the GTH pseudopotential
1039 ! and the transformation matrices "pgf" -> "prj_ppnl"
1040 CALL init_cprj_ppnl(potential)
1041
1042 ! Initialise the h(i,j) projector coefficients of the non-local part of the
1043 ! GTH pseudopotential
1044 CALL init_vprj_ppnl(potential)
1045
1046 END IF
1047
1048 END SUBROUTINE init_gth_potential
1049
1050! **************************************************************************************************
1051!> \brief Initialise the h(i,j) projector coefficients of the non-local part
1052!> of the GTH pseudopotential (and k(i,j) for SOC, see Hartwigsen, Goedecker, Hutter, PRB 1998).
1053!> \param potential ...
1054!> \date 24.10.2000
1055!> \author MK
1056!> \version 1.0
1057! **************************************************************************************************
1058 ELEMENTAL SUBROUTINE init_vprj_ppnl(potential)
1059
1060 TYPE(gth_potential_type), INTENT(INOUT) :: potential
1061
1062 INTEGER :: i, ico, iprj, iprj_ppnl, iso, j, jco, &
1063 jprj, jprj_ppnl, l, nprj
1064
1065 nprj = 0
1066
1067 DO l = 0, potential%lppnl
1068 DO iprj_ppnl = 1, potential%nprj_ppnl(l)
1069 iprj = nprj + (iprj_ppnl - 1)*nco(l)
1070 DO jprj_ppnl = 1, potential%nprj_ppnl(l)
1071 jprj = nprj + (jprj_ppnl - 1)*nco(l)
1072 DO ico = 1, nco(l)
1073 i = iprj + ico
1074 DO jco = 1, nco(l)
1075 j = jprj + jco
1076 DO iso = 1, nso(l)
1077 potential%vprj_ppnl(i, j) = potential%vprj_ppnl(i, j) + &
1078 orbtramat(l)%slm(iso, ico)* &
1079 potential%hprj_ppnl(iprj_ppnl, &
1080 jprj_ppnl, l)* &
1081 orbtramat(l)%slm(iso, jco)
1082 IF (potential%soc) THEN
1083 ! Transform spin-orbit part
1084 potential%wprj_ppnl(i, j) = potential%wprj_ppnl(i, j) + &
1085 orbtramat(l)%slm(iso, ico)* &
1086 potential%kprj_ppnl(iprj_ppnl, &
1087 jprj_ppnl, l)* &
1088 orbtramat(l)%slm(iso, jco)
1089 END IF
1090 END DO
1091 END DO
1092 END DO
1093 END DO
1094 END DO
1095 nprj = nprj + potential%nprj_ppnl(l)*nco(l)
1096 END DO
1097
1098 END SUBROUTINE init_vprj_ppnl
1099
1100! **************************************************************************************************
1101!> \brief ...
1102!> \param potential ...
1103!> \param itype ...
1104!> \param zeff ...
1105!> \param zeff_correction ...
1106! **************************************************************************************************
1107 PURE SUBROUTINE init_all_potential(potential, itype, zeff, zeff_correction)
1108
1109 TYPE(all_potential_type), INTENT(INOUT), POINTER :: potential
1110 CHARACTER(LEN=*), INTENT(IN), OPTIONAL :: itype
1111 REAL(kind=dp), INTENT(IN), OPTIONAL :: zeff, zeff_correction
1112
1113 INTEGER :: dz
1114
1115 IF (.NOT. ASSOCIATED(potential)) RETURN
1116
1117 IF (PRESENT(zeff)) potential%zeff = zeff
1118 IF (PRESENT(zeff_correction)) potential%zeff_correction = zeff_correction
1119 dz = potential%z - int(potential%zeff - potential%zeff_correction)
1120 SELECT CASE (dz)
1121 CASE DEFAULT
1122 CASE (2)
1123 potential%elec_conf(0) = potential%elec_conf(0) - 2
1124 CASE (10)
1125 potential%elec_conf(0) = potential%elec_conf(0) - 4
1126 potential%elec_conf(1) = potential%elec_conf(1) - 6
1127 CASE (18)
1128 potential%elec_conf(0) = potential%elec_conf(0) - 6
1129 potential%elec_conf(1) = potential%elec_conf(1) - 12
1130 CASE (28)
1131 potential%elec_conf(0) = potential%elec_conf(0) - 6
1132 potential%elec_conf(1) = potential%elec_conf(1) - 12
1133 potential%elec_conf(2) = potential%elec_conf(2) - 10
1134 CASE (30)
1135 potential%elec_conf(0) = potential%elec_conf(0) - 8
1136 potential%elec_conf(1) = potential%elec_conf(1) - 12
1137 potential%elec_conf(2) = potential%elec_conf(2) - 10
1138 CASE (36)
1139 potential%elec_conf(0) = potential%elec_conf(0) - 8
1140 potential%elec_conf(1) = potential%elec_conf(1) - 18
1141 potential%elec_conf(2) = potential%elec_conf(2) - 10
1142 CASE (46)
1143 potential%elec_conf(0) = potential%elec_conf(0) - 8
1144 potential%elec_conf(1) = potential%elec_conf(1) - 18
1145 potential%elec_conf(2) = potential%elec_conf(2) - 20
1146 CASE (48)
1147 potential%elec_conf(0) = potential%elec_conf(0) - 10
1148 potential%elec_conf(1) = potential%elec_conf(1) - 18
1149 potential%elec_conf(2) = potential%elec_conf(2) - 20
1150 CASE (54)
1151 potential%elec_conf(0) = potential%elec_conf(0) - 10
1152 potential%elec_conf(1) = potential%elec_conf(1) - 24
1153 potential%elec_conf(2) = potential%elec_conf(2) - 20
1154 CASE (68)
1155 potential%elec_conf(0) = potential%elec_conf(0) - 10
1156 potential%elec_conf(1) = potential%elec_conf(1) - 24
1157 potential%elec_conf(2) = potential%elec_conf(2) - 20
1158 potential%elec_conf(3) = potential%elec_conf(3) - 14
1159 CASE (78)
1160 potential%elec_conf(0) = potential%elec_conf(0) - 10
1161 potential%elec_conf(1) = potential%elec_conf(1) - 24
1162 potential%elec_conf(2) = potential%elec_conf(2) - 30
1163 potential%elec_conf(3) = potential%elec_conf(3) - 14
1164 CASE (80)
1165 potential%elec_conf(0) = potential%elec_conf(0) - 12
1166 potential%elec_conf(1) = potential%elec_conf(1) - 24
1167 potential%elec_conf(2) = potential%elec_conf(2) - 30
1168 potential%elec_conf(3) = potential%elec_conf(3) - 14
1169 CASE (86)
1170 potential%elec_conf(0) = potential%elec_conf(0) - 12
1171 potential%elec_conf(1) = potential%elec_conf(1) - 30
1172 potential%elec_conf(2) = potential%elec_conf(2) - 30
1173 potential%elec_conf(3) = potential%elec_conf(3) - 14
1174 CASE (100)
1175 potential%elec_conf(0) = potential%elec_conf(0) - 12
1176 potential%elec_conf(1) = potential%elec_conf(1) - 30
1177 potential%elec_conf(2) = potential%elec_conf(2) - 30
1178 potential%elec_conf(3) = potential%elec_conf(3) - 28
1179 END SELECT
1180
1181 IF (PRESENT(itype)) THEN
1182 IF (itype == "BARE") THEN
1183 potential%description(1) = "Bare Coulomb Potential"
1184 IF (dz > 0) THEN
1185 potential%description(2) = "Valence charge only"
1186 ELSE
1187 potential%description(2) = "Full atomic charge"
1188 END IF
1189 END IF
1190 END IF
1191
1192 END SUBROUTINE init_all_potential
1193! **************************************************************************************************
1194!> \brief Initialise a SGP potential data set structure.
1195!> \param potential ...
1196!> \version 1.0
1197! **************************************************************************************************
1198 SUBROUTINE init_sgp_potential(potential)
1199 TYPE(sgp_potential_type), INTENT(IN), POINTER :: potential
1200
1201 INTEGER :: i1, i2, j1, j2, l, la, lb, n1, n2, nnl, &
1202 nprj
1203 INTEGER, ALLOCATABLE, DIMENSION(:, :) :: ind1, ind2
1204 REAL(KIND=dp), DIMENSION(:, :), POINTER :: cprj, hnl
1205 REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: cn
1206
1207 IF (ASSOCIATED(potential)) THEN
1208 IF (potential%nppnl > 0) THEN
1209 !
1210 IF (ASSOCIATED(potential%cprj_ppnl)) THEN
1211 DEALLOCATE (potential%cprj_ppnl)
1212 END IF
1213 nnl = potential%n_nonlocal
1214 nprj = 0
1215 DO l = 0, potential%lmax
1216 nprj = nprj + nnl*nso(l)
1217 END DO
1218 ALLOCATE (potential%cprj_ppnl(potential%nppnl, nprj))
1219 cprj => potential%cprj_ppnl
1220 cprj = 0.0_dp
1221 cn => potential%c_nonlocal
1222 !
1223 ALLOCATE (ind1(potential%nppnl, 3))
1224 n1 = 0
1225 DO i1 = 1, nnl
1226 DO la = 0, potential%lmax
1227 DO j1 = 1, nco(la)
1228 n1 = n1 + 1
1229 ind1(n1, 1) = la
1230 ind1(n1, 2) = j1
1231 ind1(n1, 3) = i1
1232 END DO
1233 END DO
1234 END DO
1235 !
1236 ALLOCATE (ind2(nprj, 3))
1237 n2 = 0
1238 DO i2 = 1, nnl
1239 DO lb = 0, potential%lmax
1240 DO j2 = 1, nso(lb)
1241 n2 = n2 + 1
1242 ind2(n2, 1) = lb
1243 ind2(n2, 2) = j2
1244 ind2(n2, 3) = i2
1245 END DO
1246 END DO
1247 END DO
1248 !
1249 DO n1 = 1, SIZE(ind1, 1)
1250 la = ind1(n1, 1)
1251 j1 = ind1(n1, 2)
1252 i1 = ind1(n1, 3)
1253 DO n2 = 1, SIZE(ind2, 1)
1254 lb = ind2(n2, 1)
1255 IF (la /= lb) cycle
1256 j2 = ind2(n2, 2)
1257 i2 = ind2(n2, 3)
1258 cprj(n1, n2) = orbtramat(la)%c2s(j2, j1)*cn(i1, i2, la)
1259 END DO
1260 END DO
1261 !
1262 hnl => potential%h_nonlocal
1263 IF (ASSOCIATED(potential%vprj_ppnl)) THEN
1264 DEALLOCATE (potential%vprj_ppnl)
1265 END IF
1266 ALLOCATE (potential%vprj_ppnl(nprj))
1267 potential%vprj_ppnl = 0.0_dp
1268 DO n2 = 1, SIZE(ind2, 1)
1269 lb = ind2(n2, 1)
1270 i2 = ind2(n2, 3)
1271 potential%vprj_ppnl(n2) = hnl(i2, lb)
1272 END DO
1273 !
1274 DEALLOCATE (ind1, ind2)
1275 END IF
1276 END IF
1277
1278 END SUBROUTINE init_sgp_potential
1279
1280! **************************************************************************************************
1281!> \brief Read an atomic all-electron potential data set.
1282!> \param element_symbol ...
1283!> \param potential_name ...
1284!> \param potential ...
1285!> \param zeff_correction ...
1286!> \param para_env ...
1287!> \param potential_file_name ...
1288!> \param potential_section ...
1289!> \param update_input ...
1290!> \date 14.05.2000
1291!> \author MK
1292!> \version 1.0
1293! **************************************************************************************************
1294 SUBROUTINE read_all_potential(element_symbol, potential_name, potential, zeff_correction, &
1295 para_env, potential_file_name, potential_section, update_input)
1296
1297 CHARACTER(LEN=*), INTENT(IN) :: element_symbol, potential_name
1298 TYPE(all_potential_type), INTENT(INOUT) :: potential
1299 REAL(KIND=dp), INTENT(IN) :: zeff_correction
1300 TYPE(mp_para_env_type), INTENT(IN), POINTER :: para_env
1301 CHARACTER(len=default_path_length), INTENT(IN) :: potential_file_name
1302 TYPE(section_vals_type), INTENT(IN), POINTER :: potential_section
1303 LOGICAL, INTENT(IN) :: update_input
1304
1305 CHARACTER(LEN=240) :: line
1306 CHARACTER(LEN=242) :: line2
1307 CHARACTER(len=5*default_string_length) :: line_att
1308 CHARACTER(LEN=LEN(element_symbol)) :: symbol
1309 CHARACTER(LEN=LEN(element_symbol)+2) :: symbol2
1310 CHARACTER(LEN=LEN(potential_name)) :: apname
1311 CHARACTER(LEN=LEN(potential_name)+2) :: apname2
1312 INTEGER :: irep, l, strlen1, strlen2
1313 INTEGER, DIMENSION(:), POINTER :: elec_conf
1314 LOGICAL :: found, is_ok, match, read_from_input
1315 REAL(KIND=dp) :: alpha, r
1316 TYPE(cp_parser_type), POINTER :: parser
1317 TYPE(cp_sll_val_type), POINTER :: list
1318 TYPE(val_type), POINTER :: val
1319
1320 line2 = ""
1321 symbol2 = ""
1322 apname2 = ""
1323 NULLIFY (parser)
1324 CALL cite_reference(krack2000)
1325
1326 potential%name = potential_name
1327 read_from_input = .false.
1328 CALL section_vals_get(potential_section, explicit=read_from_input)
1329 IF (.NOT. read_from_input) THEN
1330 ALLOCATE (parser)
1331 CALL parser_create(parser, potential_file_name, para_env=para_env)
1332 END IF
1333
1334 ! Search for the requested potential in the potential file
1335 ! until the potential is found or the end of file is reached
1336
1337 apname = potential_name
1338 symbol = element_symbol
1339 irep = 0
1340 search_loop: DO
1341 IF (read_from_input) THEN
1342 NULLIFY (list, val)
1343 found = .true.
1344 CALL section_vals_list_get(potential_section, "_DEFAULT_KEYWORD_", list=list)
1345 ELSE
1346 CALL parser_search_string(parser, trim(apname), .true., found, line)
1347 END IF
1348 IF (found) THEN
1349 CALL uppercase(symbol)
1350 CALL uppercase(apname)
1351
1352 IF (read_from_input) THEN
1353 match = .true.
1354 ELSE
1355 ! Check both the element symbol and the atomic potential name
1356 match = .false.
1357 CALL uppercase(line)
1358 line2 = " "//line//" "
1359 symbol2 = " "//trim(symbol)//" "
1360 apname2 = " "//trim(apname)//" "
1361 strlen1 = len_trim(symbol2) + 1
1362 strlen2 = len_trim(apname2) + 1
1363
1364 IF ((index(line2, symbol2(:strlen1)) > 0) .AND. &
1365 (index(line2, apname2(:strlen2)) > 0)) match = .true.
1366 END IF
1367 IF (match) THEN
1368 ! Read the electronic configuration
1369 NULLIFY (elec_conf)
1370 l = 0
1371 CALL reallocate(elec_conf, 0, l)
1372 IF (read_from_input) THEN
1373 is_ok = cp_sll_val_next(list, val)
1374 IF (.NOT. is_ok) &
1375 CALL cp_abort(__location__, &
1376 "Error reading the Potential from input file!")
1377 CALL val_get(val, c_val=line_att)
1378 READ (line_att, *) elec_conf(l)
1379 CALL remove_word(line_att)
1380 DO WHILE (len_trim(line_att) /= 0)
1381 l = l + 1
1382 CALL reallocate(elec_conf, 0, l)
1383 READ (line_att, *) elec_conf(l)
1384 CALL remove_word(line_att)
1385 END DO
1386 ELSE
1387 CALL parser_get_object(parser, elec_conf(l), newline=.true.)
1388 DO WHILE (parser_test_next_token(parser) == "INT")
1389 l = l + 1
1390 CALL reallocate(elec_conf, 0, l)
1391 CALL parser_get_object(parser, elec_conf(l))
1392 END DO
1393 irep = irep + 1
1394 IF (update_input) THEN
1395 WRITE (unit=line_att, fmt="(T8,*(1X,I0))") elec_conf(:)
1396 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
1397 c_val=trim(line_att))
1398 END IF
1399 END IF
1400
1401 CALL reallocate(potential%elec_conf, 0, l)
1402 potential%elec_conf(:) = elec_conf(:)
1403
1404 potential%zeff_correction = zeff_correction
1405 potential%zeff = real(sum(elec_conf), dp) + zeff_correction
1406
1407 DEALLOCATE (elec_conf)
1408
1409 ! Read r(loc) to define the exponent of the core charge
1410 ! distribution and calculate the corresponding coefficient
1411
1412 IF (read_from_input) THEN
1413 is_ok = cp_sll_val_next(list, val)
1414 IF (.NOT. is_ok) &
1415 CALL cp_abort(__location__, &
1416 "Error reading the Potential from input file!")
1417 CALL val_get(val, c_val=line_att)
1418 READ (line_att, *) r
1419 ELSE
1420 CALL parser_get_object(parser, r, newline=.true.)
1421 irep = irep + 1
1422 IF (update_input) THEN
1423 WRITE (unit=line_att, fmt="(T9,ES25.16E3)") r
1424 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
1425 c_val=trim(line_att))
1426 END IF
1427 END IF
1428 alpha = 1.0_dp/(2.0_dp*r**2)
1429
1430 potential%alpha_core_charge = alpha
1431 potential%ccore_charge = potential%zeff*sqrt((alpha/pi)**3)
1432
1433 EXIT search_loop
1434 END IF
1435 ELSE
1436 ! Stop program, if the end of file is reached
1437 CALL cp_abort(__location__, &
1438 "The requested atomic potential <"// &
1439 trim(potential_name)// &
1440 "> for element <"// &
1441 trim(symbol)// &
1442 "> was not found in the potential file <"// &
1443 trim(potential_file_name)//">")
1444 END IF
1445 END DO search_loop
1446
1447 IF (.NOT. read_from_input) THEN
1448 ! Dump the potential info in the potential section
1449 IF (match .AND. update_input) THEN
1450 irep = irep + 1
1451 WRITE (unit=line_att, fmt="(T9,A)") &
1452 "# Potential name: "//trim(adjustl(apname2(:strlen2)))// &
1453 " for element symbol: "//trim(adjustl(symbol2(:strlen1)))
1454 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
1455 c_val=trim(line_att))
1456 irep = irep + 1
1457 WRITE (unit=line_att, fmt="(T9,A)") &
1458 "# Potential read from the potential filename: "//trim(adjustl(potential_file_name))
1459 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
1460 c_val=trim(line_att))
1461 END IF
1462 CALL parser_release(parser)
1463 DEALLOCATE (parser)
1464 END IF
1465
1466 END SUBROUTINE read_all_potential
1467
1468! **************************************************************************************************
1469!> \brief Read an atomic local potential data set.
1470!> \param element_symbol ...
1471!> \param potential_name ...
1472!> \param potential ...
1473!> \param para_env ...
1474!> \param potential_file_name ...
1475!> \param potential_section ...
1476!> \param update_input ...
1477!> \date 24.12.2014
1478!> \author JGH
1479!> \version 1.0
1480! **************************************************************************************************
1481 SUBROUTINE read_local_potential(element_symbol, potential_name, potential, &
1482 para_env, potential_file_name, potential_section, update_input)
1483
1484 CHARACTER(LEN=*), INTENT(IN) :: element_symbol, potential_name
1485 TYPE(local_potential_type), INTENT(INOUT) :: potential
1486 TYPE(mp_para_env_type), INTENT(IN), POINTER :: para_env
1487 CHARACTER(len=default_path_length), INTENT(IN) :: potential_file_name
1488 TYPE(section_vals_type), INTENT(IN), POINTER :: potential_section
1489 LOGICAL, INTENT(IN) :: update_input
1490
1491 REAL(KIND=dp), PARAMETER :: eps_tpot = 1.0e-10_dp
1492
1493 CHARACTER(LEN=240) :: line
1494 CHARACTER(LEN=242) :: line2
1495 CHARACTER(len=5*default_string_length) :: line_att
1496 CHARACTER(LEN=LEN(element_symbol)) :: symbol
1497 CHARACTER(LEN=LEN(element_symbol)+2) :: symbol2
1498 CHARACTER(LEN=LEN(potential_name)) :: apname
1499 CHARACTER(LEN=LEN(potential_name)+2) :: apname2
1500 INTEGER :: igau, ipol, irep, l, ngau, npol, &
1501 strlen1, strlen2
1502 LOGICAL :: found, is_ok, match, read_from_input
1503 REAL(KIND=dp), DIMENSION(:), POINTER :: alpha
1504 REAL(KIND=dp), DIMENSION(:, :), POINTER :: cval
1505 TYPE(cp_parser_type), POINTER :: parser
1506 TYPE(cp_sll_val_type), POINTER :: list
1507 TYPE(val_type), POINTER :: val
1508
1509 line2 = ""
1510 symbol2 = ""
1511 apname2 = ""
1512 NULLIFY (parser, alpha, cval)
1513
1514 potential%name = potential_name
1515 read_from_input = .false.
1516 CALL section_vals_get(potential_section, explicit=read_from_input)
1517 IF (.NOT. read_from_input) THEN
1518 ALLOCATE (parser)
1519 CALL parser_create(parser, potential_file_name, para_env=para_env)
1520 END IF
1521
1522 ! Search for the requested potential in the potential file
1523 ! until the potential is found or the end of file is reached
1524
1525 apname = potential_name
1526 symbol = element_symbol
1527 irep = 0
1528 search_loop: DO
1529 IF (read_from_input) THEN
1530 NULLIFY (list, val)
1531 found = .true.
1532 CALL section_vals_list_get(potential_section, "_DEFAULT_KEYWORD_", list=list)
1533 ELSE
1534 CALL parser_search_string(parser, trim(apname), .true., found, line)
1535 END IF
1536 IF (found) THEN
1537 CALL uppercase(symbol)
1538 CALL uppercase(apname)
1539
1540 IF (read_from_input) THEN
1541 match = .true.
1542 ELSE
1543 ! Check both the element symbol and the atomic potential name
1544 match = .false.
1545 CALL uppercase(line)
1546 line2 = " "//line//" "
1547 symbol2 = " "//trim(symbol)//" "
1548 apname2 = " "//trim(apname)//" "
1549 strlen1 = len_trim(symbol2) + 1
1550 strlen2 = len_trim(apname2) + 1
1551
1552 IF ((index(line2, symbol2(:strlen1)) > 0) .AND. &
1553 (index(line2, apname2(:strlen2)) > 0)) match = .true.
1554 END IF
1555 IF (match) THEN
1556
1557 ! Read ngau and npol
1558 IF (read_from_input) THEN
1559 is_ok = cp_sll_val_next(list, val)
1560 IF (.NOT. is_ok) &
1561 CALL cp_abort(__location__, &
1562 "Error reading the Potential from input file!")
1563 CALL val_get(val, c_val=line_att)
1564 READ (line_att, *) ngau, npol
1565 CALL remove_word(line_att)
1566 ELSE
1567 CALL parser_get_object(parser, ngau, newline=.true.)
1568 CALL parser_get_object(parser, npol)
1569 irep = irep + 1
1570 IF (update_input) THEN
1571 WRITE (unit=line_att, fmt="(2(1X,I0))") ngau, npol
1572 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
1573 c_val=trim(line_att))
1574 END IF
1575 END IF
1576
1577 CALL reallocate(alpha, 1, ngau)
1578 CALL reallocate(cval, 1, ngau, 1, npol)
1579 DO igau = 1, ngau
1580 IF (read_from_input) THEN
1581 is_ok = cp_sll_val_next(list, val)
1582 IF (.NOT. is_ok) &
1583 CALL cp_abort(__location__, &
1584 "Error reading the Potential from input file!")
1585 CALL val_get(val, c_val=line_att)
1586 READ (line_att, *) alpha(igau), (cval(igau, ipol), ipol=1, npol)
1587 ELSE
1588 CALL parser_get_object(parser, alpha(igau), newline=.true.)
1589 DO ipol = 1, npol
1590 CALL parser_get_object(parser, cval(igau, ipol), newline=.false.)
1591 END DO
1592 irep = irep + 1
1593 IF (update_input) THEN
1594 WRITE (unit=line_att, fmt="(*(ES25.16E3))") alpha(igau), (cval(igau, ipol), ipol=1, npol)
1595 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
1596 c_val=trim(line_att))
1597 END IF
1598 END IF
1599 END DO
1600 alpha = 1.0_dp/(2.0_dp*alpha**2)
1601
1602 potential%ngau = ngau
1603 potential%npol = npol
1604
1605 potential%alpha => alpha
1606 potential%cval => cval
1607
1608 potential%radius = 0.0_dp
1609 DO igau = 1, ngau
1610 DO ipol = 1, npol
1611 l = 2*(ipol - 1)
1612 potential%radius = max(potential%radius, &
1613 exp_radius(l, alpha(igau), eps_tpot, cval(igau, ipol), &
1614 rlow=potential%radius))
1615 END DO
1616 END DO
1617
1618 EXIT search_loop
1619 END IF
1620 ELSE
1621 ! Stop program, if the end of file is reached
1622 CALL cp_abort(__location__, &
1623 "The requested local atomic potential <"// &
1624 trim(potential_name)// &
1625 "> for element <"// &
1626 trim(symbol)// &
1627 "> was not found in the potential file <"// &
1628 trim(potential_file_name)//">")
1629 END IF
1630 END DO search_loop
1631
1632 IF (.NOT. read_from_input) THEN
1633 ! Dump the potential info in the potential section
1634 IF (match .AND. update_input) THEN
1635 irep = irep + 1
1636 WRITE (unit=line_att, fmt="(A)") &
1637 "# Potential name: "//trim(adjustl(apname2(:strlen2)))// &
1638 " for element symbol: "//trim(adjustl(symbol2(:strlen1)))
1639 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
1640 c_val=trim(line_att))
1641 irep = irep + 1
1642 WRITE (unit=line_att, fmt="(A)") &
1643 "# Potential read from the potential filename: "//trim(adjustl(potential_file_name))
1644 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
1645 c_val=trim(line_att))
1646 END IF
1647 CALL parser_release(parser)
1648 DEALLOCATE (parser)
1649 END IF
1650
1651 END SUBROUTINE read_local_potential
1652
1653! **************************************************************************************************
1654!> \brief Read an atomic GTH potential data set.
1655!> \param element_symbol ...
1656!> \param potential_name ...
1657!> \param potential ...
1658!> \param zeff_correction ...
1659!> \param para_env ...
1660!> \param potential_file_name ...
1661!> \param potential_section ...
1662!> \param update_input ...
1663!> \param monovalent ...
1664!> \date 14.05.2000
1665!> \par Literature
1666!> - S. Goedecker, M. Teter and J. Hutter,
1667!> Phys. Rev. B 54, 1703 (1996)
1668!> - C. Hartwigsen, S. Goedecker and J. Hutter,
1669!> Phys. Rev. B 58, 3641 (1998)
1670!> \par History
1671!> - Add SOC key (27.06.2023, MK)
1672!> \author MK
1673!> \version 1.0
1674! **************************************************************************************************
1675 SUBROUTINE read_gth_potential(element_symbol, potential_name, potential, zeff_correction, &
1676 para_env, potential_file_name, potential_section, update_input, &
1677 monovalent)
1678
1679 CHARACTER(LEN=*), INTENT(IN) :: element_symbol, potential_name
1680 TYPE(gth_potential_type), INTENT(INOUT) :: potential
1681 REAL(KIND=dp), INTENT(IN) :: zeff_correction
1682 TYPE(mp_para_env_type), INTENT(IN), POINTER :: para_env
1683 CHARACTER(len=default_path_length), INTENT(IN) :: potential_file_name
1684 TYPE(section_vals_type), INTENT(IN), POINTER :: potential_section
1685 LOGICAL, INTENT(IN) :: update_input
1686 LOGICAL, INTENT(IN), OPTIONAL :: monovalent
1687
1688 CHARACTER(LEN=240) :: line
1689 CHARACTER(LEN=242) :: line2
1690 CHARACTER(len=5*default_string_length) :: line_att
1691 CHARACTER(LEN=LEN(element_symbol)) :: symbol
1692 CHARACTER(LEN=LEN(element_symbol)+2) :: symbol2
1693 CHARACTER(LEN=LEN(potential_name)) :: apname
1694 CHARACTER(LEN=LEN(potential_name)+2) :: apname2
1695 INTEGER :: i, ic, ipot, irep, istr, j, l, lppnl, &
1696 lprj_ppnl_max, maxlppl, n, nppnl, &
1697 nprj_ppnl, nprj_ppnl_max, strlen1, &
1698 strlen2
1699 INTEGER, DIMENSION(:), POINTER :: elec_conf
1700 LOGICAL :: found, is_ok, match, read_from_input
1701 REAL(KIND=dp) :: alpha, ci, r, rc2
1702 REAL(KIND=dp), DIMENSION(:), POINTER :: tmp_vals
1703 REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: hprj_ppnl, kprj_ppnl
1704 TYPE(cp_parser_type), POINTER :: parser
1705 TYPE(cp_sll_val_type), POINTER :: list
1706 TYPE(val_type), POINTER :: val
1707
1708 line2 = ""
1709 symbol2 = ""
1710 apname2 = ""
1711 NULLIFY (parser, tmp_vals)
1712 CALL cite_reference(goedecker1996)
1713 CALL cite_reference(hartwigsen1998)
1714 CALL cite_reference(krack2005)
1715
1716 potential%monovalent = .false.
1717 IF (PRESENT(monovalent)) potential%monovalent = monovalent
1718
1719 potential%name = potential_name
1720 potential%aliases = potential_name
1721 read_from_input = .false.
1722 CALL section_vals_get(potential_section, explicit=read_from_input)
1723 IF (.NOT. read_from_input) THEN
1724 ALLOCATE (parser)
1725 CALL parser_create(parser, potential_file_name, para_env=para_env)
1726 END IF
1727
1728 ! Initialize extended form
1729 potential%lpotextended = .false.
1730 potential%nexp_lpot = 0
1731 potential%lsdpot = .false.
1732 potential%nexp_lsd = 0
1733 potential%nlcc = .false.
1734 potential%nexp_nlcc = 0
1735
1736 ! Search for the requested potential in the potential file
1737 ! until the potential is found or the end of file is reached
1738 apname = potential_name
1739 symbol = element_symbol
1740 irep = 0
1741 search_loop: DO
1742 IF (read_from_input) THEN
1743 NULLIFY (list, val)
1744 found = .true.
1745 CALL section_vals_list_get(potential_section, "_DEFAULT_KEYWORD_", list=list)
1746 ELSE
1747 CALL parser_search_string(parser, trim(apname), .true., found, line)
1748 END IF
1749 IF (found) THEN
1750 CALL uppercase(symbol)
1751 CALL uppercase(apname)
1752 IF (read_from_input) THEN
1753 match = .true.
1754 ELSE
1755 ! Check both the element symbol and the atomic potential name
1756 match = .false.
1757 CALL uppercase(line)
1758 line2 = " "//line//" "
1759 symbol2 = " "//trim(symbol)//" "
1760 apname2 = " "//trim(apname)//" "
1761 strlen1 = len_trim(symbol2) + 1
1762 strlen2 = len_trim(apname2) + 1
1763 i = index(line2, symbol2(:strlen1))
1764 j = index(line2, apname2(:strlen2))
1765 IF (i > 0 .AND. j > 0) THEN
1766 match = .true.
1767 i = i + 1 + index(line2(i + 1:), " ")
1768 potential%aliases = line2(i:) ! copy all names into aliases field
1769 END IF
1770 END IF
1771 IF (match) THEN
1772 ! Read the electronic configuration
1773 NULLIFY (elec_conf)
1774 l = 0
1775 CALL reallocate(elec_conf, 0, l)
1776 IF (read_from_input) THEN
1777 is_ok = cp_sll_val_next(list, val)
1778 IF (.NOT. is_ok) &
1779 CALL cp_abort(__location__, &
1780 "Error while reading GTH potential from input file")
1781 CALL val_get(val, c_val=line_att)
1782 READ (line_att, *) elec_conf(l)
1783 CALL remove_word(line_att)
1784 DO WHILE (len_trim(line_att) /= 0)
1785 l = l + 1
1786 CALL reallocate(elec_conf, 0, l)
1787 READ (line_att, *) elec_conf(l)
1788 CALL remove_word(line_att)
1789 END DO
1790 ELSE
1791 CALL parser_get_object(parser, elec_conf(l), newline=.true.)
1792 DO WHILE (parser_test_next_token(parser) == "INT")
1793 l = l + 1
1794 CALL reallocate(elec_conf, 0, l)
1795 CALL parser_get_object(parser, elec_conf(l))
1796 END DO
1797 irep = irep + 1
1798 IF (update_input) THEN
1799 WRITE (unit=line_att, fmt="(T8,*(1X,I0))") elec_conf(:)
1800 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
1801 c_val=trim(line_att))
1802 END IF
1803 END IF
1804
1805 CALL reallocate(potential%elec_conf, 0, l)
1806 IF (potential%monovalent) THEN
1807 potential%elec_conf(0) = 1
1808 ELSE
1809 potential%elec_conf(:) = elec_conf(:)
1810 END IF
1811
1812 potential%zeff_correction = zeff_correction
1813 potential%zeff = real(sum(potential%elec_conf), dp) + zeff_correction
1814
1815 DEALLOCATE (elec_conf)
1816
1817 ! Read r(loc) to define the exponent of the core charge
1818 ! distribution and calculate the corresponding coefficient
1819 IF (read_from_input) THEN
1820 is_ok = cp_sll_val_next(list, val)
1821 IF (.NOT. is_ok) &
1822 CALL cp_abort(__location__, &
1823 "Error while reading GTH potential from input file")
1824 CALL val_get(val, c_val=line_att)
1825 READ (line_att, *) r
1826 CALL remove_word(line_att)
1827 ELSE
1828 line_att = ""
1829 CALL parser_get_object(parser, r, newline=.true.)
1830 istr = len_trim(line_att) + 1
1831 WRITE (unit=line_att(istr:), fmt="(T9,ES25.16E3)") r
1832 END IF
1833 alpha = 1.0_dp/(2.0_dp*r**2)
1834
1835 potential%alpha_core_charge = alpha
1836 potential%ccore_charge = potential%zeff*sqrt((alpha/pi)**3)
1837
1838 potential%alpha_ppl = alpha
1839 potential%cerf_ppl = potential%zeff*sqrt((alpha/pi)**3)
1840
1841 ! Read the parameters for the local part of the GTH pseudopotential (ppl)
1842 IF (read_from_input) THEN
1843 READ (line_att, *) n
1844 CALL remove_word(line_att)
1845 ELSE
1846 CALL parser_get_object(parser, n)
1847 istr = len_trim(line_att) + 1
1848 WRITE (unit=line_att(istr:), fmt="(1X,I0)") n
1849 END IF
1850 potential%nexp_ppl = n
1851 CALL reallocate(potential%cexp_ppl, 1, n)
1852
1853 DO i = 1, n
1854 IF (read_from_input) THEN
1855 READ (line_att, *) ci
1856 CALL remove_word(line_att)
1857 ELSE
1858 CALL parser_get_object(parser, ci)
1859 istr = len_trim(line_att) + 1
1860 WRITE (unit=line_att(istr:), fmt="(ES25.16E3)") ci
1861 END IF
1862 rc2 = (2.0_dp*potential%alpha_ppl)
1863 potential%cexp_ppl(i) = rc2**(i - 1)*ci
1864 END DO
1865
1866 IF (.NOT. read_from_input) THEN
1867 irep = irep + 1
1868 IF (update_input) THEN
1869 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
1870 c_val=trim(line_att))
1871 END IF
1872 line_att = ""
1873 ELSE
1874 IF (len_trim(line_att) /= 0) THEN
1875 CALL cp_abort(__location__, &
1876 "Error while reading GTH potential from input file")
1877 END IF
1878 END IF
1879 maxlppl = 2*(n - 1)
1880
1881 IF (maxlppl > -1) CALL init_orbital_pointers(maxlppl)
1882
1883 ! Read extended form of GTH pseudopotential
1884 ! local potential, NLCC, LSD potential, spin-orbit coupling (SOC)
1885 IF (read_from_input) THEN
1886 read_keywords_from_input: DO
1887 is_ok = cp_sll_val_next(list, val)
1888 cpassert(is_ok)
1889 CALL val_get(val, c_val=line_att)
1890 IF (index(line_att, "LPOT") /= 0) THEN
1891 potential%lpotextended = .true.
1892 CALL remove_word(line_att)
1893 READ (line_att, *) potential%nexp_lpot
1894 n = potential%nexp_lpot
1895 maxlppl = 2*(n - 1)
1896 IF (maxlppl > -1) CALL init_orbital_pointers(maxlppl)
1897 NULLIFY (potential%alpha_lpot, potential%nct_lpot, potential%cval_lpot)
1898 CALL reallocate(potential%alpha_lpot, 1, n)
1899 CALL reallocate(potential%nct_lpot, 1, n)
1900 CALL reallocate(potential%cval_lpot, 1, 4, 1, n)
1901 DO ipot = 1, potential%nexp_lpot
1902 is_ok = cp_sll_val_next(list, val)
1903 cpassert(is_ok)
1904 CALL val_get(val, c_val=line_att)
1905 READ (line_att, *) r
1906 potential%alpha_lpot(ipot) = 0.5_dp/(r*r)
1907 CALL remove_word(line_att)
1908 READ (line_att, *) potential%nct_lpot(ipot)
1909 CALL remove_word(line_att)
1910 DO ic = 1, potential%nct_lpot(ipot)
1911 READ (line_att, *) ci
1912 rc2 = (2._dp*potential%alpha_lpot(ipot))**(ic - 1)
1913 potential%cval_lpot(ic, ipot) = ci*rc2
1914 CALL remove_word(line_att)
1915 END DO
1916 END DO
1917 ELSE IF (index(line_att, "NLCC") /= 0) THEN
1918 potential%nlcc = .true.
1919 CALL remove_word(line_att)
1920 READ (line_att, *) potential%nexp_nlcc
1921 n = potential%nexp_nlcc
1922 NULLIFY (potential%alpha_nlcc, potential%nct_nlcc, potential%cval_nlcc)
1923 CALL reallocate(potential%alpha_nlcc, 1, n)
1924 CALL reallocate(potential%nct_nlcc, 1, n)
1925 CALL reallocate(potential%cval_nlcc, 1, 4, 1, n)
1926 DO ipot = 1, potential%nexp_nlcc
1927 is_ok = cp_sll_val_next(list, val)
1928 cpassert(is_ok)
1929 CALL val_get(val, c_val=line_att)
1930 READ (line_att, *) potential%alpha_nlcc(ipot)
1931 CALL remove_word(line_att)
1932 READ (line_att, *) potential%nct_nlcc(ipot)
1933 CALL remove_word(line_att)
1934 DO ic = 1, potential%nct_nlcc(ipot)
1935 READ (line_att, *) potential%cval_nlcc(ic, ipot)
1936 ! Make it compatible with BigDFT style
1937 potential%cval_nlcc(ic, ipot) = potential%cval_nlcc(ic, ipot)/(4.0_dp*pi)
1938 CALL remove_word(line_att)
1939 END DO
1940 END DO
1941 ELSE IF (index(line_att, "LSD") /= 0) THEN
1942 potential%lsdpot = .true.
1943 CALL remove_word(line_att)
1944 READ (line_att, *) potential%nexp_lsd
1945 n = potential%nexp_lsd
1946 NULLIFY (potential%alpha_lsd, potential%nct_lsd, potential%cval_lsd)
1947 CALL reallocate(potential%alpha_lsd, 1, n)
1948 CALL reallocate(potential%nct_lsd, 1, n)
1949 CALL reallocate(potential%cval_lsd, 1, 4, 1, n)
1950 DO ipot = 1, potential%nexp_lsd
1951 is_ok = cp_sll_val_next(list, val)
1952 cpassert(is_ok)
1953 CALL val_get(val, c_val=line_att)
1954 READ (line_att, *) r
1955 potential%alpha_lsd(ipot) = 0.5_dp/(r*r)
1956 CALL remove_word(line_att)
1957 READ (line_att, *) potential%nct_lsd(ipot)
1958 CALL remove_word(line_att)
1959 DO ic = 1, potential%nct_lsd(ipot)
1960 READ (line_att, *) ci
1961 rc2 = (2._dp*potential%alpha_lsd(ipot))**(ic - 1)
1962 potential%cval_lsd(ic, ipot) = ci*rc2
1963 CALL remove_word(line_att)
1964 END DO
1965 END DO
1966 ELSE
1967 EXIT read_keywords_from_input
1968 END IF
1969 END DO read_keywords_from_input
1970 ELSE
1971 read_keywords: DO
1972 CALL parser_get_next_line(parser, 1)
1973 IF (parser_test_next_token(parser) == "INT") THEN
1974 EXIT read_keywords
1975 ELSE IF (parser_test_next_token(parser) == "STR") THEN
1976 CALL parser_get_object(parser, line)
1977 IF (index(line, "LPOT") /= 0) THEN
1978 ! Local potential
1979 potential%lpotextended = .true.
1980 CALL parser_get_object(parser, potential%nexp_lpot)
1981 n = potential%nexp_lpot
1982 NULLIFY (potential%alpha_lpot, potential%nct_lpot, potential%cval_lpot)
1983 CALL reallocate(potential%alpha_lpot, 1, n)
1984 CALL reallocate(potential%nct_lpot, 1, n)
1985 CALL reallocate(potential%cval_lpot, 1, 4, 1, n)
1986 ! Add to input section
1987 irep = irep + 1
1988 IF (update_input) THEN
1989 WRITE (unit=line_att, fmt="(T9,A,1X,I0)") "LPOT", n
1990 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
1991 c_val=trim(line_att))
1992 END IF
1993 DO ipot = 1, potential%nexp_lpot
1994 CALL parser_get_object(parser, r, newline=.true.)
1995 potential%alpha_lpot(ipot) = 0.5_dp/(r*r)
1996 CALL parser_get_object(parser, potential%nct_lpot(ipot))
1997 CALL reallocate(tmp_vals, 1, potential%nct_lpot(ipot))
1998 DO ic = 1, potential%nct_lpot(ipot)
1999 CALL parser_get_object(parser, ci)
2000 tmp_vals(ic) = ci
2001 rc2 = (2._dp*potential%alpha_lpot(ipot))**(ic - 1)
2002 potential%cval_lpot(ic, ipot) = ci*rc2
2003 END DO
2004 ! Add to input section
2005 irep = irep + 1
2006 IF (update_input) THEN
2007 WRITE (unit=line_att, fmt="(T9,ES25.16E3,1X,I0,*(ES25.16E3))") &
2008 r, potential%nct_lpot(ipot), tmp_vals(1:potential%nct_lpot(ipot))
2009 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
2010 c_val=trim(line_att))
2011 END IF
2012 END DO
2013 ELSE IF (index(line, "NLCC") /= 0) THEN
2014 ! NLCC
2015 potential%nlcc = .true.
2016 CALL parser_get_object(parser, potential%nexp_nlcc)
2017 n = potential%nexp_nlcc
2018 NULLIFY (potential%alpha_nlcc, potential%nct_nlcc, potential%cval_nlcc)
2019 CALL reallocate(potential%alpha_nlcc, 1, n)
2020 CALL reallocate(potential%nct_nlcc, 1, n)
2021 CALL reallocate(potential%cval_nlcc, 1, 4, 1, n)
2022 ! Add to input section
2023 WRITE (unit=line_att, fmt="(T9,A,1X,I0)") "NLCC", n
2024 irep = irep + 1
2025 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
2026 c_val=trim(line_att))
2027 DO ipot = 1, potential%nexp_nlcc
2028 CALL parser_get_object(parser, potential%alpha_nlcc(ipot), newline=.true.)
2029 CALL parser_get_object(parser, potential%nct_nlcc(ipot))
2030 CALL reallocate(tmp_vals, 1, potential%nct_nlcc(ipot))
2031 DO ic = 1, potential%nct_nlcc(ipot)
2032 CALL parser_get_object(parser, potential%cval_nlcc(ic, ipot))
2033 tmp_vals(ic) = potential%cval_nlcc(ic, ipot)
2034 ! Make it compatible with BigDFT style
2035 potential%cval_nlcc(ic, ipot) = potential%cval_nlcc(ic, ipot)/(4.0_dp*pi)
2036 END DO
2037 ! Add to input section
2038 irep = irep + 1
2039 IF (update_input) THEN
2040 WRITE (unit=line_att, fmt="(T9,ES25.16E3,1X,I0,*(ES25.16E3))") &
2041 potential%alpha_nlcc(ipot), potential%nct_nlcc(ipot), &
2042 tmp_vals(1:potential%nct_nlcc(ipot))
2043 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
2044 c_val=trim(line_att))
2045 END IF
2046 END DO
2047 ELSE IF (index(line, "LSD") /= 0) THEN
2048 ! LSD potential
2049 potential%lsdpot = .true.
2050 CALL parser_get_object(parser, potential%nexp_lsd)
2051 n = potential%nexp_lsd
2052 NULLIFY (potential%alpha_lsd, potential%nct_lsd, potential%cval_lsd)
2053 CALL reallocate(potential%alpha_lsd, 1, n)
2054 CALL reallocate(potential%nct_lsd, 1, n)
2055 CALL reallocate(potential%cval_lsd, 1, 4, 1, n)
2056 ! Add to input section
2057 irep = irep + 1
2058 IF (update_input) THEN
2059 WRITE (unit=line_att, fmt="(T9,A,1X,I0)") "LSD", n
2060 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
2061 c_val=trim(line_att))
2062 END IF
2063 DO ipot = 1, potential%nexp_lsd
2064 CALL parser_get_object(parser, r, newline=.true.)
2065 potential%alpha_lsd(ipot) = 0.5_dp/(r*r)
2066 CALL parser_get_object(parser, potential%nct_lsd(ipot))
2067 CALL reallocate(tmp_vals, 1, potential%nct_lsd(ipot))
2068 DO ic = 1, potential%nct_lsd(ipot)
2069 CALL parser_get_object(parser, ci)
2070 tmp_vals(ic) = ci
2071 rc2 = (2._dp*potential%alpha_lsd(ipot))**(ic - 1)
2072 potential%cval_lsd(ic, ipot) = ci*rc2
2073 END DO
2074 ! Add to input section
2075 irep = irep + 1
2076 IF (update_input) THEN
2077 WRITE (unit=line_att, fmt="(T9,ES25.16E3,1X,I0,*(ES25.16E3))") r, potential%nct_lsd(ipot), &
2078 tmp_vals(1:potential%nct_lsd(ipot))
2079 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
2080 c_val=trim(line_att))
2081 END IF
2082 END DO
2083 ELSE
2084 CALL cp_abort(__location__, &
2085 "Syntax error for <"// &
2086 trim(element_symbol)// &
2087 "> in the atomic potential <"// &
2088 trim(potential_name)// &
2089 "> potential file <"// &
2090 trim(potential_file_name)//">: "// &
2091 "Expected LPOT/NLCC/LSD keyword, got: <"// &
2092 trim(line)//">")
2093 END IF
2094 ELSE
2095 CALL parser_get_object(parser, line)
2096 CALL cp_abort(__location__, &
2097 "Syntax error for <"// &
2098 trim(element_symbol)// &
2099 "> in the atomic potential <"// &
2100 trim(potential_name)// &
2101 "> potential file <"// &
2102 trim(potential_file_name)//">: "// &
2103 "Expected LPOT/NLCC/LSD keyword or INTEGER, got: <"// &
2104 trim(line)//">")
2105 END IF
2106 END DO read_keywords
2107 END IF
2108
2109 ! Read the parameters for the non-local part of the GTH pseudopotential (ppnl)
2110 IF (read_from_input) THEN
2111 READ (line_att, *) n
2112 CALL remove_word(line_att)
2113 IF (index(line_att, "SOC") /= 0) THEN
2114 potential%soc = .true.
2115 CALL remove_word(line_att)
2116 END IF
2117 ELSE
2118 CALL parser_get_object(parser, n)
2119 IF (parser_test_next_token(parser) == "STR") THEN
2120 CALL parser_get_object(parser, line)
2121 IF (index(line, "SOC") /= 0) potential%soc = .true.
2122 END IF
2123 irep = irep + 1
2124 IF (update_input) THEN
2125 IF (potential%soc) THEN
2126 WRITE (unit=line_att, fmt="(T9,I0,2X,A)") n, "SOC"
2127 ELSE
2128 WRITE (unit=line_att, fmt="(T9,I0)") n
2129 END IF
2130 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
2131 c_val=trim(line_att))
2132 END IF
2133 END IF
2134 potential%lppnl = n - 1
2135 potential%nppnl = 0
2136
2137 potential%lprj_ppnl_max = n - 1
2138 potential%nprj_ppnl_max = 0
2139
2140 IF (n > 0) THEN
2141
2142 lppnl = potential%lppnl
2143 nppnl = potential%nppnl
2144
2145 CALL init_orbital_pointers(lppnl)
2146
2147 NULLIFY (hprj_ppnl, kprj_ppnl)
2148
2149 ! Load the parameter for n non-local projectors
2150
2151 CALL reallocate(potential%alpha_ppnl, 0, lppnl)
2152 CALL reallocate(potential%nprj_ppnl, 0, lppnl)
2153
2154 lprj_ppnl_max = -1
2155 nprj_ppnl_max = 0
2156
2157 DO l = 0, lppnl
2158 IF (read_from_input) THEN
2159 is_ok = cp_sll_val_next(list, val)
2160 IF (.NOT. is_ok) &
2161 CALL cp_abort(__location__, &
2162 "Error while reading GTH potential from input file")
2163 CALL val_get(val, c_val=line_att)
2164 READ (line_att, *) r
2165 CALL remove_word(line_att)
2166 READ (line_att, *) nprj_ppnl
2167 CALL remove_word(line_att)
2168 ELSE
2169 line_att = ""
2170 CALL parser_get_object(parser, r, newline=.true.)
2171 CALL parser_get_object(parser, nprj_ppnl)
2172 istr = len_trim(line_att) + 1
2173 WRITE (unit=line_att(istr:), fmt="(T9,ES25.16E3,1X,I0)") r, nprj_ppnl
2174 END IF
2175 IF (r == 0.0_dp .AND. nprj_ppnl /= 0) THEN
2176 CALL cp_abort(__location__, &
2177 "An error was detected in the atomic potential <"// &
2178 trim(potential_name)// &
2179 "> potential file <"// &
2180 trim(potential_file_name)//">")
2181 END IF
2182 potential%alpha_ppnl(l) = 0.0_dp
2183 IF (r /= 0.0_dp .AND. n /= 0) potential%alpha_ppnl(l) = 1.0_dp/(2.0_dp*r**2)
2184 potential%nprj_ppnl(l) = nprj_ppnl
2185 nppnl = nppnl + nprj_ppnl*nco(l)
2186 IF (nprj_ppnl > nprj_ppnl_max) THEN
2187 nprj_ppnl_max = nprj_ppnl
2188 CALL reallocate(hprj_ppnl, 1, nprj_ppnl_max, &
2189 1, nprj_ppnl_max, &
2190 0, lppnl)
2191 CALL reallocate(kprj_ppnl, 1, nprj_ppnl_max, &
2192 1, nprj_ppnl_max, &
2193 0, lppnl)
2194 END IF
2195 DO i = 1, nprj_ppnl
2196 IF (i == 1) THEN
2197 IF (read_from_input) THEN
2198 READ (line_att, *) hprj_ppnl(i, i, l)
2199 CALL remove_word(line_att)
2200 ELSE
2201 CALL parser_get_object(parser, hprj_ppnl(i, i, l))
2202 istr = len_trim(line_att) + 1
2203 WRITE (unit=line_att(istr:), fmt="(ES25.16E3)") hprj_ppnl(i, i, l)
2204 END IF
2205 ELSE
2206 IF (read_from_input) THEN
2207 IF (len_trim(line_att) /= 0) &
2208 CALL cp_abort(__location__, &
2209 "Error while reading GTH potential from input file")
2210 is_ok = cp_sll_val_next(list, val)
2211 IF (.NOT. is_ok) &
2212 CALL cp_abort(__location__, &
2213 "Error while reading GTH potential from input file")
2214 CALL val_get(val, c_val=line_att)
2215 READ (line_att, *) hprj_ppnl(i, i, l)
2216 CALL remove_word(line_att)
2217 ELSE
2218 IF (update_input) THEN
2219 irep = irep + 1
2220 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
2221 c_val=trim(line_att))
2222 END IF
2223 line_att = ""
2224 CALL parser_get_object(parser, hprj_ppnl(i, i, l), newline=.true.)
2225 istr = len_trim(line_att) + 1
2226 WRITE (unit=line_att(istr:), fmt="(T36,A,ES25.16E3)") &
2227 repeat(" ", 25*(i - 1)), hprj_ppnl(i, i, l)
2228 END IF
2229 END IF
2230 DO j = i + 1, nprj_ppnl
2231 IF (read_from_input) THEN
2232 READ (line_att, *) hprj_ppnl(i, j, l)
2233 CALL remove_word(line_att)
2234 ELSE
2235 CALL parser_get_object(parser, hprj_ppnl(i, j, l))
2236 istr = len_trim(line_att) + 1
2237 WRITE (unit=line_att(istr:), fmt="(ES25.16E3)") hprj_ppnl(i, j, l)
2238 END IF
2239 END DO
2240 END DO
2241 IF (.NOT. read_from_input) THEN
2242 IF (update_input) THEN
2243 irep = irep + 1
2244 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
2245 c_val=trim(line_att))
2246 END IF
2247 line_att = ""
2248 ELSE
2249 IF (len_trim(line_att) /= 0) THEN
2250 CALL cp_abort(__location__, &
2251 "Error while reading GTH potential from input file")
2252 END IF
2253 END IF
2254 IF (nprj_ppnl > 1) THEN
2255 CALL symmetrize_matrix(hprj_ppnl(:, :, l), "upper_to_lower")
2256 END IF
2257 IF (potential%soc .AND. (l > 0)) THEN
2258 ! Read non-local parameters for spin-orbit coupling
2259 DO i = 1, nprj_ppnl
2260 IF (read_from_input) THEN
2261 IF (len_trim(line_att) /= 0) &
2262 CALL cp_abort(__location__, &
2263 "Error while reading GTH potential from input file")
2264 is_ok = cp_sll_val_next(list, val)
2265 IF (.NOT. is_ok) &
2266 CALL cp_abort(__location__, &
2267 "Error while reading GTH potential from input file")
2268 CALL val_get(val, c_val=line_att)
2269 READ (line_att, *) kprj_ppnl(i, i, l)
2270 CALL remove_word(line_att)
2271 ELSE
2272 IF (i > 1 .AND. update_input) THEN
2273 irep = irep + 1
2274 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
2275 c_val=trim(line_att))
2276 END IF
2277 line_att = ""
2278 CALL parser_get_object(parser, kprj_ppnl(i, i, l), newline=.true.)
2279 istr = len_trim(line_att) + 1
2280 WRITE (unit=line_att(istr:), fmt="(T36,A,ES25.16E3)") &
2281 repeat(" ", 25*(i - 1)), kprj_ppnl(i, i, l)
2282 END IF
2283 DO j = i + 1, nprj_ppnl
2284 IF (read_from_input) THEN
2285 READ (line_att, *) kprj_ppnl(i, j, l)
2286 CALL remove_word(line_att)
2287 ELSE
2288 CALL parser_get_object(parser, kprj_ppnl(i, j, l))
2289 istr = len_trim(line_att) + 1
2290 WRITE (unit=line_att(istr:), fmt="(ES25.16E3)") kprj_ppnl(i, j, l)
2291 END IF
2292 END DO
2293 END DO
2294 IF (read_from_input) THEN
2295 IF (len_trim(line_att) /= 0) THEN
2296 CALL cp_abort(__location__, &
2297 "Error while reading GTH potential from input file")
2298 END IF
2299 ELSE
2300 IF (update_input) THEN
2301 irep = irep + 1
2302 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
2303 c_val=trim(line_att))
2304 END IF
2305 line_att = ""
2306 END IF
2307 IF (nprj_ppnl > 1) THEN
2308 CALL symmetrize_matrix(kprj_ppnl(:, :, l), "upper_to_lower")
2309 END IF
2310 END IF ! SOC
2311 lprj_ppnl_max = max(lprj_ppnl_max, l + 2*(nprj_ppnl - 1))
2312 END DO ! lppnl
2313
2314 potential%nppnl = nppnl
2315 CALL init_orbital_pointers(lprj_ppnl_max)
2316
2317 potential%lprj_ppnl_max = lprj_ppnl_max
2318 potential%nprj_ppnl_max = nprj_ppnl_max
2319 CALL reallocate(potential%hprj_ppnl, 1, nprj_ppnl_max, &
2320 1, nprj_ppnl_max, &
2321 0, lppnl)
2322 potential%hprj_ppnl(:, :, :) = hprj_ppnl(:, :, :)
2323 CALL reallocate(potential%kprj_ppnl, 1, nprj_ppnl_max, &
2324 1, nprj_ppnl_max, &
2325 0, lppnl)
2326 potential%kprj_ppnl(:, :, :) = kprj_ppnl(:, :, :)
2327
2328 CALL reallocate(potential%cprj, 1, ncoset(lprj_ppnl_max), 1, nppnl)
2329 CALL reallocate(potential%cprj_ppnl, 1, nprj_ppnl_max, 0, lppnl)
2330 CALL reallocate(potential%vprj_ppnl, 1, nppnl, 1, nppnl)
2331 CALL reallocate(potential%wprj_ppnl, 1, nppnl, 1, nppnl)
2332
2333 DEALLOCATE (hprj_ppnl, kprj_ppnl)
2334 END IF
2335 EXIT search_loop
2336 END IF
2337 ELSE
2338 ! Stop program, if the end of file is reached
2339 CALL cp_abort(__location__, &
2340 "The requested atomic potential <"// &
2341 trim(potential_name)// &
2342 "> for element <"// &
2343 trim(symbol)// &
2344 "> was not found in the potential file <"// &
2345 trim(potential_file_name)//">")
2346 END IF
2347 END DO search_loop
2348
2349 IF (.NOT. read_from_input) THEN
2350 ! Dump the potential info in the potential section
2351 IF (match .AND. update_input) THEN
2352 irep = irep + 1
2353 WRITE (unit=line_att, fmt="(T9,A)") &
2354 "# Potential name: "//trim(adjustl(apname2(:strlen2)))// &
2355 " for element symbol: "//trim(adjustl(symbol2(:strlen1)))
2356 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
2357 c_val=trim(line_att))
2358 irep = irep + 1
2359 WRITE (unit=line_att, fmt="(T9,A)") &
2360 "# Potential read from the potential filename: "//trim(adjustl(potential_file_name))
2361 CALL section_vals_val_set(potential_section, "_DEFAULT_KEYWORD_", i_rep_val=irep, &
2362 c_val=trim(line_att))
2363 END IF
2364 CALL parser_release(parser)
2365 DEALLOCATE (parser)
2366 END IF
2367
2368 IF (ASSOCIATED(tmp_vals)) DEALLOCATE (tmp_vals)
2369
2370 END SUBROUTINE read_gth_potential
2371
2372! **************************************************************************************************
2373!> \brief ...
2374!> \param potential ...
2375!> \param z ...
2376!> \param zeff_correction ...
2377! **************************************************************************************************
2378 SUBROUTINE set_default_all_potential(potential, z, zeff_correction)
2379
2380 TYPE(all_potential_type), INTENT(INOUT) :: potential
2381 INTEGER, INTENT(IN) :: z
2382 REAL(kind=dp), INTENT(IN) :: zeff_correction
2383
2384 CHARACTER(LEN=default_string_length) :: name
2385 INTEGER, DIMENSION(:), POINTER :: elec_conf
2386 REAL(kind=dp) :: alpha, alpha_core_charge, ccore_charge, &
2387 core_charge_radius, r, zeff
2388
2389 ALLOCATE (elec_conf(0:3))
2390 elec_conf(0:3) = ptable(z)%e_conv(0:3)
2391 zeff = real(sum(elec_conf), dp) + zeff_correction
2392 name = ptable(z)%name
2393
2394 r = ptable(z)%covalent_radius*0.5_dp
2395 r = max(r, 0.2_dp)
2396 r = min(r, 1.0_dp)
2397 alpha = 1.0_dp/(2.0_dp*r**2)
2398
2399 core_charge_radius = r
2400 alpha_core_charge = alpha
2401 ccore_charge = zeff*sqrt((alpha/pi)**3)
2402
2403 CALL set_all_potential(potential, &
2404 name=name, &
2405 alpha_core_charge=alpha_core_charge, &
2406 ccore_charge=ccore_charge, &
2407 core_charge_radius=core_charge_radius, &
2408 z=z, &
2409 zeff=zeff, &
2410 zeff_correction=zeff_correction, &
2411 elec_conf=elec_conf)
2412
2413 DEALLOCATE (elec_conf)
2414
2415 END SUBROUTINE set_default_all_potential
2416
2417! **************************************************************************************************
2418!> \brief Set the attributes of an all-electron potential data set.
2419!> \param potential ...
2420!> \param name ...
2421!> \param alpha_core_charge ...
2422!> \param ccore_charge ...
2423!> \param core_charge_radius ...
2424!> \param z ...
2425!> \param zeff ...
2426!> \param zeff_correction ...
2427!> \param elec_conf ...
2428!> \date 11.01.2002
2429!> \author MK
2430!> \version 1.0
2431! **************************************************************************************************
2432 SUBROUTINE set_all_potential(potential, name, alpha_core_charge, &
2433 ccore_charge, core_charge_radius, z, zeff, &
2434 zeff_correction, elec_conf)
2435
2436 TYPE(all_potential_type), INTENT(INOUT) :: potential
2437 CHARACTER(LEN=default_string_length), INTENT(IN), &
2438 OPTIONAL :: name
2439 REAL(kind=dp), INTENT(IN), OPTIONAL :: alpha_core_charge, ccore_charge, &
2440 core_charge_radius
2441 INTEGER, INTENT(IN), OPTIONAL :: z
2442 REAL(kind=dp), INTENT(IN), OPTIONAL :: zeff, zeff_correction
2443 INTEGER, DIMENSION(:), OPTIONAL, POINTER :: elec_conf
2444
2445 IF (PRESENT(name)) potential%name = name
2446 IF (PRESENT(alpha_core_charge)) &
2447 potential%alpha_core_charge = alpha_core_charge
2448 IF (PRESENT(ccore_charge)) potential%ccore_charge = ccore_charge
2449 IF (PRESENT(core_charge_radius)) &
2450 potential%core_charge_radius = core_charge_radius
2451 IF (PRESENT(z)) potential%z = z
2452 IF (PRESENT(zeff)) potential%zeff = zeff
2453 IF (PRESENT(zeff_correction)) potential%zeff_correction = zeff_correction
2454 IF (PRESENT(elec_conf)) THEN
2455 IF (.NOT. ASSOCIATED(potential%elec_conf)) THEN
2456 CALL reallocate(potential%elec_conf, 0, SIZE(elec_conf) - 1)
2457 END IF
2458 potential%elec_conf(:) = elec_conf(:)
2459 END IF
2460
2461 END SUBROUTINE set_all_potential
2462
2463! **************************************************************************************************
2464!> \brief Set the attributes of an atomic local potential data set.
2465!> \param potential ...
2466!> \param name ...
2467!> \param alpha ...
2468!> \param cval ...
2469!> \param radius ...
2470!> \date 24.01.2014
2471!> \author JGH
2472!> \version 1.0
2473! **************************************************************************************************
2474 SUBROUTINE set_local_potential(potential, name, alpha, cval, radius)
2475
2476 TYPE(local_potential_type), INTENT(INOUT) :: potential
2477 CHARACTER(LEN=default_string_length), INTENT(IN), &
2478 OPTIONAL :: name
2479 REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: alpha
2480 REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER :: cval
2481 REAL(KIND=dp), INTENT(IN), OPTIONAL :: radius
2482
2483 IF (PRESENT(name)) potential%name = name
2484 IF (PRESENT(alpha)) potential%alpha => alpha
2485 IF (PRESENT(cval)) potential%cval => cval
2486 IF (PRESENT(radius)) potential%radius = radius
2487
2488 END SUBROUTINE set_local_potential
2489
2490! **************************************************************************************************
2491!> \brief Set the attributes of an effective charge and inducible point
2492!> dipole potential data set.
2493!> \param potential ...
2494!> \param apol ...
2495!> \param cpol ...
2496!> \param qeff ...
2497!> \param mm_radius ...
2498!> \param qmmm_corr_radius ...
2499!> \param qmmm_radius ...
2500!> \date 05.03.2010
2501!> \author Toon.Verstraelen@gmail.com
2502! **************************************************************************************************
2503 SUBROUTINE set_fist_potential(potential, apol, cpol, qeff, mm_radius, &
2504 qmmm_corr_radius, qmmm_radius)
2505
2506 TYPE(fist_potential_type), INTENT(INOUT) :: potential
2507 REAL(kind=dp), INTENT(IN), OPTIONAL :: apol, cpol, qeff, mm_radius, &
2508 qmmm_corr_radius, qmmm_radius
2509
2510 IF (PRESENT(apol)) potential%apol = apol
2511 IF (PRESENT(cpol)) potential%cpol = cpol
2512 IF (PRESENT(mm_radius)) potential%mm_radius = mm_radius
2513 IF (PRESENT(qeff)) potential%qeff = qeff
2514 IF (PRESENT(qmmm_corr_radius)) potential%qmmm_corr_radius = qmmm_corr_radius
2515 IF (PRESENT(qmmm_radius)) potential%qmmm_radius = qmmm_radius
2516
2517 END SUBROUTINE set_fist_potential
2518
2519! **************************************************************************************************
2520!> \brief Set the attributes of a GTH potential data set.
2521!> \param potential ...
2522!> \param name ...
2523!> \param alpha_core_charge ...
2524!> \param alpha_ppl ...
2525!> \param ccore_charge ...
2526!> \param cerf_ppl ...
2527!> \param core_charge_radius ...
2528!> \param ppl_radius ...
2529!> \param ppnl_radius ...
2530!> \param lppnl ...
2531!> \param lprj_ppnl_max ...
2532!> \param nexp_ppl ...
2533!> \param nppnl ...
2534!> \param nprj_ppnl_max ...
2535!> \param z ...
2536!> \param zeff ...
2537!> \param zeff_correction ...
2538!> \param alpha_ppnl ...
2539!> \param cexp_ppl ...
2540!> \param elec_conf ...
2541!> \param nprj_ppnl ...
2542!> \param cprj ...
2543!> \param cprj_ppnl ...
2544!> \param vprj_ppnl ...
2545!> \param wprj_ppnl ...
2546!> \param hprj_ppnl ...
2547!> \param kprj_ppnl ...
2548!> \date 11.01.2002
2549!> \author MK
2550!> \version 1.0
2551! **************************************************************************************************
2552 SUBROUTINE set_gth_potential(potential, name, alpha_core_charge, alpha_ppl, &
2553 ccore_charge, cerf_ppl, core_charge_radius, &
2554 ppl_radius, ppnl_radius, lppnl, lprj_ppnl_max, &
2555 nexp_ppl, nppnl, nprj_ppnl_max, z, zeff, zeff_correction, &
2556 alpha_ppnl, cexp_ppl, elec_conf, nprj_ppnl, cprj, cprj_ppnl, &
2557 vprj_ppnl, wprj_ppnl, hprj_ppnl, kprj_ppnl)
2558
2559 TYPE(gth_potential_type), INTENT(INOUT) :: potential
2560 CHARACTER(LEN=default_string_length), INTENT(IN), &
2561 OPTIONAL :: name
2562 REAL(kind=dp), INTENT(IN), OPTIONAL :: alpha_core_charge, alpha_ppl, &
2563 ccore_charge, cerf_ppl, &
2564 core_charge_radius, ppl_radius, &
2565 ppnl_radius
2566 INTEGER, INTENT(IN), OPTIONAL :: lppnl, lprj_ppnl_max, nexp_ppl, nppnl, &
2567 nprj_ppnl_max, z
2568 REAL(kind=dp), INTENT(IN), OPTIONAL :: zeff, zeff_correction
2569 REAL(kind=dp), DIMENSION(:), OPTIONAL, POINTER :: alpha_ppnl, cexp_ppl
2570 INTEGER, DIMENSION(:), OPTIONAL, POINTER :: elec_conf, nprj_ppnl
2571 REAL(kind=dp), DIMENSION(:, :), OPTIONAL, POINTER :: cprj, cprj_ppnl, vprj_ppnl, wprj_ppnl
2572 REAL(kind=dp), DIMENSION(:, :, :), OPTIONAL, &
2573 POINTER :: hprj_ppnl, kprj_ppnl
2574
2575 IF (PRESENT(name)) potential%name = name
2576 IF (PRESENT(alpha_core_charge)) &
2577 potential%alpha_core_charge = alpha_core_charge
2578 IF (PRESENT(alpha_ppl)) potential%alpha_ppl = alpha_ppl
2579 IF (PRESENT(ccore_charge)) potential%ccore_charge = ccore_charge
2580 IF (PRESENT(cerf_ppl)) potential%cerf_ppl = cerf_ppl
2581 IF (PRESENT(core_charge_radius)) &
2582 potential%core_charge_radius = core_charge_radius
2583 IF (PRESENT(ppl_radius)) potential%ppl_radius = ppl_radius
2584 IF (PRESENT(ppnl_radius)) potential%ppnl_radius = ppnl_radius
2585 IF (PRESENT(lppnl)) potential%lppnl = lppnl
2586 IF (PRESENT(lprj_ppnl_max)) potential%lprj_ppnl_max = lprj_ppnl_max
2587 IF (PRESENT(nexp_ppl)) potential%nexp_ppl = nexp_ppl
2588 IF (PRESENT(nppnl)) potential%nppnl = nppnl
2589 IF (PRESENT(nprj_ppnl_max)) potential%nprj_ppnl_max = nprj_ppnl_max
2590 IF (PRESENT(z)) potential%z = z
2591 IF (PRESENT(zeff)) potential%zeff = zeff
2592 IF (PRESENT(zeff_correction)) potential%zeff_correction = zeff_correction
2593 IF (PRESENT(alpha_ppnl)) potential%alpha_ppnl => alpha_ppnl
2594 IF (PRESENT(cexp_ppl)) potential%cexp_ppl => cexp_ppl
2595 IF (PRESENT(elec_conf)) THEN
2596 IF (ASSOCIATED(potential%elec_conf)) THEN
2597 DEALLOCATE (potential%elec_conf)
2598 END IF
2599 ALLOCATE (potential%elec_conf(0:SIZE(elec_conf) - 1))
2600 potential%elec_conf(:) = elec_conf(:)
2601 END IF
2602 IF (PRESENT(nprj_ppnl)) potential%nprj_ppnl => nprj_ppnl
2603 IF (PRESENT(cprj)) potential%cprj => cprj
2604 IF (PRESENT(cprj_ppnl)) potential%cprj_ppnl => cprj_ppnl
2605 IF (PRESENT(hprj_ppnl)) potential%hprj_ppnl => hprj_ppnl
2606 IF (PRESENT(kprj_ppnl)) potential%kprj_ppnl => kprj_ppnl
2607 IF (PRESENT(vprj_ppnl)) potential%vprj_ppnl => vprj_ppnl
2608 IF (PRESENT(wprj_ppnl)) potential%wprj_ppnl => wprj_ppnl
2609
2610 END SUBROUTINE set_gth_potential
2611
2612! **************************************************************************************************
2613!> \brief ...
2614!> \param potential ...
2615!> \param name ...
2616!> \param description ...
2617!> \param aliases ...
2618!> \param elec_conf ...
2619!> \param z ...
2620!> \param zeff ...
2621!> \param zeff_correction ...
2622!> \param alpha_core_charge ...
2623!> \param ccore_charge ...
2624!> \param core_charge_radius ...
2625!> \param ppl_radius ...
2626!> \param ppnl_radius ...
2627!> \param ecp_local ...
2628!> \param n_local ...
2629!> \param a_local ...
2630!> \param c_local ...
2631!> \param nloc ...
2632!> \param nrloc ...
2633!> \param aloc ...
2634!> \param bloc ...
2635!> \param ecp_semi_local ...
2636!> \param sl_lmax ...
2637!> \param npot ...
2638!> \param nrpot ...
2639!> \param apot ...
2640!> \param bpot ...
2641!> \param n_nonlocal ...
2642!> \param nppnl ...
2643!> \param lmax ...
2644!> \param is_nonlocal ...
2645!> \param a_nonlocal ...
2646!> \param h_nonlocal ...
2647!> \param c_nonlocal ...
2648!> \param has_nlcc ...
2649!> \param n_nlcc ...
2650!> \param a_nlcc ...
2651!> \param c_nlcc ...
2652! **************************************************************************************************
2653 SUBROUTINE set_sgp_potential(potential, name, description, aliases, elec_conf, &
2654 z, zeff, zeff_correction, alpha_core_charge, &
2655 ccore_charge, core_charge_radius, &
2656 ppl_radius, ppnl_radius, &
2657 ecp_local, n_local, a_local, c_local, &
2658 nloc, nrloc, aloc, bloc, &
2659 ecp_semi_local, sl_lmax, npot, nrpot, apot, bpot, &
2660 n_nonlocal, nppnl, lmax, is_nonlocal, a_nonlocal, h_nonlocal, c_nonlocal, &
2661 has_nlcc, n_nlcc, a_nlcc, c_nlcc)
2662
2663 TYPE(sgp_potential_type), INTENT(INOUT) :: potential
2664 CHARACTER(LEN=default_string_length), INTENT(IN), &
2665 OPTIONAL :: name
2666 CHARACTER(LEN=default_string_length), &
2667 DIMENSION(4), INTENT(IN), OPTIONAL :: description
2668 CHARACTER(LEN=default_string_length), INTENT(IN), &
2669 OPTIONAL :: aliases
2670 INTEGER, DIMENSION(:), OPTIONAL, POINTER :: elec_conf
2671 INTEGER, INTENT(IN), OPTIONAL :: z
2672 REAL(kind=dp), INTENT(IN), OPTIONAL :: zeff, zeff_correction, &
2673 alpha_core_charge, ccore_charge, &
2674 core_charge_radius, ppl_radius, &
2675 ppnl_radius
2676 LOGICAL, INTENT(IN), OPTIONAL :: ecp_local
2677 INTEGER, INTENT(IN), OPTIONAL :: n_local
2678 REAL(kind=dp), DIMENSION(:), OPTIONAL, POINTER :: a_local, c_local
2679 INTEGER, INTENT(IN), OPTIONAL :: nloc
2680 INTEGER, DIMENSION(1:10), INTENT(IN), OPTIONAL :: nrloc
2681 REAL(dp), DIMENSION(1:10), INTENT(IN), OPTIONAL :: aloc, bloc
2682 LOGICAL, INTENT(IN), OPTIONAL :: ecp_semi_local
2683 INTEGER, INTENT(IN), OPTIONAL :: sl_lmax
2684 INTEGER, DIMENSION(0:10), OPTIONAL :: npot
2685 INTEGER, DIMENSION(1:15, 0:10), OPTIONAL :: nrpot
2686 REAL(dp), DIMENSION(1:15, 0:10), OPTIONAL :: apot, bpot
2687 INTEGER, INTENT(IN), OPTIONAL :: n_nonlocal, nppnl, lmax
2688 LOGICAL, DIMENSION(0:5), INTENT(IN), OPTIONAL :: is_nonlocal
2689 REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: a_nonlocal
2690 REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER :: h_nonlocal
2691 REAL(KIND=dp), DIMENSION(:, :, :), OPTIONAL, &
2692 POINTER :: c_nonlocal
2693 LOGICAL, INTENT(IN), OPTIONAL :: has_nlcc
2694 INTEGER, INTENT(IN), OPTIONAL :: n_nlcc
2695 REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: a_nlcc, c_nlcc
2696
2697 IF (PRESENT(name)) potential%name = name
2698 IF (PRESENT(aliases)) potential%aliases = aliases
2699 IF (PRESENT(description)) potential%description = description
2700
2701 IF (PRESENT(elec_conf)) THEN
2702 IF (ASSOCIATED(potential%elec_conf)) THEN
2703 DEALLOCATE (potential%elec_conf)
2704 END IF
2705 ALLOCATE (potential%elec_conf(0:SIZE(elec_conf) - 1))
2706 potential%elec_conf(:) = elec_conf(:)
2707 END IF
2708
2709 IF (PRESENT(z)) potential%z = z
2710 IF (PRESENT(zeff)) potential%zeff = zeff
2711 IF (PRESENT(zeff_correction)) potential%zeff_correction = zeff_correction
2712 IF (PRESENT(alpha_core_charge)) potential%alpha_core_charge = alpha_core_charge
2713 IF (PRESENT(ccore_charge)) potential%ccore_charge = ccore_charge
2714 IF (PRESENT(core_charge_radius)) potential%core_charge_radius = core_charge_radius
2715
2716 IF (PRESENT(ppl_radius)) potential%ppl_radius = ppl_radius
2717 IF (PRESENT(ppnl_radius)) potential%ppnl_radius = ppnl_radius
2718
2719 IF (PRESENT(ecp_local)) potential%ecp_local = ecp_local
2720 IF (PRESENT(n_local)) potential%n_local = n_local
2721 IF (PRESENT(a_local)) potential%a_local => a_local
2722 IF (PRESENT(c_local)) potential%c_local => c_local
2723
2724 IF (PRESENT(nloc)) potential%nloc = nloc
2725 IF (PRESENT(nrloc)) potential%nrloc = nrloc
2726 IF (PRESENT(aloc)) potential%aloc = aloc
2727 IF (PRESENT(bloc)) potential%bloc = bloc
2728
2729 IF (PRESENT(ecp_semi_local)) potential%ecp_semi_local = ecp_semi_local
2730 IF (PRESENT(sl_lmax)) potential%sl_lmax = sl_lmax
2731 IF (PRESENT(npot)) potential%npot = npot
2732 IF (PRESENT(nrpot)) potential%nrpot = nrpot
2733 IF (PRESENT(apot)) potential%apot = apot
2734 IF (PRESENT(bpot)) potential%bpot = bpot
2735
2736 IF (PRESENT(n_nonlocal)) potential%n_nonlocal = n_nonlocal
2737 IF (PRESENT(nppnl)) potential%nppnl = nppnl
2738 IF (PRESENT(lmax)) potential%lmax = lmax
2739 IF (PRESENT(is_nonlocal)) potential%is_nonlocal(:) = is_nonlocal(:)
2740 IF (PRESENT(a_nonlocal)) potential%a_nonlocal => a_nonlocal
2741 IF (PRESENT(c_nonlocal)) potential%c_nonlocal => c_nonlocal
2742 IF (PRESENT(h_nonlocal)) potential%h_nonlocal => h_nonlocal
2743
2744 IF (PRESENT(has_nlcc)) potential%has_nlcc = has_nlcc
2745 IF (PRESENT(n_nlcc)) potential%n_nlcc = n_nlcc
2746 IF (PRESENT(a_nlcc)) potential%a_nlcc => a_nlcc
2747 IF (PRESENT(c_nlcc)) potential%c_nlcc => c_nlcc
2748
2749 END SUBROUTINE set_sgp_potential
2750
2751! **************************************************************************************************
2752!> \brief Write an atomic all-electron potential data set to the output unit
2753!> \param potential ...
2754!> \param output_unit ...
2755!> \par History
2756!> - Creation (09.02.2002, MK)
2757! **************************************************************************************************
2758 SUBROUTINE write_all_potential(potential, output_unit)
2759
2760 TYPE(all_potential_type), INTENT(IN) :: potential
2761 INTEGER, INTENT(in) :: output_unit
2762
2763 CHARACTER(LEN=20) :: string
2764
2765 IF (output_unit > 0) THEN
2766 WRITE (unit=output_unit, fmt="(/,T6,A,T41,A40,/)") &
2767 "AE Potential information for", adjustr(trim(potential%name))
2768 WRITE (unit=output_unit, fmt="(T8,A,T41,A40)") &
2769 "Description: ", trim(potential%description(1)), &
2770 " ", trim(potential%description(2))
2771 WRITE (unit=output_unit, fmt="(/,T8,A,T69,F12.6)") &
2772 "Gaussian exponent of the core charge distribution: ", &
2773 potential%alpha_core_charge
2774 WRITE (unit=string, fmt="(5I4)") potential%elec_conf
2775 WRITE (unit=output_unit, fmt="(T8,A,T61,A20)") &
2776 "Electronic configuration (s p d ...):", &
2777 adjustr(trim(string))
2778 END IF
2779
2780 END SUBROUTINE write_all_potential
2781
2782! **************************************************************************************************
2783!> \brief Write an atomic local potential data set to the output unit
2784!> \param potential ...
2785!> \param output_unit ...
2786!> \par History
2787!> - Creation (24.01.2014, JGH)
2788! **************************************************************************************************
2789 SUBROUTINE write_local_potential(potential, output_unit)
2790
2791 TYPE(local_potential_type), INTENT(IN) :: potential
2792 INTEGER, INTENT(in) :: output_unit
2793
2794 INTEGER :: igau, ipol
2795
2796 IF (output_unit > 0) THEN
2797 WRITE (unit=output_unit, fmt="(/,T6,A,T41,A40)") &
2798 "Local Potential information for", adjustr(trim(potential%name))
2799 WRITE (unit=output_unit, fmt="(T8,A,T41,A40)") &
2800 "Description: ", trim(potential%description(1))
2801 DO igau = 1, potential%ngau
2802 WRITE (unit=output_unit, fmt="(T8,A,F12.6,T50,A,4(T68,I2,F10.4))") &
2803 "Exponent: ", potential%alpha(igau), &
2804 "Coefficients: ", (2*ipol - 2, potential%cval(igau, ipol), ipol=1, potential%npol)
2805 END DO
2806 END IF
2807
2808 END SUBROUTINE write_local_potential
2809
2810! **************************************************************************************************
2811!> \brief Write an atomic GTH potential data set to the output unit
2812!> \param potential ...
2813!> \param output_unit ...
2814!> \par History
2815!> - Creation (09.02.2002, MK)
2816! **************************************************************************************************
2817 SUBROUTINE write_gth_potential(potential, output_unit)
2818
2819 TYPE(gth_potential_type), INTENT(IN) :: potential
2820 INTEGER, INTENT(in) :: output_unit
2821
2822 CHARACTER(LEN=20) :: string
2823 INTEGER :: i, j, l
2824 REAL(KIND=dp) :: r
2825
2826 IF (output_unit > 0) THEN
2827 WRITE (unit=output_unit, fmt="(/,T6,A,T41,A40,/)") &
2828 "GTH Potential information for", adjustr(trim(potential%name))
2829 WRITE (unit=output_unit, fmt="(T8,A,T41,A40)") &
2830 "Description: ", adjustr(trim(potential%description(1))), &
2831 " ", adjustr(trim(potential%description(2))), &
2832 " ", adjustr(trim(potential%description(3))), &
2833 " ", adjustr(trim(potential%description(4)))
2834 WRITE (unit=output_unit, fmt="(/,T8,A,T69,F12.6)") &
2835 "Gaussian exponent of the core charge distribution: ", &
2836 potential%alpha_core_charge
2837 WRITE (unit=string, fmt="(5I4)") potential%elec_conf
2838 WRITE (unit=output_unit, fmt="(T8,A,T61,A20)") &
2839 "Electronic configuration (s p d ...):", &
2840 adjustr(trim(string))
2841
2842 r = 1.0_dp/sqrt(2.0_dp*potential%alpha_ppl)
2843
2844 WRITE (unit=output_unit, fmt="(/,T8,A,/,/,T27,A,/,T21,5F12.6)") &
2845 "Parameters of the local part of the GTH pseudopotential:", &
2846 "rloc C1 C2 C3 C4", &
2847 r, (potential%cexp_ppl(i)*r**(2*(i - 1)), i=1, potential%nexp_ppl)
2848
2849 IF (potential%lppnl > -1) THEN
2850 IF (potential%soc) THEN
2851 WRITE (unit=output_unit, fmt="(/,T8,A,/,/,(T20,A))") &
2852 "Parameters of the non-local part of the GTH (SOC) pseudopotential:", &
2853 "l r(l) h(i,j,l)", &
2854 " k(i,j,l)"
2855 ELSE
2856 WRITE (unit=output_unit, fmt="(/,T8,A,/,/,T20,A,/)") &
2857 "Parameters of the non-local part of the GTH pseudopotential:", &
2858 "l r(l) h(i,j,l)"
2859 END IF
2860 DO l = 0, potential%lppnl
2861 r = sqrt(0.5_dp/potential%alpha_ppnl(l))
2862 WRITE (unit=output_unit, fmt="(T19,I2,5F12.6)") &
2863 l, r, (potential%hprj_ppnl(1, j, l), j=1, potential%nprj_ppnl(l))
2864 DO i = 2, potential%nprj_ppnl(l)
2865 WRITE (unit=output_unit, fmt="(T33,4F12.6)") &
2866 (potential%hprj_ppnl(i, j, l), j=1, potential%nprj_ppnl(l))
2867 END DO
2868 IF (potential%soc .AND. (l > 0)) THEN
2869 DO i = 1, potential%nprj_ppnl(l)
2870 WRITE (unit=output_unit, fmt="(T33,4F12.6)") &
2871 (potential%kprj_ppnl(i, j, l), j=1, potential%nprj_ppnl(l))
2872 END DO
2873 END IF
2874 END DO
2875 END IF
2876 END IF
2877
2878 END SUBROUTINE write_gth_potential
2879
2880! **************************************************************************************************
2881!> \brief ...
2882!> \param potential ...
2883!> \param output_unit ...
2884! **************************************************************************************************
2885 SUBROUTINE write_sgp_potential(potential, output_unit)
2886
2887 TYPE(sgp_potential_type), INTENT(IN) :: potential
2888 INTEGER, INTENT(in) :: output_unit
2889
2890 CHARACTER(LEN=40) :: string
2891 INTEGER :: i, l
2892 CHARACTER(LEN=1), DIMENSION(0:10), PARAMETER :: &
2893 slqval = ["s", "p", "d", "f", "g", "h", "j", "k", "l", "m", "n"]
2894
2895 IF (output_unit > 0) THEN
2896 WRITE (unit=output_unit, fmt="(/,T6,A,T41,A40,/)") &
2897 "SGP Potential information for", adjustr(trim(potential%name))
2898 WRITE (unit=output_unit, fmt="(T8,A,T25,A56)") &
2899 "Description: ", adjustr(trim(potential%description(1))), &
2900 " ", adjustr(trim(potential%description(2))), &
2901 " ", adjustr(trim(potential%description(3))), &
2902 " ", adjustr(trim(potential%description(4)))
2903 WRITE (unit=output_unit, fmt="(/,T8,A,T69,F12.6)") &
2904 "Gaussian exponent of the core charge distribution: ", &
2905 potential%alpha_core_charge
2906 WRITE (unit=string, fmt="(10I4)") potential%elec_conf
2907 WRITE (unit=output_unit, fmt="(T8,A,T61,A20)") &
2908 "Electronic configuration (s p d ...):", &
2909 adjustr(trim(string))
2910 IF (potential%ecp_local) THEN
2911 IF (potential%nloc > 0) THEN
2912 WRITE (unit=output_unit, fmt="(/,T8,'Local pseudopotential')")
2913 WRITE (unit=output_unit, fmt="(T20,'r**(n-2)',T50,'Coefficient',T73,'Exponent')")
2914 DO i = 1, potential%nloc
2915 WRITE (unit=output_unit, fmt="(T20,I5,T47,F14.8,T69,F12.6)") &
2916 potential%nrloc(i), potential%aloc(i), potential%bloc(i)
2917 END DO
2918 END IF
2919 ELSE
2920 IF (potential%n_local > 0) THEN
2921 WRITE (unit=output_unit, fmt="(/,T8,'Local pseudopotential')")
2922 WRITE (unit=output_unit, fmt="(T8,A,10(T21,6F10.4,/))") &
2923 'Exponents:', potential%a_local(1:potential%n_local)
2924 WRITE (unit=output_unit, fmt="(T8,A,10(T21,6F10.4,/))") &
2925 'Coefficients:', potential%c_local(1:potential%n_local)
2926 END IF
2927 END IF
2928 IF (potential%ecp_semi_local) THEN
2929 WRITE (unit=output_unit, fmt="(/,T8,'Semi-local pseudopotential')")
2930 DO l = 0, potential%sl_lmax
2931 WRITE (unit=output_unit, fmt="(T8,A,A)") 'l-value: ', slqval(l)
2932 DO i = 1, potential%npot(l)
2933 WRITE (unit=output_unit, fmt="(T21,I5,2F20.8)") &
2934 potential%nrpot(i, l), potential%bpot(i, l), potential%apot(i, l)
2935 END DO
2936 END DO
2937 END IF
2938 ! nonlocal PP
2939 IF (potential%n_nonlocal > 0) THEN
2940 WRITE (unit=output_unit, fmt="(/,T8,'Nonlocal pseudopotential')")
2941 WRITE (unit=output_unit, fmt="(T8,A,T71,I10)") 'Total number of projectors:', potential%nppnl
2942 WRITE (unit=output_unit, fmt="(T8,A,10(T21,6F10.4,/))") &
2943 'Exponents:', potential%a_nonlocal(1:potential%n_nonlocal)
2944 DO l = 0, potential%lmax
2945 WRITE (unit=output_unit, fmt="(T8,'Coupling for l=',I4)") l
2946 WRITE (unit=output_unit, fmt="(10(T21,6F10.4,/))") &
2947 potential%h_nonlocal(1:potential%n_nonlocal, l)
2948 END DO
2949 END IF
2950 !
2951 IF (potential%has_nlcc) THEN
2952 WRITE (unit=output_unit, fmt="(/,T8,'Nonlinear Core Correction')")
2953 WRITE (unit=output_unit, fmt="(T8,A,10(T21,6F10.4,/))") &
2954 'Exponents:', potential%a_nlcc(1:potential%n_nlcc)
2955 WRITE (unit=output_unit, fmt="(T8,A,10(T21,6F10.4,/))") &
2956 'Coefficients:', potential%c_nlcc(1:potential%n_nlcc)
2957 END IF
2958 END IF
2959
2960 END SUBROUTINE write_sgp_potential
2961
2962! **************************************************************************************************
2963!> \brief Copy an all_potential_type to a new, unallocated variable
2964!> \param pot_in the input potential to copy
2965!> \param pot_out the newly copied and allocated potential
2966!> \par History
2967!> - Creation (12.2019, A. Bussy)
2968! **************************************************************************************************
2969 SUBROUTINE copy_all_potential(pot_in, pot_out)
2970
2971 TYPE(all_potential_type), INTENT(IN) :: pot_in
2972 TYPE(all_potential_type), INTENT(INOUT), POINTER :: pot_out
2973
2974 CALL allocate_all_potential(pot_out)
2975
2976 pot_out%name = pot_in%name
2977 pot_out%alpha_core_charge = pot_in%alpha_core_charge
2978 pot_out%ccore_charge = pot_in%ccore_charge
2979 pot_out%core_charge_radius = pot_in%core_charge_radius
2980 pot_out%zeff = pot_in%zeff
2981 pot_out%zeff_correction = pot_in%zeff_correction
2982 pot_out%z = pot_in%z
2983
2984 IF (ASSOCIATED(pot_in%elec_conf)) THEN
2985 ALLOCATE (pot_out%elec_conf(lbound(pot_in%elec_conf, 1):ubound(pot_in%elec_conf, 1)))
2986 pot_out%elec_conf(:) = pot_in%elec_conf(:)
2987 END IF
2988
2989 END SUBROUTINE copy_all_potential
2990
2991! **************************************************************************************************
2992!> \brief Copy a gth_potential_type to a new, unallocated variable
2993!> \param pot_in the input potential to copy
2994!> \param pot_out the newly copied and allocated potential
2995!> \par History
2996!> - Creation (12.2019, A. Bussy)
2997! **************************************************************************************************
2998 SUBROUTINE copy_gth_potential(pot_in, pot_out)
2999
3000 TYPE(gth_potential_type), INTENT(IN) :: pot_in
3001 TYPE(gth_potential_type), INTENT(INOUT), POINTER :: pot_out
3002
3003 CALL allocate_gth_potential(pot_out)
3004
3005 pot_out%name = pot_in%name
3006 pot_out%aliases = pot_in%aliases
3007 pot_out%alpha_core_charge = pot_in%alpha_core_charge
3008 pot_out%alpha_ppl = pot_in%alpha_ppl
3009 pot_out%ccore_charge = pot_in%ccore_charge
3010 pot_out%cerf_ppl = pot_in%cerf_ppl
3011 pot_out%zeff = pot_in%zeff
3012 pot_out%core_charge_radius = pot_in%core_charge_radius
3013 pot_out%ppl_radius = pot_in%ppl_radius
3014 pot_out%ppnl_radius = pot_in%ppnl_radius
3015 pot_out%zeff_correction = pot_in%zeff_correction
3016 pot_out%lppnl = pot_in%lppnl
3017 pot_out%lprj_ppnl_max = pot_in%lprj_ppnl_max
3018 pot_out%nexp_ppl = pot_in%nexp_ppl
3019 pot_out%nppnl = pot_in%nppnl
3020 pot_out%nprj_ppnl_max = pot_in%nprj_ppnl_max
3021 pot_out%z = pot_in%z
3022 pot_out%nlcc = pot_in%nlcc
3023 pot_out%nexp_nlcc = pot_in%nexp_nlcc
3024 pot_out%lsdpot = pot_in%lsdpot
3025 pot_out%nexp_lsd = pot_in%nexp_lsd
3026 pot_out%lpotextended = pot_in%lpotextended
3027 pot_out%nexp_lpot = pot_in%nexp_lpot
3028
3029 IF (ASSOCIATED(pot_in%alpha_ppnl)) THEN
3030 ALLOCATE (pot_out%alpha_ppnl(lbound(pot_in%alpha_ppnl, 1):ubound(pot_in%alpha_ppnl, 1)))
3031 pot_out%alpha_ppnl(:) = pot_in%alpha_ppnl(:)
3032 END IF
3033 IF (ASSOCIATED(pot_in%cexp_ppl)) THEN
3034 ALLOCATE (pot_out%cexp_ppl(lbound(pot_in%cexp_ppl, 1):ubound(pot_in%cexp_ppl, 1)))
3035 pot_out%cexp_ppl(:) = pot_in%cexp_ppl(:)
3036 END IF
3037 IF (ASSOCIATED(pot_in%elec_conf)) THEN
3038 ALLOCATE (pot_out%elec_conf(lbound(pot_in%elec_conf, 1):ubound(pot_in%elec_conf, 1)))
3039 pot_out%elec_conf(:) = pot_in%elec_conf(:)
3040 END IF
3041 IF (ASSOCIATED(pot_in%nprj_ppnl)) THEN
3042 ALLOCATE (pot_out%nprj_ppnl(lbound(pot_in%nprj_ppnl, 1):ubound(pot_in%nprj_ppnl, 1)))
3043 pot_out%nprj_ppnl(:) = pot_in%nprj_ppnl(:)
3044 END IF
3045 IF (ASSOCIATED(pot_in%cprj)) THEN
3046 ALLOCATE (pot_out%cprj(lbound(pot_in%cprj, 1):ubound(pot_in%cprj, 1), &
3047 lbound(pot_in%cprj, 2):ubound(pot_in%cprj, 2)))
3048 pot_out%cprj(:, :) = pot_in%cprj(:, :)
3049 END IF
3050 IF (ASSOCIATED(pot_in%cprj_ppnl)) THEN
3051 ALLOCATE (pot_out%cprj_ppnl(lbound(pot_in%cprj_ppnl, 1):ubound(pot_in%cprj_ppnl, 1), &
3052 lbound(pot_in%cprj_ppnl, 2):ubound(pot_in%cprj_ppnl, 2)))
3053 pot_out%cprj_ppnl(:, :) = pot_in%cprj_ppnl(:, :)
3054 END IF
3055 IF (ASSOCIATED(pot_in%hprj_ppnl)) THEN
3056 ALLOCATE (pot_out%hprj_ppnl(lbound(pot_in%hprj_ppnl, 1):ubound(pot_in%hprj_ppnl, 1), &
3057 lbound(pot_in%hprj_ppnl, 2):ubound(pot_in%hprj_ppnl, 2), &
3058 lbound(pot_in%hprj_ppnl, 3):ubound(pot_in%hprj_ppnl, 3)))
3059 pot_out%hprj_ppnl(:, :, :) = pot_in%hprj_ppnl(:, :, :)
3060 END IF
3061 IF (ASSOCIATED(pot_in%kprj_ppnl)) THEN
3062 ALLOCATE (pot_out%kprj_ppnl(lbound(pot_in%kprj_ppnl, 1):ubound(pot_in%kprj_ppnl, 1), &
3063 lbound(pot_in%kprj_ppnl, 2):ubound(pot_in%kprj_ppnl, 2), &
3064 lbound(pot_in%kprj_ppnl, 3):ubound(pot_in%kprj_ppnl, 3)))
3065 pot_out%kprj_ppnl(:, :, :) = pot_in%kprj_ppnl(:, :, :)
3066 END IF
3067 IF (ASSOCIATED(pot_in%vprj_ppnl)) THEN
3068 ALLOCATE (pot_out%vprj_ppnl(lbound(pot_in%vprj_ppnl, 1):ubound(pot_in%vprj_ppnl, 1), &
3069 lbound(pot_in%vprj_ppnl, 2):ubound(pot_in%vprj_ppnl, 2)))
3070 pot_out%vprj_ppnl(:, :) = pot_in%vprj_ppnl(:, :)
3071 END IF
3072 IF (ASSOCIATED(pot_in%wprj_ppnl)) THEN
3073 ALLOCATE (pot_out%wprj_ppnl(lbound(pot_in%wprj_ppnl, 1):ubound(pot_in%wprj_ppnl, 1), &
3074 lbound(pot_in%wprj_ppnl, 2):ubound(pot_in%wprj_ppnl, 2)))
3075 pot_out%wprj_ppnl(:, :) = pot_in%wprj_ppnl(:, :)
3076 END IF
3077 IF (ASSOCIATED(pot_in%alpha_nlcc)) THEN
3078 ALLOCATE (pot_out%alpha_nlcc(lbound(pot_in%alpha_nlcc, 1):ubound(pot_in%alpha_nlcc, 1)))
3079 pot_out%alpha_nlcc(:) = pot_in%alpha_nlcc(:)
3080 END IF
3081 IF (ASSOCIATED(pot_in%nct_nlcc)) THEN
3082 ALLOCATE (pot_out%nct_nlcc(lbound(pot_in%nct_nlcc, 1):ubound(pot_in%nct_nlcc, 1)))
3083 pot_out%nct_nlcc(:) = pot_in%nct_nlcc(:)
3084 END IF
3085 IF (ASSOCIATED(pot_in%cval_nlcc)) THEN
3086 ALLOCATE (pot_out%cval_nlcc(lbound(pot_in%cval_nlcc, 1):ubound(pot_in%cval_nlcc, 1), &
3087 lbound(pot_in%cval_nlcc, 2):ubound(pot_in%cval_nlcc, 2)))
3088 pot_out%cval_nlcc(:, :) = pot_in%cval_nlcc(:, :)
3089 END IF
3090 IF (ASSOCIATED(pot_in%alpha_lsd)) THEN
3091 ALLOCATE (pot_out%alpha_lsd(lbound(pot_in%alpha_lsd, 1):ubound(pot_in%alpha_lsd, 1)))
3092 pot_out%alpha_lsd(:) = pot_in%alpha_lsd(:)
3093 END IF
3094 IF (ASSOCIATED(pot_in%nct_lsd)) THEN
3095 ALLOCATE (pot_out%nct_lsd(lbound(pot_in%nct_lsd, 1):ubound(pot_in%nct_lsd, 1)))
3096 pot_out%nct_lsd(:) = pot_in%nct_lsd(:)
3097 END IF
3098 IF (ASSOCIATED(pot_in%cval_lsd)) THEN
3099 ALLOCATE (pot_out%cval_lsd(lbound(pot_in%cval_lsd, 1):ubound(pot_in%cval_lsd, 1), &
3100 lbound(pot_in%cval_lsd, 2):ubound(pot_in%cval_lsd, 2)))
3101 pot_out%cval_lsd(:, :) = pot_in%cval_lsd(:, :)
3102 END IF
3103 IF (ASSOCIATED(pot_in%alpha_lpot)) THEN
3104 ALLOCATE (pot_out%alpha_lpot(lbound(pot_in%alpha_lpot, 1):ubound(pot_in%alpha_lpot, 1)))
3105 pot_out%alpha_lpot(:) = pot_in%alpha_lpot(:)
3106 END IF
3107 IF (ASSOCIATED(pot_in%nct_lpot)) THEN
3108 ALLOCATE (pot_out%nct_lpot(lbound(pot_in%nct_lpot, 1):ubound(pot_in%nct_lpot, 1)))
3109 pot_out%nct_lpot(:) = pot_in%nct_lpot(:)
3110 END IF
3111 IF (ASSOCIATED(pot_in%cval_lpot)) THEN
3112 ALLOCATE (pot_out%cval_lpot(lbound(pot_in%cval_lpot, 1):ubound(pot_in%cval_lpot, 1), &
3113 lbound(pot_in%cval_lpot, 2):ubound(pot_in%cval_lpot, 2)))
3114 pot_out%cval_lpot(:, :) = pot_in%cval_lpot(:, :)
3115 END IF
3116
3117 END SUBROUTINE copy_gth_potential
3118
3119! **************************************************************************************************
3120!> \brief Copy a sgp_potential_type to a new, unallocated variable
3121!> \param pot_in the input potential to copy
3122!> \param pot_out the newly copied and allocated potential
3123!> \par History
3124!> - Creation (12.2019, A. Bussy)
3125! **************************************************************************************************
3126 SUBROUTINE copy_sgp_potential(pot_in, pot_out)
3127
3128 TYPE(sgp_potential_type), INTENT(IN) :: pot_in
3129 TYPE(sgp_potential_type), INTENT(INOUT), POINTER :: pot_out
3130
3131 CALL allocate_sgp_potential(pot_out)
3132
3133 pot_out%name = pot_in%name
3134 pot_out%aliases = pot_in%aliases
3135 pot_out%z = pot_in%z
3136 pot_out%zeff = pot_in%zeff
3137 pot_out%zeff_correction = pot_in%zeff_correction
3138 pot_out%alpha_core_charge = pot_in%alpha_core_charge
3139 pot_out%ccore_charge = pot_in%ccore_charge
3140 pot_out%core_charge_radius = pot_in%core_charge_radius
3141 pot_out%ppl_radius = pot_in%ppl_radius
3142 pot_out%ppnl_radius = pot_in%ppnl_radius
3143 pot_out%ecp_local = pot_in%ecp_local
3144 pot_out%n_local = pot_in%n_local
3145 pot_out%nloc = pot_in%nloc
3146 pot_out%nrloc = pot_in%nrloc
3147 pot_out%aloc = pot_in%aloc
3148 pot_out%bloc = pot_in%bloc
3149 pot_out%ecp_semi_local = pot_in%ecp_semi_local
3150 pot_out%sl_lmax = pot_in%sl_lmax
3151 pot_out%npot = pot_in%npot
3152 pot_out%nrpot = pot_in%nrpot
3153 pot_out%apot = pot_in%apot
3154 pot_out%bpot = pot_in%bpot
3155 pot_out%n_nonlocal = pot_in%n_nonlocal
3156 pot_out%nppnl = pot_in%nppnl
3157 pot_out%lmax = pot_in%lmax
3158 pot_out%is_nonlocal = pot_in%is_nonlocal
3159 pot_out%has_nlcc = pot_in%has_nlcc
3160 pot_out%n_nlcc = pot_in%n_nlcc
3161
3162 IF (ASSOCIATED(pot_in%elec_conf)) THEN
3163 ALLOCATE (pot_out%elec_conf(lbound(pot_in%elec_conf, 1):ubound(pot_in%elec_conf, 1)))
3164 pot_out%elec_conf(:) = pot_in%elec_conf(:)
3165 END IF
3166 IF (ASSOCIATED(pot_in%a_local)) THEN
3167 ALLOCATE (pot_out%a_local(lbound(pot_in%a_local, 1):ubound(pot_in%a_local, 1)))
3168 pot_out%a_local(:) = pot_in%a_local(:)
3169 END IF
3170 IF (ASSOCIATED(pot_in%c_local)) THEN
3171 ALLOCATE (pot_out%c_local(lbound(pot_in%c_local, 1):ubound(pot_in%c_local, 1)))
3172 pot_out%c_local(:) = pot_in%c_local(:)
3173 END IF
3174 IF (ASSOCIATED(pot_in%a_nonlocal)) THEN
3175 ALLOCATE (pot_out%a_nonlocal(lbound(pot_in%a_nonlocal, 1):ubound(pot_in%a_nonlocal, 1)))
3176 pot_out%a_nonlocal(:) = pot_in%a_nonlocal(:)
3177 END IF
3178 IF (ASSOCIATED(pot_in%h_nonlocal)) THEN
3179 ALLOCATE (pot_out%h_nonlocal(lbound(pot_in%h_nonlocal, 1):ubound(pot_in%h_nonlocal, 1), &
3180 lbound(pot_in%h_nonlocal, 2):ubound(pot_in%h_nonlocal, 2)))
3181 pot_out%h_nonlocal(:, :) = pot_in%h_nonlocal(:, :)
3182 END IF
3183 IF (ASSOCIATED(pot_in%c_nonlocal)) THEN
3184 ALLOCATE (pot_out%c_nonlocal(lbound(pot_in%c_nonlocal, 1):ubound(pot_in%c_nonlocal, 1), &
3185 lbound(pot_in%c_nonlocal, 2):ubound(pot_in%c_nonlocal, 2), &
3186 lbound(pot_in%c_nonlocal, 3):ubound(pot_in%c_nonlocal, 3)))
3187 pot_out%c_nonlocal(:, :, :) = pot_in%c_nonlocal(:, :, :)
3188 END IF
3189 IF (ASSOCIATED(pot_in%cprj_ppnl)) THEN
3190 ALLOCATE (pot_out%cprj_ppnl(lbound(pot_in%cprj_ppnl, 1):ubound(pot_in%cprj_ppnl, 1), &
3191 lbound(pot_in%cprj_ppnl, 2):ubound(pot_in%cprj_ppnl, 2)))
3192 pot_out%cprj_ppnl(:, :) = pot_in%cprj_ppnl(:, :)
3193 END IF
3194 IF (ASSOCIATED(pot_in%vprj_ppnl)) THEN
3195 ALLOCATE (pot_out%vprj_ppnl(lbound(pot_in%vprj_ppnl, 1):ubound(pot_in%vprj_ppnl, 1)))
3196 pot_out%vprj_ppnl(:) = pot_in%vprj_ppnl(:)
3197 END IF
3198 IF (ASSOCIATED(pot_in%a_nlcc)) THEN
3199 ALLOCATE (pot_out%a_nlcc(lbound(pot_in%a_nlcc, 1):ubound(pot_in%a_nlcc, 1)))
3200 pot_out%a_nlcc(:) = pot_in%a_nlcc(:)
3201 END IF
3202 IF (ASSOCIATED(pot_in%c_nlcc)) THEN
3203 ALLOCATE (pot_out%c_nlcc(lbound(pot_in%c_nlcc, 1):ubound(pot_in%c_nlcc, 1)))
3204 pot_out%c_nlcc(:) = pot_in%c_nlcc(:)
3205 END IF
3206
3207 END SUBROUTINE copy_sgp_potential
3208
3209END MODULE external_potential_types
All kind of helpful little routines.
Definition ao_util.F:14
real(kind=dp) function, public exp_radius(l, alpha, threshold, prefactor, epsabs, epsrel, rlow)
The radius of a primitive Gaussian function for a given threshold is calculated. g(r) = prefactor*r**...
Definition ao_util.F:96
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public goedecker1996
integer, save, public hartwigsen1998
integer, save, public krack2000
integer, save, public krack2005
logical function, public cp_sll_val_next(iterator, el_att)
returns true if the actual element is valid (i.e. iterator ont at end) moves the iterator to the next...
Utility routines to read data from files. Kept as close as possible to the old parser because.
subroutine, public parser_get_next_line(parser, nline, at_end)
Read the next input line and broadcast the input information. Skip (nline-1) lines and skip also all ...
character(len=3) function, public parser_test_next_token(parser, string_length)
Test next input object.
subroutine, public parser_search_string(parser, string, ignore_case, found, line, begin_line, search_from_begin_of_file)
Search a string pattern in a file defined by its logical unit number "unit". A case sensitive search ...
Utility routines to read data from files. Kept as close as possible to the old parser because.
subroutine, public parser_release(parser)
releases the parser
subroutine, public parser_create(parser, file_name, unit_nr, para_env, end_section_label, separator_chars, comment_char, continuation_char, quote_char, section_char, parse_white_lines, initial_variables, apply_preprocessing)
Start a parser run. Initial variables allow to @SET stuff before opening the file.
Definition of the atomic potential types.
subroutine, public set_default_all_potential(potential, z, zeff_correction)
...
objects that represent the structure of input sections and the data contained in an input section
subroutine, public section_vals_val_set(section_vals, keyword_name, i_rep_section, i_rep_val, val, l_val, i_val, r_val, c_val, l_vals_ptr, i_vals_ptr, r_vals_ptr, c_vals_ptr)
sets the requested value
subroutine, public section_vals_list_get(section_vals, keyword_name, i_rep_section, list)
returns the requested list
subroutine, public section_vals_get(section_vals, ref_count, n_repetition, n_subs_vals_rep, section, explicit)
returns various attributes about the section_vals
a wrapper for basic fortran types.
subroutine, public val_get(val, has_l, has_i, has_r, has_lc, has_c, l_val, l_vals, i_val, i_vals, r_val, r_vals, c_val, c_vals, len_c, type_of_var, enum)
returns the stored values
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
An array-based list which grows on demand. When the internal array is full, a new array of twice the ...
Definition list.F:24
Definition of mathematical constants and functions.
real(kind=dp), parameter, public pi
real(kind=dp), dimension(-1:2 *maxfac+1), parameter, public dfac
real(kind=dp), parameter, public rootpi
real(kind=dp), dimension(0:maxfac), parameter, public fac
Collection of simple mathematical functions and subroutines.
Definition mathlib.F:15
subroutine, public symmetrize_matrix(a, option)
Symmetrize the matrix a.
Definition mathlib.F:1204
Utility routines for the memory handling.
Interface to the message passing library MPI.
Provides Cartesian and spherical orbital pointers and indices.
subroutine, public init_orbital_pointers(maxl)
Initialize or update the orbital pointers.
integer, dimension(:, :, :), allocatable, public co
integer, dimension(:), allocatable, public nco
integer, dimension(:), allocatable, public ncoset
integer, dimension(:, :, :), allocatable, public coset
integer, dimension(:), allocatable, public nso
Calculation of the spherical harmonics and the corresponding orbital transformation matrices.
type(orbtramat_type), dimension(:), pointer, public orbtramat
Periodic Table related data definitions.
type(atom), dimension(0:nelem), public ptable
Utilities for string manipulations.
character(len=1), parameter, public newline
subroutine, public remove_word(string)
remove a word from a string (words are separated by white spaces)
elemental subroutine, public uppercase(string)
Convert all lower case characters in a string to upper case.
represent a single linked list that stores pointers to the elements
a type to have a wrapper that stores any basic fortran type
stores all the informations relevant to an mpi environment