forked from cp2k/cp2k
-
Notifications
You must be signed in to change notification settings - Fork 0
/
atom_xc.F
1338 lines (1150 loc) · 59.3 KB
/
atom_xc.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 routines that build the integrals of the Vxc potential calculated
!> for the atomic code
! **************************************************************************************************
MODULE atom_xc
USE atom_types, ONLY: atom_basis_type,&
atom_type,&
gth_pseudo,&
opmat_type,&
sgp_pseudo,&
upf_pseudo
USE atom_utils, ONLY: atom_core_density,&
atom_density,&
coulomb_potential_numeric,&
integrate_grid,&
numpot_matrix
USE cp_files, ONLY: close_file,&
get_unit_number,&
open_file
USE input_constants, ONLY: xc_none
USE input_section_types, ONLY: section_vals_get_subs_vals,&
section_vals_type,&
section_vals_val_get
USE kinds, ONLY: dp
USE mathconstants, ONLY: fourpi,&
pi
USE xc_atom, ONLY: xc_rho_set_atom_update
USE xc_derivative_desc, ONLY: &
deriv_norm_drho, deriv_norm_drhoa, deriv_norm_drhob, deriv_rho, deriv_rhoa, deriv_rhob, &
deriv_tau, deriv_tau_a, deriv_tau_b
USE xc_derivative_set_types, ONLY: xc_derivative_set_type,&
xc_dset_create,&
xc_dset_get_derivative,&
xc_dset_release,&
xc_dset_zero_all
USE xc_derivative_types, ONLY: xc_derivative_get,&
xc_derivative_type
USE xc_derivatives, ONLY: xc_functionals_eval,&
xc_functionals_get_needs
USE xc_rho_cflags_types, ONLY: xc_rho_cflags_type
USE xc_rho_set_types, ONLY: xc_rho_set_create,&
xc_rho_set_release,&
xc_rho_set_type
#include "./base/base_uses.f90"
IMPLICIT NONE
PRIVATE
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'atom_xc'
! ZMP added coulomb_potential_numeric
! ZMP added public subroutines calculate_atom_zmp, calculate_atom_ext_vxc
PUBLIC :: calculate_atom_vxc_lda, calculate_atom_vxc_lsd, &
calculate_atom_zmp, calculate_atom_ext_vxc
PUBLIC :: atom_vxc_lda, atom_vxc_lsd, atom_dpot_lda
CONTAINS
! **************************************************************************************************
!> \brief ZMP subroutine for the calculation of the effective constraint
!> potential in the atomic code.
!> \param ext_density external electron density
!> \param atom information about the atomic kind
!> \param lprint save exchange-correlation potential to the file 'linear_potentials.dat'
!> \param xcmat exchange-correlation potential matrix
!> \author D. Varsano [[email protected]]
! **************************************************************************************************
SUBROUTINE calculate_atom_zmp(ext_density, atom, lprint, xcmat)
REAL(KIND=dp), DIMENSION(:), INTENT(INOUT) :: ext_density
TYPE(atom_type), INTENT(INOUT) :: atom
LOGICAL :: lprint
TYPE(opmat_type), OPTIONAL, POINTER :: xcmat
CHARACTER(LEN=*), PARAMETER :: routineN = 'calculate_atom_zmp'
INTEGER :: extunit, handle, ir, nr, z
REAL(KIND=dp) :: int1, int2
REAL(KIND=dp), DIMENSION(:), POINTER :: deltarho, rho_dum, vxc, vxc1, vxc2
REAL(KIND=dp), DIMENSION(:, :), POINTER :: rho
CALL timeset(routineN, handle)
nr = atom%basis%grid%nr
z = atom%z
ALLOCATE (rho(nr, 1), vxc(nr), vxc1(nr), vxc2(nr), rho_dum(nr), deltarho(nr))
CALL atom_density(rho(:, 1), atom%orbitals%pmat, atom%basis, atom%state%maxl_occ, typ="RHO")
!Vxc1
int1 = integrate_grid(ext_density, atom%basis%grid)
int1 = fourpi*int1
CALL coulomb_potential_numeric(vxc1, rho(:, 1), atom%basis%grid)
vxc1 = -vxc1/z
!Vxc2
rho_dum = rho(:, 1)*int1/z
deltarho = rho_dum - ext_density
int2 = integrate_grid(deltarho, atom%basis%grid)
CALL coulomb_potential_numeric(vxc2, deltarho, atom%basis%grid)
vxc2 = vxc2*atom%lambda
!Vxc
vxc = vxc1 + vxc2
atom%energy%exc = fourpi*integrate_grid(vxc, rho(:, 1), atom%basis%grid)
atom%rho_diff_integral = fourpi*int2
IF (lprint) THEN
extunit = get_unit_number()
CALL open_file(file_name="linear_potentials.dat", file_status="UNKNOWN", &
file_form="FORMATTED", file_action="WRITE", &
unit_number=extunit)
WRITE (extunit, fmt='("#",T10,"R[bohr]",T36,"Rho[au]",T61,"DRho[au]",T86,"V_XC[au]",T111,"V_XC1[au]",T136,"V_XC2[au]")')
DO ir = 1, nr
WRITE (extunit, fmt='(T1,E24.15,T26,E24.15,T51,E24.15,T76,E24.15,T101,E24.15,T126,E24.15)') &
atom%basis%grid%rad(ir), rho(ir, 1), deltarho(ir), vxc(ir), vxc1(ir), vxc2(ir)
END DO
CALL close_file(unit_number=extunit)
END IF
IF (PRESENT(xcmat)) CALL numpot_matrix(xcmat%op, vxc, atom%basis, 0)
DEALLOCATE (rho, vxc, vxc1, vxc2, rho_dum, deltarho)
CALL timestop(handle)
END SUBROUTINE calculate_atom_zmp
! **************************************************************************************************
!> \brief ZMP subroutine for the scf external density from the external v_xc.
!> \param vxc exchange-correlation potential
!> \param atom information about the atomic kind
!> \param lprint save exchange-correlation potential to the file 'linear_potentials.dat'
!> \param xcmat exchange-correlation potential matrix
!> \author D. Varsano [[email protected]]
! **************************************************************************************************
SUBROUTINE calculate_atom_ext_vxc(vxc, atom, lprint, xcmat)
REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: vxc
TYPE(atom_type), INTENT(INOUT) :: atom
LOGICAL, INTENT(in) :: lprint
TYPE(opmat_type), INTENT(inout), OPTIONAL, POINTER :: xcmat
CHARACTER(LEN=*), PARAMETER :: routineN = 'calculate_atom_ext_vxc'
INTEGER :: extunit, handle, ir, nr
REAL(KIND=dp), DIMENSION(:, :), POINTER :: rho
CALL timeset(routineN, handle)
nr = atom%basis%grid%nr
ALLOCATE (rho(nr, 1))
CALL atom_density(rho(:, 1), atom%orbitals%pmat, atom%basis, atom%state%maxl_occ, typ="RHO")
IF (lprint) THEN
extunit = get_unit_number()
CALL open_file(file_name="linear_potentials.dat", file_status="UNKNOWN", &
file_form="FORMATTED", file_action="WRITE", &
unit_number=extunit)
WRITE (extunit, fmt='("#",T10,"R[bohr]",T36,"Rho[au]",T61,"V_XC[au]")')
DO ir = 1, nr
WRITE (extunit, fmt='(T1,E24.15,T26,E24.15,T51,E24.15)') &
atom%basis%grid%rad(ir), rho(ir, 1), vxc(ir)
END DO
CALL close_file(unit_number=extunit)
END IF
atom%energy%exc = fourpi*integrate_grid(vxc, rho(:, 1), atom%basis%grid)
IF (PRESENT(xcmat)) CALL numpot_matrix(xcmat%op, vxc, atom%basis, 0)
DEALLOCATE (rho)
CALL timestop(handle)
END SUBROUTINE calculate_atom_ext_vxc
! **************************************************************************************************
!> \brief Calculate atomic exchange-correlation potential in spin-restricted case.
!> \param xcmat exchange-correlation potential expressed in the atomic basis set
!> \param atom information about the atomic kind
!> \param xc_section XC input section
!> \par History
!> * 03.2010 extension of GTH pseudo-potential definition [Juerg Hutter]
!> * 02.2010 renamed to calculate_atom_vxc_lda() [Juerg Hutter]
!> * 08.2008 created [Juerg Hutter]
! **************************************************************************************************
SUBROUTINE calculate_atom_vxc_lda(xcmat, atom, xc_section)
TYPE(opmat_type), POINTER :: xcmat
TYPE(atom_type), INTENT(INOUT) :: atom
TYPE(section_vals_type), POINTER :: xc_section
CHARACTER(LEN=*), PARAMETER :: routineN = 'calculate_atom_vxc_lda'
INTEGER :: deriv_order, handle, i, l, myfun, n1, &
n2, n3, nr, nspins, unit_nr
INTEGER, DIMENSION(2, 3) :: bounds
LOGICAL :: lsd, nlcc
REAL(KIND=dp) :: density_cut, gradient_cut, tau_cut
REAL(KIND=dp), DIMENSION(:), POINTER :: exc, vxc
REAL(KIND=dp), DIMENSION(:, :), POINTER :: drho, lap, rho, tau
REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: taumat, xcpot
TYPE(section_vals_type), POINTER :: xc_fun_section
TYPE(xc_derivative_set_type) :: deriv_set
TYPE(xc_derivative_type), POINTER :: deriv
TYPE(xc_rho_cflags_type) :: needs
TYPE(xc_rho_set_type) :: rho_set
CALL timeset(routineN, handle)
nlcc = .FALSE.
IF (atom%potential%ppot_type == gth_pseudo) THEN
nlcc = atom%potential%gth_pot%nlcc
ELSE IF (atom%potential%ppot_type == upf_pseudo) THEN
nlcc = atom%potential%upf_pot%core_correction
ELSE IF (atom%potential%ppot_type == sgp_pseudo) THEN
nlcc = atom%potential%sgp_pot%has_nlcc
END IF
IF (ASSOCIATED(xc_section)) THEN
xc_fun_section => section_vals_get_subs_vals(xc_section, "XC_FUNCTIONAL")
CALL section_vals_val_get(xc_fun_section, "_SECTION_PARAMETERS_", i_val=myfun)
IF (myfun == xc_none) THEN
atom%energy%exc = 0._dp
ELSE
CALL section_vals_val_get(xc_section, "DENSITY_CUTOFF", r_val=density_cut)
CALL section_vals_val_get(xc_section, "GRADIENT_CUTOFF", r_val=gradient_cut)
CALL section_vals_val_get(xc_section, "TAU_CUTOFF", r_val=tau_cut)
lsd = .FALSE.
nspins = 1
needs = xc_functionals_get_needs(xc_fun_section, lsd=lsd, calc_potential=.FALSE.)
! Prepare the structures needed to calculate and store the xc derivatives
! Array dimension: here anly one dimensional arrays are used,
! i.e. only the first column of deriv_data is read.
! The other to dimensions are set to size equal 1
nr = atom%basis%grid%nr
bounds(1:2, 1:3) = 1
bounds(2, 1) = nr
! create a place where to put the derivatives
CALL xc_dset_create(deriv_set, local_bounds=bounds)
! create the place where to store the argument for the functionals
CALL xc_rho_set_create(rho_set, bounds, rho_cutoff=density_cut, &
drho_cutoff=gradient_cut, tau_cutoff=tau_cut)
! allocate the required 3d arrays where to store rho and drho
CALL xc_rho_set_atom_update(rho_set, needs, nspins, bounds)
NULLIFY (rho, drho, tau)
IF (needs%rho) THEN
ALLOCATE (rho(nr, 1))
CALL atom_density(rho(:, 1), atom%orbitals%pmat, atom%basis, atom%state%maxl_occ, typ="RHO")
IF (nlcc) THEN
CALL atom_core_density(rho(:, 1), atom%potential, typ="RHO", rr=atom%basis%grid%rad)
END IF
END IF
IF (needs%norm_drho) THEN
ALLOCATE (drho(nr, 1))
CALL atom_density(drho(:, 1), atom%orbitals%pmat, atom%basis, atom%state%maxl_occ, typ="DER")
IF (nlcc) THEN
CALL atom_core_density(drho(:, 1), atom%potential, typ="DER", rr=atom%basis%grid%rad)
END IF
END IF
IF (needs%tau) THEN
ALLOCATE (tau(nr, 1))
CALL atom_density(tau(:, 1), atom%orbitals%pmat, atom%basis, atom%state%maxl_occ, &
typ="KIN", rr=atom%basis%grid%rad2)
END IF
IF (needs%laplace_rho) THEN
ALLOCATE (lap(nr, 1))
CALL atom_density(lap(:, 1), atom%orbitals%pmat, atom%basis, atom%state%maxl_occ, &
typ="LAP", rr=atom%basis%grid%rad2)
IF (nlcc) THEN
CALL atom_core_density(lap(:, 1), atom%potential, typ="LAP", rr=atom%basis%grid%rad)
END IF
END IF
CALL fill_rho_set(rho_set, nspins, needs, rho, drho, tau, lap, nr)
CALL xc_dset_zero_all(deriv_set)
deriv_order = 1
CALL xc_functionals_eval(xc_fun_section, lsd=lsd, rho_set=rho_set, deriv_set=deriv_set, &
deriv_order=deriv_order)
! Integration to get the matrix elements and energy
deriv => xc_dset_get_derivative(deriv_set, [INTEGER::], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
atom%energy%exc = fourpi*integrate_grid(xcpot(:, 1, 1), atom%basis%grid)
! dump grid density and xcpot (xc energy?)
IF (.FALSE.) THEN
CALL open_file(unit_number=unit_nr, file_status="UNKNOWN", file_action="WRITE", file_name="atom.dat")
DO i = 1, SIZE(xcpot(:, 1, 1))
WRITE (unit_nr, *) atom%basis%grid%rad(i), rho(i, 1), xcpot(i, 1, 1)
END DO
CALL close_file(unit_nr)
END IF
IF (needs%rho) THEN
deriv => xc_dset_get_derivative(deriv_set, [deriv_rho], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
CALL numpot_matrix(xcmat%op, xcpot(:, 1, 1), atom%basis, 0)
IF (.FALSE.) THEN
CALL open_file(unit_number=unit_nr, file_status="UNKNOWN", file_action="WRITE", file_name="atompot.dat")
DO i = 1, SIZE(xcpot(:, 1, 1))
WRITE (unit_nr, *) atom%basis%grid%rad(i), rho(i, 1), xcpot(i, 1, 1)
END DO
CALL close_file(unit_nr)
END IF
DEALLOCATE (rho)
END IF
IF (needs%norm_drho) THEN
deriv => xc_dset_get_derivative(deriv_set, [deriv_norm_drho], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
CALL numpot_matrix(xcmat%op, xcpot(:, 1, 1), atom%basis, 1)
IF (.FALSE.) THEN
CALL open_file(unit_number=unit_nr, file_status="UNKNOWN", file_action="WRITE", file_name="atomdpot.dat")
DO i = 1, SIZE(xcpot(:, 1, 1))
WRITE (unit_nr, *) atom%basis%grid%rad(i), drho(i, 1), xcpot(i, 1, 1)
END DO
CALL close_file(unit_nr)
END IF
DEALLOCATE (drho)
END IF
IF (needs%tau) THEN
deriv => xc_dset_get_derivative(deriv_set, [deriv_tau], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
n1 = SIZE(xcmat%op, 1)
n2 = SIZE(xcmat%op, 2)
n3 = SIZE(xcmat%op, 3)
ALLOCATE (taumat(n1, n2, 0:n3 - 1))
taumat = 0._dp
xcpot(:, 1, 1) = 0.5_dp*xcpot(:, 1, 1)
CALL numpot_matrix(xcmat%op, xcpot(:, 1, 1), atom%basis, 2)
xcpot(:, 1, 1) = xcpot(:, 1, 1)/atom%basis%grid%rad2(:)
CALL numpot_matrix(taumat, xcpot(:, 1, 1), atom%basis, 0)
DO l = 0, 3
xcmat%op(:, :, l) = xcmat%op(:, :, l) + REAL(l*(l + 1), dp)*taumat(:, :, l)
END DO
DEALLOCATE (tau)
DEALLOCATE (taumat)
END IF
! assume lap is not really needed
IF (needs%laplace_rho) THEN
DEALLOCATE (lap)
END IF
! Release the xc structure used to store the xc derivatives
CALL xc_dset_release(deriv_set)
CALL xc_rho_set_release(rho_set)
END IF !xc_none
ELSE
! we don't have an xc_section, use a default setup
nr = atom%basis%grid%nr
ALLOCATE (rho(nr, 1), exc(nr), vxc(nr))
CALL atom_density(rho(:, 1), atom%orbitals%pmat, atom%basis, atom%state%maxl_occ, typ="RHO")
IF (nlcc) THEN
CALL atom_core_density(rho(:, 1), atom%potential, typ="RHO", rr=atom%basis%grid%rad)
END IF
CALL lda_pade(rho(:, 1), exc, vxc)
atom%energy%exc = fourpi*integrate_grid(exc, atom%basis%grid)
CALL numpot_matrix(xcmat%op, vxc, atom%basis, 0)
DEALLOCATE (rho, exc, vxc)
END IF
CALL timestop(handle)
END SUBROUTINE calculate_atom_vxc_lda
! **************************************************************************************************
!> \brief Calculate atomic exchange-correlation potential in spin-polarized case.
!> \param xcmata spin-alpha exchange-correlation potential matrix
!> \param xcmatb spin-beta exchange-correlation potential matrix
!> \param atom information about the atomic kind
!> \param xc_section XC input section
!> \par History
!> * 02.2010 created [Juerg Hutter]
! **************************************************************************************************
SUBROUTINE calculate_atom_vxc_lsd(xcmata, xcmatb, atom, xc_section)
TYPE(opmat_type), POINTER :: xcmata, xcmatb
TYPE(atom_type), INTENT(INOUT) :: atom
TYPE(section_vals_type), POINTER :: xc_section
CHARACTER(LEN=*), PARAMETER :: routineN = 'calculate_atom_vxc_lsd'
INTEGER :: deriv_order, handle, l, myfun, n1, n2, &
n3, nr, nspins
INTEGER, DIMENSION(2, 3) :: bounds
LOGICAL :: lsd, nlcc
REAL(KIND=dp) :: density_cut, gradient_cut, tau_cut
REAL(KIND=dp), DIMENSION(:), POINTER :: exc, vxca, vxcb, xfun
REAL(KIND=dp), DIMENSION(:, :), POINTER :: drho, lap, rho, tau
REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: taumat, xcpot
TYPE(section_vals_type), POINTER :: xc_fun_section
TYPE(xc_derivative_set_type) :: deriv_set
TYPE(xc_derivative_type), POINTER :: deriv
TYPE(xc_rho_cflags_type) :: needs
TYPE(xc_rho_set_type) :: rho_set
CALL timeset(routineN, handle)
nlcc = .FALSE.
IF (atom%potential%ppot_type == gth_pseudo) THEN
nlcc = atom%potential%gth_pot%nlcc
ELSE IF (atom%potential%ppot_type == upf_pseudo) THEN
nlcc = atom%potential%upf_pot%core_correction
CPASSERT(.NOT. nlcc)
ELSE IF (atom%potential%ppot_type == sgp_pseudo) THEN
nlcc = atom%potential%sgp_pot%has_nlcc
CPASSERT(.NOT. nlcc)
END IF
IF (ASSOCIATED(xc_section)) THEN
IF (nlcc) THEN
nr = atom%basis%grid%nr
ALLOCATE (xfun(nr))
END IF
xc_fun_section => section_vals_get_subs_vals(xc_section, "XC_FUNCTIONAL")
CALL section_vals_val_get(xc_fun_section, "_SECTION_PARAMETERS_", i_val=myfun)
IF (myfun == xc_none) THEN
atom%energy%exc = 0._dp
ELSE
CALL section_vals_val_get(xc_section, "DENSITY_CUTOFF", r_val=density_cut)
CALL section_vals_val_get(xc_section, "GRADIENT_CUTOFF", r_val=gradient_cut)
CALL section_vals_val_get(xc_section, "TAU_CUTOFF", r_val=tau_cut)
lsd = .TRUE.
nspins = 2
needs = xc_functionals_get_needs(xc_fun_section, lsd=lsd, calc_potential=.FALSE.)
! Prepare the structures needed to calculate and store the xc derivatives
! Array dimension: here anly one dimensional arrays are used,
! i.e. only the first column of deriv_data is read.
! The other to dimensions are set to size equal 1
nr = atom%basis%grid%nr
bounds(1:2, 1:3) = 1
bounds(2, 1) = nr
! create a place where to put the derivatives
CALL xc_dset_create(deriv_set, local_bounds=bounds)
! create the place where to store the argument for the functionals
CALL xc_rho_set_create(rho_set, bounds, rho_cutoff=density_cut, &
drho_cutoff=gradient_cut, tau_cutoff=tau_cut)
! allocate the required 3d arrays where to store rho and drho
CALL xc_rho_set_atom_update(rho_set, needs, nspins, bounds)
NULLIFY (rho, drho, tau)
IF (needs%rho_spin) THEN
ALLOCATE (rho(nr, 2))
CALL atom_density(rho(:, 1), atom%orbitals%pmata, atom%basis, atom%state%maxl_occ, typ="RHO")
CALL atom_density(rho(:, 2), atom%orbitals%pmatb, atom%basis, atom%state%maxl_occ, typ="RHO")
IF (nlcc) THEN
xfun = 0.0_dp
CALL atom_core_density(xfun(:), atom%potential, typ="RHO", rr=atom%basis%grid%rad)
rho(:, 1) = rho(:, 1) + 0.5_dp*xfun(:)
rho(:, 2) = rho(:, 2) + 0.5_dp*xfun(:)
END IF
END IF
IF (needs%norm_drho_spin) THEN
ALLOCATE (drho(nr, 2))
CALL atom_density(drho(:, 1), atom%orbitals%pmata, atom%basis, atom%state%maxl_occ, typ="DER")
CALL atom_density(drho(:, 2), atom%orbitals%pmatb, atom%basis, atom%state%maxl_occ, typ="DER")
IF (nlcc) THEN
xfun = 0.0_dp
CALL atom_core_density(xfun(:), atom%potential, typ="DER", rr=atom%basis%grid%rad)
drho(:, 1) = drho(:, 1) + 0.5_dp*xfun(:)
drho(:, 2) = drho(:, 2) + 0.5_dp*xfun(:)
END IF
END IF
IF (needs%tau_spin) THEN
ALLOCATE (tau(nr, 2))
CALL atom_density(tau(:, 1), atom%orbitals%pmata, atom%basis, atom%state%maxl_occ, &
typ="KIN", rr=atom%basis%grid%rad2)
CALL atom_density(tau(:, 2), atom%orbitals%pmatb, atom%basis, atom%state%maxl_occ, &
typ="KIN", rr=atom%basis%grid%rad2)
END IF
IF (needs%laplace_rho_spin) THEN
ALLOCATE (lap(nr, 2))
CALL atom_density(lap(:, 1), atom%orbitals%pmata, atom%basis, atom%state%maxl_occ, &
typ="LAP", rr=atom%basis%grid%rad2)
CALL atom_density(lap(:, 2), atom%orbitals%pmatb, atom%basis, atom%state%maxl_occ, &
typ="LAP", rr=atom%basis%grid%rad2)
IF (nlcc) THEN
xfun = 0.0_dp
CALL atom_core_density(xfun(:), atom%potential, typ="LAP", rr=atom%basis%grid%rad)
lap(:, 1) = lap(:, 1) + 0.5_dp*xfun(:)
lap(:, 2) = lap(:, 2) + 0.5_dp*xfun(:)
END IF
END IF
CALL fill_rho_set(rho_set, nspins, needs, rho, drho, tau, lap, nr)
CALL xc_dset_zero_all(deriv_set)
deriv_order = 1
CALL xc_functionals_eval(xc_fun_section, lsd=lsd, rho_set=rho_set, deriv_set=deriv_set, &
deriv_order=deriv_order)
! Integration to get the matrix elements and energy
deriv => xc_dset_get_derivative(deriv_set, [INTEGER::], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
atom%energy%exc = fourpi*integrate_grid(xcpot(:, 1, 1), atom%basis%grid)
IF (needs%rho_spin) THEN
deriv => xc_dset_get_derivative(deriv_set, [deriv_rhoa], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
CALL numpot_matrix(xcmata%op, xcpot(:, 1, 1), atom%basis, 0)
deriv => xc_dset_get_derivative(deriv_set, [deriv_rhob], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
CALL numpot_matrix(xcmatb%op, xcpot(:, 1, 1), atom%basis, 0)
DEALLOCATE (rho)
END IF
IF (needs%norm_drho_spin) THEN
! drhoa
NULLIFY (deriv)
deriv => xc_dset_get_derivative(deriv_set, [deriv_norm_drhoa], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
CALL numpot_matrix(xcmata%op, xcpot(:, 1, 1), atom%basis, 1)
! drhob
NULLIFY (deriv)
deriv => xc_dset_get_derivative(deriv_set, [deriv_norm_drhob], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
CALL numpot_matrix(xcmatb%op, xcpot(:, 1, 1), atom%basis, 1)
! Cross Terms
NULLIFY (deriv)
deriv => xc_dset_get_derivative(deriv_set, [deriv_norm_drho])
IF (ASSOCIATED(deriv)) THEN
CALL xc_derivative_get(deriv, deriv_data=xcpot)
CALL numpot_matrix(xcmata%op, xcpot(:, 1, 1), atom%basis, 1)
CALL numpot_matrix(xcmatb%op, xcpot(:, 1, 1), atom%basis, 1)
END IF
DEALLOCATE (drho)
END IF
IF (needs%tau_spin) THEN
n1 = SIZE(xcmata%op, 1)
n2 = SIZE(xcmata%op, 2)
n3 = SIZE(xcmata%op, 3)
ALLOCATE (taumat(n1, n2, 0:n3 - 1))
deriv => xc_dset_get_derivative(deriv_set, [deriv_tau_a], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
taumat = 0._dp
xcpot(:, 1, 1) = 0.5_dp*xcpot(:, 1, 1)
CALL numpot_matrix(xcmata%op, xcpot(:, 1, 1), atom%basis, 2)
xcpot(:, 1, 1) = xcpot(:, 1, 1)/atom%basis%grid%rad2(:)
CALL numpot_matrix(taumat, xcpot(:, 1, 1), atom%basis, 0)
DO l = 0, 3
xcmata%op(:, :, l) = xcmata%op(:, :, l) + REAL(l*(l + 1), dp)*taumat(:, :, l)
END DO
deriv => xc_dset_get_derivative(deriv_set, [deriv_tau_b], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
taumat = 0._dp
xcpot(:, 1, 1) = 0.5_dp*xcpot(:, 1, 1)
CALL numpot_matrix(xcmatb%op, xcpot(:, 1, 1), atom%basis, 2)
xcpot(:, 1, 1) = xcpot(:, 1, 1)/atom%basis%grid%rad2(:)
CALL numpot_matrix(taumat, xcpot(:, 1, 1), atom%basis, 0)
DO l = 0, 3
xcmatb%op(:, :, l) = xcmatb%op(:, :, l) + REAL(l*(l + 1), dp)*taumat(:, :, l)
END DO
DEALLOCATE (tau)
DEALLOCATE (taumat)
END IF
! assume lap is not really needed
IF (needs%laplace_rho_spin) THEN
DEALLOCATE (lap)
END IF
! Release the xc structure used to store the xc derivatives
CALL xc_dset_release(deriv_set)
CALL xc_rho_set_release(rho_set)
END IF !xc_none
IF (nlcc) DEALLOCATE (xfun)
ELSE
! we don't have an xc_section, use a default setup
nr = atom%basis%grid%nr
ALLOCATE (rho(nr, 2), exc(nr), vxca(nr), vxcb(nr))
IF (nlcc) ALLOCATE (xfun(nr))
CALL atom_density(rho(:, 1), atom%orbitals%pmata, atom%basis, atom%state%maxl_occ, typ="RHO")
CALL atom_density(rho(:, 2), atom%orbitals%pmatb, atom%basis, atom%state%maxl_occ, typ="RHO")
IF (nlcc) THEN
xfun(:) = 0.0_dp
CALL atom_core_density(xfun(:), atom%potential, typ="RHO", rr=atom%basis%grid%rad)
rho(:, 1) = rho(:, 1) + 0.5_dp*xfun(:)
rho(:, 2) = rho(:, 2) + 0.5_dp*xfun(:)
END IF
CALL lsd_pade(rho(:, 1), rho(:, 2), exc, vxca, vxcb)
atom%energy%exc = fourpi*integrate_grid(exc, atom%basis%grid)
CALL numpot_matrix(xcmata%op, vxca, atom%basis, 0)
CALL numpot_matrix(xcmatb%op, vxcb, atom%basis, 0)
DEALLOCATE (rho, exc, vxca, vxcb)
IF (nlcc) DEALLOCATE (xfun)
END IF
CALL timestop(handle)
END SUBROUTINE calculate_atom_vxc_lsd
! **************************************************************************************************
!> \brief Alternative subroutine that calculates atomic exchange-correlation potential
!> in spin-restricted case.
!> \param basis atomic basis set
!> \param pmat density matrix
!> \param maxl maximum angular momentum
!> \param xc_section XC input section
!> \param fexc exchange-correlation energy
!> \param xcmat exchange-correlation potential matrix
!> \par History
!> * 07.2016 ADMM analysis [Juerg Hutter]
! **************************************************************************************************
SUBROUTINE atom_vxc_lda(basis, pmat, maxl, xc_section, fexc, xcmat)
TYPE(atom_basis_type), POINTER :: basis
REAL(KIND=dp), DIMENSION(:, :, :), INTENT(IN) :: pmat
INTEGER, INTENT(IN) :: maxl
TYPE(section_vals_type), POINTER :: xc_section
REAL(KIND=dp), INTENT(OUT) :: fexc
REAL(KIND=dp), DIMENSION(:, :, :), INTENT(INOUT) :: xcmat
CHARACTER(LEN=*), PARAMETER :: routineN = 'atom_vxc_lda'
INTEGER :: deriv_order, handle, l, myfun, n1, n2, &
n3, nr, nspins
INTEGER, DIMENSION(2, 3) :: bounds
LOGICAL :: lsd
REAL(KIND=dp) :: density_cut, gradient_cut, tau_cut
REAL(KIND=dp), DIMENSION(:, :), POINTER :: drho, lap, rho, tau
REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: taumat, xcpot
TYPE(section_vals_type), POINTER :: xc_fun_section
TYPE(xc_derivative_set_type) :: deriv_set
TYPE(xc_derivative_type), POINTER :: deriv
TYPE(xc_rho_cflags_type) :: needs
TYPE(xc_rho_set_type) :: rho_set
CALL timeset(routineN, handle)
CPASSERT(ASSOCIATED(xc_section))
xc_fun_section => section_vals_get_subs_vals(xc_section, "XC_FUNCTIONAL")
CALL section_vals_val_get(xc_fun_section, "_SECTION_PARAMETERS_", i_val=myfun)
IF (myfun == xc_none) THEN
fexc = 0._dp
ELSE
CALL section_vals_val_get(xc_section, "DENSITY_CUTOFF", r_val=density_cut)
CALL section_vals_val_get(xc_section, "GRADIENT_CUTOFF", r_val=gradient_cut)
CALL section_vals_val_get(xc_section, "TAU_CUTOFF", r_val=tau_cut)
lsd = .FALSE.
nspins = 1
needs = xc_functionals_get_needs(xc_fun_section, lsd=lsd, calc_potential=.FALSE.)
! Prepare the structures needed to calculate and store the xc derivatives
! Array dimension: here anly one dimensional arrays are used,
! i.e. only the first column of deriv_data is read.
! The other to dimensions are set to size equal 1
nr = basis%grid%nr
bounds(1:2, 1:3) = 1
bounds(2, 1) = nr
! create a place where to put the derivatives
CALL xc_dset_create(deriv_set, local_bounds=bounds)
! create the place where to store the argument for the functionals
CALL xc_rho_set_create(rho_set, bounds, rho_cutoff=density_cut, &
drho_cutoff=gradient_cut, tau_cutoff=tau_cut)
! allocate the required 3d arrays where to store rho and drho
CALL xc_rho_set_atom_update(rho_set, needs, nspins, bounds)
NULLIFY (rho, drho, tau)
IF (needs%rho) THEN
ALLOCATE (rho(nr, 1))
CALL atom_density(rho(:, 1), pmat, basis, maxl, typ="RHO")
END IF
IF (needs%norm_drho) THEN
ALLOCATE (drho(nr, 1))
CALL atom_density(drho(:, 1), pmat, basis, maxl, typ="DER")
END IF
IF (needs%tau) THEN
ALLOCATE (tau(nr, 1))
CALL atom_density(tau(:, 1), pmat, basis, maxl, typ="KIN", rr=basis%grid%rad2)
END IF
IF (needs%laplace_rho) THEN
ALLOCATE (lap(nr, 1))
CALL atom_density(lap(:, 1), pmat, basis, maxl, typ="LAP", rr=basis%grid%rad2)
END IF
CALL fill_rho_set(rho_set, nspins, needs, rho, drho, tau, lap, nr)
CALL xc_dset_zero_all(deriv_set)
deriv_order = 1
CALL xc_functionals_eval(xc_fun_section, lsd=lsd, rho_set=rho_set, deriv_set=deriv_set, &
deriv_order=deriv_order)
! Integration to get the matrix elements and energy
deriv => xc_dset_get_derivative(deriv_set, [INTEGER::], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
fexc = fourpi*integrate_grid(xcpot(:, 1, 1), basis%grid)
IF (needs%rho) THEN
deriv => xc_dset_get_derivative(deriv_set, [deriv_rho], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
CALL numpot_matrix(xcmat, xcpot(:, 1, 1), basis, 0)
DEALLOCATE (rho)
END IF
IF (needs%norm_drho) THEN
deriv => xc_dset_get_derivative(deriv_set, [deriv_norm_drho], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
CALL numpot_matrix(xcmat, xcpot(:, 1, 1), basis, 1)
DEALLOCATE (drho)
END IF
IF (needs%tau) THEN
deriv => xc_dset_get_derivative(deriv_set, [deriv_tau], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
n1 = SIZE(xcmat, 1)
n2 = SIZE(xcmat, 2)
n3 = SIZE(xcmat, 3)
ALLOCATE (taumat(n1, n2, 0:n3 - 1))
taumat = 0._dp
xcpot(:, 1, 1) = 0.5_dp*xcpot(:, 1, 1)
CALL numpot_matrix(xcmat, xcpot(:, 1, 1), basis, 2)
xcpot(:, 1, 1) = xcpot(:, 1, 1)/basis%grid%rad2(:)
CALL numpot_matrix(taumat, xcpot(:, 1, 1), basis, 0)
DO l = 0, 3
xcmat(:, :, l) = xcmat(:, :, l) + REAL(l*(l + 1), dp)*taumat(:, :, l)
END DO
DEALLOCATE (tau)
DEALLOCATE (taumat)
END IF
! we assume lap is not really needed
IF (needs%laplace_rho) THEN
DEALLOCATE (lap)
END IF
! Release the xc structure used to store the xc derivatives
CALL xc_dset_release(deriv_set)
CALL xc_rho_set_release(rho_set)
END IF !xc_none
CALL timestop(handle)
END SUBROUTINE atom_vxc_lda
! **************************************************************************************************
!> \brief Alternative subroutine that calculates atomic exchange-correlation potential
!> in spin-polarized case.
!> \param basis atomic basis set
!> \param pmata spin-alpha density matrix
!> \param pmatb spin-beta density matrix
!> \param maxl maximum angular momentum
!> \param xc_section XC input section
!> \param fexc exchange-correlation energy
!> \param xcmata spin-alpha exchange-correlation potential matrix
!> \param xcmatb spin-beta exchange-correlation potential matrix
!> \par History
!> * 07.2016 ADMM analysis [Juerg Hutter]
! **************************************************************************************************
SUBROUTINE atom_vxc_lsd(basis, pmata, pmatb, maxl, xc_section, fexc, xcmata, xcmatb)
TYPE(atom_basis_type), POINTER :: basis
REAL(KIND=dp), DIMENSION(:, :, :), INTENT(IN) :: pmata, pmatb
INTEGER, INTENT(IN) :: maxl
TYPE(section_vals_type), POINTER :: xc_section
REAL(KIND=dp), INTENT(OUT) :: fexc
REAL(KIND=dp), DIMENSION(:, :, :), INTENT(INOUT) :: xcmata, xcmatb
CHARACTER(LEN=*), PARAMETER :: routineN = 'atom_vxc_lsd'
INTEGER :: deriv_order, handle, l, myfun, n1, n2, &
n3, nr, nspins
INTEGER, DIMENSION(2, 3) :: bounds
LOGICAL :: lsd
REAL(KIND=dp) :: density_cut, gradient_cut, tau_cut
REAL(KIND=dp), DIMENSION(:, :), POINTER :: drho, lap, rho, tau
REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: taumat, xcpot
TYPE(section_vals_type), POINTER :: xc_fun_section
TYPE(xc_derivative_set_type) :: deriv_set
TYPE(xc_derivative_type), POINTER :: deriv
TYPE(xc_rho_cflags_type) :: needs
TYPE(xc_rho_set_type) :: rho_set
CALL timeset(routineN, handle)
CPASSERT(ASSOCIATED(xc_section))
xc_fun_section => section_vals_get_subs_vals(xc_section, "XC_FUNCTIONAL")
CALL section_vals_val_get(xc_fun_section, "_SECTION_PARAMETERS_", i_val=myfun)
IF (myfun == xc_none) THEN
fexc = 0._dp
ELSE
CALL section_vals_val_get(xc_section, "DENSITY_CUTOFF", r_val=density_cut)
CALL section_vals_val_get(xc_section, "GRADIENT_CUTOFF", r_val=gradient_cut)
CALL section_vals_val_get(xc_section, "TAU_CUTOFF", r_val=tau_cut)
lsd = .TRUE.
nspins = 2
needs = xc_functionals_get_needs(xc_fun_section, lsd=lsd, calc_potential=.FALSE.)
! Prepare the structures needed to calculate and store the xc derivatives
! Array dimension: here anly one dimensional arrays are used,
! i.e. only the first column of deriv_data is read.
! The other to dimensions are set to size equal 1
nr = basis%grid%nr
bounds(1:2, 1:3) = 1
bounds(2, 1) = nr
! create a place where to put the derivatives
CALL xc_dset_create(deriv_set, local_bounds=bounds)
! create the place where to store the argument for the functionals
CALL xc_rho_set_create(rho_set, bounds, rho_cutoff=density_cut, &
drho_cutoff=gradient_cut, tau_cutoff=tau_cut)
! allocate the required 3d arrays where to store rho and drho
CALL xc_rho_set_atom_update(rho_set, needs, nspins, bounds)
NULLIFY (rho, drho, tau)
IF (needs%rho_spin) THEN
ALLOCATE (rho(nr, 2))
CALL atom_density(rho(:, 1), pmata, basis, maxl, typ="RHO")
CALL atom_density(rho(:, 2), pmatb, basis, maxl, typ="RHO")
END IF
IF (needs%norm_drho_spin) THEN
ALLOCATE (drho(nr, 2))
CALL atom_density(drho(:, 1), pmata, basis, maxl, typ="DER")
CALL atom_density(drho(:, 2), pmatb, basis, maxl, typ="DER")
END IF
IF (needs%tau_spin) THEN
ALLOCATE (tau(nr, 2))
CALL atom_density(tau(:, 1), pmata, basis, maxl, typ="KIN", rr=basis%grid%rad2)
CALL atom_density(tau(:, 2), pmatb, basis, maxl, typ="KIN", rr=basis%grid%rad2)
END IF
IF (needs%laplace_rho_spin) THEN
ALLOCATE (lap(nr, 2))
CALL atom_density(lap(:, 1), pmata, basis, maxl, typ="LAP", rr=basis%grid%rad2)
CALL atom_density(lap(:, 2), pmatb, basis, maxl, typ="LAP", rr=basis%grid%rad2)
END IF
CALL fill_rho_set(rho_set, nspins, needs, rho, drho, tau, lap, nr)
CALL xc_dset_zero_all(deriv_set)
deriv_order = 1
CALL xc_functionals_eval(xc_fun_section, lsd=lsd, rho_set=rho_set, deriv_set=deriv_set, &
deriv_order=deriv_order)
! Integration to get the matrix elements and energy
deriv => xc_dset_get_derivative(deriv_set, [INTEGER::], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
fexc = fourpi*integrate_grid(xcpot(:, 1, 1), basis%grid)
IF (needs%rho_spin) THEN
deriv => xc_dset_get_derivative(deriv_set, [deriv_rhoa], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
CALL numpot_matrix(xcmata, xcpot(:, 1, 1), basis, 0)
deriv => xc_dset_get_derivative(deriv_set, [deriv_rhob], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
CALL numpot_matrix(xcmatb, xcpot(:, 1, 1), basis, 0)
DEALLOCATE (rho)
END IF
IF (needs%norm_drho_spin) THEN
! drhoa
NULLIFY (deriv)
deriv => xc_dset_get_derivative(deriv_set, [deriv_norm_drhoa], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
CALL numpot_matrix(xcmata, xcpot(:, 1, 1), basis, 1)
! drhob
NULLIFY (deriv)
deriv => xc_dset_get_derivative(deriv_set, [deriv_norm_drhob], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
CALL numpot_matrix(xcmatb, xcpot(:, 1, 1), basis, 1)
! Cross Terms
NULLIFY (deriv)
deriv => xc_dset_get_derivative(deriv_set, [deriv_norm_drho])
IF (ASSOCIATED(deriv)) THEN
CALL xc_derivative_get(deriv, deriv_data=xcpot)
CALL numpot_matrix(xcmata, xcpot(:, 1, 1), basis, 1)
CALL numpot_matrix(xcmatb, xcpot(:, 1, 1), basis, 1)
END IF
DEALLOCATE (drho)
END IF
IF (needs%tau_spin) THEN
n1 = SIZE(xcmata, 1)
n2 = SIZE(xcmata, 2)
n3 = SIZE(xcmata, 3)
ALLOCATE (taumat(n1, n2, 0:n3 - 1))
deriv => xc_dset_get_derivative(deriv_set, [deriv_tau_a], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
taumat = 0._dp
xcpot(:, 1, 1) = 0.5_dp*xcpot(:, 1, 1)
CALL numpot_matrix(xcmata, xcpot(:, 1, 1), basis, 2)
xcpot(:, 1, 1) = xcpot(:, 1, 1)/basis%grid%rad2(:)
CALL numpot_matrix(taumat, xcpot(:, 1, 1), basis, 0)
DO l = 0, 3
xcmata(:, :, l) = xcmata(:, :, l) + REAL(l*(l + 1), dp)*taumat(:, :, l)
END DO
deriv => xc_dset_get_derivative(deriv_set, [deriv_tau_b], allocate_deriv=.FALSE.)
CALL xc_derivative_get(deriv, deriv_data=xcpot)
taumat = 0._dp
xcpot(:, 1, 1) = 0.5_dp*xcpot(:, 1, 1)
CALL numpot_matrix(xcmatb, xcpot(:, 1, 1), basis, 2)
xcpot(:, 1, 1) = xcpot(:, 1, 1)/basis%grid%rad2(:)
CALL numpot_matrix(taumat, xcpot(:, 1, 1), basis, 0)
DO l = 0, 3
xcmatb(:, :, l) = xcmatb(:, :, l) + REAL(l*(l + 1), dp)*taumat(:, :, l)
END DO
DEALLOCATE (tau)
DEALLOCATE (taumat)
END IF
! Release the xc structure used to store the xc derivatives
CALL xc_dset_release(deriv_set)
CALL xc_rho_set_release(rho_set)
END IF !xc_none
CALL timestop(handle)
END SUBROUTINE atom_vxc_lsd
! **************************************************************************************************
!> \brief Estimate the ADMM exchange energy correction term using a model exchange functional.
!> \param basis0 reference atomic basis set
!> \param pmat0 reference density matrix
!> \param basis1 ADMM basis set
!> \param pmat1 ADMM density matrix
!> \param maxl maxumum angular momentum
!> \param functional name of the model exchange functional:
!> "LINX" -- linear extrapolation of the Slater exchange potential [?]
!> \param dfexc exchange-correlation energy difference
!> \param linxpar LINX parameters
!> \par History
!> * 07.2016 ADMM analysis [Juerg Hutter]
! **************************************************************************************************
SUBROUTINE atom_dpot_lda(basis0, pmat0, basis1, pmat1, maxl, functional, dfexc, linxpar)
TYPE(atom_basis_type), POINTER :: basis0
REAL(KIND=dp), DIMENSION(:, :, :), INTENT(IN) :: pmat0
TYPE(atom_basis_type), POINTER :: basis1
REAL(KIND=dp), DIMENSION(:, :, :), INTENT(IN) :: pmat1
INTEGER, INTENT(IN) :: maxl
CHARACTER(LEN=*), INTENT(IN) :: functional
REAL(KIND=dp), INTENT(OUT) :: dfexc
REAL(KIND=dp), DIMENSION(:), INTENT(IN), OPTIONAL :: linxpar
CHARACTER(LEN=*), PARAMETER :: routineN = 'atom_dpot_lda'
REAL(KIND=dp), PARAMETER :: alx = 0.20_dp, blx = -0.06_dp, &
fs = -0.75_dp*(3._dp/pi)**(1._dp/3._dp)
INTEGER :: handle, ir, nr
REAL(KIND=dp) :: a, b, fx
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: delta, drho0, drho1, pot0, pot1, rho0, &
rho1
CALL timeset(routineN, handle)
nr = basis0%grid%nr