-
Notifications
You must be signed in to change notification settings - Fork 0
/
m_pawrad.F90
1758 lines (1524 loc) · 55.8 KB
/
m_pawrad.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_pawrad
!! NAME
!! m_pawrad
!!
!! FUNCTION
!! Module containing all the functions related to the PAW radial meshes
!!
!! COPYRIGHT
!! Copyright (C) 2013-2022 ABINIT group (MT,FJ,MG)
!! 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
!! * Routines tagged with "@type_name" are strongly connected to the definition of the data type.
!! Strongly connected means that the proper functioning of the implementation relies on the
!! assumption that the tagged procedure is consistent with the type declaration.
!! Every time a developer changes the structure "type_name" adding new entries, he/she has to make sure
!! that all the strongly connected routines are changed accordingly to accommodate the modification of the data type
!! Typical examples of strongly connected routines are creation, destruction or reset methods.
!!
!! * 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_pawrad
USE_DEFS
USE_MSG_HANDLING
USE_MPI_WRAPPERS
USE_MEMORY_PROFILING
use m_paw_numeric, only : paw_derfc
implicit none
private
!public procedures.
public :: pawrad_init ! Main creation method
public :: pawrad_free ! Free the allocated memory
public :: pawrad_print ! Printout of the basic info
public :: pawrad_isame ! Checks whether two meshes are equivalent or have the same equation.
public :: pawrad_copy ! Returns a copy of the mesh.
public :: pawrad_ifromr ! Retrieve the Index FROM a given R value in a radial grid.
public :: pawrad_deducer0 ! Extrapolate r=0 value of a function from values near r=0.
public :: pawrad_bcast ! Broadcast pawrad datastructure over a given MPI communicator
public :: simp_gen ! Performs integral on a given (generalized) grid using Simpson rule.
public :: nderiv_gen ! Do corrected first (and 2nd) derivation on a given (generalized) grid.
public :: nderiv_lin ! Do corrected first (and 2nd) derivation on a given linear grid.
public :: bound_deriv ! Computes derivatives at boundaries of the mesh
public :: poisson ! Solves Poisson eq. for angularly dependent charge distribution of angular momentum l
public :: screened_coul_kernel ! Kernel used to compute short-range screened Coulomb integrals
public :: calc_slatradl ! Calculates the radial part of Slater integrals.
interface pawrad_free
module procedure pawrad_free_0D
module procedure pawrad_free_1D
end interface pawrad_free
! TODO: Might use bit flags, but all radmesh stuff should be encapsulated here!
integer,private,parameter :: RMESH_LINEAR = 1
integer,private,parameter :: RMESH_LOG1 = 2
integer,private,parameter :: RMESH_LOG2 = 3
integer,private,parameter :: RMESH_LOG3 = 4
integer,private,parameter :: RMESH_NL = 5
!!***
!-------------------------------------------------------------------------
!!****t* m_pawrad/pawrad_type
!! NAME
!! pawrad_type
!!
!! FUNCTION
!! For PAW, RADial mesh discretization and related data
!!
!! SOURCE
type, public :: pawrad_type
!Integer scalars
integer :: int_meshsz=0
! Mesh size used in integrals computation
! Integrals will be computed up to r(int_meshsz)
integer :: mesh_size=0
! Dimension of radial mesh
integer :: mesh_type=-1
! Type of mesh
! 1=regular grid: r(i)=(i-1)*AA
! 2=logarithmic grid: r(i)=AA*(exp[BB*(i-1)]-1)
! 3=logarithmic grid: r(i>1)=AA*exp[BB*(i-1)] and r(1)=0
! 4=logarithmic grid: r(i)=-AA*ln[1-BB*(i-1)] with BB=1/n
!Real (real(dp)) scalars
real(dp) :: lstep=zero
! Exponential step of the mesh (BB parameter above)
! Defined only if mesh type is logarithmic
real(dp) :: rmax=zero
! Max. value of r = rad(mesh_size)
real(dp) :: rstep=zero
! Radial step of the mesh (AA parameter above)
real(dp) :: stepint=zero
! Radial step used to convert any function from the
! present grid onto a regular grid in order to
! integrate it using trapeze method
!Real (real(dp)) arrays
real(dp), allocatable :: rad(:)
! rad(mesh_size)
! Coordinates of all the points of the mesh
real(dp), allocatable :: radfact(:)
! radfact(mesh_size)
! Factor used to compute radial integrals
! Before being integrated on the present mesh,
! any function is multiplied by this factor
real(dp), allocatable :: simfact(:)
! simfact(mesh_size)
! Factor used to compute radial integrals by the a Simpson scheme
! Integral[f] = Sum_i [simfact(i)*f(i)]
end type pawrad_type
!!***
CONTAINS
!===========================================================
!!***
!!****f* m_pawrad/pawrad_init
!! NAME
!! pawrad_init
!!
!! FUNCTION
!! Creation method for radial meshes.
!! Compute all points (and related weights) of a radial mesh.
!! Grid can be regular or logarithimc.
!!
!! INPUTS
!! [mesh_size]=Dimension of the radial mesh
!! If not present, take mesh%mesh_size
!! [mesh_type]=Type of mesh
!! If not present, take mesh%mesh_type
!! [rstep]=Radial step of the mesh (AA parameter above)
!! If not present, take mesh%rstep
!! [lstep]=Exponential step of the mesh (BB parameter above)
!! If not present, take mesh%lstep
!! Needed only if mesh type is logarithmic.
!! [r_for_intg]=Mesh size used in integrals computation
!! If not present, take mesh%r_for_intg
!! Integrals will be computed up to rr(r_for_intg)
!! (can be negative for an integration over the whole grid)
!!
!! OUTPUT
!!
!! SIDE EFFECTS
!! mesh<pawrad_type>=The object completely initialized (containing radial grid information).
!! The following quantities are calculated inside the routine:
!! %stepint = Radial step used to convert any function from the
!! present grid onto a regular grid in order to integrate it using trapeze method
!! %rad(mesh_size) = Coordinates of all the points of the mesh.
!! %radfact(mesh_size) = Factors used to compute radial integrals.
!! %int_meshsz = Integrals will be computed up to r(int_meshsz)
!! %simfact(mesh_size) = Factor used to compute radial integrals by the a Simpson scheme
!! Integral[f] = Sum_i [simfact(i)*f(i)]
!! %rmax = Max. value of r = rad(mesh_size)
!!
!! NOTES
!! Possible mesh types (mesh%mesh_type)
!! 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
!! mesh_type=5 ( grid): rad(i)=AA*i/(n-i)
!!
!! SOURCE
subroutine pawrad_init(mesh,mesh_size,mesh_type,rstep,lstep,r_for_intg)
!Arguments ------------------------------------
!scalars
integer,intent(in),optional :: mesh_size,mesh_type
real(dp),intent(in),optional :: rstep,lstep
real(dp),intent(in),optional :: r_for_intg
type(pawrad_type),intent(inout) :: mesh
!Local variables-------------------------------
!scalars
integer :: ir,ir_last,isim,mesh_size_,mesh_type_
real(dp) :: hh,lstep_,rstep_,r_for_intg_
character(len=500) :: msg
! *************************************************************************
!@pawrad_type
!Retrieve mesh data
mesh_size_ = mesh%mesh_size ; if (present(mesh_size)) mesh_size_ = mesh_size
mesh_type_ = mesh%mesh_type ; if (present(mesh_type)) mesh_type_ = mesh_type
rstep_ = mesh%rstep ; if (present(rstep)) rstep_ = rstep
lstep_ = mesh%lstep ; if (present(lstep)) lstep_ = lstep
r_for_intg_=-1._dp;if (present(r_for_intg)) r_for_intg_=r_for_intg
mesh%mesh_size = mesh_size_
mesh%mesh_type = mesh_type_
mesh%rstep = rstep_
mesh%lstep = lstep_
LIBPAW_ALLOCATE(mesh%rad ,(mesh%mesh_size))
LIBPAW_ALLOCATE(mesh%radfact,(mesh%mesh_size))
LIBPAW_ALLOCATE(mesh%simfact,(mesh%mesh_size))
mesh%simfact=zero
if (mesh%mesh_type==1) then
isim=3
mesh%stepint=mesh%rstep
mesh%rad(1)=zero;mesh%radfact(1)=one
do ir=2,mesh%mesh_size
mesh%rad(ir) =mesh%rstep*dble(ir-1)
mesh%radfact(ir)=one
end do
else if (mesh%mesh_type==2) then
isim=3
mesh%stepint=mesh%lstep
mesh%rad(1)=zero;mesh%radfact(1)=mesh%rstep
do ir=2,mesh%mesh_size
mesh%rad(ir) =mesh%rstep*(exp(mesh%lstep*dble(ir-1))-one)
mesh%radfact(ir)=mesh%rad(ir)+mesh%rstep
end do
else if (mesh%mesh_type==3) then
isim=4
mesh%stepint=mesh%lstep
mesh%rad(1)=zero;mesh%radfact(1)=zero
do ir=2,mesh%mesh_size
mesh%rad(ir) =mesh%rstep*exp(mesh%lstep*dble(ir-2))
mesh%radfact(ir)=mesh%rad(ir)
end do
else if (mesh%mesh_type==4) then
isim=3
mesh%lstep=one/dble(mesh%mesh_size)
mesh%stepint=mesh%lstep
mesh%rad(1)=zero;mesh%radfact(1)=mesh%rstep
do ir=2,mesh%mesh_size
mesh%rad(ir) =-mesh%rstep*log(one-mesh%lstep*dble(ir-1))
mesh%radfact(ir)=mesh%rstep/(one-mesh%lstep*dble(ir-1))
end do
else if (mesh%mesh_type==5) then
isim=3
mesh%stepint=mesh%rstep
mesh%rad(1)=zero;mesh%radfact(1)=1/mesh%lstep
do ir=2,mesh%mesh_size
mesh%rad(ir) =mesh%rstep*dble(ir-1)/(mesh%lstep-dble(ir-1))
mesh%radfact(ir)=(mesh%rstep+mesh%rad(ir))/(mesh%lstep-dble(ir-1))/mesh%rstep
end do
else ! Other values of mesh_type are not allowed (see psp7in.F90)
write(msg,'(a,i0)')" Unknown value of mesh_type: ",mesh%mesh_type
LIBPAW_ERROR(msg)
end if
mesh%int_meshsz=mesh%mesh_size
if (r_for_intg_>0.d0) then
ir=min(pawrad_ifromr(mesh,r_for_intg_),mesh%mesh_size)
if (ir<mesh%mesh_size) then
if (abs(mesh%rad(ir+1)-r_for_intg_)<abs(mesh%rad(ir)-r_for_intg_)) ir=ir+1
end if
if (ir>1) then
if (abs(mesh%rad(ir-1)-r_for_intg_)<abs(mesh%rad(ir)-r_for_intg_)) ir=ir-1
end if
mesh%int_meshsz=ir
end if
hh=mesh%stepint/3.d0
mesh%simfact(mesh%int_meshsz)=hh*mesh%radfact(mesh%int_meshsz)
mesh%simfact(1:isim-2)=zero
ir_last=1
do ir=mesh%int_meshsz,isim,-2
mesh%simfact(ir-1)=4.d0*hh*mesh%radfact(ir-1)
mesh%simfact(ir-2)=2.d0*hh*mesh%radfact(ir-2)
ir_last=ir-2
end do
mesh%simfact(ir_last)=half*mesh%simfact(ir_last)
if (mesh%int_meshsz<mesh%mesh_size) mesh%simfact(mesh%int_meshsz+1:mesh%mesh_size)=zero
mesh%rmax=mesh%rad(mesh%mesh_size)
end subroutine pawrad_init
!!***
!----------------------------------------------------------------------
!!****f* m_pawrad/pawrad_free_0D
!! NAME
!! pawrad_free_0D
!!
!! FUNCTION
!! Frees all memory allocated in the object
!!
!! SOURCE
subroutine pawrad_free_0D(Rmesh)
!Arguments ------------------------------------
!arrays
type(pawrad_type),intent(inout) :: Rmesh
!Local variables-------------------------------
! *************************************************************************
!@Pawrad_type
if (allocated(Rmesh%rad )) then
LIBPAW_DEALLOCATE(Rmesh%rad)
end if
if (allocated(Rmesh%radfact)) then
LIBPAW_DEALLOCATE(Rmesh%radfact)
end if
if (allocated(Rmesh%simfact)) then
LIBPAW_DEALLOCATE(Rmesh%simfact)
end if
Rmesh%int_meshsz=0
Rmesh%mesh_size=0
Rmesh%mesh_type=-1
end subroutine pawrad_free_0D
!!***
!----------------------------------------------------------------------
!!****f* m_pawrad/pawrad_free_1D
!! NAME
!! pawrad_free_1D
!!
!! FUNCTION
!! Destroy all objects in an array of pawrad data structures
!!
!! SOURCE
subroutine pawrad_free_1D(Rmesh)
!Arguments ------------------------------------
type(pawrad_type),intent(inout) :: Rmesh(:)
!Local variables-------------------------------
integer :: ii,nn
! *************************************************************************
!@pawrad_type
nn=size(Rmesh)
if (nn==0) return
do ii=1,nn
call pawrad_free_0D(Rmesh(ii))
end do
end subroutine pawrad_free_1D
!!***
!----------------------------------------------------------------------
!!****f* m_pawrad/pawrad_print
!! NAME
!! pawrad_print
!!
!! FUNCTION
!! Reports basic info on the object.
!!
!! INPUTS
!! Rmesh<pawrad_type>=Object defining the radial mesh
!! header=String for the header provided by the user.
!! [unit]=Unit number for output, defaults to std_out
!! [prtvol]=Verbosity level, minimal if not specified.
!! [mode_paral]=Either "COLL" or "PERS". Passed to wrtout. Defaults to "COLL"
!!
!! OUTPUT
!! Only writing.
!!
!! SOURCE
subroutine pawrad_print(Rmesh,header,unit,prtvol,mode_paral)
!Arguments ------------------------------------
integer,intent(in),optional :: prtvol,unit
character(len=4),intent(in),optional :: mode_paral
character(len=*),intent(in),optional :: header
type(pawrad_type),intent(in) :: Rmesh
!Local variables-------------------------------
!scalars
integer :: my_unt,my_prtvol
character(len=4) :: my_mode
character(len=500) :: msg
! *************************************************************************
!@pawrad_type
my_unt =std_out; if (PRESENT(unit )) my_unt =unit
my_prtvol=0 ; if (PRESENT(prtvol )) my_prtvol=prtvol
my_mode ='COLL' ; if (PRESENT(mode_paral)) my_mode =mode_paral
msg=ch10//' ==== Info on the Radial Mesh ==== '
if (PRESENT(header)) msg=ch10//' ==== '//TRIM(ADJUSTL(header))//' ==== '
call wrtout(my_unt,msg,my_mode)
SELECT CASE (Rmesh%mesh_type)
CASE (RMESH_LINEAR)
write(msg,'(a,i4,a,g12.5)')&
& ' - Linear mesh: r(i)=step*(i-1), size=',Rmesh%mesh_size,', step=',Rmesh%rstep
CASE (RMESH_LOG1)
write(msg,'(a,i4,2(a,g12.5))')&
& ' - Logarithimc mesh: r(i)=AA*[exp(BB*(i-1))-1], size=',Rmesh%mesh_size,', AA=',Rmesh%rstep,' BB=',Rmesh%lstep
CASE (RMESH_LOG2)
write(msg,'(a,i4,2(a,g12.5))')&
& ' - Logarithimc mesh: r(i)=AA*exp(BB*(i-2)), size=',Rmesh%mesh_size,', AA=',Rmesh%rstep,' BB=',Rmesh%lstep
CASE (RMESH_LOG3)
write(msg,'(a,i1,a,i4,a,g12.5)')&
& ' - Logarithimc mesh: r(i)=-AA*ln(1-(i-1)/n), n=size=',Rmesh%mesh_size,', AA=',Rmesh%rstep
CASE (RMESH_NL)
write(msg,'(a,i1,a,i4,a,g12.5)')&
& ' - Non-linear mesh: r(i)=-AA*i/(n-i), n=size=',Rmesh%mesh_size,', AA=',Rmesh%rstep
CASE DEFAULT
msg = ' Unknown mesh type! Action : check your pseudopotential or input file.'
LIBPAW_ERROR(msg)
END SELECT
call wrtout(my_unt,msg,my_mode)
if (my_prtvol>1) then
write(msg,'(a,i4)')' Mesh size for integrals = ',Rmesh%int_meshsz
call wrtout(my_unt,msg,my_mode)
write(msg,'(a,g12.5)')' rmax=rad(mesh_size) = ',Rmesh%rmax
call wrtout(my_unt,msg,my_mode)
write(msg,'(a,g12.5)')' Value of stepint = ',Rmesh%stepint
call wrtout(my_unt,msg,my_mode)
end if
end subroutine pawrad_print
!!***
!----------------------------------------------------------------------
!!****f* m_pawrad/pawrad_isame
!! NAME
!! pawrad_isame
!!
!! FUNCTION
!! Check two radial meshes, returns a logical flag defining
!! whether the meshes have the same equation and an integer
!! flag
!!
!! INPUTS
!! Rmesh1,Rmesh2<pawrad_type>=The two radial meshes.
!!
!! OUTPUT
!! hasameq=.true. if the two meshes are defined by the same equation.
!! whichdenser=
!! * 0 if meshes are not compatible
!! * 1 if Rmesh1 is denser than Rmesh2
!! * 2 if Rmesh2 is denser than Rmesh1
!!
!! SOURCE
subroutine pawrad_isame(Rmesh1,Rmesh2,hasameq,whichdenser)
!Arguments ------------------------------------
integer,intent(out) :: whichdenser
logical,intent(out) :: hasameq
type(pawrad_type),intent(in) :: Rmesh1
type(pawrad_type),intent(in) :: Rmesh2
!Local variables-------------------------------
character(len=50) :: msg
! *************************************************************************
!@pawrad_type
whichdenser =0 ; hasameq=.FALSE.
if (Rmesh1%mesh_type /= Rmesh2%mesh_type) RETURN
SELECT CASE (Rmesh1%mesh_type)
CASE (RMESH_LINEAR) !check for linear meshes
hasameq = (Rmesh1%rstep == Rmesh2%rstep)
CASE (RMESH_LOG1,& !check for logarithmic meshes
& RMESH_LOG2,&
& RMESH_LOG3)
hasameq = ( Rmesh1%rstep == Rmesh2%rstep &
& .and.Rmesh1%lstep == Rmesh2%lstep )
CASE (RMESH_NL) !check for linear meshes
hasameq = (Rmesh1%rstep == Rmesh2%rstep)
CASE DEFAULT
msg='Unknown mesh type'
LIBPAW_ERROR(msg)
END SELECT
! === If meshes have same equation, check whether they are equal ===
! * Note that also int_meshsz must be equal
if (hasameq) then
whichdenser= 1
if (Rmesh2%mesh_size > Rmesh1%mesh_size ) whichdenser = 2
!if (Rmesh1%int_meshsz == Rmesh2%int_meshsz) whichdenser = 2
end if
end subroutine pawrad_isame
!!***
!----------------------------------------------------------------------
!!****f* m_pawrad/pawrad_copy
!! NAME
!! pawrad_copy
!!
!! FUNCTION
!! Copy one radial mesh (in a generalized format) to another
!!
!! INPUTS
!! mesh1 <type(pawrad_type)>=data containing radial grid information of input mesh
!!
!! OUTPUT
!! mesh2 <type(pawrad_type)>=data containing radial grid information of output mesh
!!
!! NOTES
!! Possible mesh types (mesh%mesh_type)
!! 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
!! mesh_type=5 ( grid): rad(i)=AA*i/(n-i)
!!
!! SOURCE
subroutine pawrad_copy(mesh1,mesh2)
!Arguments ------------------------------------
!scalars
type(pawrad_type),intent(in) :: mesh1
type(pawrad_type),intent(out) :: mesh2
!Local variables-------------------------------
!scalars
integer :: ir
! *************************************************************************
mesh2%mesh_type =mesh1%mesh_type
mesh2%mesh_size =mesh1%mesh_size
mesh2%int_meshsz=mesh1%int_meshsz
mesh2%lstep =mesh1%lstep
mesh2%rstep =mesh1%rstep
mesh2%stepint =mesh1%stepint
mesh2%rmax =mesh1%rmax
LIBPAW_ALLOCATE(mesh2%rad,(mesh1%mesh_size))
LIBPAW_ALLOCATE(mesh2%radfact,(mesh1%mesh_size))
LIBPAW_ALLOCATE(mesh2%simfact,(mesh1%mesh_size))
do ir=1,mesh1%mesh_size
mesh2%rad(ir) =mesh1%rad(ir)
mesh2%radfact(ir)=mesh1%radfact(ir)
mesh2%simfact(ir)=mesh1%simfact(ir)
end do
end subroutine pawrad_copy
!!***
!----------------------------------------------------------------------
!!****f* m_pawrad/pawrad_deducer0
!! NAME
!! pawrad_deducer0
!!
!! FUNCTION
!! Extrapolate r=0 value of a function from values near r=0
!! using a 3 points formula
!!
!! INPUTS
!! funcsz=size of array func
!! radmesh <type(pawrad_type)>=data containing radial grid information
!!
!! SIDE EFFECTS
!! func(funcsz)=array containing values of function to extrapolate
!!
!! SOURCE
subroutine pawrad_deducer0(func,funcsz,radmesh)
!Arguments ------------------------------------
!scalars
integer,intent(in) :: funcsz
type(pawrad_type),intent(in) :: radmesh
!arrays
real(dp),intent(inout) :: func(funcsz)
! *************************************************************************
if (radmesh%mesh_type==1.or.radmesh%mesh_type==2.or.radmesh%mesh_type==4.or.radmesh%mesh_type==5) then
func(1)=func(4)+3*(func(2)-func(3))
else if (radmesh%mesh_type==3) then
func(1)=func(4)+exp(two*radmesh%lstep)/(exp(radmesh%lstep)-one)*(func(2)-func(3))
end if
end subroutine pawrad_deducer0
!!***
!----------------------------------------------------------------------
!!****f* m_pawrad/pawrad_bcast
!! NAME
!! pawrad_bcast
!!
!! FUNCTION
!! Communicate pawrad data over a given MPI communicator
!!
!! INPUTS
!! comm_mpi= communicator used to broadcast data
!!
!! SIDE EFFECTS
!! pawrad=<type pawrad_type>= a radial mesh datastructure for PAW
!!
!! SOURCE
subroutine pawrad_bcast(pawrad,comm_mpi)
!Arguments ------------------------------------
!scalars
integer,intent(in) :: comm_mpi
type(pawrad_type),intent(inout) :: pawrad
!Local variables-------------------------------
!scalars
integer :: ierr,indx,me,nn,isz1
integer :: if_rad,if_radfact,if_simfact !flags used to communicate
character(len=500) :: msg
!arrays
integer,allocatable :: list_int(:)
real(dp),allocatable :: list_dpr(:)
!*************************************************************************
me=xpaw_mpi_comm_rank(comm_mpi)
!Initializations
if_rad=0; if_radfact=0; if_simfact=0
!calculate the size of the reals
if(me==0) then
if (allocated(pawrad%rad)) then
if_rad=1 !communicate rad
isz1=size(pawrad%rad)
if(isz1/=pawrad%mesh_size) then
msg='rad: sz1 /= pawrad%mesh_size (1)'
LIBPAW_BUG(msg)
end if
end if
if (allocated(pawrad%radfact)) then
if_radfact=1 !communicate radfact
isz1=size(pawrad%radfact)
if(isz1/=pawrad%mesh_size) then
msg='radfact: sz1 /= pawrad%mesh_size (2)'
LIBPAW_BUG(msg)
end if
end if
if (allocated(pawrad%simfact)) then
if_simfact=1 !communicate simfact
isz1=size(pawrad%simfact)
if(isz1/=pawrad%mesh_size) then
msg='simfact: sz1 /= pawrad%mesh_size (3)'
LIBPAW_BUG(msg)
end if
end if
end if
!Brodcast the integers
LIBPAW_ALLOCATE(list_int,(6))
if(me==0) then
list_int(1)=pawrad%int_meshsz
list_int(2)=pawrad%mesh_size
list_int(3)=pawrad%mesh_type
list_int(4)=if_rad
list_int(5)=if_radfact
list_int(6)=if_simfact
end if
call xpaw_mpi_bcast(list_int,0,comm_mpi,ierr)
if(me/=0) then
pawrad%int_meshsz =list_int(1)
pawrad%mesh_size =list_int(2)
pawrad%mesh_type =list_int(3)
if_rad=list_int(4)
if_radfact=list_int(5)
if_simfact=list_int(6)
end if
LIBPAW_DEALLOCATE(list_int)
!Broadcast the reals
nn=4+pawrad%mesh_size*(if_rad+if_radfact+if_simfact)
LIBPAW_ALLOCATE(list_dpr,(nn))
if(me==0) then
list_dpr(1)=pawrad%lstep
list_dpr(2)=pawrad%rmax
list_dpr(3)=pawrad%rstep
list_dpr(4)=pawrad%stepint
indx=5
if (if_rad==1) then
isz1=pawrad%mesh_size
list_dpr(indx:indx+isz1-1)=pawrad%rad(1:isz1)
indx=indx+isz1
end if
if (if_radfact==1) then
isz1=pawrad%mesh_size
list_dpr(indx:indx+isz1-1)=pawrad%radfact(1:isz1)
indx=indx+isz1
end if
if (if_simfact==1) then
isz1=pawrad%mesh_size
list_dpr(indx:indx+isz1-1)=pawrad%simfact(1:isz1)
indx=indx+isz1
end if
end if
call xpaw_mpi_bcast(list_dpr,0,comm_mpi,ierr)
if(me/=0) then
pawrad%lstep=list_dpr(1)
pawrad%rmax=list_dpr(2)
pawrad%rstep=list_dpr(3)
pawrad%stepint=list_dpr(4)
indx=5
! Deallocate all arrays:
if (allocated(pawrad%rad)) then
LIBPAW_DEALLOCATE(pawrad%rad)
end if
if (allocated(pawrad%radfact)) then
LIBPAW_DEALLOCATE(pawrad%radfact)
end if
if (allocated(pawrad%simfact)) then
LIBPAW_DEALLOCATE(pawrad%simfact)
end if
! Communicate if flag is set to 1:
if(if_rad==1) then
isz1=pawrad%mesh_size
LIBPAW_ALLOCATE(pawrad%rad,(isz1))
pawrad%rad(1:isz1)=list_dpr(indx:indx+isz1-1)
indx=indx+isz1
end if
if(if_radfact==1) then
isz1=pawrad%mesh_size
LIBPAW_ALLOCATE(pawrad%radfact,(isz1))
pawrad%radfact(1:isz1)=list_dpr(indx:indx+isz1-1)
indx=indx+isz1
end if
if(if_simfact==1) then
isz1=pawrad%mesh_size
LIBPAW_ALLOCATE(pawrad%simfact,(isz1))
pawrad%simfact(1:isz1)=list_dpr(indx:indx+isz1-1)
indx=indx+isz1
end if
end if
LIBPAW_DEALLOCATE(list_dpr)
end subroutine pawrad_bcast
!!***
!----------------------------------------------------------------------
!!****f* m_pawrad/simp_gen
!! NAME
!! simp_gen
!!
!! FUNCTION
!! Do integral on a given (generalized) grid using Simpson rule
!!
!! INPUTS
!! radmesh <type(pawrad_type)>=data containing radial grid information
!! func(:)=integrand values
!! r_for_intg=upper boundary for (future) integration over the radial grid
!! (can be negative for an integration over the whole grid)
!!
!! OUTPUT
!! intg=resulting integral by Simpson rule
!!
!! NOTES
!! Possible mesh types (radmesh%mesh_type)
!! 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
!! mesh_type=5 ( grid): rad(i)=AA*i/(n-i)
!!
!! SOURCE
subroutine simp_gen(intg,func,radmesh,r_for_intg)
#if defined HAVE_AVX_SAFE_MODE
!DEC$ NOOPTIMIZE
#endif
!Arguments ------------------------------------
!scalars
real(dp),intent(out) :: intg
real(dp),intent(in),optional :: r_for_intg
type(pawrad_type),intent(in) :: radmesh
!arrays
real(dp),intent(in) :: func(:)
!Local variables-------------------------------
!scalars
integer :: ii,int_meshsz,ir,ir_last,isim,nn
real(dp) :: hh,resid,simp
real(dp),allocatable :: simfact(:)
character(len=500) :: msg
! *************************************************************************
if (present(r_for_intg)) then
if (r_for_intg>0.d0) then
ir=min(pawrad_ifromr(radmesh,r_for_intg),radmesh%mesh_size)
if (ir<radmesh%mesh_size) then
if (abs(radmesh%rad(ir+1)-r_for_intg)<abs(radmesh%rad(ir)-r_for_intg)) ir=ir+1
end if
if (ir>1) then
if (abs(radmesh%rad(ir-1)-r_for_intg)<abs(radmesh%rad(ir)-r_for_intg)) ir=ir-1
end if
int_meshsz=ir
else
int_meshsz=radmesh%mesh_size
end if
if (int_meshsz>radmesh%mesh_size.or.int_meshsz>size(func)) then
write(msg,'(3(a,i4))')"int_meshsz= ",int_meshsz," > mesh_size=",radmesh%mesh_size,&
& ", size(func)=",size(func)
LIBPAW_BUG(msg)
end if
isim=3; if (radmesh%mesh_type==3)isim=4
LIBPAW_ALLOCATE(simfact,(radmesh%mesh_size))
hh=radmesh%stepint/3.d0
simfact(int_meshsz)=hh*radmesh%radfact(int_meshsz)
simfact(1:isim-2)=zero
ir_last=1
do ir=int_meshsz,isim,-2
simfact(ir-1)=4.d0*hh*radmesh%radfact(ir-1)
simfact(ir-2)=2.d0*hh*radmesh%radfact(ir-2)
ir_last=ir-2
end do
simfact(ir_last)=half*simfact(ir_last)
if (int_meshsz<radmesh%mesh_size) simfact(int_meshsz+1:radmesh%mesh_size)=zero
nn=int_meshsz
simp=zero
do ii=1,nn
simp=simp+func(ii)*simfact(ii)
end do
LIBPAW_DEALLOCATE(simfact)
else
if (radmesh%int_meshsz>size(func)) then
write(msg,'(2(a,i4))')"int_meshsz= ",int_meshsz," > size(func)=",size(func)
LIBPAW_BUG(msg)
end if
nn=radmesh%int_meshsz
simp=zero
do ii=1,nn
simp=simp+func(ii)*radmesh%simfact(ii)
end do
end if
resid=zero
if (radmesh%mesh_type==3) then
resid=half*(func(2)+func(1))*(radmesh%rad(2)-radmesh%rad(1))
if (mod(nn,2)==1) resid=resid+radmesh%stepint/3.d0*(1.25d0*func(2)*radmesh%radfact(2) &
& +2.d0*func(3)*radmesh%radfact(3)-0.25d0*func(4)*radmesh%radfact(4))
else if (mod(nn,2)==0) then
resid=radmesh%stepint/3.d0*(1.25d0*func(1)*radmesh%radfact(1)+2.d0*func(2)*radmesh%radfact(2) &
& -0.25d0*func(3)*radmesh%radfact(3))
end if
intg=simp+resid
end subroutine simp_gen
!!***
!----------------------------------------------------------------------
!!****f* m_pawrad/nderiv_gen
!! NAME
!! nderiv_gen
!!
!! FUNCTION
!! Do corrected first (and -if requested- second) derivation on a given (generalized) grid.
!! This routine interfaces nderiv_lin (derivation on a linear grid).
!!
!! INPUTS
!! func(:)=input function
!! radmesh <type(pawrad_type)>=data containing radial grid information
!!
!! OUTPUT
!! der(:)= 1st derivative of input function
!! [der2(:)]= -- optional -- 2nd derivative of input function
!!
!! NOTES
!! Possible mesh types (radmesh%mesh_type)
!! 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
!!
!! SOURCE
subroutine nderiv_gen(der,func,radmesh,der2)
!Arguments ------------------------------------
!scalars
type(pawrad_type),intent(in) :: radmesh
!arrays
real(dp),intent(in) :: func(:)
real(dp),intent(out) :: der(:)
real(dp),optional,intent(out) :: der2(:)
!Local variables-------------------------------
!scalars
integer :: msz
logical :: compute_2der
character(len=500) :: msg
! *************************************************************************
msz=size(func)
if (size(der)/=msz.or.msz>radmesh%mesh_size) then
msg='wrong sizes for in/out arrays!'
LIBPAW_BUG(msg)
end if
compute_2der=(present(der2))
if (radmesh%mesh_type==1) then
call nderiv_lin(radmesh%rstep,func,der,msz,1)
if (compute_2der) then
call nderiv_lin(radmesh%rstep,func,der2,msz,2)
end if
else if (radmesh%mesh_type==2) then
call nderiv_lin(radmesh%lstep,func,der,msz,1)
der(1:msz)=der(1:msz)/radmesh%radfact(1:msz)
if (compute_2der)then
call nderiv_lin(radmesh%lstep,func,der2,msz,2)
der2(1:msz)=(der2(1:msz)/radmesh%radfact(1:msz)-der(1:msz))/radmesh%radfact(1:msz)
end if
else if (radmesh%mesh_type==3) then
call nderiv_lin(radmesh%lstep,func(2:msz),der(2:msz),msz-1,1)
der(2:msz)=der(2:msz)/radmesh%radfact(2:msz)
call pawrad_deducer0(der,msz,radmesh)
if (compute_2der)then
call nderiv_lin(radmesh%lstep,func(2:msz),der2(2:msz),msz-1,2)
der2(2:msz)=(der2(2:msz)/radmesh%radfact(2:msz)-der(2:msz))/radmesh%radfact(2:msz)
call pawrad_deducer0(der2,msz,radmesh)
end if
else if (radmesh%mesh_type==4) then
call nderiv_lin(radmesh%lstep,func,der,msz,1)
der(1:msz)=der(1:msz)/radmesh%radfact(1:msz)
if (compute_2der)then
call nderiv_lin(radmesh%lstep,func,der2,msz,2)
der2(1:msz)=der2(1:msz)/radmesh%radfact(1:msz)**2-der(1:msz)/radmesh%rstep
end if
else if (radmesh%mesh_type==5) then
call nderiv_lin(one,func,der,msz,1)
der(1:msz)=der(1:msz)/(radmesh%radfact(1:msz)*radmesh%rstep)
if (compute_2der)then
call nderiv_lin(one,func,der2,msz,2)
der2(1:msz)=der2(1:msz)/(radmesh%radfact(1:msz)*radmesh%rstep)**2-two*der(1:msz)/&
& (radmesh%rstep+radmesh%rad(1:msz))