-
Notifications
You must be signed in to change notification settings - Fork 145
/
Copy pathparameters.F90
6230 lines (5520 loc) · 272 KB
/
parameters.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
!-*- mode: F90 -*-!
!------------------------------------------------------------!
! This file is distributed as part of the Wannier90 code and !
! under the terms of the GNU General Public License. See the !
! file `LICENSE' in the root directory of the Wannier90 !
! distribution, or http://www.gnu.org/copyleft/gpl.txt !
! !
! The webpage of the Wannier90 code is www.wannier.org !
! !
! The Wannier90 code is hosted on GitHub: !
! !
! https://github.com/wannier-developers/wannier90 !
!------------------------------------------------------------!
module w90_parameters
!! This module contains parameters to control the actions of wannier90.
!! Also routines to read the parameters and write them out again.
use w90_constants, only: dp
use w90_io, only: stdout, maxlen
implicit none
private
!Input
integer, public, save :: iprint
!! Controls the verbosity of the output
character(len=20), public, save :: energy_unit
!! Units for energy
character(len=20), public, save :: length_unit
!! Units for length
logical, public, save :: wvfn_formatted
!! Read the wvfn from fortran formatted file
logical, public, save :: spn_formatted
!! Read the spin from fortran formatted file
logical, public, save :: uHu_formatted
logical, public, save :: berry_uHu_formatted
!! Read the uHu from fortran formatted file
integer, public, save :: spin
!! Spin up=1 down=2
integer, public, save :: num_bands
!! Number of bands
integer, public, save :: num_dump_cycles
!! Number of steps before writing checkpoint
integer, public, save :: num_print_cycles
!! Number of steps between writing output
integer, public, save :: slwf_num
!! Number of objective Wannier functions (others excluded from spread functional)
logical, public, save :: selective_loc
!! Selective localization
logical, public, save :: slwf_constrain
!! Constrained centres
real(kind=dp), allocatable, public, save :: ccentres_frac(:, :)
real(kind=dp), allocatable, public, save :: ccentres_cart(:, :)
real(kind=dp), public, save :: slwf_lambda
!! Centre constraints for each Wannier function. Co-ordinates of centre constraint defaults
!! to centre of trial orbital. Individual Lagrange multipliers, lambdas, default to global Lagrange multiplier.
character(len=50), public, save :: devel_flag
! Adaptive vs. fixed smearing stuff [GP, Jul 12, 2012]
! Only internal, always use the local variables defined by each module
! that take this value as default
logical :: adpt_smr
real(kind=dp) :: adpt_smr_fac
real(kind=dp) :: adpt_smr_max
real(kind=dp) :: smr_fixed_en_width
! GP: added a flag to check if this is the first run of param_read in library mode or not
logical, public, save :: library_param_read_first_pass
!IVO
logical, public, save :: spin_moment
real(kind=dp), public, save :: spin_axis_polar
real(kind=dp), public, save :: spin_axis_azimuth
logical, public, save :: use_degen_pert
real(kind=dp), public, save :: degen_thr
logical, public, save :: spin_decomp
integer, public, save :: num_valence_bands
logical :: found_fermi_energy
real(kind=dp), public, save :: scissors_shift
!IVO_END
! [gp-begin, Apr 20, 2012] Smearing type
! The prefactor is given with the above parameters smr_...
! This is an internal variable, obtained from the input string smr_type
! Only internal, always use the local variables defined by each module
! that take this value as default
integer :: smr_index
! [gp-end]
integer, allocatable, public, save :: exclude_bands(:)
integer, public, save :: num_wann
!! number of wannier functions
integer, public, save :: mp_grid(3)
!! Dimensions of the Monkhorst-Pack grid
! logical, public, save :: automatic_mp_grid
logical, public, save :: gamma_only
!! Use the special Gamma-point routines
real(kind=dp), public, save :: dis_win_min
!! lower bound of the disentanglement outer window
real(kind=dp), public, save :: dis_win_max
!! upper bound of the disentanglement outer window
real(kind=dp), public, save :: dis_froz_min
!! lower bound of the disentanglement inner (frozen) window
real(kind=dp), public, save :: dis_froz_max
!! upper bound of the disentanglement inner (frozen) window
integer, public, save :: dis_num_iter
!! number of disentanglement iteration steps
real(kind=dp), public, save :: dis_mix_ratio
!! Mixing ratio for the disentanglement routine
real(kind=dp), public, save :: dis_conv_tol
!! Convergence tolerance for the disentanglement
integer, public, save :: dis_conv_window
!! Size of the convergence window for disentanglement
! GS-start
integer, public, save :: dis_spheres_first_wann
integer, public, save :: dis_spheres_num
real(kind=dp), allocatable, public, save :: dis_spheres(:, :)
! GS-end
integer, public, save :: num_iter
!! Number of wannierisation iterations
integer, public, save :: num_cg_steps
!! Number of Conjugate Gradient steps
real(kind=dp), public, save :: conv_tol
integer, public, save :: conv_window
logical, public, save :: wannier_plot
integer, allocatable, public, save :: wannier_plot_list(:)
integer, public, save :: wannier_plot_supercell(3)
character(len=20), public, save :: wannier_plot_format
character(len=20), public, save :: wannier_plot_mode
character(len=20), public, save :: wannier_plot_spinor_mode
logical, public, save :: wannier_plot_spinor_phase
logical, public, save :: write_u_matrices
logical, public, save :: bands_plot
logical, public, save :: write_bvec
integer, public, save :: bands_num_points
character(len=20), public, save :: bands_plot_format
character(len=20), public, save :: bands_plot_mode
integer, allocatable, public, save :: bands_plot_project(:)
integer, public, save :: bands_plot_dim
logical, public, save :: write_hr
logical, public, save :: write_rmn
logical, public, save :: write_tb
real(kind=dp), public, save :: hr_cutoff
real(kind=dp), public, save :: dist_cutoff
character(len=20), public, save :: dist_cutoff_mode
real(kind=dp), public, save :: dist_cutoff_hc
character(len=20), public, save :: one_dim_axis
logical, public, save :: use_ws_distance
real(kind=dp), public, save :: ws_distance_tol
!! absolute tolerance for the distance to equivalent positions
logical, public, save :: fermi_surface_plot
integer, public, save :: fermi_surface_num_points
character(len=20), public, save :: fermi_surface_plot_format
real(kind=dp), save :: fermi_energy
! module k p a t h
logical, public, save :: kpath
character(len=20), public, save :: kpath_task
integer, public, save :: kpath_num_points
character(len=20), public, save :: kpath_bands_colour
! module k s l i c e
logical, public, save :: kslice
character(len=20), public, save :: kslice_task
real(kind=dp), public, save :: kslice_corner(3)
real(kind=dp), public, save :: kslice_b1(3)
real(kind=dp), public, save :: kslice_b2(3)
integer, public, save :: kslice_2dkmesh(2)
character(len=20), public, save :: kslice_fermi_lines_colour
! module d o s
logical, public, save :: dos
! No need to save 'dos_plot', only used here (introduced 'dos_task')
logical, public :: dos_plot
character(len=20), public, save :: dos_task
logical, public, save :: dos_adpt_smr
real(kind=dp), public, save :: dos_adpt_smr_fac
integer, public, save :: dos_smr_index
real(kind=dp), public, save :: dos_smr_fixed_en_width
real(kind=dp), public, save :: dos_adpt_smr_max
real(kind=dp), public, save :: dos_energy_max
real(kind=dp), public, save :: dos_energy_min
real(kind=dp), public, save :: dos_energy_step
integer, public, save :: num_dos_project
integer, allocatable, public, save :: dos_project(:)
! character(len=20), public, save :: dos_plot_format
real(kind=dp), public, save :: dos_kmesh_spacing
integer, public, save :: dos_kmesh(3)
! real(kind=dp), public, save :: dos_gaussian_width
! Module b e r r y
logical, public, save :: berry
character(len=120), public, save :: berry_task
real(kind=dp), public, save :: berry_kmesh_spacing
integer, public, save :: berry_kmesh(3)
! --------------remove eventually----------------
! integer, public, save :: alpha
! integer, public, save :: beta
! integer, public, save :: gamma
! --------------remove eventually----------------
integer, public, save :: berry_curv_adpt_kmesh
real(kind=dp), public, save :: berry_curv_adpt_kmesh_thresh
character(len=20), public, save :: berry_curv_unit
logical, public, save :: kubo_adpt_smr
real(kind=dp), public, save :: kubo_adpt_smr_fac
integer, public, save :: kubo_smr_index
real(kind=dp), public, save :: kubo_smr_fixed_en_width
real(kind=dp), public, save :: kubo_adpt_smr_max
integer, public, save :: sc_phase_conv
real(kind=dp), public, save :: sc_eta
real(kind=dp), public, save :: sc_w_thr
logical, public, save :: wanint_kpoint_file
! logical, public, save :: sigma_abc_onlyorb
logical, public, save :: transl_inv
logical, public, save :: gyrotropic
character(len=120), public, save :: gyrotropic_task
integer, public, save :: gyrotropic_kmesh(3)
real(kind=dp), public, save :: gyrotropic_kmesh_spacing
integer, public, save :: gyrotropic_smr_index
real(kind=dp), public, save :: gyrotropic_smr_fixed_en_width
real(kind=dp) :: gyrotropic_freq_min
real(kind=dp) :: gyrotropic_freq_max
real(kind=dp) :: gyrotropic_freq_step
integer, public, save :: gyrotropic_nfreq
complex(kind=dp), allocatable, public, save :: gyrotropic_freq_list(:)
real(kind=dp), public, save :: gyrotropic_box_corner(3), gyrotropic_box(3, 3)
real(kind=dp) :: gyrotropic_box_tmp(3)
real(kind=dp), public, save :: gyrotropic_degen_thresh
integer, allocatable, public, save :: gyrotropic_band_list(:)
integer, public, save :: gyrotropic_num_bands
real(kind=dp) :: smr_max_arg
real(kind=dp), public, save :: gyrotropic_smr_max_arg
real(kind=dp), public, save :: gyrotropic_eigval_max
logical :: fermi_energy_scan
real(kind=dp) :: fermi_energy_min
real(kind=dp) :: fermi_energy_max
real(kind=dp) :: fermi_energy_step
integer, public, save :: nfermi
real(kind=dp), allocatable, public, save :: fermi_energy_list(:)
real(kind=dp) :: kubo_freq_min
real(kind=dp) :: kubo_freq_max
real(kind=dp) :: kubo_freq_step
integer, public, save :: kubo_nfreq
complex(kind=dp), allocatable, public, save :: kubo_freq_list(:)
real(kind=dp), public, save :: kubo_eigval_max
! Module s p i n
real(kind=dp), public, save :: spin_kmesh_spacing
integer, public, save :: spin_kmesh(3)
! [gp-begin, Apr 13, 2012]
! Global interpolation k mesh variables
! These don't need to be public, since their values are copied in the variables of the
! local interpolation meshes. JRY: added save attribute
real(kind=dp), save :: kmesh_spacing
integer, save :: kmesh(3)
logical, save :: global_kmesh_set
! [gp-end]
! [gp-begin, Jun 1, 2012]
! GeneralInterpolator variables
logical, public, save :: geninterp
logical, public, save :: geninterp_alsofirstder
logical, public, save :: geninterp_single_file
! [gp-end, Jun 1, 2012]
! [gp-begin, Apr 12, 2012]
! BoltzWann variables
logical, public, save :: boltzwann
logical, public, save :: boltz_calc_also_dos
integer, public, save :: boltz_2d_dir_num
character(len=4), save :: boltz_2d_dir
real(kind=dp), public, save :: boltz_dos_energy_step
real(kind=dp), public, save :: boltz_dos_energy_min
real(kind=dp), public, save :: boltz_dos_energy_max
logical, public, save :: boltz_dos_adpt_smr
real(kind=dp), public, save :: boltz_dos_smr_fixed_en_width
real(kind=dp), public, save :: boltz_dos_adpt_smr_fac
real(kind=dp), public, save :: boltz_dos_adpt_smr_max
real(kind=dp), public, save :: boltz_mu_min
real(kind=dp), public, save :: boltz_mu_max
real(kind=dp), public, save :: boltz_mu_step
real(kind=dp), public, save :: boltz_temp_min
real(kind=dp), public, save :: boltz_temp_max
real(kind=dp), public, save :: boltz_temp_step
real(kind=dp), public, save :: boltz_kmesh_spacing
integer, public, save :: boltz_kmesh(3)
real(kind=dp), public, save :: boltz_tdf_energy_step
integer, public, save :: boltz_TDF_smr_index
integer, public, save :: boltz_dos_smr_index
real(kind=dp), public, save :: boltz_relax_time
real(kind=dp), public, save :: boltz_TDF_smr_fixed_en_width
logical, public, save :: boltz_bandshift
integer, public, save :: boltz_bandshift_firstband
real(kind=dp), public, save :: boltz_bandshift_energyshift
! [gp-end, Apr 12, 2012]
logical, public, save :: transport
logical, public, save :: tran_easy_fix
character(len=20), public, save :: transport_mode
real(kind=dp), public, save :: tran_win_min
real(kind=dp), public, save :: tran_win_max
real(kind=dp), public, save :: tran_energy_step
integer, public, save :: tran_num_bb
integer, public, save :: tran_num_ll
integer, public, save :: tran_num_rr
integer, public, save :: tran_num_cc
integer, public, save :: tran_num_lc
integer, public, save :: tran_num_cr
integer, public, save :: tran_num_bandc
logical, public, save :: tran_write_ht
logical, public, save :: tran_read_ht
logical, public, save :: tran_use_same_lead
integer, public, save :: tran_num_cell_ll
integer, public, save :: tran_num_cell_rr
real(kind=dp), public, save :: tran_group_threshold
real(kind=dp), public, save :: translation_centre_frac(3)
integer, public, save :: num_shells
!! no longer an input keyword
logical, public, save :: skip_B1_tests
!! do not check the B1 condition
logical, public, save :: explicit_nnkpts
!! nnkpts block is in the input file (allowed only for post-proc setup)
integer, allocatable, public, save :: shell_list(:)
real(kind=dp), allocatable, public, save :: kpt_latt(:, :)
!! kpoints in lattice vecs
real(kind=dp), public, save :: real_lattice(3, 3)
logical, public, save :: postproc_setup
logical, public, save :: cp_pp
!! Car-Parinello post-proc flag/transport
logical, public, save :: calc_only_A
logical, public, save :: use_bloch_phases
character(len=20), public, save :: restart
logical, public, save :: write_r2mn
logical, public, save :: guiding_centres
integer, public, save :: num_guide_cycles
integer, public, save :: num_no_guide_iter
real(kind=dp), public, save :: fixed_step
real(kind=dp), public, save :: trial_step
logical, public, save :: precond
logical, public, save :: write_proj
integer, public, save :: timing_level
logical, public, save :: spinors !are our WF spinors?
integer, public, save :: num_elec_per_state
logical, public, save :: translate_home_cell
logical, public, save :: write_xyz
logical, public, save :: write_hr_diag
real(kind=dp), public, save :: conv_noise_amp
integer, public, save :: conv_noise_num
real(kind=dp), public, save :: wannier_plot_radius
real(kind=dp), public, save :: wannier_plot_scale
integer, public, save :: search_shells !for kmesh
real(kind=dp), public, save :: kmesh_tol
integer, public, save :: optimisation
! aam: for WF-based calculation of vdW C6 coefficients
logical, public, save :: write_vdw_data
! Restarts
real(kind=dp), public, save :: omega_invariant
character(len=20), public, save :: checkpoint
logical, public, save :: have_disentangled
! Atom sites
real(kind=dp), allocatable, public, save :: atoms_pos_frac(:, :, :)
real(kind=dp), allocatable, public, save :: atoms_pos_cart(:, :, :)
integer, allocatable, public, save :: atoms_species_num(:)
character(len=maxlen), allocatable, public, save :: atoms_label(:)
character(len=2), allocatable, public, save :: atoms_symbol(:)
integer, public, save :: num_atoms
integer, public, save :: num_species
! Projections
real(kind=dp), allocatable, public, save :: proj_site(:, :)
integer, allocatable, public, save :: proj_l(:)
integer, allocatable, public, save :: proj_m(:)
integer, allocatable, public, save :: proj_s(:)
real(kind=dp), allocatable, public, save :: proj_s_qaxis(:, :)
real(kind=dp), allocatable, public, save :: proj_z(:, :)
real(kind=dp), allocatable, public, save :: proj_x(:, :)
integer, allocatable, public, save :: proj_radial(:)
real(kind=dp), allocatable, public, save :: proj_zona(:)
integer, public, save :: num_proj
!parameters dervied from input
integer, public, save :: num_kpts
real(kind=dp), public, save :: recip_lattice(3, 3)
real(kind=dp), public, save :: cell_volume
real(kind=dp), public, save :: real_metric(3, 3)
real(kind=dp), public, save :: recip_metric(3, 3)
integer, public, save :: bands_num_spec_points
character(len=20), allocatable, public, save ::bands_label(:)
real(kind=dp), allocatable, public, save ::bands_spec_points(:, :)
real(kind=dp), allocatable, public, save ::kpt_cart(:, :) !kpoints in cartesians
logical, public, save :: disentanglement
real(kind=dp), public, save :: lenconfac
integer, public, save :: num_wannier_plot
integer, public, save :: num_bands_project
integer, public, save :: num_exclude_bands
logical, public, save :: lfixstep
! kmesh parameters (set in kmesh)
integer, public, save :: nnh ! the number of b-directions (bka)
integer, public, save :: nntot ! total number of neighbours for each k-point
integer, public, save, allocatable :: nnlist(:, :) ! list of neighbours for each k-point
integer, public, save, allocatable :: neigh(:, :)
integer, public, save, allocatable :: nncell(:, :, :) ! gives BZ of each neighbour of each k-point
real(kind=dp), public, save :: wbtot
real(kind=dp), public, save, allocatable :: wb(:) ! weights associated with neighbours of each k-point
real(kind=dp), public, save, allocatable :: bk(:, :, :) ! the b-vectors that go from each k-point to its neighbours
real(kind=dp), public, save, allocatable :: bka(:, :) ! the b-directions from 1st k-point to its neighbours
! disentangle parameters
integer, public, save, allocatable :: ndimwin(:)
logical, public, save, allocatable :: lwindow(:, :)
logical, public, save :: frozen_states
! a_matrix and m_matrix_orig can be calculated internally from bloch states
! or read in from an ab-initio grid
! a_matrix = projection of trial orbitals on bloch states
! m_matrix_orig = overlap of bloch states
complex(kind=dp), allocatable, save, public :: a_matrix(:, :, :)
complex(kind=dp), allocatable, save, public :: m_matrix_orig(:, :, :, :)
complex(kind=dp), allocatable, save, public :: m_matrix_orig_local(:, :, :, :)
real(kind=dp), allocatable, save, public :: eigval(:, :)
logical, save, public :: eig_found
! $![ysl-b]
! $ ! ph_g = phase factor of Bloch functions at Gamma
! $ ! assuming that Bloch functions at Gamma are real except this phase factor
! $ complex(kind=dp), allocatable, save, public :: ph_g(:)
! $![ysl-e]
! u_matrix_opt gives the num_wann dimension optimal subspace from the
! original bloch states
complex(kind=dp), allocatable, save, public :: u_matrix_opt(:, :, :)
! u_matrix gives the unitary rotations from the optimal subspace to the
! optimally smooth states.
! m_matrix we store here, becuase it is needed for restart of wannierise
complex(kind=dp), allocatable, save, public :: u_matrix(:, :, :)
complex(kind=dp), allocatable, save, public :: m_matrix(:, :, :, :)
complex(kind=dp), allocatable, save, public :: m_matrix_local(:, :, :, :)
! RS: symmetry-adapted Wannier functions
logical, public, save :: lsitesymmetry = .false.
real(kind=dp), public, save :: symmetrize_eps = 1.d-3
! The maximum number of shells we need to satisfy B1 condition in kmesh
integer, parameter, public :: max_shells = 6
integer, parameter, public :: num_nnmax = 12
! Are we running as a library
logical, save, public :: library = .false.
! Are we running postw90?
logical, save, public :: ispostw90 = .false.
! IVO
! Are we running postw90 starting from an effective model?
logical, save, public :: effective_model = .false.
! Wannier centres and spreads
real(kind=dp), public, save, allocatable :: wannier_centres(:, :)
real(kind=dp), public, save, allocatable :: wannier_spreads(:)
real(kind=dp), public, save :: omega_total
real(kind=dp), public, save :: omega_tilde
! [ omega_invariant is declared above ]
! For Hamiltonian matrix in WF representation
logical, public, save :: automatic_translation
integer, public, save :: one_dim_dir
! vv: SCDM method
logical, public, save :: scdm_proj
integer, public, save :: scdm_entanglement
real(kind=dp), public, save :: scdm_mu
real(kind=dp), public, save :: scdm_sigma
! Private data
integer :: num_lines
character(len=maxlen), allocatable :: in_data(:)
character(len=maxlen) :: ctmp
logical :: ltmp
! AAM_2016-09-15: hr_plot is a deprecated input parameter. Replaced by write_hr.
logical :: hr_plot
public :: param_read
public :: param_write
public :: param_postw90_write
public :: param_dealloc
public :: param_write_header
public :: param_write_chkpt
public :: param_read_chkpt
public :: param_lib_set_atoms
public :: param_memory_estimate
public :: param_get_smearing_type
public :: param_get_convention_type
public :: param_dist
public :: param_chkpt_dist
contains
!==================================================================!
subroutine param_read()
!==================================================================!
! !
!! Read parameters and calculate derived values
!!
!! Note on parallelization: this function should be called
!! from the root node only!
!!
! !
!===================================================================
use w90_constants, only: bohr, eps6, cmplx_i
use w90_utility, only: utility_recip_lattice, utility_metric
use w90_io, only: io_error, io_file_unit, seedname, post_proc_flag
implicit none
!local variables
real(kind=dp) :: real_lattice_tmp(3, 3)
integer :: nkp, i, j, n, k, itmp, i_temp, i_temp2, eig_unit, loop, ierr, iv_temp(3), rows
logical :: found, found2, lunits, chk_found
character(len=6) :: spin_str
real(kind=dp) :: cosa(3), rv_temp(3)
integer, allocatable, dimension(:, :) :: nnkpts_block
integer, allocatable, dimension(:) :: nnkpts_idx
call param_in_file
!%%%%%%%%%%%%%%%%
! Site symmetry
!%%%%%%%%%%%%%%%%
! default value is lsitesymmetry=.false.
call param_get_keyword('site_symmetry', found, l_value=lsitesymmetry)!YN:
! default value is symmetrize_eps=0.001
call param_get_keyword('symmetrize_eps', found, r_value=symmetrize_eps)!YN:
!%%%%%%%%%%%%%%%%
! Transport
!%%%%%%%%%%%%%%%%
transport = .false.
call param_get_keyword('transport', found, l_value=transport)
tran_read_ht = .false.
call param_get_keyword('tran_read_ht', found, l_value=tran_read_ht)
tran_easy_fix = .false.
call param_get_keyword('tran_easy_fix', found, l_value=tran_easy_fix)
if (transport .and. tran_read_ht) restart = ' '
!%%%%%%%%%%%%%%%%
!System variables
!%%%%%%%%%%%%%%%%
timing_level = 1 ! Verbosity of timing output info
call param_get_keyword('timing_level', found, i_value=timing_level)
iprint = 1 ! Verbosity
call param_get_keyword('iprint', found, i_value=iprint)
optimisation = 3 ! Verbosity
call param_get_keyword('optimisation', found, i_value=optimisation)
if (transport .and. tran_read_ht) goto 301
!ivo
call param_get_keyword('effective_model', found, l_value=effective_model)
energy_unit = 'ev' !
call param_get_keyword('energy_unit', found, c_value=energy_unit)
length_unit = 'ang' !
lenconfac = 1.0_dp
call param_get_keyword('length_unit', found, c_value=length_unit)
if (length_unit .ne. 'ang' .and. length_unit .ne. 'bohr') &
call io_error('Error: value of length_unit not recognised in param_read')
if (length_unit .eq. 'bohr') lenconfac = 1.0_dp/bohr
wvfn_formatted = .false. ! formatted or "binary" file
call param_get_keyword('wvfn_formatted', found, l_value=wvfn_formatted)
spn_formatted = .false. ! formatted or "binary" file
call param_get_keyword('spn_formatted', found, l_value=spn_formatted)
uHu_formatted = .false. ! formatted or "binary" file
call param_get_keyword('uhu_formatted', found, l_value=uHu_formatted)
spin = 1
call param_get_keyword('spin', found, c_value=spin_str)
if (found) then
if (index(spin_str, 'up') > 0) then
spin = 1
elseif (index(spin_str, 'down') > 0) then
spin = 2
else
call io_error('Error: unrecognised value of spin found: '//trim(spin_str))
end if
end if
num_wann = -99
call param_get_keyword('num_wann', found, i_value=num_wann)
if (.not. found) call io_error('Error: You must specify num_wann')
if (num_wann <= 0) call io_error('Error: num_wann must be greater than zero')
num_exclude_bands = 0
call param_get_range_vector('exclude_bands', found, num_exclude_bands, lcount=.true.)
if (found) then
if (num_exclude_bands < 1) call io_error('Error: problem reading exclude_bands')
if (allocated(exclude_bands)) deallocate (exclude_bands)
allocate (exclude_bands(num_exclude_bands), stat=ierr)
if (ierr /= 0) call io_error('Error allocating exclude_bands in param_read')
call param_get_range_vector('exclude_bands', found, num_exclude_bands, .false., exclude_bands)
if (any(exclude_bands < 1)) &
call io_error('Error: exclude_bands must contain positive numbers')
end if
! AAM_2016-09-16: some changes to logic to patch a problem with uninitialised num_bands in library mode
! num_bands = -1
call param_get_keyword('num_bands', found, i_value=i_temp)
if (found .and. library) write (stdout, '(/a)') ' Ignoring <num_bands> in input file'
if (.not. library .and. .not. effective_model) then
if (found) num_bands = i_temp
if (.not. found) num_bands = num_wann
end if
! GP: I subtract it here, but only the first time when I pass the total number of bands
! In later calls, I need to pass instead num_bands already subtracted.
if (library .and. library_param_read_first_pass) num_bands = num_bands - num_exclude_bands
if (.not. effective_model) then
if (found .and. num_bands < num_wann) then
write (stdout, *) 'num_bands', num_bands
write (stdout, *) 'num_wann', num_wann
call io_error('Error: num_bands must be greater than or equal to num_wann')
endif
endif
num_dump_cycles = 100 ! frequency to write backups at
call param_get_keyword('num_dump_cycles', found, i_value=num_dump_cycles)
if (num_dump_cycles < 0) call io_error('Error: num_dump_cycles must be positive')
num_print_cycles = 1 ! frequency to write at
call param_get_keyword('num_print_cycles', found, i_value=num_print_cycles)
if (num_print_cycles < 0) call io_error('Error: num_print_cycles must be positive')
devel_flag = ' ' !
call param_get_keyword('devel_flag', found, c_value=devel_flag)
! mp_grid=-99
call param_get_keyword_vector('mp_grid', found, 3, i_value=iv_temp)
if (found .and. library) write (stdout, '(a)') ' Ignoring <mp_grid> in input file'
if (.not. library .and. .not. effective_model) then
if (found) mp_grid = iv_temp
if (.not. found) then
call io_error('Error: You must specify dimensions of the Monkhorst-Pack grid by setting mp_grid')
elseif (any(mp_grid < 1)) then
call io_error('Error: mp_grid must be greater than zero')
end if
num_kpts = mp_grid(1)*mp_grid(2)*mp_grid(3)
end if
![ysl-b]
ltmp = .false.
call param_get_keyword('gamma_only', found, l_value=ltmp)
if (.not. library) then
gamma_only = ltmp
if (gamma_only .and. (num_kpts .ne. 1)) &
call io_error('Error: gamma_only is true, but num_kpts > 1')
else
if (found) write (stdout, '(a)') ' Ignoring <gamma_only> in input file'
endif
![ysl-e]
! aam: automatic_mp_grid no longer used
! automatic_mp_grid = .false.
! call param_get_keyword('automatic_mp_grid',found,l_value=automatic_mp_grid)
postproc_setup = .false. ! set to true to write .nnkp file and exit
call param_get_keyword('postproc_setup', found, l_value=postproc_setup)
! We allow this keyword to be overriden by a command line arg -pp
if (post_proc_flag) postproc_setup = .true.
cp_pp = .false. ! set to true if doing CP post-processing
call param_get_keyword('cp_pp', found, l_value=cp_pp)
calc_only_A = .false.
call param_get_keyword('calc_only_A', found, l_value=calc_only_A)
restart = ' '
call param_get_keyword('restart', found, c_value=restart)
if (found) then
if ((restart .ne. 'default') .and. (restart .ne. 'wannierise') &
.and. (restart .ne. 'plot') .and. (restart .ne. 'transport')) then
call io_error('Error in input file: value of restart not recognised')
else
inquire (file=trim(seedname)//'.chk', exist=chk_found)
if (.not. chk_found) &
call io_error('Error: restart requested but '//trim(seedname)//'.chk file not found')
endif
endif
!post processing takes priority (user is not warned of this)
if (postproc_setup) restart = ' '
write_r2mn = .false.
call param_get_keyword('write_r2mn', found, l_value=write_r2mn)
write_proj = .false.
call param_get_keyword('write_proj', found, l_value=write_proj)
ltmp = .false. ! by default our WF are not spinors
call param_get_keyword('spinors', found, l_value=ltmp)
if (.not. library) then
spinors = ltmp
else
if (found) write (stdout, '(a)') ' Ignoring <spinors> in input file'
endif
! if(spinors .and. (2*(num_wann/2))/=num_wann) &
! call io_error('Error: For spinor WF num_wann must be even')
! We need to know if the bands are double degenerate due to spin, e.g. when
! calculating the DOS
if (spinors) then
num_elec_per_state = 1
else
num_elec_per_state = 2
endif
call param_get_keyword('num_elec_per_state', found, i_value=num_elec_per_state)
if ((num_elec_per_state /= 1) .and. (num_elec_per_state /= 2)) &
call io_error('Error: num_elec_per_state can be only 1 or 2')
if (spinors .and. num_elec_per_state /= 1) &
call io_error('Error: when spinors = T num_elec_per_state must be 1')
translate_home_cell = .false.
call param_get_keyword('translate_home_cell', found, l_value=translate_home_cell)
write_xyz = .false.
call param_get_keyword('write_xyz', found, l_value=write_xyz)
write_hr_diag = .false.
call param_get_keyword('write_hr_diag', found, l_value=write_hr_diag)
!%%%%%%%%%%%
! Wannierise
!%%%%%%%%%%%
num_iter = 100
call param_get_keyword('num_iter', found, i_value=num_iter)
if (num_iter < 0) call io_error('Error: num_iter must be positive')
num_cg_steps = 5
call param_get_keyword('num_cg_steps', found, i_value=num_cg_steps)
if (num_cg_steps < 0) call io_error('Error: num_cg_steps must be positive')
conv_tol = 1.0e-10_dp
call param_get_keyword('conv_tol', found, r_value=conv_tol)
if (conv_tol < 0.0_dp) call io_error('Error: conv_tol must be positive')
conv_noise_amp = -1.0_dp
call param_get_keyword('conv_noise_amp', found, r_value=conv_noise_amp)
conv_window = -1
if (conv_noise_amp > 0.0_dp) conv_window = 5
call param_get_keyword('conv_window', found, i_value=conv_window)
conv_noise_num = 3
call param_get_keyword('conv_noise_num', found, i_value=conv_noise_num)
if (conv_noise_num < 0) call io_error('Error: conv_noise_num must be positive')
guiding_centres = .false.
call param_get_keyword('guiding_centres', found, l_value=guiding_centres)
num_guide_cycles = 1
call param_get_keyword('num_guide_cycles', found, i_value=num_guide_cycles)
if (num_guide_cycles < 0) call io_error('Error: num_guide_cycles must be >= 0')
num_no_guide_iter = 0
call param_get_keyword('num_no_guide_iter', found, i_value=num_no_guide_iter)
if (num_no_guide_iter < 0) call io_error('Error: num_no_guide_iter must be >= 0')
fixed_step = -999.0_dp; lfixstep = .false.
call param_get_keyword('fixed_step', found, r_value=fixed_step)
if (found .and. (fixed_step < 0.0_dp)) call io_error('Error: fixed_step must be > 0')
if (fixed_step > 0.0_dp) lfixstep = .true.
trial_step = 2.0_dp
call param_get_keyword('trial_step', found, r_value=trial_step)
if (found .and. lfixstep) call io_error('Error: cannot specify both fixed_step and trial_step')
precond = .false.
call param_get_keyword('precond', found, l_value=precond)
slwf_num = num_wann
selective_loc = .false.
call param_get_keyword('slwf_num', found, i_value=slwf_num)
if (found) then
if (slwf_num .gt. num_wann .or. slwf_num .lt. 1) then
call io_error('Error: slwf_num must be an integer between 1 and num_wann')
end if
if (slwf_num .lt. num_wann) selective_loc = .true.
end if
slwf_constrain = .false.
call param_get_keyword('slwf_constrain', found, l_value=slwf_constrain)
if (found .and. slwf_constrain) then
if (selective_loc) then
allocate (ccentres_frac(num_wann, 3), stat=ierr)
if (ierr /= 0) call io_error('Error allocating ccentres_frac in param_get_centre_constraints')
allocate (ccentres_cart(num_wann, 3), stat=ierr)
if (ierr /= 0) call io_error('Error allocating ccentres_cart in param_get_centre_constraints')
else
write (stdout, *) ' No selective localisation requested. Ignoring constraints on centres'
slwf_constrain = .false.
end if
end if
slwf_lambda = 1.0_dp
call param_get_keyword('slwf_lambda', found, r_value=slwf_lambda)
if (found) then
if (slwf_lambda < 0.0_dp) call io_error('Error: slwf_lambda must be positive.')
endif
!%%%%%%%%%
! Plotting
!%%%%%%%%%
wannier_plot = .false.
call param_get_keyword('wannier_plot', found, l_value=wannier_plot)
wannier_plot_supercell = 2
call param_get_vector_length('wannier_plot_supercell', found, length=i)
if (found) then
if (i .eq. 1) then
call param_get_keyword_vector('wannier_plot_supercell', found, 1, &
i_value=wannier_plot_supercell)
wannier_plot_supercell(2) = wannier_plot_supercell(1)
wannier_plot_supercell(3) = wannier_plot_supercell(1)
elseif (i .eq. 3) then
call param_get_keyword_vector('wannier_plot_supercell', found, 3, &
i_value=wannier_plot_supercell)
else
call io_error('Error: wannier_plot_supercell must be provided as either one integer or a vector of three integers')
end if
if (any(wannier_plot_supercell <= 0)) &
call io_error('Error: wannier_plot_supercell elements must be greater than zero')
end if
wannier_plot_format = 'xcrysden'
call param_get_keyword('wannier_plot_format', found, c_value=wannier_plot_format)
wannier_plot_mode = 'crystal'
call param_get_keyword('wannier_plot_mode', found, c_value=wannier_plot_mode)
wannier_plot_spinor_mode = 'total'
call param_get_keyword('wannier_plot_spinor_mode', found, c_value=wannier_plot_spinor_mode)
wannier_plot_spinor_phase = .true.
call param_get_keyword('wannier_plot_spinor_phase', found, l_value=wannier_plot_spinor_phase)
call param_get_range_vector('wannier_plot_list', found, num_wannier_plot, lcount=.true.)
if (found) then
if (num_wannier_plot < 1) call io_error('Error: problem reading wannier_plot_list')
if (allocated(wannier_plot_list)) deallocate (wannier_plot_list)
allocate (wannier_plot_list(num_wannier_plot), stat=ierr)
if (ierr /= 0) call io_error('Error allocating wannier_plot_list in param_read')
call param_get_range_vector('wannier_plot_list', found, num_wannier_plot, .false., wannier_plot_list)
if (any(wannier_plot_list < 1) .or. any(wannier_plot_list > num_wann)) &
call io_error('Error: wannier_plot_list asks for a non-valid wannier function to be plotted')
else
! we plot all wannier functions
num_wannier_plot = num_wann
if (allocated(wannier_plot_list)) deallocate (wannier_plot_list)
allocate (wannier_plot_list(num_wannier_plot), stat=ierr)
if (ierr /= 0) call io_error('Error allocating wannier_plot_list in param_read')
do loop = 1, num_wann
wannier_plot_list(loop) = loop
end do
end if
wannier_plot_radius = 3.5_dp
call param_get_keyword('wannier_plot_radius', found, r_value=wannier_plot_radius)
wannier_plot_scale = 1.0_dp
call param_get_keyword('wannier_plot_scale', found, r_value=wannier_plot_scale)
! checks
if (wannier_plot) then
if ((index(wannier_plot_format, 'xcrys') .eq. 0) .and. (index(wannier_plot_format, 'cub') .eq. 0)) &
call io_error('Error: wannier_plot_format not recognised')
if ((index(wannier_plot_mode, 'crys') .eq. 0) .and. (index(wannier_plot_mode, 'mol') .eq. 0)) &
call io_error('Error: wannier_plot_mode not recognised')
if ((index(wannier_plot_spinor_mode, 'total') .eq. 0) .and. (index(wannier_plot_spinor_mode, 'up') .eq. 0) &
.and. (index(wannier_plot_spinor_mode, 'down') .eq. 0)) &
call io_error('Error: wannier_plot_spinor_mode not recognised')
if (wannier_plot_radius < 0.0_dp) call io_error('Error: wannier_plot_radius must be positive')
if (wannier_plot_scale < 0.0_dp) call io_error('Error: wannier_plot_scale must be positive')
endif
write_u_matrices = .false.
call param_get_keyword('write_u_matrices', found, l_value=write_u_matrices)
bands_plot = .false.
call param_get_keyword('bands_plot', found, l_value=bands_plot)
write_bvec = .false.
call param_get_keyword('write_bvec', found, l_value=write_bvec)
bands_num_points = 100
call param_get_keyword('bands_num_points', found, i_value=bands_num_points)
bands_plot_format = 'gnuplot'
call param_get_keyword('bands_plot_format', found, c_value=bands_plot_format)
bands_plot_mode = 's-k'
call param_get_keyword('bands_plot_mode', found, c_value=bands_plot_mode)
bands_plot_dim = 3
call param_get_keyword('bands_plot_dim', found, i_value=bands_plot_dim)
num_bands_project = 0
call param_get_range_vector('bands_plot_project', found, num_bands_project, lcount=.true.)
if (found) then
if (num_bands_project < 1) call io_error('Error: problem reading bands_plot_project')
if (allocated(bands_plot_project)) deallocate (bands_plot_project)
allocate (bands_plot_project(num_bands_project), stat=ierr)
if (ierr /= 0) call io_error('Error allocating bands_plot_project in param_read')
call param_get_range_vector('bands_plot_project', found, num_bands_project, .false., bands_plot_project)
if (any(bands_plot_project < 1) .or. any(bands_plot_project > num_wann)) &
call io_error('Error: bands_plot_project asks for a non-valid wannier function to be projected')
endif
bands_num_spec_points = 0
call param_get_block_length('kpoint_path', found, i_temp)
if (found) then
bands_num_spec_points = i_temp*2
if (allocated(bands_label)) deallocate (bands_label)
allocate (bands_label(bands_num_spec_points), stat=ierr)
if (ierr /= 0) call io_error('Error allocating bands_label in param_read')
if (allocated(bands_spec_points)) deallocate (bands_spec_points)
allocate (bands_spec_points(3, bands_num_spec_points), stat=ierr)
if (ierr /= 0) call io_error('Error allocating bands_spec_points in param_read')
call param_get_keyword_kpath
end if
if (.not. found .and. bands_plot) &
call io_error('A bandstructure plot has been requested but there is no kpoint_path block')
! checks
if (bands_plot) then
if ((index(bands_plot_format, 'gnu') .eq. 0) .and. (index(bands_plot_format, 'xmgr') .eq. 0)) &
call io_error('Error: bands_plot_format not recognised')
if ((index(bands_plot_mode, 's-k') .eq. 0) .and. (index(bands_plot_mode, 'cut') .eq. 0)) &
call io_error('Error: bands_plot_mode not recognised')
if (bands_num_points < 0) call io_error('Error: bands_num_points must be positive')
endif
fermi_surface_plot = .false.
call param_get_keyword('fermi_surface_plot', found, l_value=fermi_surface_plot)
fermi_surface_num_points = 50
call param_get_keyword('fermi_surface_num_points', found, i_value=fermi_surface_num_points)
fermi_surface_plot_format = 'xcrysden'
call param_get_keyword('fermi_surface_plot_format', &
found, c_value=fermi_surface_plot_format)
nfermi = 0
found_fermi_energy = .false.
call param_get_keyword('fermi_energy', found, r_value=fermi_energy)
if (found) then
found_fermi_energy = .true.
nfermi = 1
endif
!
fermi_energy_scan = .false.
call param_get_keyword('fermi_energy_min', found, r_value=fermi_energy_min)
if (found) then
if (found_fermi_energy) call io_error( &
'Error: Cannot specify both fermi_energy and fermi_energy_min')
fermi_energy_scan = .true.
fermi_energy_max = fermi_energy_min + 1.0_dp
call param_get_keyword('fermi_energy_max', found, &
r_value=fermi_energy_max)
if (found .and. fermi_energy_max <= fermi_energy_min) call io_error( &
'Error: fermi_energy_max must be larger than fermi_energy_min')
fermi_energy_step = 0.01_dp
call param_get_keyword('fermi_energy_step', found, &
r_value=fermi_energy_step)
if (found .and. fermi_energy_step <= 0.0_dp) call io_error( &
'Error: fermi_energy_step must be positive')
nfermi = nint(abs((fermi_energy_max - fermi_energy_min)/fermi_energy_step)) + 1
endif
!
if (found_fermi_energy) then
if (allocated(fermi_energy_list)) deallocate (fermi_energy_list)