-
Notifications
You must be signed in to change notification settings - Fork 8
/
me_numer.f90
2246 lines (2171 loc) · 67.4 KB
/
me_numer.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
module me_numer
!
! Numerical matrix elements based on the Numerov procedure
!
use accuracy
use molecules
use moltype
use lapack
use timer
implicit none
private
public ik, rk, out
public ME_numerov,simpsonintegral,simpsonintegral_ark,integral_rect_ark
!
real(ark), parameter :: enermax= 700000.0_rk ! Upper limit for the energy search (1/cm)
real(ark), parameter :: enerdev= 0.1_rk ! In the i-th solution has been already defined at
! previouse steps, we solve the eigen-problem anyway
! but within +/- "enerdev" (cm-1) interval from
! the found solution
real(ark), parameter :: enerstep = 100.0_ark ! Energy interval needed for the step-wise search of the
! solution, when the iterative procedure does not help
!
real(ark), parameter :: enersmall= 0.000001_ark ! Small energy value to check whether the search interval not collapsed
!
real(ark), parameter :: thrsh_int = 10.0_ark**(-ark) ! Threshold in Numerov-integration
real(ark), parameter :: thrsh_dif = 0.1e-2_ark ! Threshold in numerical differentiation
real(ark), parameter :: thrsh_upper = 1.0e20_ark ! Upper limit for fncts in out-numerov integr.
!
integer(ik), parameter :: itermax= 500000 ! Maximal possible number of iterations (tries) for Numerov
integer(ik), parameter :: maxslots = 500 ! Maximal possible number of slots for found energies
integer(ik), parameter :: epoints = 40 ! N of points used for extrapolation in lsq_fit
!
real(ark) :: rhostep ! mesh size
!
real(ark),parameter :: wave_init=1.0_ark ! initial arbitrary value at rho = rho_ref
!
!
integer(ik), parameter :: iterlimit=5000 ! Iteration limit in numerov integration
!
integer(ik) :: npoints ! number of grid points
real(ark) :: rho_b(2) ! rhomin..rhomax
integer(ik) :: isingular ! point with singularity, none if <0
integer(ik) :: iperiod ! the periodicity (can be negative for the reflecation)
!
integer(ik) :: verbose = 6 ! Verbosity level
!
logical :: periodic ! periodic boundary conditions: true or false
!
integer(ik) :: vmax,maxorder,imin
!
integer(ik) :: io_slot ! unit numeber to store the numerov eigenvectors and their derivatives
!
integer(ik) :: imode ! the current mode
!
integer(ik) :: Nr = 4 ! 2*Nr+1 is the number of interpolation points
!
integer(ik) :: iparity !
!
character(len=cl),parameter :: deriv_method = 'ML_diffs' ! We use this method to estimate d pvi_v / d rho
! where phi_v is a numerical eigenfunction from numerov
! deriv_method can be either 'd04aaf' of '5 points'
! '5 points'
real(ark) :: rho_switch = .0174532925199432957692369_ark ! the value of abcisse rho of the switch between regions (1 deg)
integer(ik) :: iswitch ! the grid point of switch
!
contains
!
! Matrix elements with Numerov-eigenfunctions
!
subroutine ME_numerov(vmax_,maxorder_,rho_b_,isingular_,npoints_,numerpoints_,drho_,xton_,poten_,mu_rr_,icoord,iperiod_,&
verbose_,g_numerov,energy)
!
integer(ik),intent(in) :: vmax_,maxorder_,npoints_,isingular_,numerpoints_,iperiod_
real(ark),intent(out) :: g_numerov(-1:3,0:maxorder_,0:vmax_,0:vmax_)
real(ark),intent(out) :: energy(0:vmax_)
!
real(ark),intent(in) :: rho_b_(2)
real(ark),intent(in) :: poten_(0:npoints_),mu_rr_(0:npoints_),drho_(0:npoints_,3),xton_(0:npoints_,0:maxorder_)
integer(ik),intent(in) :: icoord ! coordinate number for which the numerov is employed
integer(ik),intent(in) :: verbose_ ! Verbosity level
!
real(ark) :: rho
real(ark) :: h_t,sigma,sigma_t,rms,psipsi_t,characvalue,rhostep_,step_scale,fval,df_t
!
integer(ik) :: vl,vr,lambda,alloc,i,rec_len,k,i_,i1,i2
!
real(ark),allocatable :: phil(:),phir(:),dphil(:),dphir(:),phivphi(:),rho_kinet(:),rho_poten(:),rho_extF(:)
real(ark),allocatable :: f(:),poten(:),mu_rr(:),d2fdr2(:),dfdr(:),rho_(:),xton(:,:)
character(len=cl) :: unitfname
real(ark),allocatable :: enerslot(:),enerslot_(:)
!
if (verbose>=1) write (out,"(/'Numerov matrix elements calculations')")
!
! global variables
!
vmax = vmax_
maxorder = maxorder_
npoints = numerpoints_
rho_b = rho_b_
isingular = isingular_
iperiod = iperiod_
verbose = verbose_
imode=icoord
rho_switch = molec%specparam(icoord)
!
periodic = .false.
if (iperiod>0) periodic = .true.
!
allocate(phil(0:npoints_),phir(0:npoints_),dphil(0:npoints_),dphir(0:npoints_), &
phivphi(0:npoints_),rho_kinet(0:npoints_),rho_poten(0:npoints_),rho_extF(0:npoints_),enerslot(0:maxslots), &
f(0:npoints),dfdr(0:npoints),d2fdr2(0:npoints),poten(0:npoints),mu_rr(0:npoints),&
xton(0:npoints,0:maxorder_),stat=alloc)
if (alloc/=0) then
write (out,"('phi - out of memory')")
stop 'phi - out of memory'
endif
!
! numerov step size
rhostep = (rho_b(2)-rho_b(1))/real(npoints,kind=ark)
!
! real step size
rhostep_ = (rho_b(2)-rho_b(1))/real(npoints_,kind=ark)
!
! define the rho-type coordinate
!
rho_kinet(:) = drho_(:,1)
rho_poten(:) = drho_(:,2)
rho_extF(:) = drho_(:,3)
!
step_scale = rhostep/rhostep_
!
! interpolation
!
if (npoints==npoints_) then
!
poten = poten_
mu_rr = mu_rr_
!
do lambda = 0,maxorder_
!
xton(:,lambda) = xton_(:,lambda)
!
enddo
!
else
!
allocate(rho_(0:npoints_),stat=alloc)
if (alloc/=0) stop 'rho_ - out of memory'
!
forall(i_ = 0:npoints_) rho_(i_) = rho_b(1)+real(i_,kind=ark)*rhostep_
!
do i = 0,npoints
!
rho = rho_b(1)+real(i,kind=ark)*rhostep
!
i_ = int( real(i,ark)*step_scale )
!
i1 = max(0,i_-Nr) ; i2 = min(npoints_,i_+Nr)
!
call polintark(rho_(i1:i2),poten_(i1:i2),rho,fval,df_t)
poten(i) = fval
!
call polintark(rho_(i1:i2),mu_rr_(i1:i2),rho,fval,df_t)
mu_rr(i) = fval
!
do lambda = 0,maxorder
!
call polintark(rho_(i1:i2),xton_(i1:i2,lambda),rho,fval,df_t)
xton(i,lambda) = fval
!
enddo
!
enddo
!
deallocate(rho_)
!
endif
!
! Do some reporting
!
if (verbose>=3) then
write (out,"('vmax = ',i8)") vmax
write (out,"('maxorder = ',i8)") maxorder
write (out,"('icoord = ',i4)") icoord
write (out,"('rho_b (x) = ',2f12.4)") rho_b(1:2) !*180.0_rk/pi
write (out,"('rhostep (x) = ',2f12.4)") rhostep !*180.0_rk/pi
endif
!
call diff_2d_4points_ark(npoints,rho_b,mu_rr,periodic,0_ik,dfdr,d2fdr2)
!
f(:) = (d2fdr2(:)/mu_rr(:)-dfdr(:)**2/mu_rr(:)**2*0.5_ark)*0.5_ark
!
! Print out
!
if (verbose>=3) then
write(out,"('grid values (i,rho,rho_kinet,rho_poten,poten, mu_rr, f): ')")
do i_=0,npoints_,2
i = int( real(i_,ark)/step_scale )
rho = rho_b(1)+real(i_,kind=ark)*rhostep_
write(out,"(i8,3f14.6,3g14.6)") i_,rho,rho_kinet(i_),rho_poten(i_),poten(i),mu_rr(i),f(i)
enddo
endif
!
inquire(iolength=rec_len) phil(:),dphil(:)
!
write(unitfname,"('Numerov basis set # ',i6)") icoord
call IOStart(trim(unitfname),io_slot)
!
open(unit=io_slot,status='scratch',access='direct',recl=rec_len)
!
! solve 1D Schroedinger equation with Numerov algorithm
!
iparity = 0
!
iswitch = max(1,int( rho_switch/rhostep_))
!
if (isingular>-1.and.(iswitch<=0.or.iswitch>=npoints_)) then
write(out,"('me_numerov: illegal iswitch or npoints',2i8)") iswitch,npoints_
stop 'iswitch inconsistent with npoints'
endif
!
if (iperiod/=0) vmax = vmax/2
!
call numerov(npoints,npoints_,step_scale,poten,mu_rr,f,enerslot)
!
if (iperiod/=0) then
allocate(enerslot_(0:maxslots),stat=alloc)
if (alloc/=0) then
write (out,"('phi - out of memory')")
stop 'phi - out of memory'
endif
!
iparity = 1
!
call numerov(npoints,npoints_,step_scale,poten,mu_rr,f,enerslot_)
!
do vl = vmax,0,-1
!
enerslot(2*vl ) = enerslot(vl)
if (2*vl+1<=vmax_) enerslot(2*vl+1) = enerslot_(vl)
!
enddo
!
vmax = vmax_
!
deallocate(enerslot_)
!
endif
!
! Matrix elements
!
sigma = 0.0_ark
rms = 0.0_ark
characvalue = maxval(enerslot(0:vmax))
energy(0:vmax) = enerslot(0:vmax)-enerslot(0)
!
do vl = 0,vmax
!
read (io_slot,rec=vl+1) (phil(i),i=0,npoints_),(dphil(i),i=0,npoints_)
!
do vr = vl,vmax
!
if (iperiod/=0.and.mod(abs(vl-vr),2)==1) cycle
!
if (vl==vr) then
phir = phil
dphir = dphil
else
read (io_slot,rec=vr+1) (phir(i),i=0,npoints_),(dphir(i),i=0,npoints_)
endif
!
! Here we prepare integrals of the potential
! <vl|poten|vr> and use to check the solution of the Schroedinger eq-n
! obtained above by the Numerov
!
phivphi(:) = phil(:)*poten_(:)*phir(:)
!
h_t = simpsonintegral_ark(npoints_,rho_b(2)-rho_b(1),phivphi)
!
! momenta-quadratic part
!
phivphi(:) =-dphil(:)*mu_rr_(:)*dphir(:)
!
psipsi_t = simpsonintegral_ark(npoints_,rho_b(2)-rho_b(1),phivphi)
!
! Add the diagonal kinetic part to the tested mat. elem-s
!
h_t = h_t - 0.5_ark*psipsi_t
!
psipsi_t = 0
!
do lambda = 0,maxorder
!
! momenta-free part in potential part
!
if (lambda==0) then
phivphi(:) = phil(:)*phir(:)
else
phivphi(:) = phil(:)*rho_poten(:)**lambda*phir(:)
endif
!
g_numerov(0,lambda,vl,vr) = simpsonintegral_ark(npoints_,rho_b(2)-rho_b(1),phivphi)
!
! external field expansion
!
if (lambda==0) then
phivphi(:) = phil(:)*phir(:)
else
phivphi(:) = phil(:)*rho_extF(:)**lambda*phir(:)
endif
!
g_numerov(3,lambda,vl,vr) = simpsonintegral_ark(npoints_,rho_b(2)-rho_b(1),phivphi)
if (vl/=vr) g_numerov(3,lambda,vr,vl) = g_numerov(3,lambda,vl,vr)
!
! momenta-free in kinetic part
!
if (lambda==0) then
phivphi(:) = phil(:)*phir(:)
else
phivphi(:) = phil(:)*rho_kinet(:)**lambda*phir(:)
endif
!
phivphi(:) = phil(:)*xton(:,lambda)*phir(:)
!
g_numerov(-1,lambda,vl,vr) = simpsonintegral_ark(npoints_,rho_b(2)-rho_b(1),phivphi)
!
! We also control the orthogonality of the basis set
!
if (lambda==0) psipsi_t = g_numerov(0,lambda,vl,vr)
!
if (vl/=vr) g_numerov(-1:0,lambda,vr,vl) = g_numerov(-1:0,lambda,vl,vr)
!
! momenta-quadratic part
!
if (lambda==0) then
phivphi(:) =-dphil(:)*dphir(:)
else
phivphi(:) =-dphil(:)*rho_kinet(:)**lambda*dphir(:)
endif
!
phivphi(:) =-dphil(:)*xton(:,lambda)*dphir(:)
!
g_numerov(2,lambda,vl,vr) = simpsonintegral_ark(npoints_,rho_b(2)-rho_b(1),phivphi)
!
if (vl/=vr) g_numerov(2,lambda,vr,vl) = g_numerov(2,lambda,vl,vr)
!
! momenta-linear part:
! < vl | d/dx g(x) | vr > = - < vr | g(x) d/dx | vl >
!
!
if (lambda==0) then
phivphi(:) = phil(:)*dphir(:)
else
phivphi(:) = phil(:)*rho_kinet(:)**lambda*dphir(:)
endif
!
phivphi(:) = phil(:)*xton(:,lambda)*dphir(:)
!
g_numerov(1,lambda,vl,vr) = simpsonintegral_ark(npoints_,rho_b(2)-rho_b(1),phivphi)
!
if (vl/=vr) then
!
if (lambda==0) then
phivphi(:) = dphil(:)*phir(:)
else
phivphi(:) = dphil(:)*rho_kinet(:)**lambda*phir(:)
endif
!
phivphi(:) = dphil(:)*xton(:,lambda)*phir(:)
!
g_numerov(1,lambda,vr,vl) = simpsonintegral_ark(npoints_,rho_b(2)-rho_b(1),phivphi)
!
endif
!
!
if (verbose>=7) then
write(out,"('g_numerov(0,',i4,i4,i4,') = ',f18.8)") lambda,vl,vr,g_numerov(0,lambda,vl,vr)
write(out,"('g_numerov(1,',i4,i4,i4,') = ',f18.8)") lambda,vl,vr,g_numerov(1,lambda,vl,vr)
write(out,"('g_numerov(2,',i4,i4,i4,') = ',f18.8)") lambda,vl,vr,g_numerov(2,lambda,vl,vr)
write(out,"('g_numerov(3,',i4,i4,i4,') = ',f18.8)") lambda,vl,vr,g_numerov(3,lambda,vl,vr)
if (vl/=vr) then
write(out,"('g_numerov(0,',i4,i4,i4,') = ',f18.8)") lambda,vr,vl,g_numerov(0,lambda,vr,vl)
write(out,"('g_numerov(1,',i4,i4,i4,') = ',f18.8)") lambda,vr,vl,g_numerov(1,lambda,vr,vl)
write(out,"('g_numerov(2,',i4,i4,i4,') = ',f18.8)") lambda,vr,vl,g_numerov(2,lambda,vr,vl)
write(out,"('g_numerov(3,',i4,i4,i4,') = ',f18.8)") lambda,vr,vl,g_numerov(3,lambda,vr,vl)
endif
endif
!
enddo
!
! Count the error, as a maximal deviation sigma = | <i|H|j>-E delta_ij |
!
sigma_t = abs(h_t)
if (vl==vr) sigma_t = abs(h_t-enerslot(vl))
!
sigma = max(sigma,sigma_t)
rms = rms + sigma_t**2
!
! Now we test the h_t = <vl|h|vr> matrix elements and check if Numerov cracked
! the Schroedinger all right
if (vl/=vr.and.abs(h_t)>sqrt(small_)*abs(characvalue)*1e4) then
write(out,"('ME_numerov: wrong Numerovs solution for <',i4,'|H|',i4,'> = ',f20.10)") vl,vr,h_t
stop 'ME_numerov: bad Numerov solution'
endif
!
if (vl==vr.and.abs(h_t-enerslot(vl))>sqrt(small_)*abs(characvalue)*1e4) then
write(out,"('ME_numerov: wrong <',i4,'|H|',i4,'> (',f16.6,') =/= energy (',f16.6,')')") vl,vr,h_t,enerslot(vl)
stop 'ME_numerov: bad Numerov solution'
endif
!
! Reporting the quality of the matrix elemenst
!
if (verbose>=3) then
if (vl/=vr) then
write(out,"('<',i4,'|H|',i4,'> = ',e16.2,'<-',8x,'0.0',5x,'; <',i4,'|',i4,'> = ',e16.2,'<-',8x,'0.0')") &
vl,vr,h_t,vl,vr,psipsi_t
else
write(out,"('<',i4,'|H|',i4,'> = ',f16.6,'<-',f16.6,'; <',i4,'|',i4,'> = ',f16.6)")&
vl,vr,h_t,enerslot(vl),vl,vr,psipsi_t
endif
endif
!
enddo
enddo
!
rms = sqrt(rms/real((vmax+1)*(vmax+2)/2,kind=ark))
!
if (verbose>=1) then
write(out,"('Maximal deviation sigma = | <i|H|j>-E delta_ij | is ',f18.8)") sigma
write(out,"('rms deviation is ',f18.8)") sqrt(rms)
endif
!
deallocate(phil,phir,dphil,dphir,phivphi,rho_kinet,rho_poten,rho_extF,enerslot,f,poten,mu_rr,d2fdr2,dfdr,xton)
!
!
end subroutine ME_numerov
!
!********************************************************************!
! !
! numerov calculates bending functions and stores them !
! !
! !
!********************************************************************!
subroutine numerov(npoints,npoints_,step_scale,poten,mu_rr,f,enerslot)
!
integer(ik),intent(in) :: npoints,npoints_
real(ark),intent(in) :: step_scale
real(ark),intent(in) :: poten(0:npoints),mu_rr(0:npoints),f(0:npoints)
!
real(ark),intent(out) :: enerslot(0:maxslots)
integer(ik) :: v,i_
integer(ik) :: ierr,numnod,ntries,vrec,k0,ipoint,mtries
!
real(ark) :: eold,tsum
!
real(ark) :: pcout
!
real(ark),allocatable :: pot_eff(:),i0(:),phi_f(:),phi_t(:),psi_t(:),phi_f_(:),phi_d_(:)
!
real(ark) :: eguess,enerlow,enerupp,enerdelta,enershift,enermid
!
!real(rk) :: simpsonintegral
!
logical :: notfound
!
real(ark) :: potmin,oldphi,newphi,v_t(-2:2),der(4),erest(4),v_tt(-4:4)
!
real(ark) :: r_t(1:4),r_tt,f_t(1:4),f_tt,df_t1,df_t2
!
integer(ik) :: alloc,i,j,jp2,iter,nm2,ic,imin_l,imin_r
integer(ik) :: icslots(0:maxslots),iright,ileft,ireflect
allocate(pot_eff(0:npoints),i0(0:npoints),phi_f(0:npoints),phi_t(0:npoints),psi_t(-30:npoints+30),phi_f_(0:npoints),&
phi_d_(0:npoints),stat=alloc)
if (alloc/=0) then
write (6,"('numerov: allocation is faild - out of memory')")
stop 'numerov: allocation is faild - out of memory'
endif
!fcoef = planck*avogno*1.0d+16/( 4.0d+00*pi*pi*vellgt )
!
! Check for the potential "1/0" probelm in mu_rr
!
do i=0,npoints
if ( abs(mu_rr(i))<small_*100.0_rk ) then
write(out,"('numerov: mu_rr is zero at i = ',i8)") i
stop 'numerov: mu_rr has a zero element'
endif
enddo
!
! before commencing the numerov-cooley numerical integration we
! store the following function in the f1 array:
!
! 0
! w = f (p ) + 2 i v (p )
! i 1 i 0 i
!
! in the renner-teller version of morbid, this function is
! extended with a term originating in the non-zero
! electronic angular momentum.
!
! -1 0 -1
! v is given in cm , and i is given in cm
! 0
!
! calculate effective potential w_i
!
pot_eff(:) = f(:) + 2.0_ark*( poten(:) )/mu_rr(:)
!
! Effective inertia moment
!
i0(:) = 2.0_ark/mu_rr(:)
!
! find minimum : position
!
potmin = safe_max
!
do i=0,npoints
!
if (pot_eff(i)<potmin) then
imin = i
potmin = pot_eff(i)
endif
!
enddo
! and value
!
potmin = pot_eff(imin)
!
potmin = max(0.0_ark,potmin)
!
!if (imin==0) imin = npoints/2
!
!imin = minloc(pot_eff(:),dim=1)
!
if (imin<0.or.imin>npoints) then
write(out,"('numerov: pot_eff has no minimum',i8)")
stop 'numerov: pot_eff has no minimum'
endif
!
imin_l = minloc(pot_eff(0:imin-1),dim=1)
imin_r = minloc(pot_eff(imin+1:npoints),dim=1)+imin
!
!
ileft = npoints
iright = 0
!
! Print out
!
if (verbose>=5) then
write(out,"('grid values (poten, mu_rr, poten_eff): ')")
do i=0,npoints,2
write(out,"(i8,3f18.8)") i,poten(i),mu_rr(i),pot_eff(i)
enddo
write(out,"('potmin(eff) = ',f18.8,' at i = ',i8)") potmin,imin
endif
!
! assign initial values
!
enerslot(0:maxslots)=poten(imin)-sqrt(safe_max)
icslots(0:maxslots) = npoints
!
! Parameters for the Simson rule integration
!
! facodd=2.0_rk*rhostep/3.0_rk
! faceve=2.0_rk*facodd
!
! Here we start solving Schroedinger equation for the v-th eigenvalue
!
enerlow=potmin/2.0_ark*mu_rr(imin)
enerupp=enermax
!
do v=0,vmax+1
!
if (v/=0) enerlow = enerslot(v-1)
!
! Determine the highest undefined energy slot
i = v
do while (enerslot(i)<poten(imin).and.i<maxslots)
i = i + 1
enddo
!
if (enerslot(i)>=poten(imin)) then
enerupp=enerslot(i)
else
enerupp=enermax
endif
!
if (enerslot(v)>poten(imin)) then
!
enerlow=max(enerslot(v)-enerdev,poten(imin))
enerupp=min(enerslot(v)+enerdev,enermax)
ic = icslots(v)
!
endif
!
enerdelta= enerstep
enershift = enerstep
ntries=0
mtries = 0
!
! start itererative procedure for searching the current solution
!
iter = 0
!
notfound = .true.
!
search_loop : do while(notfound.and.iter<itermax)
!
iter = iter + 1
!
enermid=(enerlow+enerupp)/2.0
!
!enerdelta = enerdelta*0.5_ark
!
eguess=enermid
!
eold=eguess
!
! Numerov procedure
!
ic = icslots(v)
!
if (ic/=npoints.and.enerslot(v)>=poten(imin).and.iter==1) eguess = enerslot(v)
!
if (verbose>=5) write(out,"('eguess = ',e14.7)") eguess
!
! special case of a singular point i = 0
!
if (isingular==0) then
!
forall(i = 0:npoints) phi_f(i) = sqrt(rhostep*real(i,ark))!*sqrt(mu_rr(i))
!
endif
!
call numcoo ( v, pot_eff, i0, eguess, enerlow, ic, phi_f, pcout, ierr)
!
if (ierr>1) then
!
write(out,"('numerov: no solution found in numcoo, ierr = ',i8)") ierr
!
do i = 0,maxslots
if (enerslot(i)>poten(imin)) write (out,"(' v,ener = ',i8,f20.10)") i,enerslot(i)
enddo
!
stop 'numerov: no solution found in numcoo'
!
endif
!
!
if (ierr == 0) then
!if (trim(molec%coords_transform)=='R-EXPRHO'.and.imode==molec%Nmodes) then
! numnod=0
! nm2=npoints-2
! oldphi=phi_f(0)*exp((rho_b(1)+rhostep*real(0,rk)))
! newphi=phi_f(2)*exp((rho_b(1)+rhostep*real(2,rk)))
! do i=2,nm2
!
! if (oldphi*phi_f(i)*exp((rho_b(1)+rhostep*real(i,rk)))<0.0_ark .or. &
! ( abs(phi_f(i)*exp((rho_b(1)+rhostep*real(i,rk))))==0.0_ark .and. newphi*oldphi<0.0_ark) ) then
! if (exp( (rho_b(1)+rhostep*real(i,rk)) )>1e-4) numnod=numnod+1
! endif
! !
! newphi=phi_f(i+2)*exp((rho_b(1)+rhostep*real(i+2,rk)))
! oldphi=phi_f(i)*exp((rho_b(1)+rhostep*real(i,rk)))
! enddo
!else
numnod=0
nm2=npoints-2
oldphi=phi_f(0)
newphi=phi_f(2)
do i=2,nm2
if ( oldphi*phi_f(i)<0.0 .or. &
( phi_f(i)==0.0_ark .and. newphi*oldphi<0.0_ark) ) numnod=numnod+1
!
newphi=phi_f(i+2)
oldphi=phi_f(i)
enddo
!
endif
!
if (ierr==0.and.periodic) then
!
phi_t(:) = phi_f(:)/sqrt(mu_rr(:))
!
tsum = simpsonintegral_ark(npoints,rho_b(2)-rho_b(1),phi_t(:)**2)
!
phi_t(:)=phi_t(:)/sqrt(tsum)
!
ireflect = 0
!
if (iparity/=0) ireflect = -1
!
call diff_2d_4points_ark(Npoints,rho_b,phi_t,periodic,ireflect,psi_t(0:Npoints))
!
do i=2,5
r_t(i-1) = rho_b(1)+real(i,kind=ark)*rhostep
enddo
r_tt = rho_b(1)
f_t(1:4) = psi_t(2:5)
call polintark(r_t(1:4),f_t,r_tt,df_t1,f_tt)
!
do i=1,4
r_t(i) = rho_b(1)+real(npoints-6+i,kind=ark)*rhostep
enddo
r_tt = rho_b(2)
f_t(1:4) = psi_t(npoints-5:npoints-2)
call polintark(r_t(1:4),f_t,r_tt,df_t2,f_tt)
!
if (abs(phi_t(0)-(-1.0_ark)**iparity*phi_t(npoints ))>sqrt(thrsh_int).or.&
abs(df_t1-(-1.0_ark)**iparity*df_t2)>1000.*sqrt(thrsh_int)) then
!
if (verbose>=5) write(out,"(/'f(0),f(1),.,f(N-1),f(N),f(0)-tau*f(N),energy: ',5g18.8,', energy= ',f12.4,',n=',i9)")&
phi_t(0),phi_t(npoints),phi_t(1),phi_t(npoints-1),&
phi_t(0)-(-1.0_rk)**iparity*phi_t(npoints ),eguess,numnod
if (verbose>=5) write(out,"(/'df(0)-df(N): ',3g18.8)") df_t1,df_t2,df_t1-(-1.0_ark)**iparity*df_t2
!
if ( numnod+1<maxslots.and.numnod>v) then
enerslot(numnod)=eguess
icslots(numnod) = ic
endif
!
numnod = maxslots+1
!
endif
endif
!
!
if (ierr==0.and.iperiod<0) then
!
phi_t(:) = phi_f(:)/sqrt(mu_rr(:))
!
tsum = simpsonintegral_ark(npoints,rho_b(2)-rho_b(1),phi_t(:)**2)
!
phi_t(:)=phi_t(:)/sqrt(tsum)
!
ireflect = 1
if (iparity/=0) ireflect = -1
!
call diff_2d_4points_ark(Npoints,rho_b,phi_t,periodic,ireflect,psi_t(0:Npoints))
!
if ( (iparity==0.and.abs(psi_t(npoints ))>sqrt(thrsh_int)).or.&
(iparity==1.and.abs(phi_t(npoints ))>sqrt(thrsh_int))) then
!
if (verbose>=5) write(out,"(/'phi(N) and psi(N) = ',2g18.8,', energy = ',f12.4,', n = ',9i7)") &
phi_t(npoints),psi_t(npoints ),eguess,numnod
!
if ( numnod+1<maxslots.and.numnod>v) then
enerslot(numnod)=eguess
icslots(numnod) = ic
endif
!
numnod = maxslots+1
!
endif
!
endif
!
!if (iperiod<0) then
! !
! if (abs(phi_f(npoints))<small_.and.mod(v,2)/=0) then
! !
! numnod = numnod*2+1
! !
! else
! !
! numnod = numnod*2
! !
! endif
! !
!endif
!
if ( ierr==0.and.numnod+1<maxslots.and.numnod>=v ) then
enerslot(numnod)=eguess
icslots(numnod) = ic
if (verbose>=6) then
write (out,"('v,numnod,ener = ',2i8,f20.10)") v,numnod,enerslot(numnod)
!
do i=0,npoints
!
write(out,"(i8,2f18.8)") i,phi_f(i),exp((rho_b(1)+rhostep*real(i,rk)))
!
enddo
!
endif
endif
!
if (verbose>=5) write(out,"('v = ',i5,'; numnod = ',i8,', efound = ',f16.7,', ierr = ',i9)") v,numnod,eguess,ierr
!
if ( ierr == 0.and.v==numnod ) notfound = .false.
!
! if it cannot find a solution by Numerov procedure and ierr = 1, we can still try to
! devide the searching interval and repeate
! the Numerov procedure with a closer initial guess
! We do the same thing if the found solution is not what we need (/=v)
!
if (abs(enerupp-enerlow) < enersmall.and.(ierr == 1 .or. numnod /= v)) then
!
ntries=0
!
!enerdelta = enerdelta*0.5_ark
!
! Are we still at the begining?
!
if (v.eq.0) then
!
enerlow=poten(imin)
!
!elseif (iperiod<0) then
! !
! enerlow=poten(imin)
! !
else
!
enerlow=enerslot(v-1)-0.001
!
!enerlow=poten(imin)
endif
!
! Determine the highest undefined energy slot
i = v+1
do while (enerslot(i)<poten(imin).and.i<maxslots)
i = i + 1
enddo
!
if (enerslot(i)>=poten(imin)) then
enerupp=enerslot(i)
else
enerupp=enermax
endif
!
if (enermid>=enerupp) then
enerdelta = enerdelta*0.5_ark
enershift = 0.0_ark
endif
!
enerlow=enerlow+enershift
enershift=enershift+enerdelta
enerupp=enerlow
!
cycle search_loop
endif
!
if (ierr <=1 .and. numnod /= v) then
!
! We need to do something if the seacrh interval is collapsed
! this must have happend at ierr=1, and no solution found
!
if (mtries==0) then
!
mtries = 1
!
! Are we still at the begining?
!
if (v.eq.0) then
enerlow=poten(imin)
else
enerlow=enerslot(v-1)
!enerlow=poten(imin)
endif
!
! Determine the highest undefined energy slot
i = v+1
do while (enerslot(i)<poten(imin).and.i<maxslots)
i = i + 1
enddo
!
if (enerslot(i)>=poten(imin)) then
enerupp=enerslot(i)
else
enerupp=enermax
endif
!
endif
!
if (numnod<v.or.ierr==1) then
enerlow=enermid
else
enerupp=enermid
endif
!
endif
!
enddo search_loop
!
if( notfound ) then
!
write(out,"('numerov: no solution found after ',i8,' iterations')") itermax
write (out,"(' v = ',i8,', eguess= ',e20.10)") v,eguess
!
do i = 0,maxslots
if (enerslot(i)>poten(imin)) write (out,"(' v,ener = ',i8,f20.10)") i,enerslot(i)
enddo
!
stop 'numerov: no solution found'
!
endif
!
if (ierr==0.and.periodic) then
if (verbose>=5) write(out,"(/'phi_t(0)-phi_t(npoints ) : ',3g18.8)") phi_t(0),phi_t(npoints ),phi_t(0)-phi_t(npoints )
if (verbose>=5) write(out,"(/'dphi_t(0)-dphi_t(npoints) : ',5g18.8)") phi_t(1),phi_t(npoints-1),df_t1,df_t2,df_t1-df_t2
endif
!
! multiply the wavefunction with sqrt(irr) (eq. (6.4) of jensen)
!
phi_f(:) = phi_f(:)/sqrt(mu_rr(:))
!
phi_t(:) = phi_f(:)*phi_f(:)
!
!if (trim(molec%coords_transform)=='R-RHO'.and.imode==molec%Nmodes) then
! !
! do i = 0,Npoints
! !
! phi_t(i) = phi_t(i)*sqrt(rho_b(1)+rhostep*real(i,rk) )**1
! !
! enddo
! !
!endif
!
! normalize the wavefunction using the Simpson rule integration
!
!tsum = sum(phi_f(:)*phi_f(:):0:2)*facodd + sum(phi_f(:)*phi_f(:):0:2)*faceve
!
!do i=0,npoints-1,2
! tsum=tsum+facodd*phi_f(i)*phi_f(i)+faceve*phi_f(i+1)*phi_f(i+1)
!enddo ! --- i
!
! numerical intagration with simpson's rule #2
!
tsum = simpsonintegral_ark(npoints,rho_b(2)-rho_b(1),phi_t)
!
phi_f(:)=phi_f(:)/sqrt(tsum)
!
!if (iperiod<0) phi_f(:)=phi_f(:)/sqrt(2.0_ark)
!
!
if (verbose>=5) then
write (out,"('v,ener = ',i8,f20.10)") v,enerslot(numnod)
endif
!
if (verbose>=6) then
write(out,"(f18.8)") phi_f
endif
!
! direvatives of the wavefunctions method I
!
if(deriv_method == '5 points') then
!
do i=0,npoints ! --- np
!
do k0 = -2, 2
!
if (periodic) then
ipoint = min(max(i+k0,0),npoints)
v_t(k0) = phi_f(ipoint)
else
if (i+k0<0) then
v_t(k0) = phi_f(-(i+k0))
elseif(i+k0>npoints) then
v_t(k0) = phi_f(npoints)
else
v_t(k0) = phi_f(i+k0)
endif
endif
!
enddo
!
! 5-points expression
!
phi_t(i) = (-v_t( 2)/12.0_ark+2.0_ark/3.0_ark*v_t(1) &
+v_t(-2)/12.0_ark-2.0_ark/3.0_ark*v_t(-1) )/rhostep
!
!phi_t(i) = (v_t(1)-v_t(-1) )/(rhostep*2.0_ark)
!
enddo
!
! direvatives of the wavefunctions method II
!
elseif(deriv_method == 'ML_diffs') then
!
if (iperiod<0) then
!
ireflect = 1
if (iparity/=0) ireflect = -1