-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathwavefunction_tools.f90
1652 lines (1652 loc) · 69.5 KB
/
wavefunction_tools.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
!
! SCID-TDSE: Simple 1-electron atomic TDSE solver
! Copyright (C) 2015-2021 Serguei Patchkovskii, Serguei.Patchkovskii@mbi-berlin.de
!
! This program is free software: you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation, either version 3 of the License, or
! (at your option) any later version.
!
! This program is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! You should have received a copy of the GNU General Public License
! along with this program. If not, see <https://www.gnu.org/licenses/>.
!
!
! Sundry tools for manipulating wavefunction
!
! The routines in this modules are shared-memory parallelized (OpenMP). They
! are aware of a possible distributed-memory parallelization, but rely on the
! upstream routines to set up the distributed-memory calculations.
!
module wavefunction_tools
use accuracy
use constants
use hacks
use timer
use spherical_data
use spherical_data_initialize
use tridiagonal_tools
use node_tools
use sort_tools
use lapack
use math
implicit none
private
public wt_atomic_cache_prefix, wt_iterative_improvement, wt_disable_orthogonalization
public wt_max_solution_iterations, wt_enable_memory_caches
public wt_init_memory_caches
public wt_atomic_solutions, wt_one_atomic_solution
public wt_normalize, wt_energy, wt_dipole
public wt_transition_matrix_elements, wt_tme, wt_tme_file
public wt_adaptive_r_buffer
public wt_reset_lrmax, wt_update_lrmax
public wt_resize
public wt_r2r_scale
public wt_fake_left
public wt_reconstruct_left
public rcsid_wavefunction_tools
!
character(len=clen), save :: rcsid_wavefunction_tools = "$Id: wavefunction_tools.f90,v 1.65 2023/12/23 11:09:20 ps Exp $"
!
integer, parameter :: wt_adaptive_r_buffer = 8 ! Maximum number of points to examine for adaptive nradial determination
character(len=clen), save :: wt_atomic_cache_prefix = ' ' ! Atomic solution for each L will be cached
! here; these quantities are expensive to
! recompute, especially if we run a 2D simulation
! Blank means no caching
logical, save :: wt_iterative_improvement = .true.
! Set to .true. if the atomic eigenstates are to be
! iteratively refined
logical, save :: wt_disable_orthogonalization = .false.
! Set to .true. to skip explicit Gram-Schmidt orthogonalization
! for the atomic solutions
integer(ik), save :: wt_max_solution_iterations = 20_ik
! Max number of iterations in trying to find the eigenvectors
! in wt_one_atomic_solution()
logical, save :: wt_enable_memory_caches(2) = (/ .false., .false. /)
! Set to .true. to enable in-memory caching of eigensolutions (1)
! and right-to-left transformation matrices (2)
! If eigensolutions and transformation matrices exist in the disk
! cache specified by wt_atomic_cache_prefix, they'll be loaded.
logical, save :: wt_tme = .false.
! Calculate transition matrix elements between all bound states
character(len=clen), save :: wt_tme_file = "tme.dat"
! Output file for the transition matrix elements
!
! In-memory caches. These caches can be enormously large; do not activate them
! unless you have plenty of RAM!
!
complex(rk), allocatable, save :: cache_eval(:, :) ! First index: 1..sd_nradial; Second index: 0..sd_lmax
complex(rk), allocatable, save :: cache_evec(:,:,:,:) ! First & second indices: 1..sd_nradial; Third index: 1,2; Fourth index: 0..sd_lmax
complex(rk), allocatable, save :: cache_r2l (:,:, :) ! First & second indices: 1..sd_nradial; Third index: 0..sd_lmax
!
contains
!
subroutine wt_init_memory_caches(verbose)
integer(ik), intent(in) :: verbose
!
integer(ik) :: alloc, lval
real(rk) :: ram
!
if (.not.any(wt_enable_memory_caches)) return
!
call TimerStart('Atomic memory cache')
!
! Eigensolution cache
!
if (wt_enable_memory_caches(1)) then
ram = real(sd_lmax+1,rk)*real(sd_nradial+2*sd_nradial**2,rk)*2._rk*rk_bytes()
write (out,"(/'Allocating ',f0.3,' Gbytes of RAM for the atomic eigenfunctions cache.'/)") ram/1024._rk**3
call flush_wrapper(out)
!
! Cache field-free eigenvalues and eigenfunctions
!
wt_enable_memory_caches(1) = .false. ! To stop wt_atomic_solutions from going weird
allocate (cache_eval(sd_nradial,0:sd_lmax),cache_evec(sd_nradial,sd_nradial,2,0:sd_lmax),stat=alloc)
if (alloc/=0) then
write (out,"('Error ',i0,' allocating memory for eigenfunctions cache')") alloc
stop 'wavefunction_tools%wt_init_memory_caches (1)'
end if
!$omp parallel do default(shared) private(lval)
fill_eigenfunctions_cache: do lval=0,sd_lmax
call wt_atomic_solutions(verbose,lval,cache_eval(:,lval),cache_evec(:,:,:,lval))
end do fill_eigenfunctions_cache
!$omp end parallel do
wt_enable_memory_caches(1) = .true. ! We can start caching
end if
!
! R2L cache
!
if (wt_enable_memory_caches(2)) then
ram = real(sd_lmax+1,rk)*real(sd_nradial**2,rk)*2._rk*rk_bytes()
write (out,"(/'Allocating ',f0.3,' Gbytes of RAM for the R2L cache.'/)") ram/1024._rk**3
call flush_wrapper(out)
!
! Cache R2L transformation matrices
!
wt_enable_memory_caches(2) = .false. ! To stop wr_r2l_transform from going wonky
allocate (cache_r2l(sd_nradial,sd_nradial,0:sd_lmax),stat=alloc)
if (alloc/=0) then
write (out,"('Error ',i0,' allocating memory for R2L cache')") alloc
stop 'wavefunction_tools%wt_init_memory_caches (2)'
end if
!$omp parallel do default(shared) private(lval)
fill_r2l_cache: do lval=0,sd_lmax
call wt_r2l_transform(verbose,lval,cache_r2l(:,:,lval))
end do fill_r2l_cache
!$omp end parallel do
wt_enable_memory_caches(2) = .true. ! We can start caching
end if
!
call TimerStop('Atomic memory cache')
end subroutine wt_init_memory_caches
!
! Calculate a block of atomic solutions. We return the non-zero part of the
! eigenvectors only, instead of the much larger total wavefunction.
!
subroutine wt_atomic_solutions(verbose,lval,eval,evec)
integer(ik), intent(in) :: verbose ! Reporting level
integer(ik), intent(in) :: lval ! Desired angular momentum
complex(rk), intent(out) :: eval(:) ! Eigenvalues of the Hamiltonian matrix block for L=lval
complex(rk), intent(out) :: evec(:,:,:) ! Left (:,:,1) and right (:,:,2) column-eigenvectors
!
character(len=clen) :: cache
integer(ik) :: iu_temp, ios
!
call TimerStart('Atomic solutions')
if (size(eval)/=sd_nradial .or. ubound(evec,1)/=sd_nradial .or. ubound(evec,2)/=sd_nradial .or. ubound(evec,3)/=2 &
.or. lval<0 .or. lval>sd_lmax) then
stop 'wavefunction_tools%wt_atomic_solutions - bad arguments'
end if
!
! First, go to the in-memory cache if it's present
!
if (wt_enable_memory_caches(1)) then
eval = cache_eval(:, lval)
evec = cache_evec(:,:,:,lval)
call TimerStop('Atomic solutions')
return
end if
!
! Try to load from the disk cache next
!
if (wt_atomic_cache_prefix/=' ') then
!
! This function may be called from a parallel region; make sure it
! does not try to access the same file from multiple threads!
!
write (cache,"(a,'-L=',i5.5)") trim(wt_atomic_cache_prefix), lval
!*ibc !$omp critical
open (newunit=iu_temp,action='read',status='old',form='unformatted',iostat=ios,file=trim(cache))
!*ibc !$omp end critical
if (ios==0) then
read (iu_temp,iostat=ios) eval, evec
close (iu_temp)
end if
! Read errors (e.g. because another thread is creating the file right now)
! will cause us to drop through and recompute the eigensolution
if (ios==0) then
call TimerStop('Atomic solutions')
return
end if
!
! Open was unsuccessful; proceed through calculation, and try to save the
! results once we have them
!
end if
!
! No cached result available; compute it from scratch.
!
! Begin by computing the guess eigendecomposition using LAPACK.
! Unfortunately, LAPACK sometimes produces very inaccurate results for general
! complex matrices, so that these eigenvalues/eigenvectors cannot be used as is
!
call wt_atomic_solution_guess(verbose,lval,eval,evec)
!
call wt_atomic_solution_biorthogonalize(evec)
!
if (wt_iterative_improvement) then
call wt_atomic_solution_improve(verbose,lval,eval,evec)
end if
!
call wt_atomic_solution_verify(lval,evec)
!
if (wt_atomic_cache_prefix/=' ') then
!*ibc !$omp critical
open (newunit=iu_temp,action='write',status='new',form='unformatted',iostat=ios,file=trim(cache))
!*ibc !$omp end critical
if (ios/=0) then
write (out,"('Error ',i0,' creating new cache file ',a)") ios, trim(cache)
write (out,"('Skipping wavefunction save for L=',i0,' and continuing')") lval
else
write (iu_temp) eval, evec
close (iu_temp)
end if
end if
call TimerStop('Atomic solutions')
end subroutine wt_atomic_solutions
!
! Compute a single atomic solution using inverse iterations.
! Requires approximate energy of the solution as the starting guess.
!
subroutine wt_one_atomic_solution(verbose,lval,eval,evec)
integer(ik), intent(in) :: verbose ! Reporting level
integer(ik), intent(in) :: lval ! Desired angular momentum
complex(rk), intent(inout) :: eval ! In: guess eigenvalue
! Out: improved eigenvalue
complex(rk), intent(out) :: evec(:,:) ! Left (:,1) and right (:,2) eigenvectors
!
real(rk) :: guess_buf(sd_nradial)
complex(rk) :: eval_new, delta
integer(ik) :: iter
logical :: converged
!
call TimerStart('Atomic solution (single)')
if (ubound(evec,1)/=sd_nradial .or. ubound(evec,2)/=2 .or. lval<0 .or. lval>sd_lmax) then
stop 'wavefunction_tools%wt_one_atomic_solution - bad arguments'
end if
!
! Start with a random real guess. Don't forget to normalize!
!
call random_number(guess_buf)
evec(:,1) = guess_buf
call random_number(guess_buf)
evec(:,2) = guess_buf
call wt_normalize_one_atomic_solution(evec)
!
write (out,"('Solving for atomic eigenvector L=',i0,' E(guess)=',2(g26.12,1x),'using inverse iterations.')") lval, eval
converged = .false.
inverse_iteration_passes: do iter=1,wt_max_solution_iterations
eval_new = eval
call wt_improve_one_atomic_eigenstate(verbose,lval,eval_new,evec)
delta = eval_new - eval
eval = eval_new
if (verbose>=0) then
write (out,"(1x,'Iteration ',i0,' E= ',2(g35.24e3,1x),' change = ',2(g15.6e3,1x))") iter, eval, delta
end if
if (abs(delta)<=100._rk*spacing(abs(eval_new))) then
!
! We are converged, but do couple more iterations for a good measure
!
eval_new = eval
call wt_improve_one_atomic_eigenstate(verbose,lval,eval_new,evec)
call wt_improve_one_atomic_eigenstate(verbose,lval,eval_new,evec)
delta = eval_new - eval
eval = eval_new
if (verbose>=0) then
write (out,"(1x,'Final iteration E= ',2(g35.24e3,1x),' change = ',2(g15.6e3,1x))") eval, delta
end if
converged = .true.
exit inverse_iteration_passes
end if
end do inverse_iteration_passes
!
write (out,"('Final energy ',2(g35.24e3,1x))") eval
if (.not.converged) then
write (out,"('Inverse iterations failed to converge after ',i0,' passes.')") wt_max_solution_iterations
write (out,"('Continuing with possibly unconverged eigenvectors.')")
end if
!
call TimerStop('Atomic solution (single)')
end subroutine wt_one_atomic_solution
!
subroutine wt_atomic_solution_guess(verbose,lval,eval,evec)
integer(ik), intent(in) :: verbose ! Reporting level
integer(ik), intent(in) :: lval ! Desired angular momentum
complex(rk), intent(out) :: eval(:) ! Eigenvalues of the Hamiltonian matrix block for L=lval
complex(rk), intent(out) :: evec(:,:,:) ! Left (:,:,1) and right (:,:,2) eigenvectors
!
real(rk), allocatable :: tmat(:,:)
complex(lrk), allocatable :: evec_guess(:,:,:), eval_guess(:)
real(rk), allocatable :: eval_guess_real(:)
integer(ik) :: order(sd_nradial)
integer(ik) :: alloc, ipt
!
call TimerStart('Atomic solutions: Guess')
allocate (tmat(sd_nradial,sd_nradial),evec_guess(sd_nradial,sd_nradial,2), &
eval_guess(sd_nradial),eval_guess_real(sd_nradial),stat=alloc)
if (alloc/=0) stop 'wavefunction_tools%wt_atomic_solution_guess - no memory for tmat'
!
if (lval==0) then
call sd_expand_implicit_operator(sd_d2n_l0,sd_m2n_l0,sd_m2nf_l0,sd_m2np_l0,tmat)
else
call sd_expand_implicit_operator(sd_d2n_lx,sd_m2n_lx,sd_m2nf_lx,sd_m2np_lx,tmat)
end if
evec_guess(:,:,1) = cmplx((-1._rk/(2._rk*electron_mass))*tmat,kind=lrk)
!
deallocate (tmat)
!
add_potential: do ipt=1,sd_nradial
evec_guess(ipt,ipt,1) = evec_guess(ipt,ipt,1) + real(sd_pottab(ipt,1,lval),kind=lrk)
end do add_potential
if (sd_capped) then
add_cap: do ipt=sd_cap_start,sd_nradial
evec_guess(ipt,ipt,1) = evec_guess(ipt,ipt,1) + cmplx(sd_captab(ipt),kind=lrk)
end do add_cap
end if
!
if (verbose>=2) then
write (out,"(' Largest matrix element of H for L=',i0,' is ',g24.13)") &
lval, maxval(abs(evec_guess(:,:,1)))
write (out,"('Largest deviation from hermiticity for L=',i0,' is ',g24.13)") &
lval, maxval(abs(evec_guess(:,:,1)-conjg(transpose(evec_guess(:,:,1)))))
end if
!
call lapack_geev(sd_nradial,evec_guess,eval_guess)
eval_guess_real = real(eval_guess,kind=rk)
call order_keys(eval_guess_real,order)
eval = eval_guess(order)
!
! LAPACK _GEEV conjugates the left eigenvectors. Let's begin by undoing the damage.
!
evec(:,:,1) = conjg(evec_guess(:,order,1))
evec(:,:,2) = evec_guess(:,order,2)
!
deallocate (evec_guess,eval_guess,eval_guess_real)
call TimerStop('Atomic solutions: Guess')
end subroutine wt_atomic_solution_guess
!
subroutine wt_atomic_solution_biorthogonalize(evec)
complex(rk), intent(inout) :: evec(:,:,:) ! Left (:,:,1) and right (:,:,2) eigenvectors
!
integer(ik) :: imo, jmo
complex(rk) :: scl
! complex(rk) :: fact_l(sd_nradial), fact_r(sd_nradial)
!
if (wt_disable_orthogonalization) return
call TimerStart('Atomic solutions: Orthogonalize')
!
! LAPACK _GEEV routines do not give biorthogonal eigenvectors
! in the present of multiple roots. Although this should not
! occur in our case, it is better to be safe than sorry. Since
! we also do not like LAPACK's normalization convention for the
! eigenvectors (Cartesian norm = 1 separately for the left and
! right vectors), we have to do a bit of work here
!
normalize_eigenvectors: do imo=1,sd_nradial
!
! Enforce bi-orthogonality to all lower eigenvectors,
! using a variation of the Gram-Schmidt process.
! We assume that all lower eigenvectors are already normalized to
! our conventions (L.R = 1)
!
orthogonalize_lower: do jmo=1,imo-1
! Right eigenvector
scl = sum(evec(:,jmo,1)*evec(:,imo,2))
evec(:,imo,2) = evec(:,imo,2) - scl*evec(:,jmo,2)
! Left eigenvector
scl = sum(evec(:,jmo,2)*evec(:,imo,1))
evec(:,imo,1) = evec(:,imo,1) - scl*evec(:,jmo,1)
end do orthogonalize_lower
! The matrix-vector code below should be faster; in fact, it is
! actually (much) slower. Oops.
! fact_l(:imo-1) = matmul(evec(:,imo,2),evec(:,:imo-1,1))
! fact_r(:imo-1) = matmul(evec(:,imo,1),evec(:,:imo-1,2))
! evec(:,imo,1) = evec(:,imo,1) - matmul(evec(:,:imo-1,1),fact_r(:imo-1))
! evec(:,imo,2) = evec(:,imo,2) - matmul(evec(:,:imo-1,2),fact_l(:imo-1))
!
call wt_normalize_one_atomic_solution(evec(:,imo,:))
end do normalize_eigenvectors
call TimerStop('Atomic solutions: Orthogonalize')
end subroutine wt_atomic_solution_biorthogonalize
!
! Normalize right eigenvector to Euclidian norm of 1,
! simultaneously making the largest coefficient positive real.
!
subroutine wt_normalize_one_atomic_solution(vecx)
complex(rk), intent(inout) :: vecx(:,:)
!
real(rk) :: nrm
complex(rk) :: scl
integer(ik) :: imax
!
if (ubound(vecx,1)/=sd_nradial .or. ubound(vecx,2)/=2) then
stop 'wavefunction_tools%wt_normalize_one_atomic_solution - bad arguments'
end if
!
! Normalize right eigenvector. Iterative refinement tends to produce very large eigenvectors,
! with norms which can easily exceed the dynamics range of the underlying real type. We therefore
! have to be very careful here to avoid overflowing our data types.
!
! Start by making the largest element of the right eigenvector equal to 1.
!
imax = maxloc(abs(vecx(:,2)),dim=1)
scl = 1._rk / vecx(imax,2)
vecx(:,2) = scl * vecx(:,2)
!
! Do the same for the left eigenvector
!
imax = maxloc(abs(vecx(:,1)),dim=1)
scl = 1._rk / vecx(imax,1)
vecx(:,1) = scl * vecx(:,1)
!
! Now set the 2-norm of the right eigenvector to 1.
!
nrm = sum(abs(vecx(:,2))**2)
vecx(:,2) = (1._rk/sqrt(nrm)) * vecx(:,2)
!
! Finally, normalize the left eigenvector
!
scl = sum(vecx(:,1)*vecx(:,2))
vecx(:,1) = (1._rk/scl) * vecx(:,1)
!
! Debugging code: check for NaNs. Should not happen, but being paranoid costs us very little
!
! if (any(vecx/=vecx)) then
! write (out,"(/'FATAL: atomic eigenvectors contain NaNs in wt_normalize_one_atomic_solution')")
! write (out,"('Left eigenvector:')")
! write (out,"(6(1x,g24.14))") vecx(:,1)
! write (out,"('Right eigenvector:')")
! write (out,"(6(1x,g24.14))") vecx(:,2)
! stop 'wavefunction_tools%wt_normalize_one_atomic_solution - got NaN'
! end if
end subroutine wt_normalize_one_atomic_solution
!
subroutine wt_atomic_solution_improve(verbose,lval,eval,evec)
integer(ik), intent(in) :: verbose ! Reporting level
integer(ik), intent(in) :: lval ! Desired angular momentum
complex(rk), intent(inout) :: eval(:) ! Eigenvalues of the Hamiltonian matrix block for L=lval
complex(rk), intent(inout) :: evec(:,:,:) ! Left (:,:,1) and right (:,:,2) eigenvectors
!
integer(ik) :: ipass
integer(ik) :: iev
complex(rk) :: eval_initial(size(eval))
!
call TimerStart('Atomic solutions: Improve')
eval_initial = eval
improvement_passes: do ipass=1,1
scan_eigenvalues: do iev=1,sd_nradial
call wt_improve_one_atomic_eigenstate(verbose,lval,eval(iev),evec(:,iev,:))
if (verbose>=2) then
write (out,"('Iterative improvement pass ',i0,' L= ',i0,' I= ',i0,' E= ',2g26.18)") ipass, lval, iev, eval(iev)
end if
end do scan_eigenvalues
call wt_atomic_solution_biorthogonalize(evec)
end do improvement_passes
!
if (verbose>=1) then
eval_initial = eval_initial - eval
iev = maxloc(abs(eval_initial),dim=1)
write (out,"(/'Iterative update of L=',i0,' solutions complete')") lval
write (out,"(' ground-state energy is ',g34.22,1x,g34.22,' change = ',g18.7,1x,g18.7)") eval(1), eval_initial(1)
write (out,"(' most affected eigenvalue is ',g34.22,1x,g34.22,' change = ',g18.7,1x,g18.7)") eval(iev), eval_initial(iev)
end if
!
call TimerStop('Atomic solutions: Improve')
end subroutine wt_atomic_solution_improve
!
! Improve a single eigenstate using inverse iterations
!
subroutine wt_improve_one_atomic_eigenstate(verbose,lval,eval,evec)
integer(ik), intent(in) :: verbose ! Reporting level
integer(ik), intent(in) :: lval ! Desired angular momentum
complex(rk), intent(inout) :: eval ! Eigenvalue to improve
complex(rk), intent(inout) :: evec(:,:) ! Left (:,1) and right (:,2) eigenvectors to improve
!
integer(ik) :: iter_evec, iter_eval, imax
logical :: fail(2)
logical :: upd_nan(3)
real(rk) :: mn (sd_nradial,3), mt (sd_nradial,3)
complex(rk) :: leqn (sd_nradial,3), leqt (sd_nradial,3)
complex(rk) :: leqnf(sd_nradial,m3d_dc_size)
complex(rk) :: leqtf(sd_nradial,m3d_dc_size)
logical :: leqnp(sd_nradial), leqtp(sd_nradial)
complex(rk) :: tmp (sd_nradial)
complex(rk) :: vecx(sd_nradial,2)
complex(rk) :: sclx ! Largest element of vecx(:,2)
complex(rk) :: scr (sd_nradial,m3d_sc_size)
complex(rk) :: delta_m1 ! Correction to the eigenvalue^{-1}
complex(rk) :: corr ! Correction to the eigenvalue
complex(rk) :: eval0 ! Eigenvalue from the previous iteration
!
improve_eval: do iter_eval=1,5
eval0 = eval
call build_right_linear_system(eval,leqn,leqnf,leqnp,mn,fail(2))
call build_left_linear_system (eval,leqt,leqtf,leqtp,mt,fail(1))
if (any(fail)) then
if (iter_eval<=1 .or. verbose>2) then
write (out,"('WARNING: wt_improve_one_atomic_eigenstate: Update iteration ',i0,' failed, leave solutions unchanged')") &
iter_eval
end if
return
end if
!
! There is non-negligible cost to building linear system solutions;
! perform a few eigenvector updates before updating the eigenvalue
!
! There is a possibility of blowing the dynamic range of the floating-
! point representation during the inverse iterations. We need to to be
! a little careful with the results.
!
improve_evec: do iter_evec=1,3
! Inverse iteration: right vector
call m3d_multiply(mn,evec(:,2),tmp)
call m3d_solve(leqn,leqnf,leqnp,tmp,vecx(:,2),scr)
! Inverse iteration: left vector
call m3d_solve(leqt,leqtf,leqtp,evec(:,1),tmp,scr)
call m3d_multiply(mt,tmp,vecx(:,1))
! Adjust the range of vecx(:,1) to prevent overflows
! If we already had an infinity somewhere, it should become a NaN, and trigger an error exit below.
imax = maxloc(abs(vecx(:,1)),dim=1)
sclx = 1._rk / vecx(imax,1)
vecx(:,1) = sclx * vecx(:,1)
! Adjust the range of vecx(:,2) to prevent overflows
imax = maxloc(abs(vecx(:,2)),dim=1)
sclx = 1._rk / vecx(imax,2)
vecx(:,2) = sclx * vecx(:,2)
! Compute correction to the eigenvalue. Don't forget to apply the scale factor later!
delta_m1 = sum(evec(:,1)*vecx(:,2))
corr = sclx/delta_m1
!
call wt_normalize_one_atomic_solution(vecx)
!
! Make sure none of the updated quantities are NaN - this could happen if we have
! insufficient dynamic range.
!
upd_nan(1) = isnan_wrapper(corr)
upd_nan(2) = any(isnan_wrapper(vecx(:,1)))
upd_nan(3) = any(isnan_wrapper(vecx(:,2)))
!
if (verbose>2) then
write (out,"('Update ',i0,',',i0,' upd_nan = ',3l1,' eval = ',2(g32.20,1x),' corr. = ',2(g32.20,1x))") &
iter_eval, iter_evec, upd_nan, eval, corr
end if
!
if (any(upd_nan)) then
if (iter_eval<=1 .or. verbose>2) then
write (out,"('WARNING: wt_improve_one_atomic_eigenstate: Update ',i0,',',i0,' failed, taking last update.')") &
iter_eval, iter_evec
end if
return
end if
! Store the updated eigenvectors
evec(:,:) = vecx(:,:)
! Update eigenvalue, in case we'll need to do an early exit later
eval = eval0 + corr
end do improve_evec
end do improve_eval
!
contains
subroutine build_right_linear_system(eps,leqn,leqnf,leqnp,mn,fail)
complex(rk), intent(in) :: eps ! Current eigenvalue guess
complex(rk), intent(out) :: leqn (:,:) ! Linear system
complex(rk), intent(out) :: leqnf(:,:) ! Factorization of the linear system
logical, intent(out) :: leqnp( :) ! Pivot list for the linear system
real(rk), intent(out) :: mn (:,:) ! Scaling factor for the RHS
logical, intent(out) :: fail
!
complex(rk) :: vtmp( sd_nradial)
!
! Prepare linear system matrices for solving right linear system:
!
! [-(1/(2m))M^{-1}Delta+V-eps] Y = Z.
!
! This is done as:
!
! [-(1/(2m))Delta+M(V-eps)] Y = M Z
!
vtmp(:) = sd_pottab(:,1,lval) - eps
if (sd_capped) then
vtmp(sd_cap_start:) = vtmp(sd_cap_start:) + sd_captab(sd_cap_start:)
end if
if (lval==0) then
call m3d_right_scale(sd_m2n_l0,vtmp,leqn)
leqn(:,:) = leqn + (-1._rk/(2._rk*electron_mass))*sd_d2n_l0
mn (:,:) = sd_m2n_l0
else
call m3d_right_scale(sd_m2n_lx,vtmp,leqn)
leqn(:,:) = leqn + (-1._rk/(2._rk*electron_mass))*sd_d2n_lx
mn (:,:) = sd_m2n_lx
end if
call m3d_decompose(leqn,leqnf,leqnp,fail=fail)
end subroutine build_right_linear_system
!
subroutine build_left_linear_system(eps,leqt,leqtf,leqtp,mt,fail)
complex(rk), intent(in) :: eps ! Current eigenvalue guess
complex(rk), intent(out) :: leqt (:,:) ! Linear system
complex(rk), intent(out) :: leqtf(:,:) ! Factorization of the linear system
logical, intent(out) :: leqtp( :) ! Pivot list for the linear system
real(rk), intent(out) :: mt (:,:) ! Scaling factor for the solution
logical, intent(out) :: fail
!
complex(rk) :: vtmp( sd_nradial)
!
! Prepare linear system matrices for solving left linear system:
!
! [-(1/(2m))Delta^T M^{-T}+V-eps] Y = Z.
!
! This is done as:
!
! [-(1/(2m))Delta^T+(V-eps)M^T] Q = Z
! Y = M^T Q
!
vtmp(:) = sd_pottab(:,1,lval) - eps
if (sd_capped) then
vtmp(sd_cap_start:) = vtmp(sd_cap_start:) + sd_captab(sd_cap_start:)
end if
if (lval==0) then
call m3d_left_scale(vtmp,sd_m2t_l0,leqt)
leqt(:,:) = leqt(:,:) + (-1._rk/(2._rk*electron_mass))*sd_d2t_l0
mt (:,:) = sd_m2t_l0
else
call m3d_left_scale(vtmp,sd_m2t_lx,leqt)
leqt(:,:) = leqt(:,:) + (-1._rk/(2._rk*electron_mass))*sd_d2t_lx
mt (:,:) = sd_m2t_lx
end if
call m3d_decompose(leqt,leqtf,leqtp,fail=fail)
end subroutine build_left_linear_system
end subroutine wt_improve_one_atomic_eigenstate
!
subroutine wt_atomic_solution_verify(lval,evec)
integer(ik), intent(in) :: lval
complex(rk), intent(in) :: evec(:,:,:) ! Left (:,:,1) and right (:,:,2) eigenvectors
!
integer(ik) :: imo, jmo, alloc
real(rk) :: nrm
complex(rk), allocatable :: norm(:,:)
integer(ik) :: warning_count
complex(rk) :: err, err_max
integer(ik), parameter :: max_warnings = 5
!
call TimerStart('Atomic solutions: Verify')
allocate (norm(sd_nradial,sd_nradial),stat=alloc)
if (alloc/=0) stop 'wavefunction_tools%wt_atomic_solution_verify - allocation failed'
!
norm = matmul(transpose(evec(:,:,1)),evec(:,:,2))
warning_count = 0
err_max = 0
test_norm_j: do jmo=1,sd_nradial
test_norm_i: do imo=1,sd_nradial
nrm = 0
if (imo==jmo) nrm = 1
err = norm(imo,jmo)-nrm
if (abs(err)>1e5_rk*spacing(1._rk)) then
warning_count = warning_count + 1
if (abs(err)>abs(err_max)) err_max = err
if (warning_count<=max_warnings) then
write (out,"('WARNING: _GEEV eigenvectors ',i0,',',i0,' L=',i0,' are not biorthogonal. Error =',2(1x,g24.13))") &
imo, jmo, lval, err
end if
end if
end do test_norm_i
end do test_norm_j
!
if (warning_count>max_warnings) then
write (out,"('WARNING: ',i0,' _GEEV biorthogonality warnings for L=',i0,' suppressed. Max error was ',2(1x,g24.13))") &
warning_count-max_warnings, lval, err_max
end if
!
deallocate (norm)
call TimerStop('Atomic solutions: Verify')
end subroutine wt_atomic_solution_verify
!
! Choose wavefunction lmax and nradial. This is only done once; do not bother
! parallizing this routine.
!
subroutine wt_reset_lrmax(wfn)
type(sd_wfn), intent(inout) :: wfn ! Wavefunction to examine
!
integer(ik) :: lval, mval, ir
real(rk) :: max_val_l(0:sd_lmax)
real(rk) :: max_val_r(sd_nradial)
real(rk) :: eps_l, eps_r ! For the initial reset, use tigher tolerances.
! This will stop the wavefunction from growing spuriously later on
! due to the edge-reflection effects
!
wfn%lmax = sd_lmax
wfn%nradial = sd_nradial
if (.not.sd_adaptive) return
!
call TimerStart('WT reset_lrmax')
!
max_val_l = 0
max_val_r = 0
!
! Collect maximum wavefunction values across L channels
!
sense_loop_m: do mval=sd_mmin,sd_mmax
l_sense_loop_l: do lval=abs(mval),sd_lmax
max_val_l(lval) = max(max_val_l(lval),maxval(abs(wfn%wfn(:,:,lval,mval))))
end do l_sense_loop_l
sense_loop_r: do ir=1,sd_nradial
max_val_r(ir) = max(max_val_r(ir),maxval(abs(wfn%wfn(ir,:,abs(mval):sd_lmax,mval))))
end do sense_loop_r
end do sense_loop_m
!
if (sd_adaptive_l) then
eps_l = 1e-2_rk * wfn%sd_tol_l
choose_lmax: do lval=sd_lmax,0,-1
if (max_val_l(lval)>eps_l) then
wfn%lmax = min(lval+1,sd_lmax)
exit choose_lmax
end if
end do choose_lmax
end if
if (sd_adaptive_r) then
eps_r = 1e-2_rk * wfn%sd_tol_r
choose_rmax: do ir=sd_nradial,1,-1
if (max_val_r(ir)>eps_r) then
wfn%nradial = min(ir+1,sd_nradial)
exit choose_rmax
end if
end do choose_rmax
end if
call TimerStop('WT reset_lrmax')
end subroutine wt_reset_lrmax
!
! Update wavefunction lmax if necessary
!
subroutine wt_update_lrmax(verbose,wfn)
integer(ik), intent(in) :: verbose ! Verbosity level
type(sd_wfn), intent(inout) :: wfn ! Wavefunction to examine
!
integer(ik) :: mval, mmin, mmax, ir, lval, old_lmax
integer(ik) :: nrtop, old_nradial
real(rk) :: max_val_r, max_val_l(2), max_t
real(rk) :: max_r_buf(wt_adaptive_r_buffer+1)
!
if (.not.sd_adaptive) return
call TimerStart('WT update_lrmax')
!
! Determine the radial extent to examine
!
nrtop = min(sd_nradial,wfn%nradial+wt_adaptive_r_buffer)
!
! Adjust angular momentum. We may adjust both up and down!
!
if (sd_adaptive_l) then
!
! max_val_l(1)=max_l1 is the largest element for L=lmax
! max_val_l(2)=max_l2 is the largest element for L=lmax-1
!
! Older versions of Gfortran have a bug in max: reduction for
! arrays; we'll work around it by using explicit variables.
!
max_val_l = 0
!
if (nts%dist(wfn%lmax)==0) then
mmin = max(sd_mmin,-wfn%lmax)
mmax = min(sd_mmax, wfn%lmax)
max_t = 0
!$omp parallel do if(mmax>mmin) default(none) &
!$omp& shared(mmin,mmax,sd_nspin,wfn,nrtop,nts) &
!$omp& private(mval) reduction(max:max_t)
l_sense_loop_m1: do mval=mmin,mmax
if (abs(mval)>wfn%lmax) cycle l_sense_loop_m1
max_t = max(max_t,maxval(abs(wfn%wfn(:nrtop,1:sd_nspin,wfn%lmax,mval))))
end do l_sense_loop_m1
!$omp end parallel do
max_val_l(1) = max_t
end if
if (wfn%lmax>0) then
if (nts%dist(wfn%lmax-1)==0) then
mmin = max(sd_mmin,-(wfn%lmax-1))
mmax = min(sd_mmax, (wfn%lmax-1))
max_t = 0
!$omp parallel do if(mmax>mmin) default(none) &
!$omp& shared(mmin,mmax,sd_nspin,wfn,nrtop,nts) &
!$omp& private(mval) reduction(max:max_t)
l_sense_loop_m2: do mval=mmin,mmax
if (abs(mval)>wfn%lmax-1) cycle l_sense_loop_m2
max_t = max(max_t,maxval(abs(wfn%wfn(:nrtop,1:sd_nspin,wfn%lmax-1,mval))))
end do l_sense_loop_m2
!$omp end parallel do
max_val_l(2) = max_t
end if
end if
call nt_max(max_val_l)
!
if (max_val_l(1)>wfn%sd_tol_l) then
!
! Grow maximum L value by 1
!
if (verbose>0 .and. wfn%lmax<sd_lmax) then
write (out,"('wt_update_lmax: L increased from ',i0,' max_val = ',2g24.12)") wfn%lmax, max_val_l
end if
wfn%lmax = min(sd_lmax,wfn%lmax+1_ik)
wfn%lmax_top = max(wfn%lmax_top,wfn%lmax)
end if
!
if (all(max_val_l<=0.1_rk*wfn%sd_tol_l)) then
!
! Shrink maximum L value by 1
!
if (verbose>0 .and. wfn%lmax>0) then
write (out,"('wt_update_lmax: L decreased from ',i0,' max_val = ',2g24.12)") wfn%lmax, max_val_l
end if
old_lmax = wfn%lmax
wfn%lmax = max(0_ik,wfn%lmax-1_ik)
!
! If lmax decreased, clear the wavefunction array: we assume that parts of the
! wavefunction with L>wfn%lmax is rigorously zero - or we'll start picking up
! errors as we grow the array again.
!
if (old_lmax/=wfn%lmax .and. nts%dist(old_lmax)==0) then
mmin = max(sd_mmin,-old_lmax)
mmax = min(sd_mmax, old_lmax)
!$omp parallel do if(mmax>mmin) default(none) &
!$omp& private(mval) shared(mmin,mmax,wfn,sd_nspin,old_lmax)
l_clear_loop_m: do mval=mmin,mmax
if (abs(mval)>old_lmax) cycle l_clear_loop_m
wfn%wfn(:,1:sd_nspin,old_lmax,mval) = 0
end do l_clear_loop_m
!$omp end parallel do
end if
end if
end if
!
! Adjust radial extent
!
if (sd_adaptive_r .and. wfn%nradial<sd_nradial) then
mmin = max(sd_mmin,-wfn%lmax)
mmax = min(sd_mmax, wfn%lmax)
old_nradial = wfn%nradial
r_sense_loop_r: do ir=old_nradial,nrtop
max_val_r = 0
!$omp parallel do collapse(2) default(none) &
!$omp& shared(mmin,mmax,sd_nspin,wfn,nrtop,ir,nts) &
!$omp& private(lval,mval) &
!$omp& reduction(max:max_val_r)
r_sense_loop_m: do mval=mmin,mmax
r_sense_loop_l: do lval=0,wfn%lmax
if (lval<abs(mval)) cycle r_sense_loop_l
if (nts%dist(lval)>0) cycle r_sense_loop_l
max_val_r = max(max_val_r,maxval(abs(wfn%wfn(ir,1:sd_nspin,lval,mval))))
end do r_sense_loop_l
end do r_sense_loop_m
!$omp end parallel do
max_r_buf(ir-old_nradial+1) = max_val_r
end do r_sense_loop_r
call nt_max(max_r_buf)
!
r_adjust_loop: do ir=old_nradial,nrtop
max_val_r = max_r_buf(ir-old_nradial+1)
if (max_val_r>wfn%sd_tol_r) then
wfn%nradial = min(sd_nradial,ir+1)
end if
! write (out,"('wfn%nradial updated to ',i8)") wfn%nradial
end do r_adjust_loop
end if
call TimerStop('WT update_lrmax')
end subroutine wt_update_lrmax
!
! Resize wavefunction arrays. If the wavefunction is truncated, data will be discarded
! for larger values of the indices. If it is grown, new elements will be initialized to
! zero. There is no data conversion or resampling of any kind goind on.
!
! wt_resize() does not need to be aware of distributed-memory parallization at this point.
!
subroutine wt_resize(wfn,nradial)
type(sd_wfn), intent(inout) :: wfn ! Wavefunction to reallocate
integer(ik), intent(in) :: nradial ! New nradial
!
complex(rk), allocatable :: tmp_wfn(:,:,:,:)
integer(ik) :: alloc
integer(ik) :: old_nr, min_nr
integer(ik) :: lv, mv
!
if (nradial==size(wfn%wfn,dim=1)) return
!
call TimerStart('WT resize')
!
old_nr = size(wfn%wfn,dim=1)
min_nr = min(old_nr,nradial)
allocate (tmp_wfn(old_nr,sd_nspin,0:sd_lmax,sd_mmin:sd_mmax),stat=alloc)
if (alloc/=0) stop 'wavefunction_tools%wt_resize - allocate failed (1)'
!
!$omp parallel default(none) &
!$omp& shared(wfn,tmp_wfn,alloc,old_nr,min_nr,sd_nspin,sd_lmax,sd_mmin,sd_mmax,nradial) &
!$omp& private(lv,mv)
!
! Make a copy
!
!$omp do collapse(2)
copy_out_m: do mv=sd_mmin,sd_mmax
copy_out_l: do lv=0,sd_lmax
tmp_wfn(:,:,lv,mv) = wfn%wfn(:,:,lv,mv)
end do copy_out_l
end do copy_out_m
!$omp end do
!
!$omp single
deallocate (wfn%wfn,stat=alloc)
if (alloc/=0) stop 'wavefunction_tools%wt_resize - deallocate failed (1)'
allocate (wfn%wfn(nradial,sd_nspin,0:sd_lmax,sd_mmin:sd_mmax),stat=alloc)
if (alloc/=0) stop 'wavefunction_tools%wt_resize - allocate failed (2)'
!$omp end single
!
!$omp do collapse(2)
copy_in_m: do mv=sd_mmin,sd_mmax
copy_in_l: do lv=0,sd_lmax
wfn%wfn(:min_nr,:,lv,mv) = tmp_wfn(:min_nr,:,lv,mv)
wfn%wfn(min_nr+1:nradial,:,lv,mv) = 0
end do copy_in_l
end do copy_in_m
!$omp end do
!
!$omp end parallel
!
deallocate(tmp_wfn,stat=alloc)
if (alloc/=0) stop 'wavefunction_tools%wt_resize - deallocate failed (1)'
call TimerStop('WT resize')
end subroutine wt_resize
!
! Reset wavefunction norm.
! We will set Cartesian norm of the right wavefunction to 1, while maintaining its phase.
! The left wavefunction will be adjusted so that <l|r> product is also 1.
! Keep in mind that the left wavefunction does not need to be conjugated!
!
subroutine wt_normalize(wfn_l,wfn_r,norm)
type(sd_wfn), intent(inout) :: wfn_l ! Left wavefunction
type(sd_wfn), intent(inout) :: wfn_r ! Right wavefunction
! Since we are dealing with (potentially) non-Hermitian
! operators, our wavefunctions always come as a pair of
! the left and right vectors.
complex(rk), intent(out) :: norm(2) ! norm(1) - Input wavefunction norm <l|r>
! norm(2) - Cartesian norm of the input right wavefunction <r|r>
! norm(2) is guaranteed to be real.
!
integer(ik) :: lval, mval, sval, my_lmax
real(rk) :: sum_rr, scl_r
complex(rk) :: sum_lr, scl_l
!
call TimerStart('WF normalize')
sum_rr = 0
sum_lr = 0
my_lmax = max(wfn_l%lmax,wfn_r%lmax)
!$omp parallel default(none) &
!$omp& shared(sd_mmin,sd_mmax,my_lmax,sd_nspin,wfn_l,wfn_r,nts) &
!$omp& private(mval,lval,sval) &
!$omp& shared(sum_rr,sum_lr,scl_r,scl_l)
sense_loop_m: do mval=sd_mmin,sd_mmax
!$omp do reduction(+:sum_rr,sum_lr)
sense_loop_l: do lval=abs(mval),my_lmax
if (nts%dist(lval)>0) cycle sense_loop_l ! Hook for distributed-memory parallelization
sense_loop_s: do sval=1,sd_nspin
sum_rr = sum_rr + sum(abs(wfn_r%wfn(:,sval,lval,mval))**2)
sum_lr = sum_lr + sum(wfn_l%wfn(:,sval,lval,mval)*wfn_r%wfn(:,sval,lval,mval))
end do sense_loop_s
end do sense_loop_l
!$omp end do nowait
end do sense_loop_m
!$omp barrier
!$omp single
call nt_add(sum_rr)
call nt_add(sum_lr)
scl_r = 1._rk / sqrt(sum_rr)
scl_l = sqrt(sum_rr) / sum_lr
!$omp end single
!$omp barrier
scale_loop_m: do mval=sd_mmin,sd_mmax
!$omp do
scale_loop_l: do lval=abs(mval),my_lmax
if (nts%dist(lval)>0) cycle scale_loop_l ! Hook for distributed-memory parallelization
scale_loop_s: do sval=1,sd_nspin
wfn_r%wfn(:,sval,lval,mval) = scl_r*wfn_r%wfn(:,sval,lval,mval)
wfn_l%wfn(:,sval,lval,mval) = scl_l*wfn_l%wfn(:,sval,lval,mval)
end do scale_loop_s
end do scale_loop_l
!$omp end do nowait
end do scale_loop_m
!$omp end parallel
norm(1) = sum_lr
norm(2) = sum_rr
call TimerStop('WF normalize')
end subroutine wt_normalize
!
! Calculate expectation value of the Hamiltonian