-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathm_pawtab.F90
3104 lines (2744 loc) · 104 KB
/
m_pawtab.F90
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
!!****m* ABINIT/m_pawtab
!! NAME
!! m_pawtab
!!
!! FUNCTION
!! This module contains the definition of the pawtab_type structured datatype,
!! as well as related functions and methods.
!! pawtab_type variables define TABulated data for PAW (from pseudopotential)
!!
!! COPYRIGHT
!! Copyright (C) 2013-2022 ABINIT group (MT)
!! This file is distributed under the terms of the
!! GNU General Public License, see ~abinit/COPYING
!! or http://www.gnu.org/copyleft/gpl.txt .
!!
!! NOTES
!! FOR DEVELOPPERS: in order to preserve the portability of libPAW library,
!! please consult ~abinit/src/??_libpaw/libpaw-coding-rules.txt
!!
!! SOURCE
#include "libpaw.h"
MODULE m_pawtab
USE_DEFS
USE_MSG_HANDLING
USE_MPI_WRAPPERS
USE_MEMORY_PROFILING
implicit none
private
!!***
!----------------------------------------------------------------------
!!****t* m_pawtab/wvlpaw_rholoc_type
!! NAME
!! wvlpaw_rholoc_type
!!
!! FUNCTION
!! Objects for WVL+PAW
!!
!! SOURCE
type,public :: wvlpaw_rholoc_type
integer :: msz
! mesh size
real(dp),allocatable :: d(:,:)
! local rho and derivatives
real(dp),allocatable :: rad(:)
! radial mesh
end type wvlpaw_rholoc_type
public :: wvlpaw_rholoc_free ! Free memory
public :: wvlpaw_rholoc_nullify ! Nullify content
!!***
!----------------------------------------------------------------------
!!****t* m_pawtab/wvlpaw_type
!! NAME
!! wvlpaw_type
!!
!! FUNCTION
!! Objects for WVL+PAW
!!
!! SOURCE
type,public :: wvlpaw_type
! WARNING : if you modify this datatype, please check whether there might be creation/destruction/copy routines,
! declared in another part of ABINIT, that might need to take into account your modification.
!Integer scalars
integer :: npspcode_init_guess
! This is for the PAW-WVL case, only for the initial guess
integer :: ptotgau
! total number of complex gaussians
! for tproj
integer,allocatable :: pngau(:)
! number of complex gaussians per basis element
! for tproj
!Real pointers
real(dp),allocatable :: parg(:,:)
!argument of Gaussians
real(dp),allocatable :: pfac(:,:)
!factors of Gaussians
!Other scalars
type(wvlpaw_rholoc_type) :: rholoc
! local density
! d(:,1): local rho
! d(:,2): local rho 2nd-derivative
! d(:,3): local pot
! d(:,4): local pot 2nd-derivative
end type wvlpaw_type
public :: wvlpaw_allocate ! Allocate memory
public :: wvlpaw_free ! Free memory
public :: wvlpaw_nullify
!!***
!----------------------------------------------------------------------
!!****t* m_pawtab/pawtab_type
!! NAME
!! pawtab_type
!!
!! FUNCTION
!! This structured datatype contains TABulated data for PAW (from pseudopotential)
!! used in PAW calculations.
!!
!! SOURCE
type,public :: pawtab_type
! WARNING : if you modify this datatype, please check whether there might be creation/destruction/copy routines,
! declared in another part of ABINIT, that might need to take into account your modification.
!Integer scalars
integer :: basis_size
! Number of elements for the paw nl basis on the considered atom type
integer :: has_coretau
! Flag controling use of core kinetic enrgy density (AE and pseudo)
! if 1, [t]coretau() is allocated.
! if 2, [t]coretau() is computed and stored.
integer :: has_fock
! if 1, onsite matrix elements of the core-valence Fock operator are allocated
! if 2, onsite matrix elements of the core-valence Fock operator are computed and stored
integer :: has_kij
! if 1, onsite matrix elements of the kinetic operator are allocated
! if 2, onsite matrix elements of the kinetic operator are computed and stored
integer :: has_shapefncg
! if 1, the spherical Fourier transforms of the radial shape functions are allocated
! if 2, the spherical Fourier transforms of the radial shape functions are computed and stored
integer :: has_nabla
! if 1, onsite matrix elements of the nabla operator are allocated.
! if 2, onsite matrix elements of the nabla operator are computed and stored.
! if 3, onsite matrix elements of the nabla operator between valence and core states are computed and stored.
! if 4, onsite matrix elements of the nabla operator between valence and core spinor states are are computed and stored.
integer :: has_nablaphi
! if 1, nablaphi are allocated
! if 2, nablaphi are computed for MetaGGA and stored.
integer :: has_tproj
! Flag controling use of projectors in real space (0 if tnval is unknown)
! if 1, tproj() is allocated.
! if 2, tproj() is computed and stored.
integer :: has_tvale
! Flag controling use of pseudized valence density (0 if tnval is unknown)
! if 1, tvalespl() is allocated.
! if 2, tvalespl() is computed and stored.
integer :: has_vhtnzc
! if 1, space for vhtnzc is allocated
! if 2, vhtnzc has been read from PAW file and stored
integer :: has_vhnzc
! if 1, space for vhnzc is allocated
! if 2, vhnzc has been computed and stored
integer :: has_vminushalf
! has_vminushalf=0 ; vminushal is not allocated
! has_vminushalf=1 ; vminushal is not allocated and stored
integer :: has_wvl
! if 1, data for wavelets (pawwvl) are allocated
! if 2, data for wavelets (pawwvl) are computed and stored
integer :: ij_proj
! Number of (i,j) elements for the orbitals on which U acts (PAW+U only)
! on the considered atom type (ij_proj=1 (1 projector), 3 (2 projectors)...)
! Also used for local exact-exchange
integer :: ij_size
! Number of (i,j) elements for the symetric paw basis
! on the considered atom type (ij_size=basis_size*(basis_size+1)/2)
integer :: lcut_size
! Maximum value of l+1 leading to non zero Gaunt coeffs
! modified by dtset%pawlcutd
! lcut_size=min(2*l_max,dtset%pawlcutd)+1
integer :: l_size
! Maximum value of l+1 leading to non zero Gaunt coeffs
! l_size=2*l_max-1
integer :: lexexch
! lexexch gives l on which local exact-exchange is applied for a given type of atom.
integer :: lmn_size
! Number of (l,m,n) elements for the paw basis
integer :: lmn2_size
! lmn2_size=lmn_size*(lmn_size+1)/2
! where lmn_size is the number of (l,m,n) elements for the paw basis
integer :: lmnmix_sz
! lmnmix_sz=number of klmn=(lmn,lmn_prime) verifying l<=lmix and l_prime<=lmix
integer :: lpawu
! lpawu gives l on which U is applied for a given type of atom.
integer :: nproju
! nproju is the number of projectors for orbitals on which paw+u acts.
! Also used for local exact-exchange
integer :: option_interaction_pawu
! Option for interaction energy (PAW+U) in case of non-collinear magnetism:
! 1: E_int=-J/4.N.(N-2)
! 2: E_int=-J/2.(Nup.(Nup-1)+Ndn.(Ndn-1)) (Nup and Ndn are ill-defined)
! 3: E_int=-J/4.( N.(N-2) + mx^2 + my^2 + mz^2 )
integer :: mesh_size
! Dimension of radial mesh for generic arrays contained in this pawtab datastructure
! The mesh is usually defined up to the PAW augmentation region boundary
! (+ a few additional points). May be different from pawrad%mesh_size
integer :: core_mesh_size
! Dimension of radial mesh for core density
integer :: coretau_mesh_size
! Dimension of radial mesh for core density
integer :: vminus_mesh_size
! Dimension of radial mesh for vminushalf
integer :: partialwave_mesh_size
! Dimension of radial mesh for partial waves (phi, tphi)
! May be different from pawrad%mesh_size and pawtab%mesh_size
integer :: tnvale_mesh_size
! Dimension of radial mesh for tnvale
integer :: mqgrid
! Number of points in the reciprocal space grid on which
! the radial functions (tcorespl, tvalespl...) are specified
! Same as psps%mqgrid_vl
integer :: mqgrid_shp
! Number of points in the reciprocal space grid on which
! the radial shape functions (shapefncg) are given
integer :: usespnorb
! usespnorb=0 ; no spin-orbit coupling
! usespnorb=1 ; spin-orbit coupling
integer :: shape_lambda
! Lambda parameter in gaussian shapefunction (shape_type=2)
integer :: shape_type
! Radial shape function type
! shape_type=-1 ; g(r)=numeric (read from psp file)
! shape_type= 1 ; g(r)=[sin(pi*r/rshp)/(pi*r/rshp)]**2 if r<=rshp, zero if r>rshp
! shape_type= 2 ; g(r)=exp[-(r/sigma)**lambda]
! shape_type= 3 ; gl(r)=Alpha(1,l)*jl(q(1,l)*r)+Alpha(2,l)*jl(q(2,l)*r) for each l
integer :: useexexch
! useexexch=0 ; do not use local exact-exchange
! useexexch=1 ; use local exact-exchange
integer :: usepawu
! usepawu= 0 ; do not use PAW+U formalism
! usepawu= 1 ; use PAW+U formalism (Full localized limit)
! usepawu= 2 ; use PAW+U formalism (Around Mean Field)
! usepawu= 3 ; use PAW+U formalism (Around Mean Field) - Alternative
! usepawu= 4 ; use PAW+U formalism (FLL) without polarization in the XC
! usepawu=-1 ; use PAW+U formalism (FLL) - No use of the occupation matrix - Experimental
! usepawu=-2 ; use PAW+U formalism (AMF) - No use of the occupation matrix - Experimental
! usepawu=-4 ; use PAW+U formalism (FLL) without polarization in the XC - No use of the occupation matrix - Experimental
! usepawu=10 ; use PAW+U within DMFT
! usepawu=14 ; use PAW+U within DMFT without polarization in the XC
integer :: usepotzero
! usepotzero=0 if it is the Kresse-Joubert convention
! usepotzero=1 if it is the new convention
! usepotzero=2 if it is the PWscf convention
integer :: usetcore
! Flag controling use of pseudized core density (0 if tncore=zero)
integer :: usexcnhat
! 0 if compensation charge density is not included in XC terms
! 1 if compensation charge density is included in XC terms
!Real (real(dp)) scalars
real(dp) :: beta
! contains the integral of the difference between vH[nZc] and vH[tnZc]
real(dp) :: dncdq0
! Gives 1/q d(tNcore(q))/dq for q=0
! (tNcore(q) = FT of pseudo core density)
real(dp) :: d2ncdq0
! Gives contribution of d2(tNcore(q))/d2q for q=0
! \int{(16/15)*pi^5*n(r)*r^6* dr}
! (tNcore(q) = FT of pseudo core density)
real(dp) :: dnvdq0
! Gives 1/q d(tNvale(q))/dq for q=0
! (tNvale(q) = FT of pseudo valence density)
real(dp) :: dtaucdq0
! Gives 1/q d(tTAUcore(q))/dq for q=0
! (tTAUcore(q) = FT of pseudo core kinetic density)
real(dp) :: ex_cc
! Exchange energy for the core-core interaction of the Fock operator
real(dp) :: exccore
! Exchange-correlation energy for the core density
real(dp) :: exchmix
! mixing of exact exchange; default is 0.25 (PBE0)
real(dp) :: f4of2_sla
! Ratio of Slater Integrals F4 and F2
real(dp) :: f6of2_sla
! Ratio of Slater Integrals F6 and F4
real(dp) :: jpawu
! jpawu
! Value of J parameter for paw+u for a given type.
real(dp) :: rpaw
! Radius of PAW sphere
real(dp) :: rshp
! Compensation charge radius (if r>rshp, g(r)=zero)
real(dp) :: rcore
! Radius of core corrections (rcore >= rpaw)
real(dp) :: rcoretau
! Radius of kinetic core corrections (rcoretau >= rpaw)
real(dp) :: shape_sigma
! Sigma parameter in gaussian shapefunction (shape_type=2)
real(dp) :: upawu
! upawu
! Value of U parameter for paw+u for a given type.
!Objects
type(wvlpaw_type), pointer :: wvl
!variable containing objects needed
!for wvl+paw implementation
!Warning: it is a pointer; it has to be allocated before use
!Integer arrays
integer, allocatable :: indklmn(:,:)
! indklmn(8,lmn2_size)
! Array giving klm, kln, abs(il-jl), (il+jl), ilm and jlm, ilmn and jlmn for each klmn=(ilmn,jlmn)
! Note: ilmn=(il,im,in) and ilmn<=jlmn
integer, allocatable :: indlmn(:,:)
! indlmn(6,lmn_size)
! For each type of psp,
! array giving l,m,n,lm,ln,spin for i=lmn (if useylm=1)
integer, allocatable :: klmntomn(:,:)
! klmntomn(4,lmn2_size)
! Array giving im, jm ,in, and jn for each klmn=(ilmn,jlmn)
! Note: ilmn=(il,im,in) and ilmn<=jlmn
! NB: klmntomn is an application and not a bijection
integer, allocatable :: kmix(:)
! kmix(lmnmix_sz)
! Indirect array selecting the klmn=(lmn,lmn_prime) verifying l<=lmix and l_prime<=lmix
integer, allocatable :: lnproju(:)
! lnproju(nproju) gives ln (index for phi) for each projectors on which U acts (PAW+U only)
! nproju is 1 or 2 and is the number of projectors for correlated orbitals
! Also used for local exact-exchange
integer, allocatable :: orbitals(:)
! gives the l quantum number per basis element
!Real (real(dp)) arrays
real(dp), allocatable :: coredens(:)
! coredens(mesh_size)
! Gives the core density of the atom
real(dp), allocatable :: coretau(:)
! coretau(mesh_size)
! Gives the kinetic energy density of the atom
real(dp), allocatable :: dij0(:)
! dij0(lmn2_size)
! Part of the Dij term (non-local operator) completely
! calculated in the atomic data part
real(dp), allocatable :: dltij(:)
! dltij(lmn2_size)
! Factor used to compute sums over klmn=(ilmn,jlmn)
! ((ilmn,ilmn) term has to be added once)
! dltij(klmn)=1 if ilmn=jlmn, else dltij(klmn)=2
real(dp), allocatable :: dshpfunc(:,:,:)
! shapefunc(mesh_size,l_size,4)
! Gives the 4 first derivatives of radial shape function
! for each l component; used only if shape_type=-1
real(dp), allocatable :: eijkl(:,:)
! eijkl(lmn2_size,lmn2_size)
! Hartree kernel for the on-site terms (E_hartree=Sum_ijkl[rho_ij rho_kl e_ijkl])
! Used for Hartree and/or Fock contributions
real(dp), allocatable :: eijkl_sr(:,:)
! eijkl_sr(lmn2_size,lmn2_size)
! Screened Hartree kernel for the on-site terms (E_hartree=Sum_ijkl[rho_ij rho_kl e_ijkl_sr])
! Used for screened Fock contributions
real(dp), allocatable :: euijkl(:,:,:,:,:)
! euijkl(3,lmn_size,lmn_size,lmn_size,lmn_size)
! PAW+U kernel for the on-site terms ( E_PAW+U = 0.5 * Sum_ijkl Sum_s1s2 [rho_ij^s1 rho_kl^s2 euijkl(s1,s2)] )
! Contrary to eijkl and eijkl_sr, euijkl is not invariant with respect to the permutations i <--> j or k <--> l
! However, it is still invariant with respect to the permutation i,k <--> j,l, see pawpuxinit.F90
! Also, it depends on two spin indexes
! Used for PAW+U contributions
real(dp), allocatable :: euij_fll(:)
! euij_fll(lmn2_size)
! Double counting part of the PAW+U kernel in the "fully localized limit".This term is only linear with respect to rho_ij,
! while euijkl is quadratic.
! Used for PAW+U contributions
real(dp), allocatable :: ex_cvij(:)
! ex_cvij(lmn2_size))
! Onsite exact_exchange matrix elements for core-valence interactions of the Fock operator
real(dp), allocatable :: fk(:,:)
! fk(6,4)
! Slater integrals used for local exact exchange
real(dp), allocatable :: gammaij(:)
! gammaij(lmn2_size)
! background contribution from the densities
real(dp), allocatable :: gnorm(:)
! gnorm(l_size)
! Give the the normalization factor of each radial shape function
real(dp), allocatable :: kij(:)
! kij(lmn2_size))
! Onsite matrix elements <phi|\kinetic|phj>-<tphi|\kinetic|tphj>
real(dp), allocatable :: nabla_ij(:,:,:)
! nabla_ij(3,lmn_size,lmn_size)
! Onsite matrix elements <phi|\nabla|phj>-<tphi|\nabla|tphj>
real(dp), allocatable :: nabla_im_ij(:,:,:)
! nabla_im_ij(3,lmn_size,lmn_size)
! Imaginary part of onsite matrix elements <phi|\nabla|phj>-<tphi|\nabla|tphj>
! Used in case of core spinor wave functions
real(dp), allocatable :: nablaphi(:,:)
! nablaphi(partialwave_mesh_size, basis_size)
! store the results of dphi/dr-(1/r)phi
real(dp), allocatable :: phi(:,:)
! phi(partialwave_mesh_size, basis_size)
! Gives the paw electron wavefunctions on the radial grid
real(dp), allocatable :: phiphj(:,:)
! phiphj(mesh_size,ij_size)
! Useful product Phi(:,i)*Phi(:,j)
real(dp), allocatable :: phiphjint(:)
! phiphjint(ij_proj)
! Integration of Phi(:,i)*Phi(:,j) for DFT+U/local exact-exchange occupation matrix
real(dp), allocatable :: ph0phiint(:)
! ph0phjint(ij_proj)
! Integration of Phi(:,1)*Phi(:,j) for LDA+DMFT projections
real(dp), allocatable :: qgrid_shp(:)
! qgrid_shp(mqgrid_shp)
! Grid of points in reciprocal space on which the shape functions are given
real(dp), allocatable :: qijl(:,:)
! qijl(l_size**2,lmn2_size)
! The qijl are the moments of the charge density difference between
! the AE and PS partial wave for each channel (i,j). They take part
! to the building of the compensation charge
real(dp), allocatable :: rad_for_spline(:)
! rad_for_spline(mesh_size)
! Radial mesh used to spline quantities on radial mesh;
! Allocated and used only when
! shape_type=-1 (numerical shape function)
! or usedvloc=1 (use of vloc derivative)
real(dp), allocatable :: rhoij0(:)
! rhoij0(lmn2_size)
! Initial guess for rhoij
real(dp), allocatable :: shape_alpha(:,:)
! shape_alpha(2,l_size)
! Alpha_i parameters in Bessel shapefunctions (shape_type=3)
real(dp), allocatable :: shape_q(:,:)
! shape_q(2,l_size)
! Q_i parameters in Bessel shapefunctions (shape_type=3)
real(dp), allocatable :: shapefunc(:,:)
! shapefunc(mesh_size,l_size)
! Gives the normalized radial shape function for each l component
real(dp), allocatable :: shapefncg(:,:,:)
! shapefncg(mqgrid_shp,2,l_size)
! Gives the spherical Fourier transform of the radial shape function
! for each l component (for each qgrid_shp(i)) + second derivative
real(dp), allocatable :: sij(:)
! sij(lmn2_size)
! Nonlocal part of the overlap operator
real(dp), allocatable :: tcoredens(:,:)
! tcoredens(core_mesh_size,1)
! Gives the pseudo core density of the atom
! In PAW+WVL:
! tcoredens(core_mesh_size,2:6)
! are the first to the fifth derivatives of the pseudo core density.
real(dp), allocatable :: tcoretau(:)
! tcoretau(coretau_mesh_size)
! Gives the pseudo core kinetic energy density of the atom
real(dp), allocatable :: tcorespl(:,:)
! tcorespl(mqgrid,2)
! Gives the pseudo core density in reciprocal space on a regular grid
real(dp), allocatable :: tcoretauspl(:,:)
! tcoretauspl(mqgrid,2)
! Gives the pseudo kinetic core density in reciprocal space on a regular grid
real(dp), allocatable :: tnablaphi(:,:)
! tphi(partialwave_mesh_size,basis_size)
! Gives, on the radial grid, the paw atomic pseudowavefunctions
real(dp), allocatable :: tphi(:,:)
! tphi(partialwave_mesh_size,basis_size)
! Gives, on the radial grid, the paw atomic pseudowavefunctions
real(dp), allocatable :: tphitphj(:,:)
! tphitphj(mesh_size,ij_size)
! Useful product tPhi(:,i)*tPhi(:,j)
real(dp), allocatable :: tproj(:,:)
! non-local projectors
real(dp), allocatable :: tvalespl(:,:)
! tvalespl(mqgrid,2)
! Gives the pseudo valence density in reciprocal space on a regular grid
real(dp), allocatable :: Vee(:,:,:,:)
! PAW+U:
! Screened interaction matrix deduced from U and J parameters
! computed on the basis of orbitals on which U acts.
real(dp), allocatable :: Vex(:,:,:,:,:)
! Local exact-exchange:
! Screened interaction matrix deduced from calculation of Slater integrals
! computed on the basis of orbitals on which local exact exchange acts.
real(dp), allocatable :: vhtnzc(:)
! vhtnzc(mesh_size)
! Hartree potential for pseudized Zc density, v_H[\tilde{n}_{Zc}]
! read in from PAW file
real(dp), allocatable :: VHnZC(:)
! VHnZC(mesh_size)
! Hartree potential for Zc density, v_H[n_{Zc}]
! constructed from core density in PAW file (see psp7in.F90)
real(dp), allocatable :: vminushalf(:)
! vminushalf(mesh_size)
! External potential for LDA minus half calculation
! read in from PAW file
real(dp), allocatable :: zioneff(:)
! zioneff(ij_proj)
! "Effective charge"*n "seen" at r_paw, deduced from Phi at r_paw, n:
! pricipal quantum number
! good approximation to model wave function outside PAW-sphere through
end type pawtab_type
public :: pawtab_free ! Free memory
public :: pawtab_nullify ! Nullify content
public :: pawtab_get_lsize ! Get the max. l for a product of 2 partial waves
public :: pawtab_set_flags ! Set the value of the internal flags
public :: pawtab_print ! Printout of the object.
public :: pawtab_bcast ! MPI broadcast the object
interface pawtab_nullify
module procedure pawtab_nullify_0D
module procedure pawtab_nullify_1D
end interface pawtab_nullify
interface pawtab_free
module procedure pawtab_free_0D
module procedure pawtab_free_1D
end interface pawtab_free
interface pawtab_set_flags
module procedure pawtab_set_flags_0D
module procedure pawtab_set_flags_1D
end interface pawtab_set_flags
!!***
CONTAINS !===========================================================
!!***
!----------------------------------------------------------------------
!!****f* m_pawtab/pawtab_nullify_0D
!! NAME
!! pawtab_nullify_0D
!!
!! FUNCTION
!! Nullify pointers and flags in a pawtab structure
!!
!! SIDE EFFECTS
!! Pawtab<type(pawtab_type)>=PAW arrays tabulated.
!! Nullified in output
!!
!! SOURCE
subroutine pawtab_nullify_0D(Pawtab)
!Arguments ------------------------------------
!arrays
type(Pawtab_type),intent(inout) :: Pawtab
!Local variables-------------------------------
! *************************************************************************
!@Pawtab_type
nullify(Pawtab%wvl)
! === Reset all flags and sizes ===
!Flags controlling optional arrays
Pawtab%has_fock=0
Pawtab%has_kij=0
Pawtab%has_tproj=0
Pawtab%has_tvale=0
Pawtab%has_coretau=0
Pawtab%has_vhtnzc=0
Pawtab%has_vhnzc=0
Pawtab%has_vminushalf=0
Pawtab%has_nabla=0
Pawtab%has_nablaphi=0
Pawtab%has_shapefncg=0
Pawtab%has_wvl=0
Pawtab%usetcore=0
Pawtab%usexcnhat=0
Pawtab%useexexch=0
Pawtab%usepawu=0
Pawtab%usepotzero=0
Pawtab%usespnorb=0
Pawtab%mqgrid=0
Pawtab%mqgrid_shp=0
Pawtab%basis_size=0
Pawtab%ij_proj=0
Pawtab%ij_size=0
Pawtab%lcut_size=0
Pawtab%l_size=0
Pawtab%lexexch=-1
Pawtab%lmn_size=0
Pawtab%lmn2_size=0
Pawtab%lmnmix_sz=0
Pawtab%lpawu=-1
Pawtab%nproju=0
Pawtab%option_interaction_pawu=0
Pawtab%mesh_size=0
Pawtab%partialwave_mesh_size=0
Pawtab%core_mesh_size=0
Pawtab%coretau_mesh_size=0
Pawtab%vminus_mesh_size=0
Pawtab%tnvale_mesh_size=0
Pawtab%shape_type=-10
end subroutine pawtab_nullify_0D
!!***
!----------------------------------------------------------------------
!!****f* m_pawtab/pawtab_nullify_1D
!! NAME
!! pawtab_nullify_1D
!!
!! FUNCTION
!! Nullify all pointers in an array of pawtab data structures
!!
!! SOURCE
subroutine pawtab_nullify_1D(Pawtab)
!Arguments ------------------------------------
type(pawtab_type),intent(inout) :: Pawtab(:)
!Local variables-------------------------------
integer :: ii,nn
! *************************************************************************
!@pawtab_type
nn=size(Pawtab)
if (nn==0) return
do ii=1,nn
call pawtab_nullify_0D(Pawtab(ii))
end do
end subroutine pawtab_nullify_1D
!!***
!----------------------------------------------------------------------
!!****f* m_pawtab/pawtab_free_0D
!! NAME
!! pawtab_free_0D
!!
!! FUNCTION
!! Deallocate pointers and nullify flags in a pawtab structure
!!
!! SIDE EFFECTS
!! Pawtab<type(pawtab_type)>=PAW arrays tabulated.
!! All allocated arrays in Pawtab are deallocated
!!
!! SOURCE
subroutine pawtab_free_0D(Pawtab)
!Arguments ------------------------------------
!arrays
type(Pawtab_type),intent(inout) :: Pawtab
!Local variables-------------------------------
! *************************************************************************
!@Pawtab_type
if (allocated(Pawtab%indklmn)) then
LIBPAW_DEALLOCATE(Pawtab%indklmn)
end if
if (allocated(Pawtab%indlmn)) then
LIBPAW_DEALLOCATE(Pawtab%indlmn)
end if
if (allocated(Pawtab%klmntomn)) then
LIBPAW_DEALLOCATE(Pawtab%klmntomn)
end if
if (allocated(Pawtab%kmix)) then
LIBPAW_DEALLOCATE(Pawtab%kmix)
end if
if (allocated(Pawtab%lnproju)) then
LIBPAW_DEALLOCATE(Pawtab%lnproju)
end if
if (allocated(Pawtab%coredens)) then
LIBPAW_DEALLOCATE(Pawtab%coredens)
end if
if (allocated(Pawtab%coretau)) then
LIBPAW_DEALLOCATE(Pawtab%coretau)
end if
if (allocated(Pawtab%dij0)) then
LIBPAW_DEALLOCATE(Pawtab%dij0)
end if
if (allocated(Pawtab%dltij)) then
LIBPAW_DEALLOCATE(Pawtab%dltij)
end if
if (allocated(Pawtab%dshpfunc)) then
LIBPAW_DEALLOCATE(Pawtab%dshpfunc)
end if
if (allocated(Pawtab%eijkl)) then
LIBPAW_DEALLOCATE(Pawtab%eijkl)
end if
if (allocated(Pawtab%eijkl_sr)) then
LIBPAW_DEALLOCATE(Pawtab%eijkl_sr)
end if
if (allocated(Pawtab%euijkl)) then
LIBPAW_DEALLOCATE(Pawtab%euijkl)
end if
if (allocated(Pawtab%euij_fll)) then
LIBPAW_DEALLOCATE(Pawtab%euij_fll)
end if
if (allocated(Pawtab%fk)) then
LIBPAW_DEALLOCATE(Pawtab%fk)
end if
if (allocated(Pawtab%gammaij)) then
LIBPAW_DEALLOCATE(Pawtab%gammaij)
end if
if (allocated(Pawtab%gnorm)) then
LIBPAW_DEALLOCATE(Pawtab%gnorm)
end if
if (allocated(Pawtab%ex_cvij)) then
LIBPAW_DEALLOCATE(Pawtab%ex_cvij)
end if
if (allocated(Pawtab%kij)) then
LIBPAW_DEALLOCATE(Pawtab%kij)
end if
if (allocated(Pawtab%nabla_ij)) then
LIBPAW_DEALLOCATE(Pawtab%nabla_ij)
end if
if (allocated(Pawtab%nabla_im_ij)) then
LIBPAW_DEALLOCATE(Pawtab%nabla_im_ij)
end if
if (allocated(Pawtab%nablaphi)) then
LIBPAW_DEALLOCATE(Pawtab%nablaphi)
end if
if (allocated(Pawtab%orbitals)) then
LIBPAW_DEALLOCATE(Pawtab%orbitals)
end if
if (allocated(Pawtab%phi)) then
LIBPAW_DEALLOCATE(Pawtab%phi)
end if
if (allocated(Pawtab%phiphj)) then
LIBPAW_DEALLOCATE(Pawtab%phiphj)
end if
if (allocated(Pawtab%phiphjint)) then
LIBPAW_DEALLOCATE(Pawtab%phiphjint)
end if
if (allocated(Pawtab%ph0phiint)) then
LIBPAW_DEALLOCATE(Pawtab%ph0phiint)
end if
if (allocated(Pawtab%qgrid_shp)) then
LIBPAW_DEALLOCATE(Pawtab%qgrid_shp)
end if
if (allocated(Pawtab%qijl)) then
LIBPAW_DEALLOCATE(Pawtab%qijl)
end if
if (allocated(Pawtab%rad_for_spline)) then
LIBPAW_DEALLOCATE(Pawtab%rad_for_spline)
end if
if (allocated(Pawtab%rhoij0)) then
LIBPAW_DEALLOCATE(Pawtab%rhoij0)
end if
if (allocated(Pawtab%shape_alpha)) then
LIBPAW_DEALLOCATE(Pawtab%shape_alpha)
end if
if (allocated(Pawtab%shape_q)) then
LIBPAW_DEALLOCATE(Pawtab%shape_q)
end if
if (allocated(Pawtab%shapefunc)) then
LIBPAW_DEALLOCATE(Pawtab%shapefunc)
end if
if (allocated(Pawtab%shapefncg)) then
LIBPAW_DEALLOCATE(Pawtab%shapefncg)
end if
if (allocated(Pawtab%sij)) then
LIBPAW_DEALLOCATE(Pawtab%sij)
end if
if (allocated(Pawtab%tcoredens)) then
LIBPAW_DEALLOCATE(Pawtab%tcoredens)
end if
if (allocated(Pawtab%tcoretau)) then
LIBPAW_DEALLOCATE(Pawtab%tcoretau)
end if
if (allocated(Pawtab%tcorespl)) then
LIBPAW_DEALLOCATE(Pawtab%tcorespl)
end if
if (allocated(Pawtab%tcoretauspl)) then
LIBPAW_DEALLOCATE(Pawtab%tcoretauspl)
end if
if (allocated(Pawtab%tnablaphi)) then
LIBPAW_DEALLOCATE(Pawtab%tphi)
end if
if (allocated(Pawtab%tphi)) then
LIBPAW_DEALLOCATE(Pawtab%tphi)
end if
if (allocated(Pawtab%tphitphj)) then
LIBPAW_DEALLOCATE(Pawtab%tphitphj)
end if
if (allocated(Pawtab%tproj)) then
LIBPAW_DEALLOCATE(Pawtab%tproj)
end if
if (allocated(Pawtab%tvalespl)) then
LIBPAW_DEALLOCATE(Pawtab%tvalespl)
end if
if (allocated(Pawtab%vee)) then
LIBPAW_DEALLOCATE(Pawtab%vee)
end if
if (allocated(Pawtab%Vex)) then
LIBPAW_DEALLOCATE(Pawtab%Vex)
end if
if (allocated(Pawtab%vhtnzc)) then
LIBPAW_DEALLOCATE(Pawtab%vhtnzc)
end if
if (allocated(Pawtab%VHnZC)) then
LIBPAW_DEALLOCATE(Pawtab%VHnZC)
end if
if (allocated(Pawtab%vminushalf)) then
LIBPAW_DEALLOCATE(Pawtab%vminushalf)
end if
if (allocated(Pawtab%zioneff)) then
LIBPAW_DEALLOCATE(Pawtab%zioneff)
end if
call wvlpaw_free(Pawtab%wvl)
! === Reset all flags and sizes ===
!CAUTION: do not reset these flags
!They are set from input data and must be kept
!Pawtab%has_kij=0
!Pawtab%has_tproj=0
!Pawtab%has_coretau=0
!Pawtab%has_tvale=0
!Pawtab%has_vhtnzc=0
!Pawtab%has_vhnzc=0
!Pawtab%has_vminushalf=0
!Pawtab%has_nabla=0
!Pawtab%has_nablaphi=0
!Pawtab%has_shapefncg=0
!Pawtab%has_wvl=0
Pawtab%usetcore=0
Pawtab%usexcnhat=0
Pawtab%useexexch=0
Pawtab%usepawu=0
Pawtab%usepotzero=0
Pawtab%usespnorb=0
Pawtab%mqgrid=0
Pawtab%mqgrid_shp=0
Pawtab%basis_size=0
Pawtab%ij_proj=0
Pawtab%ij_size=0
Pawtab%lcut_size=0
Pawtab%l_size=0
Pawtab%lexexch=-1
Pawtab%lmn_size=0
Pawtab%lmn2_size=0
Pawtab%lmnmix_sz=0
Pawtab%lpawu=-1
Pawtab%nproju=0
Pawtab%option_interaction_pawu=0
Pawtab%mesh_size=0
Pawtab%partialwave_mesh_size=0
Pawtab%core_mesh_size=0
Pawtab%coretau_mesh_size=0
Pawtab%vminus_mesh_size=0
Pawtab%tnvale_mesh_size=0
Pawtab%shape_type=-10
end subroutine pawtab_free_0D
!!***
!----------------------------------------------------------------------
!!****f* m_pawtab/pawtab_free_1D
!! NAME
!! pawtab_free_1D
!!
!! FUNCTION
!! Destroy (deallocate) all pointers in an array of pawtab data structures
!!
!! SOURCE
subroutine pawtab_free_1D(Pawtab)
!Arguments ------------------------------------
type(pawtab_type),intent(inout) :: Pawtab(:)
!Local variables-------------------------------
integer :: ii,nn