-
Notifications
You must be signed in to change notification settings - Fork 0
/
m_paw_numeric.F90
1120 lines (989 loc) · 31.1 KB
/
m_paw_numeric.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_paw_numeric
!! NAME
!! m_paw_numeric
!!
!! FUNCTION
!! Wrappers for various numeric operations (spline, sort, ...)
!!
!! COPYRIGHT
!! Copyright (C) 2012-2022 ABINIT group (MT,TR)
!! 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_paw_numeric
USE_DEFS
USE_MSG_HANDLING
USE_MEMORY_PROFILING
implicit none
private
!public procedures
public:: paw_spline
public:: paw_splint
public:: paw_splint_der
public:: paw_uniform_splfit
public:: paw_smooth
public:: paw_sort_dp
public:: paw_jbessel
public:: paw_solvbes
public:: paw_jbessel_4spline
public:: paw_derfc
!!***
CONTAINS
!===========================================================
!!***
!----------------------------------------------------------------------
!!****f* m_paw_numeric/paw_spline
!! NAME
!! paw_spline
!!
!! FUNCTION
!! Computes the second derivatives of a cubic spline
!!
!! INPUTS
!! * Input, integer N, the number of data points; N must be at least 2.
!! In the special case where N = 2 and IBCBEG = IBCEND = 0, the
!! spline will actually be linear.
!! * Input, real(dp) T(N), the knot values, that is, the points where data
!! is specified. The knot values should be distinct, and increasing.
!! * Input, real(dp) Y(N), the data values to be interpolated.
!! * Input, real(dp) YBCBEG, YBCEND, the values to be used in the boundary
!! conditions if IBCBEG or IBCEND is equal to 1 or 2.
!!
!! OUTPUT
!! Output, real(dp) YPP(N), the second derivatives of the cubic spline.
!! Work space, real(dp) DIAG(N) - should be removed ...
!!
!! SOURCE
subroutine paw_spline(t,y,n,ybcbeg,ybcend,ypp)
!*******************************************************************************
! Discussion:
! For data interpolation, the user must call SPLINE_CUBIC_SET to
! determine the second derivative data, passing in the data to be
! interpolated, and the desired boundary conditions.
! The data to be interpolated, plus the SPLINE_CUBIC_SET output,
! defines the spline. The user may then call SPLINE_CUBIC_VAL to
! evaluate the spline at any point.
! The cubic spline is a piecewise cubic polynomial. The intervals
! are determined by the "knots" or abscissas of the data to be
! interpolated. The cubic spline has continous first and second
! derivatives over the entire interval of interpolation.
! For any point T in the interval T(IVAL), T(IVAL+1), the form of
! the spline is
! SPL(T) = A(IVAL)
! + B(IVAL) * ( T - T(IVAL) )
! + C(IVAL) * ( T - T(IVAL) )**2
! + D(IVAL) * ( T - T(IVAL) )**3
! If we assume that we know the values Y(*) and YPP(*), which represent
! the values and second derivatives of the spline at each knot, then
! the coefficients can be computed as:
! A(IVAL) = Y(IVAL)
! B(IVAL) = ( Y(IVAL+1) - Y(IVAL) ) / ( T(IVAL+1) - T(IVAL) )
! - ( YPP(IVAL+1) + 2 * YPP(IVAL) ) * ( T(IVAL+1) - T(IVAL) ) / 6
! C(IVAL) = YPP(IVAL) / 2
! D(IVAL) = ( YPP(IVAL+1) - YPP(IVAL) ) / ( 6 * ( T(IVAL+1) - T(IVAL) ) )
! Since the first derivative of the spline is
! SPL'(T) = B(IVAL)
! + 2 * C(IVAL) * ( T - T(IVAL) )
! + 3 * D(IVAL) * ( T - T(IVAL) )**2,
! the requirement that the first derivative be continuous at interior
! knot I results in a total of N-2 equations, of the form:
! B(IVAL-1) + 2 C(IVAL-1) * (T(IVAL)-T(IVAL-1))
! + 3 * D(IVAL-1) * (T(IVAL) - T(IVAL-1))**2 = B(IVAL)
! or, setting H(IVAL) = T(IVAL+1) - T(IVAL)
! ( Y(IVAL) - Y(IVAL-1) ) / H(IVAL-1)
! - ( YPP(IVAL) + 2 * YPP(IVAL-1) ) * H(IVAL-1) / 6
! + YPP(IVAL-1) * H(IVAL-1)
! + ( YPP(IVAL) - YPP(IVAL-1) ) * H(IVAL-1) / 2
! =
! ( Y(IVAL+1) - Y(IVAL) ) / H(IVAL)
! - ( YPP(IVAL+1) + 2 * YPP(IVAL) ) * H(IVAL) / 6
! or
! YPP(IVAL-1) * H(IVAL-1) + 2 * YPP(IVAL) * ( H(IVAL-1) + H(IVAL) )
! + YPP(IVAL) * H(IVAL)
! =
! 6 * ( Y(IVAL+1) - Y(IVAL) ) / H(IVAL)
! - 6 * ( Y(IVAL) - Y(IVAL-1) ) / H(IVAL-1)
! Boundary conditions must be applied at the first and last knots.
! The resulting tridiagonal system can be solved for the YPP values.
!
! Author:
! John Burkardt, modified by Xavier Gonze
!Arguments ------------------------------------
!scalars
integer,intent(in) :: n
real(dp),intent(in) :: ybcbeg,ybcend
!arrays
real(dp),intent(in) :: t(n),y(n)
real(dp),intent(out) :: ypp(n)
!Local variables-------------------------------
!scalars
integer :: ibcbeg,ibcend,i,k
real(dp) :: ratio,pinv
character(len=500) :: msg
!arrays
real(dp),allocatable :: tmp(:)
! *************************************************************************
!Check
if (n<=1) then
write(msg,'(6a,i8)') ch10, &
& 'SPLINE_CUBIC_SET - Fatal error!',ch10, &
& ' The number of knots must be at least 2.',ch10, &
& ' The input value of N = ', n
LIBPAW_ERROR(msg)
end if
LIBPAW_ALLOCATE(tmp,(n))
do i=1,n-1
if (t(i)>=t(i+1)) then
write(msg,'(6a,i8,a,es19.12,2a,i8,a,es19.12)') ch10, &
& 'SPLINE_CUBIC_SET - Fatal error!',ch10, &
& ' The knots must be strictly increasing, but',ch10, &
& ' T(', i,') = ', t(i), ch10, &
& ' T(',i+1,') = ', t(i+1)
LIBPAW_ERROR(msg)
end if
end do
ibcbeg=1;if(ybcbeg>1.0d+30)ibcbeg=0
ibcend=1;if(ybcend>1.0d+30)ibcend=0
!Set the first and last equations
if (ibcbeg==0) then
ypp(1) = 0._dp
tmp(1) = 0._dp
else if ( ibcbeg == 1 ) then
ypp(1) = -0.5_dp
tmp(1) = (3._dp/(t(2)-t(1)))*((y(2)-y(1))/(t(2)-t(1))-ybcbeg)
end if
if (ibcend==0) then
ypp(n) = 0._dp
tmp(n) = 0._dp
else if ( ibcend == 1 ) then
ypp(n) = 0.5_dp
tmp(n) = (3._dp/(t(n)-t(n-1)))*(ybcend-(y(n)-y(n-1))/(t(n)-t(n-1)))
end if
!Set the intermediate equations
do i=2,n-1
ratio=(t(i)-t(i-1))/(t(i+1)-t(i-1))
pinv = 1.0_dp/(ratio*ypp(i-1) + 2.0_dp)
ypp(i) = (ratio-1.0_dp)*pinv
tmp(i)=(6.0_dp*((y(i+1)-y(i))/(t(i+1)-t(i))-(y(i)-y(i-1)) &
& /(t(i)-t(i-1)))/(t(i+1)-t(i-1))-ratio*tmp(i-1))*pinv
if (abs(tmp(i))<1.d5*tiny(0._dp)) tmp(i)=0._dp
end do
!Solve the equations
ypp(n) = (tmp(n)-ypp(n)*tmp(n-1))/(ypp(n)*ypp(n-1)+1.0_dp)
do k=n-1,1,-1
ypp(k)=ypp(k)*ypp(k+1)+tmp(k)
end do
LIBPAW_DEALLOCATE(tmp)
end subroutine paw_spline
!!***
!----------------------------------------------------------------------
!!****f* m_paw_numeric/paw_splint
!! NAME
!! paw_splint
!!
!! FUNCTION
!! Compute spline interpolation of a tabulated function.
!! There is no hypothesis about the spacing of the input grid points.
!!
!! INPUTS
!! nspline: number of grid points of input mesh
!! xspline(nspline): input mesh
!! yspline(nspline): function on input mesh
!! ysplin2(nspline): second derivative of yspline on input mesh
!! nfit: number of points of output mesh
!! xfit(nfit): output mesh
!!
!! OUTPUT
!! yfit(nfit): function on output mesh
!! [ierr]=A non-zero value is used to signal that some points in xfit exceed xspline(nspline).
!! The input value is incremented by the number of such points.
!!
!! SOURCE
subroutine paw_splint(nspline,xspline,yspline,ysplin2,nfit,xfit,yfit,ierr)
!Arguments ------------------------------------
!scalars
integer,intent(in) :: nfit, nspline
integer,optional,intent(out) :: ierr
!arrays
real(dp),intent(in) :: xspline(nspline),yspline(nspline)
real(dp),intent(in) :: ysplin2(nspline),xfit(nfit)
real(dp),intent(out) :: yfit(nfit)
!Local variables-------------------------------
!scalars
integer :: left,i,k,right,my_err
real(dp) :: delarg,invdelarg,aa,bb
character(len=50) :: msg
!arrays
! *************************************************************************
my_err=0
left=1
do i=1,nfit
yfit(i)=0._dp ! Initialize for the unlikely event that rmax exceed r(mesh)
do k=left+1, nspline
if(xspline(k) >= xfit(i)) then
if(xspline(k-1) <= xfit(i)) then
right = k
left = k-1
else
if (k-1.eq.1 .and. i.eq.1) then
msg='xfit(1) < xspline(1)'
else
msg='xfit not properly ordered'
end if
LIBPAW_ERROR(msg)
end if
delarg= xspline(right) - xspline(left)
invdelarg= 1.0_dp/delarg
aa= (xspline(right)-xfit(i))*invdelarg
bb= (xfit(i)-xspline(left))*invdelarg
yfit(i) = aa*yspline(left)+bb*yspline(right) &
& +( (aa*aa*aa-aa)*ysplin2(left) + &
& (bb*bb*bb-bb)*ysplin2(right) ) *delarg*delarg/6.0_dp
exit
end if
end do ! k
if (k==nspline+1) my_err=my_err+1 ! xfit not found
end do ! i
if (present(ierr)) ierr=my_err
end subroutine paw_splint
!!***
!----------------------------------------------------------------------
!!****f* m_paw_numeric/paw_splint_der
!! NAME
!! paw_splint_der
!!
!! FUNCTION
!! Compute spline interpolation of the derivative of a tabulated function.
!! There is no hypothesis about the spacing of the input grid points.
!!
!! INPUTS
!! nspline: number of grid points of input mesh
!! xspline(nspline): input mesh
!! yspline(nspline): function on input mesh
!! ysplin2(nspline): second derivative of yspline on input mesh
!! nfit: number of points of output mesh
!! xfit(nfit): output mesh
!!
!! OUTPUT
!! dydxfit(nfit): 1st-derivative of function on output mesh
!! [ierr]=A non-zero value is used to signal that some points in xfit exceed xspline(nspline).
!! The input value is incremented by the number of such points.
!!
!! SOURCE
subroutine paw_splint_der(nspline,xspline,yspline,ysplin2,nfit,xfit,dydxfit,ierr)
!Arguments ------------------------------------
!scalars
integer,intent(in) :: nfit, nspline
integer,optional,intent(out) :: ierr
!arrays
real(dp),intent(in) :: xspline(nspline),yspline(nspline)
real(dp),intent(in) :: ysplin2(nspline),xfit(nfit)
real(dp),intent(out) :: dydxfit(nfit)
!Local variables-------------------------------
!scalars
integer :: left,i,k,right,my_err
real(dp) :: delarg,invdelarg,aa,bb
character(len=50) :: msg
!arrays
! *************************************************************************
my_err=0
left=1
do i=1,nfit
dydxfit(i)=0._dp ! Initialize for the unlikely event that rmax exceed r(mesh)
do k=left+1, nspline
if(xspline(k) >= xfit(i)) then
if(xspline(k-1) <= xfit(i)) then
right = k
left = k-1
else
if (k-1.eq.1 .and. i.eq.1) then
msg='xfit(1) < xspline(1)'
else
msg='xfit not properly ordered'
end if
LIBPAW_ERROR(msg)
end if
delarg= xspline(right) - xspline(left)
invdelarg= 1.0_dp/delarg
aa= (xspline(right)-xfit(i))*invdelarg
bb= (xfit(i)-xspline(left))*invdelarg
dydxfit(i) = (yspline(right)-yspline(left))*invdelarg &
& -( (3.0_dp*(aa*aa)-1.0_dp) *ysplin2(left) &
& -(3.0_dp*(bb*bb)-1.0_dp) *ysplin2(right) ) *delarg/6.0_dp
exit
end if
end do ! k
if (k==nspline+1) my_err=my_err+1 ! xfit not found
end do ! i
if (present(ierr)) ierr=my_err
end subroutine paw_splint_der
!!***
!----------------------------------------------------------------------
!!****f* m_paw_numeric/paw_uniform_splfit
!! NAME
!! paw_uniform_splfit
!!
!! FUNCTION
!! Evaluate cubic spline fit to get function values on input set
!! of ORDERED, UNIFORMLY SPACED points.
!! Optionally gives derivatives (first and second) at those points too.
!! If point lies outside the range of arg, assign the extremal
!! point values to these points, and zero derivative.
!!
!! INPUTS
!! arg(numarg)=equally spaced arguments (spacing delarg) for data
!! to which spline was fit.
!! fun(numarg,2)=function values to which spline was fit and spline
!! fit to second derivatives (from Numerical Recipes spline).
!! ider= see above
!! newarg(numnew)=new values of arguments at which function is desired.
!! numarg=number of arguments at which spline was fit.
!! numnew=number of arguments at which function values are desired.
!!
!! OUTPUT
!! derfun(numnew)=(optional) values of first or second derivative of function.
!! This is only computed for ider=1 or 2; otherwise derfun not used.
!! newfun(numnew)=values of function at newarg(numnew).
!! This is only computed for ider=0 or 1.
!!
!! NOTES
!! if ider=0, compute only the function (contained in fun)
!! if ider=1, compute the function (contained in fun) and its first derivative (in derfun)
!! if ider=2, compute only the second derivative of the function (in derfun)
!!
!! SOURCE
subroutine paw_uniform_splfit(arg,derfun,fun,ider,newarg,newfun,numarg,numnew)
!Arguments ------------------------------------
!scalars
integer, intent(in) :: ider,numarg,numnew
!arrays
real(dp), intent(in) :: arg(numarg),fun(numarg,2),newarg(numnew)
real(dp), intent(out) :: derfun(numnew)
real(dp), intent(inout) :: newfun(numnew)
!Local variables-------------------------------
!scalars
integer :: i,jspl
real(dp) :: argmin,delarg,d,aa,bb,cc,dd
character(len=500) :: msg
!arrays
! *************************************************************************
!argmin is smallest x value in spline fit; delarg is uniform spacing of spline argument
argmin=arg(1)
delarg=(arg(numarg)-argmin)/dble(numarg-1)
if(delarg<tol12)then
write(msg,'(a,es16.8)') 'delarg should be strictly positive, while delarg= ',delarg
LIBPAW_ERROR(msg)
endif
jspl=-1
!Do one loop for no grads, other for grads:
if (ider==0) then
! Spline index loop for no grads:
do i=1,numnew
if (newarg(i).ge.arg(numarg)) then
newfun(i)=fun(numarg,1)
else if (newarg(i).le.arg(1)) then
newfun(i)=fun(1,1)
else
jspl=1+int((newarg(i)-argmin)/delarg)
d=newarg(i)-arg(jspl)
bb = d/delarg
aa = 1.0d0-bb
cc = aa*(aa**2-1.0d0)*(delarg**2/6.0d0)
dd = bb*(bb**2-1.0d0)*(delarg**2/6.0d0)
newfun(i)=aa*fun(jspl,1)+bb*fun(jspl+1,1)+cc*fun(jspl,2)+dd*fun(jspl+1,2)
end if
end do
else if(ider==1)then
! Spline index loop includes grads:
do i=1,numnew
if (newarg(i).ge.arg(numarg)) then
newfun(i)=fun(numarg,1)
derfun(i)=0.0d0
else if (newarg(i).le.arg(1)) then
newfun(i)=fun(1,1)
derfun(i)=0.0d0
else
! cubic spline interpolation:
jspl=1+int((newarg(i)-arg(1))/delarg)
d=newarg(i)-arg(jspl)
bb = d/delarg
aa = 1.0d0-bb
cc = aa*(aa**2-1.0d0)*(delarg**2/6.0d0)
dd = bb*(bb**2-1.0d0)*(delarg**2/6.0d0)
newfun(i)=aa*fun(jspl,1)+bb*fun(jspl+1,1)+cc*fun(jspl,2)+dd*fun(jspl+1,2)
! spline fit to first derivative:
! note correction of Numerical Recipes sign error
derfun(i) = (fun(jspl+1,1)-fun(jspl,1))/delarg + &
& (-(3.d0*aa**2-1.d0)*fun(jspl,2)+ &
& (3.d0*bb**2-1.d0)*fun(jspl+1,2)) * delarg/6.0d0
end if
end do
else if (ider==2) then
do i=1,numnew
if (newarg(i).ge.arg(numarg)) then
derfun(i)=0.0d0
else if (newarg(i).le.arg(1)) then
derfun(i)=0.0d0
else
! cubic spline interpolation:
jspl=1+int((newarg(i)-argmin)/delarg)
d=newarg(i)-arg(jspl)
bb = d/delarg
aa = 1.0d0-bb
! second derivative of spline (piecewise linear function)
derfun(i) = aa*fun(jspl,2)+bb*fun(jspl+1,2)
end if
end do
end if
end subroutine paw_uniform_splfit
!!***
!----------------------------------------------------------------------
!!****f* m_paw_numeric/paw_smooth
!! NAME
!! paw_smooth
!!
!! FUNCTION
!! Smooth an array of given ordinates (y's) that are in order of
!! increasing abscissas (x's), but without using the abscissas themselves
!! supposed to be equally spaced.
!!
!! INPUTS
!! it=number of abscissas to treat
!! mesh=size of the array (number of abscissas)
!!
!! OUTPUT
!!
!! SIDE EFFECTS
!! a(mesh)=array to be smoothed
!!
!! SOURCE
subroutine paw_smooth(a,mesh,it)
!Arguments ------------------------------------
!scalars
integer, intent(in) :: it,mesh
!arrays
real(dp), intent(inout) :: a(mesh)
!Local variables-------------------------------
!scalars
integer :: i,k
!arrays
real(dp) :: asm(mesh)
! *************************************************************************
asm(1:4) = zero ! ?? Correct me ...
do k=1,it
asm(5)=0.2_dp*(a(3)+a(4)+a(5)+a(6)+a(7))
asm(mesh-4)=0.2_dp*(a(mesh-2)+a(mesh-3)+a(mesh-4)+&
& a(mesh-5)+a(mesh-6))
asm(mesh-3)=0.2_dp*(a(mesh-1)+a(mesh-2)+a(mesh-3)+&
& a(mesh-4)+a(mesh-5))
asm(mesh-2)=0.2_dp*(a(mesh)+a(mesh-1)+a(mesh-2)+&
& a(mesh-3)+a(mesh-4))
asm(mesh-1)=0.25_dp*(a(mesh)+a(mesh-1)+a(mesh-2)+a(mesh-3))
asm(mesh)=1.0_dp/3.0_dp*(a(mesh)+a(mesh-1)+a(mesh-2))
do i=6,mesh-5
asm(i)=0.1_dp *a(i)+0.1_dp*(a(i+1)+a(i-1))+&
& 0.1_dp *(a(i+2)+a(i-2))+&
& 0.1_dp *(a(i+3)+a(i-3))+&
& 0.1_dp *(a(i+4)+a(i-4))+&
& 0.05_dp*(a(i+5)+a(i-5))
end do
do i=1,mesh
a(i)=asm(i)
end do
end do
end subroutine paw_smooth
!!***
!----------------------------------------------------------------------
!!****f* m_paw_numeric/paw_sort_dp
!! NAME
!! paw_sort_dp
!!
!! FUNCTION
!! Sort real(dp) array list(n) into ascending numerical order using Heapsort
!! algorithm, while making corresponding rearrangement of the integer
!! array iperm. Consider that two real(dp) numbers
!! within tolerance tol are equal.
!!
!! INPUTS
!! n intent(in) dimension of the list
!! tol intent(in) numbers within tolerance are equal
!! list(n) intent(inout) list of real(dp) numbers to be sorted
!! iperm(n) intent(inout) iperm(i)=i (very important)
!!
!! OUTPUT
!! list(n) sorted list
!! iperm(n) index of permutation given the right ascending order
!!
!! SOURCE
subroutine paw_sort_dp(n,list,iperm,tol)
!Arguments ------------------------------------
!scalars
integer, intent(in) :: n
real(dp), intent(in) :: tol
!arrays
integer, intent(inout) :: iperm(n)
real(dp), intent(inout) :: list(n)
!Local variables-------------------------------
!scalars
integer :: l,ir,iap,i,j
real(dp) :: ap
character(len=500) :: msg
!arrays
! *************************************************************************
!Accomodate case of array of length 1: already sorted!
if (n==1) return
!Should not call with n<1
if (n<1) then
write(msg,'(a,i12,2a)') &
& 'paw_sort_dp has been called with array length n=',n, ch10, &
& ' having a value less than 1. This is not allowed.'
LIBPAW_ERROR(msg)
end if
!Conduct the usual sort
l=n/2+1 ; ir=n
do ! Infinite do-loop
if (l>1) then
l=l-1
ap=list(l)
iap=iperm(l)
else ! l<=1
ap=list(ir)
iap=iperm(ir)
list(ir)=list(1)
iperm(ir)=iperm(1)
ir=ir-1
if (ir==1) then
list(1)=ap
iperm(1)=iap
exit ! This is the end of this algorithm
end if
end if ! l>1
i=l
j=l+l
do while (j<=ir)
if (j<ir) then
if ( list(j)<list(j+1)-tol .or. &
& (list(j)<list(j+1)+tol.and.iperm(j)<iperm(j+1))) j=j+1
endif
if (ap<list(j)-tol.or.(ap<list(j)+tol.and.iap<iperm(j))) then
list(i)=list(j)
iperm(i)=iperm(j)
i=j
j=j+j
else
j=ir+1
end if
end do
list(i)=ap
iperm(i)=iap
end do ! End infinite do-loop
end subroutine paw_sort_dp
!!***
!----------------------------------------------------------------------
!!****f* m_paw_numeric/paw_jbessel
!! NAME
!! paw_jbessel
!!
!! FUNCTION
!! Compute spherical Bessel function j_l(x) and derivative(s)
!!
!! INPUTS
!! ll=l-order of the Bessel function
!! order=1 if first derivative is requested
!! 2 if first and second derivatives are requested
!! xx=where to compute j_l
!!
!! OUTPUT
!! bes= Bessel function j_l at xx
!! besp= first derivative of j_l at xx (only if order>=1)
!! bespp= second derivative of j_l at xx (only if order=2)
!!
!! SOURCE
subroutine paw_jbessel(bes,besp,bespp,ll,order,xx)
!Arguments ---------------------------------------------
!scalars
integer,intent(in) :: ll,order
real(dp),intent(in) :: xx
real(dp),intent(out) :: bes,besp,bespp
!Local variables ---------------------------------------
!scalars
integer,parameter :: imax=40
integer :: ii,il
real(dp),parameter :: prec=1.d-15
real(dp) :: besp1,fact,factp,factpp,jn,jnp,jnpp,jr,xx2,xxinv
character(len=200) :: msg
! *********************************************************************
if (order>2) then
msg='Wrong order in paw_jbessel!'
LIBPAW_ERROR(msg)
end if
if (abs(xx)<prec) then
bes=zero;if (ll==0) bes=one
if (order>=1) then
besp=zero;if (ll==1) besp=third
end if
if (order==2) then
bespp=zero
if (ll==0) bespp=-third
if (ll==2) bespp=2._dp/15._dp
end if
return
end if
xxinv=one/xx
if (order==0) then
factp=zero
factpp=zero
jnp=zero
jnpp=zero
end if
if (xx<one) then
xx2=0.5_dp*xx*xx
fact=one
do il=1,ll
fact=fact*xx/dble(2*il+1)
end do
jn=one;jr=one;ii=0
do while(abs(jr)>=prec.and.ii<imax)
ii=ii+1;jr=-jr*xx2/dble(ii*(2*(ll+ii)+1))
jn=jn+jr
end do
bes=jn*fact
if (abs(jr)>prec) then
msg='Bessel function did not converge!'
LIBPAW_ERROR(msg)
end if
if (order>=1) then
factp=fact*xx/dble(2*ll+3)
jnp=one;jr=one;ii=0
do while(abs(jr)>=prec.AND.ii<imax)
ii=ii+1;jr=-jr*xx2/dble(ii*(2*(ll+ii)+3))
jnp=jnp+jr
end do
besp=-jnp*factp+jn*fact*xxinv*dble(ll)
if (abs(jr)>prec) then
msg='1st der. of Bessel function did not converge!'
LIBPAW_ERROR(msg)
end if
end if
if (order==2) then
factpp=factp*xx/dble(2*ll+5)
jnpp=one;jr=one;ii=0
do while(abs(jr)>=prec.AND.ii<imax)
ii=ii+1;jr=-jr*xx2/dble(ii*(2*(ll+ii)+5))
jnpp=jnpp+jr
end do
besp1=-jnpp*factpp+jnp*factp*xxinv*dble(ll+1)
if (abs(jr)>prec) then
msg='2nd der. of Bessel function did not converge !'
LIBPAW_ERROR(msg)
end if
end if
else
jn =sin(xx)*xxinv
jnp=(-cos(xx)+jn)*xxinv
do il=2,ll+1
jr=-jn+dble(2*il-1)*jnp*xxinv
jn=jnp;jnp=jr
end do
bes=jn
if (order>=1) besp =-jnp+jn *xxinv*dble(ll)
if (order==2) besp1= jn -jnp*xxinv*dble(ll+2)
end if
if (order==2) bespp=-besp1+besp*ll*xxinv-bes*ll*xxinv*xxinv
end subroutine paw_jbessel
!!***
!----------------------------------------------------------------------
!!****f* m_paw_numeric/paw_solvbes
!! NAME
!! paw_solvbes
!!
!! FUNCTION
!! Find nq first roots of instrinsic equation:
!! alpha.jl(Q) + beta.Q.djl/dr(Q) = 0
!!
!! INPUTS
!! alpha,beta= factors in intrinsic equation
!! ll= l quantum number
!! nq= number of roots to find
!!
!! OUTPUT
!! root(nq)= roots of instrinsic equation
!!
!! SOURCE
subroutine paw_solvbes(root,alpha,beta,ll,nq)
!Arguments ------------------------------------
!scalars
integer :: ll,nq
real(dp) :: alpha,beta
!arrays
real(dp) :: root(nq)
!Local variables-------------------------------
!scalars
integer :: nroot
real(dp),parameter :: dh=0.1_dp,tol=tol14
real(dp) :: dum,hh,jbes,jbesp,qq,qx,y1,y2
! *************************************************************************
qq=dh;nroot=0
do while (nroot<nq)
call paw_jbessel(jbes,jbesp,dum,ll,1,qq)
y1=alpha*jbes+beta*qq*jbesp
qq=qq+dh
call paw_jbessel(jbes,jbesp,dum,ll,1,qq)
y2=alpha*jbes+beta*qq*jbesp
do while (y1*y2>=zero)
qq=qq+dh
call paw_jbessel(jbes,jbesp,dum,ll,1,qq)
y2=alpha*jbes+beta*qq*jbesp
end do
hh=dh;qx=qq
do while (hh>tol)
hh=half*hh
if (y1*y2<zero) then
qx=qx-hh
else
qx=qx+hh
end if
call paw_jbessel(jbes,jbesp,dum,ll,1,qx)
y2=alpha*jbes+beta*qx*jbesp
end do
nroot=nroot+1
root(nroot)=qx
end do
end subroutine paw_solvbes
!!***
!----------------------------------------------------------------------
!!****f* m_special_funcs/paw_jbessel_4spline
!! NAME
!! paw_jbessel_4spline
!!
!! FUNCTION
!! Compute spherical Bessel functions and derivatives.
!! A polynomial approximation is employed for q-->0.
!!
!! INPUTS
!! ll=l-order of the Bessel function
!! tol=tolerance below which a Polynomial approximation is employed
!! both for jl and its derivative (if required)
!! order=1 if only first derivative is requested
!! 2 if first and second derivatives are requested
!! xx=where to compute j_l
!!
!! OUTPUT
!! bes=Spherical Bessel function j_l at xx
!! besp= first derivative of j_l at xx (only if order>=1)
!!
!! TODO
!! Remove inline definitions, they are obsolete in F2003
!!
!! SOURCE
subroutine paw_jbessel_4spline(bes,besp,ll,order,xx,tol)
!Arguments ---------------------------------------------
!scalars
integer,intent(in) :: ll,order
real(dp),intent(in) :: xx,tol
real(dp),intent(out) :: bes,besp
!Local variables ---------------------------------------
!scalars
real(dp) :: bespp
!real(dp) :: arg,bes0a,bes0ap,bes0b,bes0bp,bes1a,bes1ap,bes1b,bes1bp
!real(dp) :: bes2a,bes2ap,bes2b,bes2bp,bes3a,bes3ap,bes3b,bes3bp
character(len=100) :: msg
! *********************************************************************
! === l=0,1,2 and 3 spherical Bessel functions (and derivatives) ===
! Statement functions are obsolete. Sorry ...
!bes0a(arg)=1.0_dp-arg**2/6.0_dp*(1.0_dp-arg**2/20.0_dp)
!bes0b(arg)=sin(arg)/arg
!bes1a(arg)=(10.0_dp-arg*arg)*arg/30.0_dp
!bes1b(arg)=(sin(arg)-arg*cos(arg))/arg**2
!bes2a(arg)=arg*arg/15.0_dp-arg**4/210.0_dp
!bes2b(arg)=((3.0_dp-arg**2)*sin(arg)-3.0_dp*arg*cos(arg))/arg**3
!bes3a(arg)=arg*arg*arg/105.0_dp-arg**5/1890.0_dp+arg**7/83160.0_dp
!bes3b(arg)=(15.0_dp*sin(arg)-15.0_dp*arg*cos(arg)-6.0_dp*arg**2*sin(arg)+arg**3*cos(arg))/arg**4
!bes0ap(arg)=(-10.0_dp+arg*arg)*arg/30.0_dp
!bes0bp(arg)=-(sin(arg)-arg*cos(arg))/arg**2
!bes1ap(arg)=(10.0_dp-3.0_dp*arg*arg)/30.0_dp
!bes1bp(arg)=((arg*arg-2.0_dp)*sin(arg)+2.0_dp*arg*cos(arg))/arg**3
!bes2ap(arg)=(1.0_dp-arg*arg/7.0_dp)*2.0_dp*arg/15.0_dp
!bes2bp(arg)=((4.0_dp*arg*arg-9.0_dp)*sin(arg)+(9.0_dp-arg*arg)*arg*cos(arg))/arg**4
!bes3ap(arg)=(1.0_dp/35-arg*arg/378.0_dp+arg**4/11880.0_dp)*arg*arg
!bes3bp(arg)=((-60.0_dp+27.0_dp*arg*arg-arg**4)*sin(arg)+(60.0_dp*arg-7.0_dp*arg**3)*cos(arg))/arg**5
! This is to test paw_jbessel calculation without polynomial approximation for q-->0.
! call paw_jbessel(bes,besp,bespp,ll,order,xx)
! RETURN
if (order>2) then
msg='Wrong order in paw_jbessel_4spline'
LIBPAW_ERROR(msg)
end if
select case (ll)
case (0)
if (xx<TOL) then
bes=1.0_dp-xx**2/6.0_dp*(1.0_dp-xx**2/20.0_dp)
if (order>=1) besp=(-10.0_dp+xx*xx)*xx/30.0_dp
else
bes=sin(xx)/xx
if (order>=1) besp=-(sin(xx)-xx*cos(xx))/xx**2
end if
case (1)
if (xx<TOL) then
bes=(10.0_dp-xx*xx)*xx/30.0_dp
if (order>=1) besp=(10.0_dp-3.0_dp*xx*xx)/30.0_dp
else
bes=(sin(xx)-xx*cos(xx))/xx**2
if (order>=1) besp=((xx*xx-2.0_dp)*sin(xx)+2.0_dp*xx*cos(xx))/xx**3
end if
case (2)
if (xx<TOL) then
bes=xx*xx/15.0_dp-xx**4/210.0_dp
if (order>=1) besp=(1.0_dp-xx*xx/7.0_dp)*2.0_dp*xx/15.0_dp
else
bes=((3.0_dp-xx**2)*sin(xx)-3.0_dp*xx*cos(xx))/xx**3
if (order>=1) besp=((4.0_dp*xx*xx-9.0_dp)*sin(xx)+(9.0_dp-xx*xx)*xx*cos(xx))/xx**4
end if
case (3)
if (xx<TOL) then
bes=xx*xx*xx/105.0_dp-xx**5/1890.0_dp+xx**7/83160.0_dp
if (order>=1) besp=(1.0_dp/35-xx*xx/378.0_dp+xx**4/11880.0_dp)*xx*xx
else
bes=(15.0_dp*sin(xx)-15.0_dp*xx*cos(xx)-6.0_dp*xx**2*sin(xx)+xx**3*cos(xx))/xx**4
if (order>=1) besp=((-60.0_dp+27.0_dp*xx*xx-xx**4)*sin(xx)+(60.0_dp*xx-7.0_dp*xx**3)*cos(xx))/xx**5
end if
case (4:)
call paw_jbessel(bes,besp,bespp,ll,order,xx)
case default
write(msg,'(a,i4)')' wrong value for ll = ',ll
LIBPAW_BUG(msg)
end select
end subroutine paw_jbessel_4spline
!!***
!----------------------------------------------------------------------
!!****f* m_paw_numeric/paw_derfc
!! NAME
!! paw_derfc
!!
!! FUNCTION
!! Evaluates the complementary error function in real(dp).
!!
!! INPUTS
!! yy
!!
!! OUTPUT
!! derfc_yy=complementary error function of yy
!!
!! SOURCE
elemental function paw_derfc(yy) result(derfc_yy)
!Arguments ------------------------------------
!scalars