-
Notifications
You must be signed in to change notification settings - Fork 1
/
atom_utils.F
2754 lines (2419 loc) · 110 KB
/
atom_utils.F
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!--------------------------------------------------------------------------------------------------!
! CP2K: A general program to perform molecular dynamics simulations !
! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
! !
! SPDX-License-Identifier: GPL-2.0-or-later !
!--------------------------------------------------------------------------------------------------!
! **************************************************************************************************
!> \brief Some basic routines for atomic calculations
!> \author jgh
!> \date 01.04.2008
!> \version 1.0
!>
! **************************************************************************************************
MODULE atom_utils
USE ai_onecenter, ONLY: sg_overlap,&
sto_overlap
USE ai_overlap, ONLY: overlap_ab_s,&
overlap_ab_sp
USE ao_util, ONLY: exp_radius
USE atom_types, ONLY: &
CGTO_BASIS, GTO_BASIS, NUM_BASIS, STO_BASIS, atom_basis_type, atom_gthpot_type, &
atom_hfx_type, atom_potential_type, atom_state, atom_type, ecp_pseudo, eri, gth_pseudo, &
lmat, no_pseudo, sgp_pseudo, upf_pseudo
USE basis_set_types, ONLY: srules
USE cp_files, ONLY: close_file,&
get_unit_number,&
open_file
USE input_constants, ONLY: do_rhf_atom,&
do_rks_atom,&
do_rohf_atom,&
do_uhf_atom,&
do_uks_atom
USE kahan_sum, ONLY: accurate_dot_product
USE kinds, ONLY: default_string_length,&
dp
USE lapack, ONLY: lapack_ssyev
USE mathconstants, ONLY: dfac,&
fac,&
fourpi,&
maxfac,&
pi,&
rootpi
USE mathlib, ONLY: binomial_gen,&
invmat_symm
USE orbital_pointers, ONLY: deallocate_orbital_pointers,&
init_orbital_pointers
USE orbital_transformation_matrices, ONLY: deallocate_spherical_harmonics,&
init_spherical_harmonics
USE periodic_table, ONLY: nelem,&
ptable
USE physcon, ONLY: bohr
USE qs_grid_atom, ONLY: grid_atom_type
USE splines, ONLY: spline3ders
USE string_utilities, ONLY: uppercase
#include "./base/base_uses.f90"
IMPLICIT NONE
PRIVATE
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'atom_utils'
PUBLIC :: atom_condnumber, atom_completeness, atom_basis_condnum
PUBLIC :: atom_set_occupation, get_maxl_occ, get_maxn_occ
PUBLIC :: atom_denmat, atom_density, atom_core_density
PUBLIC :: integrate_grid, atom_trace, atom_solve
PUBLIC :: coulomb_potential_numeric, coulomb_potential_analytic
PUBLIC :: exchange_numeric, exchange_semi_analytic
PUBLIC :: numpot_matrix, ceri_contract, eeri_contract, err_matrix
PUBLIC :: slater_density, wigner_slater_functional, atom_local_potential
PUBLIC :: atom_orbital_charge, atom_orbital_nodes, atom_consistent_method
PUBLIC :: atom_orbital_max, atom_wfnr0, get_rho0
PUBLIC :: contract2, contract4, contract2add
! ZMP added public subroutines
PUBLIC :: atom_read_external_density
PUBLIC :: atom_read_external_vxc
PUBLIC :: atom_read_zmp_restart
PUBLIC :: atom_write_zmp_restart
!-----------------------------------------------------------------------------!
INTERFACE integrate_grid
MODULE PROCEDURE integrate_grid_function1, &
integrate_grid_function2, &
integrate_grid_function3
END INTERFACE
! **************************************************************************************************
CONTAINS
! **************************************************************************************************
!> \brief Set occupation of atomic orbitals.
!> \param ostring list of electronic states
!> \param occupation ...
!> \param wfnocc ...
!> \param multiplicity ...
!> \par History
!> * 11.2009 unrestricted KS and HF methods [Juerg Hutter]
!> * 08.2008 created [Juerg Hutter]
! **************************************************************************************************
SUBROUTINE atom_set_occupation(ostring, occupation, wfnocc, multiplicity)
CHARACTER(LEN=default_string_length), &
DIMENSION(:), POINTER :: ostring
REAL(Kind=dp), DIMENSION(0:lmat, 10) :: occupation, wfnocc
INTEGER, INTENT(OUT), OPTIONAL :: multiplicity
CHARACTER(len=2) :: elem
CHARACTER(LEN=default_string_length) :: pstring
INTEGER :: i, i1, i2, ielem, is, jd, jf, jp, js, k, &
l, mult, n, no
REAL(Kind=dp) :: e0, el, oo
occupation = 0._dp
CPASSERT(ASSOCIATED(ostring))
CPASSERT(SIZE(ostring) > 0)
no = SIZE(ostring)
is = 1
! look for multiplicity
mult = -1 !not specified
IF (is <= no) THEN
IF (INDEX(ostring(is), "(") /= 0) THEN
i1 = INDEX(ostring(is), "(")
i2 = INDEX(ostring(is), ")")
CPASSERT((i2 - i1 - 1 > 0) .AND. (i2 - i1 - 1 < 3))
elem = ostring(is) (i1 + 1:i2 - 1)
IF (INDEX(elem, "HS") /= 0) THEN
mult = -2 !High spin
ELSE IF (INDEX(elem, "LS") /= 0) THEN
mult = -3 !Low spin
ELSE
READ (elem, *) mult
END IF
is = is + 1
END IF
END IF
IF (is <= no) THEN
IF (INDEX(ostring(is), "CORE") /= 0) is = is + 1 !Pseudopotential detected
END IF
IF (is <= no) THEN
IF (INDEX(ostring(is), "none") /= 0) is = is + 1 !no electrons, used with CORE
END IF
IF (is <= no) THEN
IF (INDEX(ostring(is), "[") /= 0) THEN
! core occupation from element [XX]
i1 = INDEX(ostring(is), "[")
i2 = INDEX(ostring(is), "]")
CPASSERT((i2 - i1 - 1 > 0) .AND. (i2 - i1 - 1 < 3))
elem = ostring(is) (i1 + 1:i2 - 1)
ielem = 0
DO k = 1, nelem
IF (elem == ptable(k)%symbol) THEN
ielem = k
EXIT
END IF
END DO
CPASSERT(ielem /= 0)
DO l = 0, MIN(lmat, UBOUND(ptable(ielem)%e_conv, 1))
el = 2._dp*(2._dp*REAL(l, dp) + 1._dp)
e0 = ptable(ielem)%e_conv(l)
DO k = 1, 10
occupation(l, k) = MIN(el, e0)
e0 = e0 - el
IF (e0 <= 0._dp) EXIT
END DO
END DO
is = is + 1
END IF
END IF
DO i = is, no
pstring = ostring(i)
CALL uppercase(pstring)
js = INDEX(pstring, "S")
jp = INDEX(pstring, "P")
jd = INDEX(pstring, "D")
jf = INDEX(pstring, "F")
CPASSERT(js + jp + jd + jf > 0)
IF (js > 0) THEN
CPASSERT(jp + jd + jf == 0)
READ (pstring(1:js - 1), *) n
READ (pstring(js + 1:), *) oo
CPASSERT(n > 0)
CPASSERT(oo >= 0._dp)
CPASSERT(occupation(0, n) == 0)
occupation(0, n) = oo
END IF
IF (jp > 0) THEN
CPASSERT(js + jd + jf == 0)
READ (pstring(1:jp - 1), *) n
READ (pstring(jp + 1:), *) oo
CPASSERT(n > 1)
CPASSERT(oo >= 0._dp)
CPASSERT(occupation(1, n - 1) == 0)
occupation(1, n - 1) = oo
END IF
IF (jd > 0) THEN
CPASSERT(js + jp + jf == 0)
READ (pstring(1:jd - 1), *) n
READ (pstring(jd + 1:), *) oo
CPASSERT(n > 2)
CPASSERT(oo >= 0._dp)
CPASSERT(occupation(2, n - 2) == 0)
occupation(2, n - 2) = oo
END IF
IF (jf > 0) THEN
CPASSERT(js + jp + jd == 0)
READ (pstring(1:jf - 1), *) n
READ (pstring(jf + 1:), *) oo
CPASSERT(n > 3)
CPASSERT(oo >= 0._dp)
CPASSERT(occupation(3, n - 3) == 0)
occupation(3, n - 3) = oo
END IF
END DO
wfnocc = 0._dp
DO l = 0, lmat
k = 0
DO i = 1, 10
IF (occupation(l, i) /= 0._dp) THEN
k = k + 1
wfnocc(l, k) = occupation(l, i)
END IF
END DO
END DO
!Check for consistency with multiplicity
IF (mult /= -1) THEN
! count open shells
js = 0
DO l = 0, lmat
k = 2*(2*l + 1)
DO i = 1, 10
IF (wfnocc(l, i) /= 0._dp .AND. wfnocc(l, i) /= REAL(k, dp)) THEN
js = js + 1
i1 = l
i2 = i
END IF
END DO
END DO
IF (js == 0 .AND. mult == -2) mult = 1
IF (js == 0 .AND. mult == -3) mult = 1
IF (js == 0) THEN
CPASSERT(mult == 1)
END IF
IF (js == 1) THEN
l = i1
i = i2
k = NINT(wfnocc(l, i))
IF (k > (2*l + 1)) k = 2*(2*l + 1) - k
IF (mult == -2) mult = k + 1
IF (mult == -3) mult = MOD(k, 2) + 1
CPASSERT(MOD(k + 1 - mult, 2) == 0)
END IF
IF (js > 1 .AND. mult /= -2) THEN
CPASSERT(mult == -2)
END IF
END IF
IF (PRESENT(multiplicity)) multiplicity = mult
END SUBROUTINE atom_set_occupation
! **************************************************************************************************
!> \brief Return the maximum orbital quantum number of occupied orbitals.
!> \param occupation ...
!> \return ...
!> \par History
!> * 08.2008 created [Juerg Hutter]
! **************************************************************************************************
PURE FUNCTION get_maxl_occ(occupation) RESULT(maxl)
REAL(Kind=dp), DIMENSION(0:lmat, 10), INTENT(IN) :: occupation
INTEGER :: maxl
INTEGER :: l
maxl = 0
DO l = 0, lmat
IF (SUM(occupation(l, :)) /= 0._dp) maxl = l
END DO
END FUNCTION get_maxl_occ
! **************************************************************************************************
!> \brief Return the maximum principal quantum number of occupied orbitals.
!> \param occupation ...
!> \return ...
!> \par History
!> * 08.2008 created [Juerg Hutter]
! **************************************************************************************************
PURE FUNCTION get_maxn_occ(occupation) RESULT(maxn)
REAL(Kind=dp), DIMENSION(0:lmat, 10), INTENT(IN) :: occupation
INTEGER, DIMENSION(0:lmat) :: maxn
INTEGER :: k, l
maxn = 0
DO l = 0, lmat
DO k = 1, 10
IF (occupation(l, k) /= 0._dp) maxn(l) = maxn(l) + 1
END DO
END DO
END FUNCTION get_maxn_occ
! **************************************************************************************************
!> \brief Calculate a density matrix using atomic orbitals.
!> \param pmat electron density matrix
!> \param wfn atomic wavefunctions
!> \param nbas number of basis functions
!> \param occ occupation numbers
!> \param maxl maximum angular momentum to consider
!> \param maxn maximum principal quantum number for each angular momentum
!> \par History
!> * 08.2008 created [Juerg Hutter]
! **************************************************************************************************
PURE SUBROUTINE atom_denmat(pmat, wfn, nbas, occ, maxl, maxn)
REAL(KIND=dp), DIMENSION(:, :, 0:), INTENT(INOUT) :: pmat
REAL(KIND=dp), DIMENSION(:, :, 0:), INTENT(IN) :: wfn
INTEGER, DIMENSION(0:lmat), INTENT(IN) :: nbas
REAL(KIND=dp), DIMENSION(0:, :), INTENT(IN) :: occ
INTEGER, INTENT(IN) :: maxl
INTEGER, DIMENSION(0:lmat), INTENT(IN) :: maxn
INTEGER :: i, j, k, l, n
pmat = 0._dp
n = SIZE(wfn, 2)
DO l = 0, maxl
DO i = 1, MIN(n, maxn(l))
DO k = 1, nbas(l)
DO j = 1, nbas(l)
pmat(j, k, l) = pmat(j, k, l) + occ(l, i)*wfn(j, i, l)*wfn(k, i, l)
END DO
END DO
END DO
END DO
END SUBROUTINE atom_denmat
! **************************************************************************************************
!> \brief Map the electron density on an atomic radial grid.
!> \param density computed electron density
!> \param pmat electron density matrix
!> \param basis atomic basis set
!> \param maxl maximum angular momentum to consider
!> \param typ type of the matrix to map:
!> RHO -- density matrix;
!> DER -- first derivatives of the electron density;
!> KIN -- kinetic energy density;
!> LAP -- second derivatives of the electron density.
!> \param rr abscissa on the radial grid (required for typ == 'KIN')
!> \par History
!> * 08.2008 created [Juerg Hutter]
! **************************************************************************************************
SUBROUTINE atom_density(density, pmat, basis, maxl, typ, rr)
REAL(KIND=dp), DIMENSION(:), INTENT(OUT) :: density
REAL(KIND=dp), DIMENSION(:, :, 0:), INTENT(IN) :: pmat
TYPE(atom_basis_type), INTENT(IN) :: basis
INTEGER, INTENT(IN) :: maxl
CHARACTER(LEN=*), OPTIONAL :: typ
REAL(KIND=dp), DIMENSION(:), INTENT(IN), OPTIONAL :: rr
CHARACTER(LEN=3) :: my_typ
INTEGER :: i, j, l, n
REAL(KIND=dp) :: ff
my_typ = "RHO"
IF (PRESENT(typ)) my_typ = typ(1:3)
IF (my_typ == "KIN") THEN
CPASSERT(PRESENT(rr))
END IF
density = 0._dp
DO l = 0, maxl
n = basis%nbas(l)
DO i = 1, n
DO j = i, n
ff = pmat(i, j, l)
IF (i /= j) ff = 2._dp*pmat(i, j, l)
IF (my_typ == "RHO") THEN
density(:) = density(:) + ff*basis%bf(:, i, l)*basis%bf(:, j, l)
ELSE IF (my_typ == "DER") THEN
density(:) = density(:) + ff*basis%dbf(:, i, l)*basis%bf(:, j, l) &
+ ff*basis%bf(:, i, l)*basis%dbf(:, j, l)
ELSE IF (my_typ == "KIN") THEN
density(:) = density(:) + 0.5_dp*ff*( &
basis%dbf(:, i, l)*basis%dbf(:, j, l) + &
REAL(l*(l + 1), dp)*basis%bf(:, i, l)*basis%bf(:, j, l)/rr(:))
ELSE IF (my_typ == "LAP") THEN
density(:) = density(:) + ff*basis%ddbf(:, i, l)*basis%bf(:, j, l) &
+ ff*basis%bf(:, i, l)*basis%ddbf(:, j, l) &
+ 2._dp*ff*basis%dbf(:, i, l)*basis%bf(:, j, l)/rr(:) &
+ 2._dp*ff*basis%bf(:, i, l)*basis%dbf(:, j, l)/rr(:)
ELSE
CPABORT("")
END IF
END DO
END DO
END DO
! this factor from the product of two spherical harmonics
density = density/fourpi
END SUBROUTINE atom_density
! **************************************************************************************************
!> \brief ZMP subroutine to write external restart file.
!> \param atom information about the atomic kind
!> \date 07.10.2013
!> \author D. Varsano [[email protected]]
!> \version 1.0
! **************************************************************************************************
SUBROUTINE atom_write_zmp_restart(atom)
TYPE(atom_type), INTENT(IN) :: atom
INTEGER :: extunit, i, k, l, n
extunit = get_unit_number()
CALL open_file(file_name=atom%zmp_restart_file, file_status="UNKNOWN", &
file_form="FORMATTED", file_action="WRITE", &
unit_number=extunit)
n = SIZE(atom%orbitals%wfn, 2)
WRITE (extunit, *) atom%basis%nbas
DO l = 0, atom%state%maxl_occ
DO i = 1, MIN(n, atom%state%maxn_occ(l))
DO k = 1, atom%basis%nbas(l)
WRITE (extunit, *) atom%orbitals%wfn(k, i, l)
END DO
END DO
END DO
CALL close_file(unit_number=extunit)
END SUBROUTINE atom_write_zmp_restart
! **************************************************************************************************
!> \brief ZMP subroutine to read external restart file.
!> \param atom information about the atomic kind
!> \param doguess flag that indicates that the restart file has not been read,
!> so the initial guess is required
!> \param iw output file unit
!> \date 07.10.2013
!> \author D. Varsano [[email protected]]
!> \version 1.0
! **************************************************************************************************
SUBROUTINE atom_read_zmp_restart(atom, doguess, iw)
TYPE(atom_type), INTENT(INOUT) :: atom
LOGICAL, INTENT(INOUT) :: doguess
INTEGER, INTENT(IN) :: iw
INTEGER :: er, extunit, i, k, l, maxl, n
INTEGER, DIMENSION(0:lmat) :: maxn, nbas
INQUIRE (file=atom%zmp_restart_file, exist=atom%doread)
IF (atom%doread) THEN
WRITE (iw, fmt="(' ZMP | Restart file found')")
extunit = get_unit_number()
CALL open_file(file_name=atom%zmp_restart_file, file_status="OLD", &
file_form="FORMATTED", file_action="READ", &
unit_number=extunit)
READ (extunit, *, IOSTAT=er) nbas
IF (er .NE. 0) THEN
WRITE (iw, fmt="(' ZMP | ERROR! Restart file unreadable')")
WRITE (iw, fmt="(' ZMP | ERROR! Starting ZMP calculation form initial atomic guess')")
doguess = .TRUE.
atom%doread = .FALSE.
ELSE IF (nbas(1) .NE. atom%basis%nbas(1)) THEN
WRITE (iw, fmt="(' ZMP | ERROR! Restart file contains a different basis set')")
WRITE (iw, fmt="(' ZMP | ERROR! Starting ZMP calculation form initial atomic guess')")
doguess = .TRUE.
atom%doread = .FALSE.
ELSE
nbas = atom%basis%nbas
maxl = atom%state%maxl_occ
maxn = atom%state%maxn_occ
n = SIZE(atom%orbitals%wfn, 2)
DO l = 0, atom%state%maxl_occ
DO i = 1, MIN(n, atom%state%maxn_occ(l))
DO k = 1, atom%basis%nbas(l)
READ (extunit, *) atom%orbitals%wfn(k, i, l)
END DO
END DO
END DO
doguess = .FALSE.
END IF
CALL close_file(unit_number=extunit)
ELSE
WRITE (iw, fmt="(' ZMP | WARNING! Restart file not found')")
WRITE (iw, fmt="(' ZMP | WARNING! Starting ZMP calculation form initial atomic guess')")
END IF
END SUBROUTINE atom_read_zmp_restart
! **************************************************************************************************
!> \brief ZMP subroutine to read external density from linear grid of density matrix.
!> \param density external density
!> \param atom information about the atomic kind
!> \param iw output file unit
!> \date 07.10.2013
!> \author D. Varsano [[email protected]]
!> \version 1.0
! **************************************************************************************************
SUBROUTINE atom_read_external_density(density, atom, iw)
REAL(KIND=dp), DIMENSION(:), INTENT(OUT) :: density
TYPE(atom_type), INTENT(INOUT) :: atom
INTEGER, INTENT(IN) :: iw
CHARACTER(LEN=default_string_length) :: filename
INTEGER :: extunit, ir, j, k, l, maxl_occ, maxnbas, &
nbas, nr
LOGICAL :: ldm
REAL(KIND=dp) :: rr
REAL(KIND=dp), ALLOCATABLE :: pmatread(:, :, :)
filename = atom%ext_file
ldm = atom%dm
extunit = get_unit_number()
CALL open_file(file_name=filename, file_status="OLD", &
file_form="FORMATTED", file_action="READ", &
unit_number=extunit)
IF (.NOT. ldm) THEN
READ (extunit, *) nr
IF (nr .NE. atom%basis%grid%nr) THEN
IF (iw > 0) WRITE (iw, fmt="(' ZMP | ERROR! External grid dimension ',I4,' differs from internal grid ',I4)") &
nr, atom%basis%grid%nr
IF (iw > 0) WRITE (iw, fmt="(' ZMP | ERROR! Stopping ZMP calculation')")
CPABORT("")
END IF
DO ir = 1, nr
READ (extunit, *) rr, density(ir)
IF (ABS(rr - atom%basis%grid%rad(ir)) .GT. atom%zmpgrid_tol) THEN
IF (iw > 0) WRITE (iw, fmt="(' ZMP | ERROR! Grid points do not coincide: ')")
IF (iw > 0) WRITE (iw, fmt='(" ZMP |",T20,"R_out[bohr]",T36,"R_in[bohr]",T61,"R_diff[bohr]")')
IF (iw > 0) WRITE (iw, fmt='(" ZMP |",T14,E24.15,T39,E24.15,T64,E24.15)') &
rr, atom%basis%grid%rad(ir), ABS(rr - atom%basis%grid%rad(ir))
CPABORT("")
END IF
END DO
CALL close_file(unit_number=extunit)
ELSE
READ (extunit, *) maxl_occ
maxnbas = MAXVAL(atom%basis%nbas)
ALLOCATE (pmatread(maxnbas, maxnbas, 0:maxl_occ))
pmatread = 0.0
DO l = 0, maxl_occ
nbas = atom%basis%nbas(l)
READ (extunit, *) ! Read empty line
DO k = 1, nbas
READ (extunit, *) (pmatread(j, k, l), j=1, k)
DO j = 1, k
pmatread(k, j, l) = pmatread(j, k, l)
END DO
END DO
END DO
CALL close_file(unit_number=extunit)
CALL atom_density(density, pmatread, atom%basis, maxl_occ, typ="RHO")
extunit = get_unit_number()
CALL open_file(file_name="rho_target.dat", file_status="UNKNOWN", &
file_form="FORMATTED", file_action="WRITE", unit_number=extunit)
IF (iw > 0) WRITE (iw, fmt="(' ZMP | Writing target density from density matrix')")
WRITE (extunit, fmt='("# Target density built from density matrix : ",A20)') filename
WRITE (extunit, fmt='("#",T10,"R[bohr]",T36,"Rho[au]")')
nr = atom%basis%grid%nr
DO ir = 1, nr
WRITE (extunit, fmt='(T1,E24.15,T26,E24.15)') &
atom%basis%grid%rad(ir), density(ir)
END DO
DEALLOCATE (pmatread)
CALL close_file(unit_number=extunit)
END IF
END SUBROUTINE atom_read_external_density
! **************************************************************************************************
!> \brief ZMP subroutine to read external v_xc in the atomic code.
!> \param vxc external exchange-correlation potential
!> \param atom information about the atomic kind
!> \param iw output file unit
!> \author D. Varsano [[email protected]]
! **************************************************************************************************
SUBROUTINE atom_read_external_vxc(vxc, atom, iw)
REAL(KIND=dp), DIMENSION(:), INTENT(OUT) :: vxc
TYPE(atom_type), INTENT(INOUT) :: atom
INTEGER, INTENT(IN) :: iw
CHARACTER(LEN=default_string_length) :: adum, filename
INTEGER :: extunit, ir, nr
REAL(KIND=dp) :: rr
filename = atom%ext_vxc_file
extunit = get_unit_number()
CALL open_file(file_name=filename, file_status="OLD", &
file_form="FORMATTED", file_action="READ", &
unit_number=extunit)
READ (extunit, *)
READ (extunit, *) adum, nr
IF (nr .NE. atom%basis%grid%nr) THEN
IF (iw > 0) WRITE (iw, fmt="(' ZMP | ERROR! External grid dimension ',I4,' differs from internal grid ',I4)") &
nr, atom%basis%grid%nr
IF (iw > 0) WRITE (iw, fmt="(' ZMP | ERROR! Stopping ZMP calculation')")
CPABORT("")
END IF
DO ir = 1, nr
READ (extunit, *) rr, vxc(ir)
IF (ABS(rr - atom%basis%grid%rad(ir)) .GT. atom%zmpvxcgrid_tol) THEN
IF (iw > 0) WRITE (iw, fmt="(' ZMP | ERROR! Grid points do not coincide: ')")
IF (iw > 0) WRITE (iw, fmt='(" ZMP |",T20,"R_out[bohr]",T36,"R_in[bohr]",T61,"R_diff[bohr]")')
IF (iw > 0) WRITE (iw, fmt='(" ZMP |",T14,E24.15,T39,E24.15,T64,E24.15)') &
rr, atom%basis%grid%rad(ir), ABS(rr - atom%basis%grid%rad(ir))
CPABORT("")
END IF
END DO
END SUBROUTINE atom_read_external_vxc
! **************************************************************************************************
!> \brief ...
!> \param charge ...
!> \param wfn ...
!> \param rcov ...
!> \param l ...
!> \param basis ...
! **************************************************************************************************
PURE SUBROUTINE atom_orbital_charge(charge, wfn, rcov, l, basis)
REAL(KIND=dp), INTENT(OUT) :: charge
REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: wfn
REAL(KIND=dp), INTENT(IN) :: rcov
INTEGER, INTENT(IN) :: l
TYPE(atom_basis_type), INTENT(IN) :: basis
INTEGER :: i, j, m, n
REAL(KIND=dp) :: ff
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: den
charge = 0._dp
m = SIZE(basis%bf, 1)
ALLOCATE (den(m))
n = basis%nbas(l)
den = 0._dp
DO i = 1, n
DO j = 1, n
ff = wfn(i)*wfn(j)
den(1:m) = den(1:m) + ff*basis%bf(1:m, i, l)*basis%bf(1:m, j, l)
END DO
END DO
DO i = 1, m
IF (basis%grid%rad(i) > rcov) den(i) = 0._dp
END DO
charge = SUM(den(1:m)*basis%grid%wr(1:m))
DEALLOCATE (den)
END SUBROUTINE atom_orbital_charge
! **************************************************************************************************
!> \brief ...
!> \param corden ...
!> \param potential ...
!> \param typ ...
!> \param rr ...
!> \par History
!> * 01.2017 rewritten [Juerg Hutter]
!> * 03.2010 extension of GTH pseudo-potential definition [Juerg Hutter]
! **************************************************************************************************
SUBROUTINE atom_core_density(corden, potential, typ, rr)
REAL(KIND=dp), DIMENSION(:), INTENT(INOUT) :: corden
TYPE(atom_potential_type), INTENT(IN) :: potential
CHARACTER(LEN=*), OPTIONAL :: typ
REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: rr
CHARACTER(LEN=3) :: my_typ
INTEGER :: i, j, m, n
LOGICAL :: reverse
REAL(KIND=dp) :: a, a2, cval, fb
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: fe, rc, rhoc, rval
my_typ = "RHO"
IF (PRESENT(typ)) my_typ = typ(1:3)
SELECT CASE (potential%ppot_type)
CASE (no_pseudo, ecp_pseudo)
! we do nothing
CASE (gth_pseudo)
IF (potential%gth_pot%nlcc) THEN
m = SIZE(corden)
ALLOCATE (fe(m), rc(m))
n = potential%gth_pot%nexp_nlcc
DO i = 1, n
a = potential%gth_pot%alpha_nlcc(i)
a2 = a*a
! note that all terms are computed with rc, not rr
rc(:) = rr(:)/a
fe(:) = EXP(-0.5_dp*rc(:)*rc(:))
DO j = 1, potential%gth_pot%nct_nlcc(i)
cval = potential%gth_pot%cval_nlcc(j, i)
IF (my_typ == "RHO") THEN
corden(:) = corden(:) + fe(:)*rc**(2*j - 2)*cval
ELSE IF (my_typ == "DER") THEN
corden(:) = corden(:) - fe(:)*rc**(2*j - 1)*cval/a
IF (j > 1) THEN
corden(:) = corden(:) + REAL(2*j - 2, dp)*fe(:)*rc**(2*j - 3)*cval/a
END IF
ELSE IF (my_typ == "LAP") THEN
fb = 2._dp*cval/a
corden(:) = corden(:) - fb*fe(:)/rr(:)*rc**(2*j - 1)
corden(:) = corden(:) + fe(:)*rc**(2*j)*cval/a2
IF (j > 1) THEN
corden(:) = corden(:) + fb*REAL(2*j - 2, dp)*fe(:)/rr(:)*rc**(2*j - 3)
corden(:) = corden(:) + REAL((2*j - 2)*(2*j - 3), dp)*fe(:)*rc**(2*j - 4)*cval/a2
corden(:) = corden(:) - REAL(2*j - 2, dp)*fe(:)*rc**(2*j - 2)*cval/a2
END IF
ELSE
CPABORT("")
END IF
END DO
END DO
DEALLOCATE (fe, rc)
END IF
CASE (upf_pseudo)
IF (potential%upf_pot%core_correction) THEN
m = SIZE(corden)
n = potential%upf_pot%mesh_size
reverse = .FALSE.
IF (rr(1) > rr(m)) reverse = .TRUE.
ALLOCATE (rhoc(m), rval(m))
IF (reverse) THEN
DO i = 1, m
rval(i) = rr(m - i + 1)
END DO
ELSE
rval(1:m) = rr(1:m)
END IF
IF (my_typ == "RHO") THEN
CALL spline3ders(potential%upf_pot%r(1:n), potential%upf_pot%rho_nlcc(1:n), rval(1:m), &
ynew=rhoc(1:m))
ELSE IF (my_typ == "DER") THEN
CALL spline3ders(potential%upf_pot%r(1:n), potential%upf_pot%rho_nlcc(1:n), rval(1:m), &
dynew=rhoc(1:m))
ELSE IF (my_typ == "LAP") THEN
CALL spline3ders(potential%upf_pot%r(1:n), potential%upf_pot%rho_nlcc(1:n), rval(1:m), &
d2ynew=rhoc(1:m))
ELSE
CPABORT("")
END IF
IF (reverse) THEN
DO i = 1, m
rval(i) = rr(m - i + 1)
corden(i) = corden(i) + rhoc(m - i + 1)
END DO
ELSE
corden(1:m) = corden(1:m) + rhoc(1:m)
END IF
DEALLOCATE (rhoc, rval)
END IF
CASE (sgp_pseudo)
IF (potential%sgp_pot%has_nlcc) THEN
CPABORT("not implemented")
END IF
CASE DEFAULT
CPABORT("Unknown PP type")
END SELECT
END SUBROUTINE atom_core_density
! **************************************************************************************************
!> \brief ...
!> \param locpot ...
!> \param gthpot ...
!> \param rr ...
! **************************************************************************************************
PURE SUBROUTINE atom_local_potential(locpot, gthpot, rr)
REAL(KIND=dp), DIMENSION(:), INTENT(INOUT) :: locpot
TYPE(atom_gthpot_type), INTENT(IN) :: gthpot
REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: rr
INTEGER :: i, j, m, n
REAL(KIND=dp) :: a
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: fe, rc
m = SIZE(locpot)
ALLOCATE (fe(m), rc(m))
rc(:) = rr(:)/gthpot%rc
DO i = 1, m
locpot(i) = -gthpot%zion*erf(rc(i)/SQRT(2._dp))/rr(i)
END DO
n = gthpot%ncl
fe(:) = EXP(-0.5_dp*rc(:)*rc(:))
DO i = 1, n
locpot(:) = locpot(:) + fe(:)*rc**(2*i - 2)*gthpot%cl(i)
END DO
IF (gthpot%lpotextended) THEN
DO j = 1, gthpot%nexp_lpot
a = gthpot%alpha_lpot(j)
rc(:) = rr(:)/a
fe(:) = EXP(-0.5_dp*rc(:)*rc(:))
n = gthpot%nct_lpot(j)
DO i = 1, n
locpot(:) = locpot(:) + fe(:)*rc**(2*i - 2)*gthpot%cval_lpot(i, j)
END DO
END DO
END IF
DEALLOCATE (fe, rc)
END SUBROUTINE atom_local_potential
! **************************************************************************************************
!> \brief ...
!> \param rmax ...
!> \param wfn ...
!> \param rcov ...
!> \param l ...
!> \param basis ...
! **************************************************************************************************
PURE SUBROUTINE atom_orbital_max(rmax, wfn, rcov, l, basis)
REAL(KIND=dp), INTENT(OUT) :: rmax
REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: wfn
REAL(KIND=dp), INTENT(IN) :: rcov
INTEGER, INTENT(IN) :: l
TYPE(atom_basis_type), INTENT(IN) :: basis
INTEGER :: i, m, n
REAL(KIND=dp) :: ff
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: dorb
m = SIZE(basis%bf, 1)
ALLOCATE (dorb(m))
n = basis%nbas(l)
dorb = 0._dp
DO i = 1, n
ff = wfn(i)
dorb(1:m) = dorb(1:m) + ff*basis%dbf(1:m, i, l)
END DO
rmax = -1._dp
DO i = 1, m - 1
IF (basis%grid%rad(i) < 2*rcov) THEN
IF (dorb(i)*dorb(i + 1) < 0._dp) THEN
rmax = MAX(rmax, basis%grid%rad(i))
END IF
END IF
END DO
DEALLOCATE (dorb)
END SUBROUTINE atom_orbital_max
! **************************************************************************************************
!> \brief ...
!> \param node ...
!> \param wfn ...
!> \param rcov ...
!> \param l ...
!> \param basis ...
! **************************************************************************************************
PURE SUBROUTINE atom_orbital_nodes(node, wfn, rcov, l, basis)
INTEGER, INTENT(OUT) :: node
REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: wfn
REAL(KIND=dp), INTENT(IN) :: rcov
INTEGER, INTENT(IN) :: l
TYPE(atom_basis_type), INTENT(IN) :: basis
INTEGER :: i, m, n
REAL(KIND=dp) :: ff
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: orb
node = 0
m = SIZE(basis%bf, 1)
ALLOCATE (orb(m))
n = basis%nbas(l)
orb = 0._dp
DO i = 1, n
ff = wfn(i)
orb(1:m) = orb(1:m) + ff*basis%bf(1:m, i, l)
END DO
DO i = 1, m - 1
IF (basis%grid%rad(i) < rcov) THEN
IF (orb(i)*orb(i + 1) < 0._dp) node = node + 1
END IF
END DO
DEALLOCATE (orb)
END SUBROUTINE atom_orbital_nodes
! **************************************************************************************************
!> \brief ...
!> \param value ...
!> \param wfn ...
!> \param basis ...
! **************************************************************************************************
PURE SUBROUTINE atom_wfnr0(value, wfn, basis)
REAL(KIND=dp), INTENT(OUT) :: value
REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: wfn
TYPE(atom_basis_type), INTENT(IN) :: basis
INTEGER :: i, m, n
value = 0._dp
m = MAXVAL(MINLOC(basis%grid%rad))
n = basis%nbas(0)
DO i = 1, n
value = value + wfn(i)*basis%bf(m, i, 0)
END DO
END SUBROUTINE atom_wfnr0
! **************************************************************************************************
!> \brief Solve the generalised eigenproblem for every angular momentum.
!> \param hmat Hamiltonian matrix
!> \param umat transformation matrix which reduces the overlap matrix to its unitary form
!> \param orb atomic wavefunctions
!> \param ener atomic orbital energies
!> \param nb number of contracted basis functions
!> \param nv number of linear-independent contracted basis functions
!> \param maxl maximum angular momentum to consider
!> \par History
!> * 08.2008 created [Juerg Hutter]
! **************************************************************************************************
SUBROUTINE atom_solve(hmat, umat, orb, ener, nb, nv, maxl)
REAL(KIND=dp), DIMENSION(:, :, 0:), INTENT(IN) :: hmat, umat
REAL(KIND=dp), DIMENSION(:, :, 0:), INTENT(INOUT) :: orb
REAL(KIND=dp), DIMENSION(:, 0:), INTENT(INOUT) :: ener
INTEGER, DIMENSION(0:), INTENT(IN) :: nb, nv
INTEGER, INTENT(IN) :: maxl
INTEGER :: info, l, lwork, m, n
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: w, work
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: a
CPASSERT(ALL(nb >= nv))
orb = 0._dp
DO l = 0, maxl
n = nb(l)
m = nv(l)
IF (n > 0 .AND. m > 0) THEN
lwork = 10*m
ALLOCATE (a(n, n), w(n), work(lwork))
a(1:m, 1:m) = MATMUL(TRANSPOSE(umat(1:n, 1:m, l)), MATMUL(hmat(1:n, 1:n, l), umat(1:n, 1:m, l)))
CALL lapack_ssyev("V", "U", m, a(1:m, 1:m), m, w(1:m), work, lwork, info)
a(1:n, 1:m) = MATMUL(umat(1:n, 1:m, l), a(1:m, 1:m))
m = MIN(m, SIZE(orb, 2))
orb(1:n, 1:m, l) = a(1:n, 1:m)
ener(1:m, l) = w(1:m)
DEALLOCATE (a, w, work)
END IF
END DO
END SUBROUTINE atom_solve
! **************************************************************************************************
!> \brief THIS FUNCTION IS NO LONGER IN USE.
!> \param fun ...
!> \param deps ...
!> \return ...
! **************************************************************************************************
FUNCTION prune_grid(fun, deps) RESULT(nc)
REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: fun
REAL(KIND=dp), INTENT(IN), OPTIONAL :: deps
INTEGER :: nc
INTEGER :: i, nr
REAL(KIND=dp) :: meps
meps = 1.e-14_dp
IF (PRESENT(deps)) meps = deps
nr = SIZE(fun)
nc = 0
DO i = nr, 1, -1
IF (ABS(fun(i)) > meps) THEN
nc = i
EXIT
END IF