-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathm_pawpsp.F90
5111 lines (4607 loc) · 172 KB
/
m_pawpsp.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_pawpsp
!! NAME
!! m_pawpsp
!!
!! FUNCTION
!! Module to read PAW atomic data
!!
!! COPYRIGHT
!! Copyright (C) 2012-2022 ABINIT group (MT, FJ,TR, GJ, FB, FrD, AF, GMR, DRH)
!! 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 DEVELOPERS: in order to preserve the portability of libPAW library,
!! please consult ~abinit/src/??_libpaw/libpaw-coding-rules.txt
!!
!! SOURCE
#include "libpaw.h"
module m_pawpsp
USE_DEFS
USE_MSG_HANDLING
USE_MPI_WRAPPERS
USE_MEMORY_PROFILING
use m_libpaw_libxc
#if defined LIBPAW_HAVE_FOX
use fox_sax
#endif
use m_libpaw_tools, only : libpaw_basename, libpaw_get_free_unit
use m_pawang, only: pawang_type
use m_pawtab, only: pawtab_type, wvlpaw_type, wvlpaw_allocate, wvlpaw_rholoc_free, &
& pawtab_free, wvlpaw_free, wvlpaw_rholoc_nullify, pawtab_bcast, &
& pawtab_set_flags, wvlpaw_allocate, wvlpaw_free, wvlpaw_rholoc_nullify, &
& wvlpaw_rholoc_free
use m_pawxmlps, only: rdpawpsxml_core, paw_setup_t, paw_setuploc, paw_setup_free
use m_pawrad, only: pawrad_type, pawrad_init, pawrad_free, pawrad_copy, &
& pawrad_bcast, pawrad_ifromr, simp_gen, nderiv_gen, bound_deriv, pawrad_deducer0, poisson
use m_paw_numeric, only: paw_splint, paw_spline, paw_smooth, paw_jbessel_4spline
use m_paw_atom, only: atompaw_shapebes, atompaw_vhnzc, atompaw_shpfun, &
& atompaw_dij0, atompaw_kij
use m_pawxc, only: pawxc, pawxcm, pawxc_get_usekden
use m_paw_gaussfit, only: gaussfit_projector
implicit none
private
public:: pawpsp_calc_d5 !calculate up to the 5th derivative
public:: pawpsp_main !main routine to read psp
public:: pawpsp_nl !make paw projector form factors f_l(q)
public:: pawpsp_read !read psp from file
public:: pawpsp_read_header !read header of psp file
public:: pawpsp_read_corewf !read core wavefunction
public:: pawpsp_read_header_2 !reads pspversion, basis_size and lmn_size
public:: pawpsp_rw_atompaw !read and writes ATOMPAW psp with gaussian |p>
public:: pawpsp_wvl !wavelet and icoulomb>0 related operations
public:: pawpsp_wvl_calc !wavelet related operations
public:: pawpsp_7in !reads non-XML atomic data
public:: pawpsp_17in !reads XML atomic data
public:: pawpsp_calc !calculates atomic quantities from psp info
public:: pawpsp_read_header_xml !read header of psp file for XML
public:: pawpsp_read_pawheader !read header variables from XML objects
public:: pawpsp_bcast ! broadcast PAW psp data
public:: pawpsp_cg !compute sin FFT transform of a density
public:: pawpsp_lo !compute sin FFT transform of local potential
! Private procedures
private:: pawpsp_wvl_sin2gauss !convert sin/cos to gaussians
!!***
!-------------------------------------------------------------------------
!!****t* m_pawpsp/pawpsp_header_type
!! NAME
!! pawpsp_header_type
!!
!! FUNCTION
!! For PAW, header related data
!!
!! SOURCE
type, public :: pawpsp_header_type
!Integer scalars
integer :: basis_size ! Number of elements of the wf basis ((l,n) quantum numbers)
integer :: l_size ! Maximum value of l+1 leading to a non zero Gaunt coefficient
integer :: lmn_size ! Number of elements of the paw basis
integer :: mesh_size ! Dimension of (main) radial mesh
integer :: pawver ! Version number of paw psp format
integer :: shape_type ! Type of shape function
real(dp) :: rpaw ! Radius for paw spheres
real(dp) :: rshp ! Cut-off radius of shape function
end type pawpsp_header_type
!!***
CONTAINS
!===========================================================
!!***
!-------------------------------------------------------------------------
!!****f* m_pawpsp/pawpsp_nl
!! NAME
!! pawpsp_nl
!!
!! FUNCTION
!! Make paw projector form factors f_l(q) for each l
!!
!! INPUTS
!! indlmn(6,lmnmax)= array giving l,m,n,lm,ln,s for i=lmn
!! lmnmax=max number of (l,m,n) components
!! lnmax=max number of (l,n) components
!! mqgrid=number of grid points for q grid
!! qgrid(mqgrid)=values at which form factors are returned
!! radmesh <type(pawrad_type)>=data containing radial grid information
!! wfll(:,lnmax)=paw projector on radial grid
!!
!! OUTPUT
!! ffspl(mqgrid,2,lnmax)= form factor f_l(q) and second derivative
!!
!! NOTES
!! u_l(r) is the paw projector (input as wfll);
!! j_l(q) is a spherical Bessel function;
!! f_l(q) = $ \int_0^{rmax}[j_l(2\pi q r) u_l(r) r dr]$
!!
!! SOURCE
subroutine pawpsp_nl(ffspl,indlmn,lmnmax,lnmax,mqgrid,qgrid,radmesh,wfll)
!Arguments ------------------------------------
!scalars
integer,intent(in) :: lmnmax,lnmax,mqgrid
type(pawrad_type),intent(in) :: radmesh
!arrays
integer,intent(in) :: indlmn(6,lmnmax)
real(dp),intent(in) :: qgrid(mqgrid)
real(dp),intent(in) :: wfll(:,:)
real(dp),intent(inout) :: ffspl(mqgrid,2,lnmax)
!Local variables-------------------------------
!scalars
integer :: ilmn,iln,iln0,iq,ir,ll,meshsz,mmax
real(dp),parameter :: eps=tol14**4,TOLJ=0.001_dp
real(dp) :: arg,argn,bes
real(dp) :: besp,qr
real(dp) :: yp1,ypn
character(len=100) :: msg
type(pawrad_type) :: tmpmesh
!arrays
real(dp),allocatable :: ff(:),gg(:),rr(:),rr2(:),rr2wf(:),rrwf(:),work(:)
!*************************************************************************
!Is mesh beginning with r=0 ?
if (radmesh%rad(1)>tol10) then
msg='Radial mesh cannot begin with r<>0!'
LIBPAW_BUG(msg)
end if
meshsz=size(wfll,1)
if (meshsz>radmesh%mesh_size) then
msg='wrong size for wfll!'
LIBPAW_BUG(msg)
end if
!Init. temporary arrays and variables
LIBPAW_ALLOCATE(ff,(meshsz))
LIBPAW_ALLOCATE(gg,(meshsz))
LIBPAW_ALLOCATE(rr,(meshsz))
LIBPAW_ALLOCATE(rr2,(meshsz))
LIBPAW_ALLOCATE(rrwf,(meshsz))
LIBPAW_ALLOCATE(rr2wf,(meshsz))
LIBPAW_ALLOCATE(work,(mqgrid))
rr(1:meshsz) =radmesh%rad(1:meshsz)
rr2(1:meshsz)=two_pi*rr(1:meshsz)*rr(1:meshsz)
argn=two_pi*qgrid(mqgrid)
mmax=meshsz
!Loop on (l,n) projectors
iln0=0
do ilmn=1,lmnmax
iln=indlmn(5,ilmn)
if(iln>iln0) then
iln0=iln;ll=indlmn(1,ilmn)
ir=meshsz
do while (abs(wfll(ir,iln))<eps)
ir=ir-1
end do
mmax=min(ir+1,meshsz)
if (mmax/=radmesh%int_meshsz) then
call pawrad_init(tmpmesh,mesh_size=meshsz,mesh_type=radmesh%mesh_type, &
& rstep=radmesh%rstep,lstep=radmesh%lstep,r_for_intg=rr(mmax))
else
call pawrad_copy(radmesh,tmpmesh)
end if
rrwf(:) =rr (:)*wfll(:,iln)
rr2wf(:)=rr2(:)*wfll(:,iln)
! 1-Compute f_l(0<q<qmax)
if (mqgrid>2) then
do iq=2,mqgrid-1
arg=two_pi*qgrid(iq)
do ir=1,mmax
qr=arg*rr(ir)
call paw_jbessel_4spline(bes,besp,ll,0,qr,TOLJ)
ff(ir)=bes*rrwf(ir)
end do
call simp_gen(ffspl(iq,1,iln),ff,tmpmesh)
end do
end if
! 2-Compute f_l(q=0) and first derivative
ffspl(1,1,iln)=zero;yp1=zero
if (ll==0) then
call simp_gen(ffspl(1,1,iln),rrwf,tmpmesh)
end if
if (ll==1) then
call simp_gen(yp1,rr2wf,tmpmesh)
yp1=yp1*third
end if
! 3-Compute f_l(q=qmax) and first derivative
if (mqgrid>1) then
! if (ll==0.or.ll==1) then
do ir=1,mmax
qr=argn*rr(ir)
call paw_jbessel_4spline(bes,besp,ll,1,qr,TOLJ)
ff(ir)=bes*rrwf(ir)
gg(ir)=besp*rr2wf(ir)
end do
call simp_gen(ffspl(mqgrid,1,iln),ff,tmpmesh)
call simp_gen(ypn,gg,tmpmesh)
else
ypn=yp1
end if
! 4-Compute second derivative of f_l(q)
call paw_spline(qgrid,ffspl(:,1,iln),mqgrid,yp1,ypn,ffspl(:,2,iln))
call pawrad_free(tmpmesh)
! End loop on (l,n) projectors
end if
end do
LIBPAW_DEALLOCATE(ff)
LIBPAW_DEALLOCATE(gg)
LIBPAW_DEALLOCATE(rr)
LIBPAW_DEALLOCATE(rr2)
LIBPAW_DEALLOCATE(rrwf)
LIBPAW_DEALLOCATE(rr2wf)
LIBPAW_DEALLOCATE(work)
end subroutine pawpsp_nl
!!***
!-------------------------------------------------------------------------
!!****f* m_pawpsp/pawpsp_lo
!! NAME
!! pawpsp_lo
!!
!! FUNCTION
!! Compute sine transform to transform from V(r) to q^2 V(q).
!! Computes integrals on (generalized) grid using corrected trapezoidal integration.
!!
!! INPUTS
!! mqgrid=number of grid points in q from 0 to qmax.
!! qgrid(mqgrid)=q grid values (bohr**-1).
!! radmesh <type(pawrad_type)>=data containing radial grid information
!! vloc(:)=V(r) on radial grid.
!! zion=nominal valence charge of atom.
!!
!! OUTPUT
!! epsatm=$ 4\pi\int[r^2 (V(r)+\frac{Zv}{r}dr]$.
!!{{\\ \begin{equation}
!! q2vq(mqgrid)
!! =q^2 V(q)
!! = -\frac{Zv}{\pi}
!! + q^2 4\pi\int[(\frac{\sin(2\pi q r)}{2\pi q r})(r^2 V(r)+r Zv)dr].
!!\end{equation} }}
!! yp1,ypn=derivatives of q^2 V(q) wrt q at q=0 and q=qmax (needed for spline fitter).
!!
!! SOURCE
subroutine pawpsp_lo(epsatm,mqgrid,qgrid,q2vq,radmesh,vloc,yp1,ypn,zion)
!Arguments----------------------------------------------------------
!scalars
integer,intent(in) :: mqgrid
real(dp),intent(in) :: zion
real(dp),intent(out) :: epsatm,yp1,ypn
type(pawrad_type),intent(in) :: radmesh
!arrays
real(dp),intent(in) :: qgrid(mqgrid)
real(dp),intent(in) :: vloc(:)
real(dp),intent(out) :: q2vq(mqgrid)
!Local variables ------------------------------
!scalars
integer :: iq,ir,irmax,mesh_size
real(dp) :: arg,r0tor1,r1torm,rmtoin
logical :: begin_r0
!arrays
real(dp),allocatable :: ff(:),rvpz(:)
!************************************************************************
mesh_size=size(vloc)
irmax=pawrad_ifromr(radmesh,min(20._dp,radmesh%rmax))
irmax=min(irmax,mesh_size)
!Particular case of a zero potential
if (maxval(abs(vloc(1:irmax)))<=1.e-20_dp) then
q2vq=zero;yp1=zero;ypn=zero;epsatm=zero
return
end if
LIBPAW_ALLOCATE(ff,(mesh_size))
LIBPAW_ALLOCATE(rvpz,(mesh_size))
ff=zero;rvpz=zero
!Is mesh beginning with r=0 ?
begin_r0=(radmesh%rad(1)<1.e-20_dp)
!Store r.V+Z
do ir=1,irmax
rvpz(ir)=radmesh%rad(ir)*vloc(ir)+zion
end do
!===========================================
!=== Compute q^2 v(q) for q=0 separately
!===========================================
!Integral from 0 to r1 (only if r1<>0)
r0tor1=zero;if (.not.begin_r0) &
& r0tor1=(zion*0.5_dp+radmesh%rad(1)*vloc(1)/3._dp)*radmesh%rad(1)**2
!Integral from r1 to rmax
do ir=1,irmax
if (abs(rvpz(ir))>1.e-20_dp) then
ff(ir)=radmesh%rad(ir)*rvpz(ir)
end if
end do
call simp_gen(r1torm,ff,radmesh)
!Integral from rmax to infinity
!This part is neglected... might be improved.
rmtoin=zero
!Some of the three parts
epsatm=four_pi*(r0tor1+r1torm+rmtoin)
q2vq(1)=-zion/pi
!===========================================
!=== Compute q^2 v(q) for other q''s
!===========================================
!Loop over q values
do iq=2,mqgrid
arg=two_pi*qgrid(iq)
! Integral from 0 to r1 (only if r1<>0)
r0tor1=zero;if (.not.begin_r0) &
& r0tor1=( vloc(1)/arg*sin(arg*radmesh%rad(1)) &
& -rvpz(1) *cos(arg*radmesh%rad(1)) +zion )/pi
! Integral from r1 to rmax
do ir=1,irmax
if (abs(rvpz(ir))>1.e-20_dp) ff(ir)=sin(arg*radmesh%rad(ir))*rvpz(ir)
end do
call simp_gen(r1torm,ff,radmesh)
! Integral from rmax to infinity
! This part is neglected... might be improved.
rmtoin=zero
! Store q^2 v(q)
q2vq(iq)=-zion/pi + two*qgrid(iq)*(r0tor1+r1torm+rmtoin)
end do
!===========================================
!=== Compute derivatives of q^2 v(q)
!=== at ends of interval
!===========================================
!yp(0)=zero
yp1=zero
!yp(qmax)=$ 2\int_0^\infty[(\sin(2\pi qmax r)+(2\pi qmax r)*\cos(2\pi qmax r)(r V(r)+Z) dr]$
arg=two_pi*qgrid(mqgrid)
!Integral from 0 to r1 (only if r1<>0)
r0tor1=zero;if (.not.begin_r0) &
& r0tor1=zion*radmesh%rad(1) *sin(arg*radmesh%rad(1)) &
& +three*radmesh%rad(1)*vloc(1)/arg *cos(arg*radmesh%rad(1)) &
& +(radmesh%rad(1)**2-one/arg**2)*vloc(1)*sin(arg*radmesh%rad(1))
!Integral from r1 to rmax
do ir=1,irmax
if (abs(rvpz(ir))>1.e-20_dp) ff(ir)=( arg*radmesh%rad(ir)*cos(arg*radmesh%rad(ir)) &
& + sin(arg*radmesh%rad(ir))) *rvpz(ir)
end do
call simp_gen(r1torm,ff,radmesh)
!Integral from rmax to infinity
!This part is neglected... might be improved.
rmtoin=zero
!Some of the three parts
ypn=two*(r0tor1+r1torm+rmtoin)
LIBPAW_DEALLOCATE(ff)
LIBPAW_DEALLOCATE(rvpz)
end subroutine pawpsp_lo
!!***
!-------------------------------------------------------------------------
!!****f* m_pawpsp/pawpsp_cg
!! NAME
!! pawpsp_cg
!!
!! FUNCTION
!! Compute sine transform to transform from n(r) to n(q).
!! Computes integrals on (generalized) grid using corrected trapezoidal integration.
!!
!! INPUTS
!! mqgrid=number of grid points in q from 0 to qmax.
!! qgrid(mqgrid)=q grid values (bohr**-1).
!! radmesh <type(pawrad_type)>=data containing radial grid information
!! nr(:)=n(r) on radial grid.
!!
!! OUTPUT
!! dnqdq0= 1/q dn(q)/dq for q=0
!! d2nqdq0 = Gives contribution of d2(tNcore(q))/d2q for q=0
!! compute \int{(16/15)*pi^5*n(r)*r^6* dr}
!!{{\\ \begin{equation}
!! nq(mqgrid)= n(q)
!! = 4\pi\int[(\frac{\sin(2\pi q r)}{2\pi q r})(r^2 n(r))dr].
!!\end{equation} }}
!! yp1,ypn=derivatives of n(q) wrt q at q=0 and q=qmax (needed for spline fitter).
!!
!! SOURCE
subroutine pawpsp_cg(dnqdq0,d2nqdq0,mqgrid,qgrid,nq,radmesh,nr,yp1,ypn)
!Arguments----------------------------------------------------------
!scalars
integer,intent(in) :: mqgrid
real(dp),intent(out) :: dnqdq0,d2nqdq0,yp1,ypn
type(pawrad_type),intent(in) :: radmesh
!arrays
real(dp),intent(in) :: nr(:)
real(dp),intent(in) :: qgrid(mqgrid)
real(dp),intent(out) :: nq(mqgrid)
!Local variables-------------------------------
!scalars
integer :: iq,ir,mesh_size
real(dp) :: aexp,arg,bexp,dn,r0tor1,r1torm,rm,rmtoin
logical :: begin_r0
!character(len=500) :: msg
!arrays
real(dp),allocatable :: ff(:),rnr(:)
! *************************************************************************
mesh_size=min(size(nr),radmesh%mesh_size)
LIBPAW_ALLOCATE(ff,(mesh_size))
LIBPAW_ALLOCATE(rnr,(mesh_size))
ff=zero;rnr=zero
do ir=1,mesh_size
rnr(ir)=radmesh%rad(ir)*nr(ir)
end do
!Is mesh beginning with r=0 ?
begin_r0=(radmesh%rad(1)<1.d-20)
!Adjustment of an exponentional at r_max (n_exp(r)=aexp*Exp[-bexp*r])
rm=radmesh%rad(mesh_size)
dn=one/(12._dp*radmesh%stepint*radmesh%radfact(mesh_size)) &
& *( 3._dp*nr(mesh_size-4) &
& -16._dp*nr(mesh_size-3) &
& +36._dp*nr(mesh_size-2) &
& -48._dp*nr(mesh_size-1) &
& +25._dp*nr(mesh_size))
if (dn<0._dp.and. &
& abs(radmesh%rad(mesh_size)*nr(mesh_size))>1.d-20) then
bexp=-dn/nr(mesh_size)
if (bexp * rm > 50._dp) then
! This solves the problem with the weird core charge used in v4[62] in which bexp x rm ~= 10^3
!write(msg,"(a,es16.8)")"Tooooo large bexp * rm: ", bexp*rm, ", setting aexp to 0"
!LIBPAW_WARNING(msg)
bexp=0.001_dp;aexp=zero
else
aexp=nr(mesh_size)*exp(bexp*rm)
if (abs(aexp)<1.d-20) then
bexp=0.001_dp;aexp=zero
end if
end if
else
bexp=0.001_dp;aexp=zero
end if
!===========================================
!=== Compute n(q) for q=0 separately
!===========================================
!Integral from 0 to r1 (only if r1<>0)
r0tor1=zero
if (.not.begin_r0) r0tor1=(rnr(1)*radmesh%rad(1)**2)/3.d0
!Integral from r1 to rmax
do ir=1,mesh_size
if (abs(rnr(ir))>1.d-20) ff(ir)=rnr(ir)*radmesh%rad(ir)
end do
call simp_gen(r1torm,ff,radmesh)
!Integral from rmax to infinity
!This part is approximated using an exponential density aexp*Exp[-bexp*r]
!(formulae obtained with mathematica)
rmtoin=aexp*exp(-bexp*rm)/bexp**3*(two+two*bexp*rm+bexp*bexp*rm*rm)
!Some of the three parts
nq(1)=four_pi*(r0tor1+r1torm+rmtoin)
!===========================================
!=== Compute n(q) for other q''s
!===========================================
!Loop over q values
do iq=2,mqgrid
arg=two_pi*qgrid(iq)
! Integral from 0 to r1 (only if r1<>0)
r0tor1=zero;if (.not.begin_r0) &
& r0tor1=nr(1)*(sin(arg*radmesh%rad(1))/arg/arg&
& -radmesh%rad(1)*cos(arg*radmesh%rad(1))/arg)
! Integral from r1 to rmax
do ir=1,mesh_size
if (abs(rnr(ir))>1.d-20) ff(ir)=sin(arg*radmesh%rad(ir))*rnr(ir)
end do
call simp_gen(r1torm,ff,radmesh)
! Integral from rmax to infinity
! This part is approximated using an exponential density aexp*Exp[-bexp*r]
! (formulae obtained with mathematica)
rmtoin=aexp*exp(-bexp*rm)/(arg**2+bexp**2)**2 &
& *(arg*(two*bexp+arg**2*rm+bexp**2*rm)*cos(arg*rm) &
& +(arg**2*(bexp*rm-one)+bexp**2*(bexp*rm+one))*sin(arg*rm))
! Store q^2 v(q)
nq(iq)=two/qgrid(iq)*(r0tor1+r1torm+rmtoin)
end do
!===========================================
!=== Compute derivatives of n(q)
!=== at ends of interval
!===========================================
!yp(0)=zero
yp1=zero
!yp(qmax)=$ 2\int_0^\infty[(-\sin(2\pi qmax r)+(2\pi qmax r)*\cos(2\pi qmax r) r n(r) dr]$
arg=two_pi*qgrid(mqgrid)
!Integral from 0 to r1 (only if r1<>0)
r0tor1=zero;if (.not.begin_r0) &
& r0tor1=two_pi*nr(1)*(3.d0*radmesh%rad(1)/arg /arg*cos(arg*radmesh%rad(1))+ &
& (radmesh%rad(1)**2/arg-3.0d0/arg**3)*sin(arg*radmesh%rad(1)))
!Integral from r1 to rmax
do ir=1,mesh_size
if (abs(rnr(ir))>1.d-20) ff(ir)=(two_pi*radmesh%rad(ir)*cos(arg*radmesh%rad(ir)) &
& - sin(arg*radmesh%rad(ir))/qgrid(mqgrid)) *rnr(ir)
end do
call simp_gen(r1torm,ff,radmesh)
!Integral from rmax to infinity
!This part is approximated using an exponential density aexp*Exp[-bexp*r]
!(formulae obtained with mathematica)
rmtoin=-one/(qgrid(mqgrid)*(arg**2+bexp**2)**3) &
& *aexp*exp(-bexp*rm) &
& *((arg**5*rm-two_pi*arg**4*qgrid(mqgrid)*rm*(bexp*rm-two) &
& +two*arg**3*bexp*(bexp*rm+one)+arg*bexp**3*(bexp*rm+two) &
& -four_pi*arg**2*bexp*qgrid(mqgrid)*(bexp**2*rm**2-three) &
& -two_pi*bexp**3*qgrid(mqgrid)*(bexp**2*rm**2+two*bexp*rm+two))*cos(arg*rm) &
& +(two*arg**2*bexp**3*rm+two_pi*arg**5*qgrid(mqgrid)*rm**2 &
& +arg**4*(bexp*rm-one)+bexp**4*(bexp*rm+one) &
& +four_pi*arg**3*qgrid(mqgrid)*(bexp**2*rm**2+two*bexp*rm-one) &
& +two_pi*arg*bexp**2*qgrid(mqgrid)*(bexp**2*rm**2+four*bexp*rm+6._dp))*sin(arg*rm))
!Some of the three parts
ypn=two/qgrid(mqgrid)*(r0tor1+r1torm+rmtoin)
!===========================================
!=== Compute 1/q dn(q)/dq at q=0
!===========================================
!Integral from 0 to r1 (only if r1<>0)
r0tor1=zero
if (.not.begin_r0) r0tor1=(rnr(1)*radmesh%rad(1)**4)/5.d0
!Integral from r1 to rmax
do ir=1,mesh_size
if (abs(rnr(ir))>1.d-20) ff(ir)=rnr(ir)*radmesh%rad(ir)**3
end do
call simp_gen(r1torm,ff,radmesh)
!Integral from rmax to infinity
!This part is approximated using an exponential density aexp*Exp[-bexp*r]
!(formulae obtained with mathematica)
rmtoin=aexp*exp(-bexp*rm)/bexp**5 &
& *(24._dp+24._dp*bexp*rm+12._dp*bexp**2*rm**2+four*bexp**3*rm**3+bexp**4*rm**4)
!Some of the three parts
dnqdq0=-(2.d0/3.d0)*two_pi**3*(r0tor1+r1torm+rmtoin)
LIBPAW_DEALLOCATE(ff)
LIBPAW_DEALLOCATE(rnr)
d2nqdq0 = 1_dp
end subroutine pawpsp_cg
!!***
!-------------------------------------------------------------------------
!!****f* m_pawpsp/pawpsp_read
!! NAME
!! pawpsp_read
!!
!! FUNCTION
!!
!! INPUTS
!!
!! OUTPUT
!!
!! SIDE EFFECTS
!!
!! NOTES
!! File format of formatted PAW psp input (the 3 first lines
!! have already been read in calling -pspatm- routine) :
!! (1) title (character) line
!! (2) psps%znuclpsp(ipsp), zion, pspdat
!! (3) pspcod, pspxc, lmax, lloc, mmax, r2well
!! (4) psp_version, creatorID
!! (5) basis_size, lmn_size
!! (6) orbitals (for l=1 to basis_size)
!! (7) number_of_meshes
!! For imsh=1 to number_of_meshes
!! (8) mesh_index, mesh_type ,mesh_size, rad_step[, log_step]
!! (9) r_cut(SPH)
!! (10) shape_type, r_shape[, shapefunction arguments]
!! For iln=1 to basis_size
!! (11) comment(character)
!! (12) radial mesh index for phi
!! (13) phi(r) (for ir=1 to phi_meshsz)
!! For iln=1 to basis_size
!! (14) comment(character)
!! (15) radial mesh index for tphi
!! (16) tphi(r) (for ir=1 to phi_mesh_size)
!! For iln=1 to basis_size
!! (17) comment(character)
!! (18) radial mesh index for tproj
!! (19) tproj(r) (for ir=1 to proj_mesh_size)
!! (20) comment(character)
!! (21) radial mesh index for core_density
!! (22) core_density (for ir=1 to core_mesh_size)
!! (23) comment(character)
!! (24) radial mesh index for pseudo_core_density
!! (25) tcore_density (for ir=1 to core_mesh_size)
!! (26) comment(character)
!! (27) Dij0 (for ij=1 to lmn_size*(lmn_size+1)/2)
!! (28) comment(character)
!! (29) Rhoij0 (for ij=1 to lmn_size*(lmn_size+1)/2)
!! (30) comment(character)
!! (31) radial mesh index for Vloc, format of Vloc (0=Vbare, 1=VH(tnzc), 2=VH(tnzc) without nhat in XC)
!! (32) Vloc(r) (for ir=1 to vloc_mesh_size)
!! ===== Following lines only if shape_type=-1 =====
!! For il=1 to 2*max(orbitals)+1
!! (33) comment(character)
!! (34) radial mesh index for shapefunc
!! (35) shapefunc(r)*gnorm(l)*r**l (for ir=1 to shape_mesh_size)
!! (36) comment(character)
!! (37) radial mesh index for pseudo_valence_density
!! (38) tvale(r) (for ir=1 to vale_mesh_size)
!!
!! Comments:
!! * psp_version= ID of PAW_psp version
!! 4 characters string of the form 'pawn' (with n varying)
!! * creatorID= ID of psp generator
!! creatorid=1xyz : psp generated from Holzwarth AtomPAW generator version x.yz
!! creatorid=2xyz : psp generated from Vanderbilt ultra-soft generator version x.yz
!! creatorid=-1: psp for tests (for developpers only)
!! * mesh_type= type of radial mesh
!! mesh_type=1 (regular grid): rad(i)=(i-1)*AA
!! mesh_type=2 (logari. grid): rad(i)=AA*(exp[BB*(i-1)]-1)
!! mesh_type=3 (logari. grid): rad(i>1)=AA*exp[BB*(i-2)] and rad(1)=0
!! mesh_type=4 (logari. grid): rad(i)=-AA*ln[1-BB*(i-1)] with BB=1/n
!! * radial shapefunction type
!! shape_type=-1 ; gl(r)=numeric (read from psp file)
!! shape_type= 1 ; gl(r)=k(r).r^l; k(r)=exp[-(r/sigma)**lambda]
!! shape_type= 2 ; gl(r)=k(r).r^l; k(r)=[sin(pi*r/rshp)/(pi*r/rshp)]**2 if r<=rshp
!! shape_type= 3 ; gl(r)=Alpha(1,l)*jl(q(1,l)*r)+Alpha(2,l)*jl(q(2,l)*r) for each l
!!
!! SOURCE
subroutine pawpsp_read(core_mesh,funit,imainmesh,lmax,&
& ncore,nmesh,pawrad,pawtab,pspversion,radmesh,save_core_msz,&
& tcoretau,tncore,tnvale,tproj,tproj_mesh,usexcnhat_in,usexcnhat_out,vale_mesh,&
& vlocopt,vlocr,vloc_mesh,znucl)
!Arguments ------------------------------------
integer,intent(in):: funit,lmax,usexcnhat_in
integer,intent(out) :: imainmesh,pspversion,usexcnhat_out,vlocopt
logical,intent(in) :: save_core_msz
real(dp),intent(in):: znucl
!arrays
real(dp),pointer :: ncore(:),tcoretau(:),tncore(:),tnvale(:),tproj(:,:),vlocr(:)
type(pawrad_type),intent(inout) :: pawrad
type(pawrad_type),intent(out)::core_mesh,tproj_mesh,vale_mesh,vloc_mesh
type(pawrad_type),pointer :: radmesh(:)
type(pawtab_type),intent(inout) :: pawtab
integer,intent(out)::nmesh
!Local variables-------------------------------
integer :: creatorid,imsh
integer :: icoremesh,ishpfmesh,ivalemesh,ivlocmesh
integer :: ib,il,ilm,ilmn,iln,iprojmesh
integer :: ii,ir,iread1,iread2,jj
integer :: msz,pngau_,ptotgau_
real(dp):: rc,rread1,rread2
real(dp) :: yp1,ypn
!arrays
integer,allocatable :: nprj(:)
real(dp),allocatable :: shpf(:,:),val(:),vhnzc(:)
real(dp),allocatable :: work1(:),work2(:),work3(:),work4(:)
character :: blank=' ',numb=' '
character(len=80) :: pspline
character(len=500) :: msg,submsg
logical :: read_gauss=.false.
type(pawrad_type)::shpf_mesh
! *************************************************************************
!==========================================================
!Read lines 4 to 11 of the header
!This is important for BigDFT in standalone mode
call pawpsp_read_header_2(funit,pspversion,pawtab%basis_size,pawtab%lmn_size)
!Check pspversion for wvl-paw
if(pspversion<4 .and. pawtab%has_wvl>0)then
write(msg, '(a,i2,a,a)' )&
& 'In reading atomic psp file, finds pspversion=',pspversion,ch10,&
& 'For WVL-PAW, pspversion >= 4 is required.'
LIBPAW_BUG(msg)
end if
!Have to maintain compatibility with Abinit v4.2.x
if (pspversion==1) then
LIBPAW_ALLOCATE(pawtab%orbitals,(pawtab%basis_size))
read(funit,*) (pawtab%orbitals(ib), ib=1,pawtab%basis_size)
pawtab%l_size=2*maxval(pawtab%orbitals)+1
nmesh=3
LIBPAW_DATATYPE_ALLOCATE(radmesh,(nmesh))
read(funit,'(a80)') pspline
radmesh(1)%lstep=zero
read(unit=pspline,fmt=*,err=10,end=10) radmesh(1)%mesh_type,&
& radmesh(1)%rstep,radmesh(1)%lstep
10 read(funit,*) pawtab%rpaw
read(funit,*) radmesh(1)%mesh_size,radmesh(2)%mesh_size,&
& radmesh(3)%mesh_size
read(funit,'(a80)') pspline
pawtab%shape_lambda=-1;pawtab%shape_sigma=1.d99
read(unit=pspline,fmt=*,err=11,end=11) pawtab%shape_type,&
& pawtab%shape_lambda,pawtab%shape_sigma
11 read(funit,*) creatorid
if (pawtab%shape_type==3) pawtab%shape_type=-1
radmesh(2)%mesh_type=radmesh(1)%mesh_type
radmesh(3)%mesh_type=radmesh(1)%mesh_type
radmesh(2)%rstep=radmesh(1)%rstep
radmesh(3)%rstep=radmesh(1)%rstep
radmesh(2)%lstep=radmesh(1)%lstep
radmesh(3)%lstep=radmesh(1)%lstep
else
! Here psp file for Abinit 4.3+
LIBPAW_ALLOCATE(pawtab%orbitals,(pawtab%basis_size))
read(funit,*) (pawtab%orbitals(ib), ib=1,pawtab%basis_size)
pawtab%l_size=2*maxval(pawtab%orbitals)+1
read(funit,*) nmesh
LIBPAW_DATATYPE_ALLOCATE(radmesh,(nmesh))
do imsh=1,nmesh
rread2=zero
read(funit,'(a80)') pspline
read(unit=pspline,fmt=*,err=20,end=20) ii,iread1,iread2,rread1,rread2
20 continue
if (ii<=nmesh) then
radmesh(ii)%mesh_type=iread1
radmesh(ii)%mesh_size=iread2
radmesh(ii)%rstep=rread1
radmesh(ii)%lstep=rread2
else
write(msg, '(3a)' )&
& 'Index of mesh out of range !',ch10,&
& 'Action : check your pseudopotential file.'
LIBPAW_ERROR(msg)
end if
end do
read(funit,*) pawtab%rpaw
read(funit,'(a80)') pspline
read(unit=pspline,fmt=*) pawtab%shape_type
pawtab%shape_lambda=-1;pawtab%shape_sigma=1.d99
end if
!Initialize radial meshes
do imsh=1,nmesh
call pawrad_init(radmesh(imsh))
end do
!==========================================================
!Initialize various dims and indexes
pawtab%l_size=2*maxval(pawtab%orbitals)+1
pawtab%lmn2_size=pawtab%lmn_size*(pawtab%lmn_size+1)/2
pawtab%ij_size=pawtab%basis_size*(pawtab%basis_size+1)/2
pawtab%usexcnhat=usexcnhat_in
!indlmn calculation (indices for (l,m,n) basis)
if (allocated(pawtab%indlmn)) then
LIBPAW_DEALLOCATE(pawtab%indlmn)
end if
LIBPAW_ALLOCATE(pawtab%indlmn,(6,pawtab%lmn_size))
LIBPAW_BOUND1_ALLOCATE(nprj,BOUNDS(0,maxval(pawtab%orbitals)))
pawtab%indlmn(:,:)=0
ilmn=0;iln=0;nprj=0
do ib=1,pawtab%basis_size
il=pawtab%orbitals(ib)
nprj(il)=nprj(il)+1
iln=iln+1
do ilm=1,2*il+1
pawtab%indlmn(1,ilmn+ilm)=il
pawtab%indlmn(2,ilmn+ilm)=ilm-(il+1)
pawtab%indlmn(3,ilmn+ilm)=nprj(il)
pawtab%indlmn(4,ilmn+ilm)=il*il+ilm
pawtab%indlmn(5,ilmn+ilm)=iln
pawtab%indlmn(6,ilmn+ilm)=1
end do
ilmn=ilmn+2*il+1
end do
LIBPAW_DEALLOCATE(nprj)
!Are ilmn (found here) and pawtab%lmn_size compatibles ?
if (ilmn/=pawtab%lmn_size) then
write(msg, '(a,a,a,a,a)' )&
& 'Calculated lmn size differs from',ch10,&
& 'lmn_size read from pseudo !',ch10,&
& 'Action: check your pseudopotential file.'
LIBPAW_ERROR(msg)
end if
!==========================================================
!Here reading shapefunction parameters
!Shapefunction parameters for Abinit 4.3...4.5
if (pspversion==2) then
if (pawtab%shape_type==1) read(unit=pspline,fmt=*) ii,pawtab%shape_lambda,pawtab%shape_sigma
if (pawtab%shape_type==3) pawtab%shape_type=-1
pawtab%rshp=zero
!Shapefunction parameters for Abinit 4.6+
else if (pspversion>=3) then
pawtab%rshp=zero
if (pawtab%shape_type==-1) read(unit=pspline,fmt=*,err=21,end=21) ii,pawtab%rshp
if (pawtab%shape_type== 1) read(unit=pspline,fmt=*,err=21,end=21) ii,pawtab%rshp, &
& pawtab%shape_lambda,pawtab%shape_sigma
if (pawtab%shape_type== 2) read(unit=pspline,fmt=*,err=21,end=21) ii,pawtab%rshp
if (pawtab%shape_type== 3) read(unit=pspline,fmt=*,err=21,end=21) ii,pawtab%rshp
end if
21 continue
!If shapefunction type is gaussian, check exponent
if (pawtab%shape_type==1) then
if (pawtab%shape_lambda<2) then
write(msg, '(3a)' )&
& 'For a gaussian shape function, exponent lambda must be >1 !',ch10,&
& 'Action: check your psp file.'
LIBPAW_ERROR(msg)
end if
end if
!If shapefunction type is Bessel, deduce here its parameters from rc
if (pawtab%shape_type==3) then
LIBPAW_ALLOCATE(pawtab%shape_alpha,(2,pawtab%l_size))
LIBPAW_ALLOCATE(pawtab%shape_q,(2,pawtab%l_size))
rc=pawtab%rshp;if (rc<1.d-8) rc=pawtab%rpaw
do il=1,pawtab%l_size
call atompaw_shapebes(pawtab%shape_alpha(1:2,il),pawtab%shape_q(1:2,il),il-1,rc)
end do
end if
!==========================================================
!Mirror pseudopotential parameters to the output and log files
write(msg,'(a,i1)')' Pseudopotential format is: paw',pspversion
call wrtout(ab_out,msg,'COLL')
call wrtout(std_out, msg,'COLL')
write(msg,'(2(a,i3),a,64i4)') &
& ' basis_size (lnmax)=',pawtab%basis_size,' (lmn_size=',&
& pawtab%lmn_size,'), orbitals=',pawtab%orbitals(1:pawtab%basis_size)
call wrtout(ab_out,msg,'COLL')
call wrtout(std_out, msg,'COLL')
write(msg,'(a,f11.8)')' Spheres core radius: rc_sph=',pawtab%rpaw
call wrtout(ab_out,msg,'COLL')
call wrtout(std_out, msg,'COLL')
write(msg,'(a,i1,a)')' ',nmesh,' radial meshes are used:'
call wrtout(ab_out,msg,'COLL')
call wrtout(std_out, msg,'COLL')
do imsh=1,nmesh
if (radmesh(imsh)%mesh_type==1) &
& write(msg,'(a,i1,a,i4,a,g12.5)') &
& ' - mesh ',imsh,': r(i)=step*(i-1), size=',radmesh(imsh)%mesh_size,&
& ' , step=',radmesh(imsh)%rstep
if (radmesh(imsh)%mesh_type==2) &
& write(msg,'(a,i1,a,i4,2(a,g12.5))') &
& ' - mesh ',imsh,': r(i)=AA*[exp(BB*(i-1))-1], size=',radmesh(imsh)%mesh_size,&
& ' , AA=',radmesh(imsh)%rstep,' BB=',radmesh(imsh)%lstep
if (radmesh(imsh)%mesh_type==3) &
& write(msg,'(a,i1,a,i4,2(a,g12.5))') &
& ' - mesh ',imsh,': r(i)=AA*exp(BB*(i-2)), size=',radmesh(imsh)%mesh_size,&
& ' , AA=',radmesh(imsh)%rstep,' BB=',radmesh(imsh)%lstep
if (radmesh(imsh)%mesh_type==4) &
& write(msg,'(a,i1,a,i4,a,g12.5)') &
& ' - mesh ',imsh,': r(i)=-AA*ln(1-(i-1)/n), n=size=',radmesh(imsh)%mesh_size,&
& ' , AA=',radmesh(imsh)%rstep
call wrtout(ab_out,msg,'COLL')
call wrtout(std_out, msg,'COLL')
end do
if (pawtab%shape_type==-1) then
write(msg,'(a)')&
' Shapefunction is NUMERIC type: directly read from atomic data file'
call wrtout(ab_out,msg,'COLL')
call wrtout(std_out, msg,'COLL')
end if
if (pawtab%shape_type==1) then
write(msg,'(2a,a,f6.3,a,i3)')&
& ' Shapefunction is EXP type: shapef(r)=exp(-(r/sigma)**lambda)',ch10,&
& ' with sigma=',pawtab%shape_sigma,' and lambda=',pawtab%shape_lambda
call wrtout(ab_out,msg,'COLL')
call wrtout(std_out, msg,'COLL')
end if
if (pawtab%shape_type==2) then
write(msg,'(a)')&
' Shapefunction is SIN type: shapef(r)=[sin(pi*r/rshp)/(pi*r/rshp)]**2'
call wrtout(ab_out,msg,'COLL')
call wrtout(std_out, msg,'COLL')
end if
if (pawtab%shape_type==3) then
write(msg,'(a)')&
& ' Shapefunction is BESSEL type: shapef(r,l)=aa(1,l)*jl(q(1,l)*r)+aa(2,l)*jl(q(2,l)*r)'
call wrtout(ab_out,msg,'COLL')
call wrtout(std_out, msg,'COLL')
end if
if (pawtab%rshp<1.d-8) then
write(msg,'(a)') ' Radius for shape functions = sphere core radius'
else
write(msg,'(a,f11.8)') ' Radius for shape functions = ',pawtab%rshp
end if
call wrtout(ab_out,msg,'COLL')
call wrtout(std_out, msg,'COLL')
!==========================================================
!Perfom tests
!Are lmax and orbitals compatibles ?
if (lmax/=maxval(pawtab%orbitals)) then
write(msg, '(a,a,a)' )&
& 'lmax /= MAX(orbitals) !',ch10,&
& 'Action: check your pseudopotential file.'
LIBPAW_ERROR(msg)
end if
!Only mesh_type=1,2, 3 or 4 allowed
do imsh=1,nmesh
if (radmesh(imsh)%mesh_type>4) then