forked from QcmPlab/CDMFT-LANC-ED
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathED_AUX_FUNX.f90
673 lines (609 loc) · 23.5 KB
/
ED_AUX_FUNX.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
MODULE ED_AUX_FUNX
USE ED_INPUT_VARS
USE ED_VARS_GLOBAL
USE SF_TIMER
USE SF_LINALG
USE SF_MISC, only: assert_shape
USE SF_IOTOOLS, only:free_unit,reg,txtfy
implicit none
private
interface lso2nnn_reshape
module procedure d_nlso2nnn
module procedure c_nlso2nnn
end interface lso2nnn_reshape
interface so2nn_reshape
module procedure d_nso2nn
module procedure c_nso2nn
end interface so2nn_reshape
interface nnn2lso_reshape
module procedure d_nnn2nlso
module procedure c_nnn2nlso
end interface nnn2lso_reshape
interface nn2so_reshape
module procedure d_nn2nso
module procedure c_nn2nso
end interface nn2so_reshape
interface print_Hloc
module procedure :: print_Hloc_nnn
module procedure :: print_Hloc_lso
end interface print_Hloc
!interface set_Hloc
!module procedure set_Hloc_lso
!module procedure set_Hloc_nnn
!end interface set_Hloc
#if __GNUC__ > 6
interface read(unformatted)
module procedure :: read_unformatted
end interface read(unformatted)
interface write(unformatted)
module procedure :: write_unformatted
end interface write(unformatted)
interface read(formatted)
module procedure :: read_formatted
end interface read(formatted)
interface write(formatted)
module procedure :: write_formatted
end interface write(formatted)
#endif
public :: index_stride_lso
!
!public :: set_Hloc
public :: print_Hloc
!
public :: save_gfprime
public :: read_gfprime
!
public :: lso2nnn_reshape
public :: so2nn_reshape
public :: nnn2lso_reshape
public :: nn2so_reshape
!
public :: ed_search_variable
!
contains
!> Get stride position in the one-particle many-body space
function index_stride_lso(ilat,ispin,iorb) result(indx)
integer :: ilat
integer :: iorb
integer :: ispin
integer :: indx
indx = iorb + (ilat-1)*Norb + (ispin-1)*Norb*Nlat
end function index_stride_lso
!##################################################################
! HLOC ROUTINES
!##################################################################
!+------------------------------------------------------------------+
!PURPOSE : Print Hloc
!+------------------------------------------------------------------+
subroutine print_Hloc_nnn(hloc,file)![Nspin][Nspin][Norb][Norb]
real(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb) :: hloc
character(len=*),optional :: file
integer :: ilat,jlat,iorb,jorb,ispin,jspin
integer :: unit
!
unit=LOGfile
if(present(file))then
open(free_unit(unit),file=reg(file))
write(LOGfile,"(A)")"print_Hloc on file :"//reg(file)
endif
do ispin=1,Nspin
do ilat=1,Nlat
do iorb=1,Norb
write(unit,"(20(F7.3,2x))")&
(((Hloc(ilat,jlat,ispin,jspin,iorb,jorb),jorb =1,Norb),jlat=1,Nlat),jspin=1,Nspin)
enddo
enddo
enddo
write(unit,*)""
if(present(file))close(unit)
end subroutine print_Hloc_nnn
subroutine print_Hloc_lso(hloc,file) ![Nlso][Nlso]
real(8),dimension(Nlat*Nspin*Norb,Nlat*Nspin*Norb) :: hloc
character(len=*),optional :: file
integer :: ilat,iorb,jorb,unit
unit=LOGfile
!
if(present(file))then
open(free_unit(unit),file=reg(file))
write(LOGfile,"(A)")"print_Hloc on file :"//reg(file)
endif
!
Nlso = Nlat*Nspin*Norb
do iorb=1,Nlso
write(unit,"(20(F7.3,2x))")(Hloc(iorb,jorb),jorb =1,Nlso)
enddo
write(unit,*)""
if(present(file))close(unit)
end subroutine print_Hloc_lso
!##################################################################
! RESHAPE ROUTINES
!##################################################################
!+-----------------------------------------------------------------------------+!
!PURPOSE:
! reshape a matrix from the [Nlso][Nlso] shape
! from/to the [Nlat][Nspin][Nspin][Norb][Norb] shape.
! _nlso2nnn : from [Nlso][Nlso] to [Nlat][Nspin][Nspin][Norb][Norb] !
! _nso2nn : from [Nso][Nso] to [Nspin][Nspin][Norb][Norb]
!+-----------------------------------------------------------------------------+!
function d_nlso2nnn(Hlso,Nlat,Nspin,Norb) result(Hnnn)
integer :: Nlat,Nspin,Norb
real(8),dimension(Nlat*Nspin*Norb,Nlat*Nspin*Norb) :: Hlso
real(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb) :: Hnnn
integer :: ilat,jlat
integer :: iorb,jorb
integer :: ispin,jspin
integer :: is,js
Hnnn=zero
do ilat=1,Nlat
do jlat=1,Nlat
do ispin=1,Nspin
do jspin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
is = index_stride_lso(ilat,ispin,iorb)
js = index_stride_lso(jlat,jspin,jorb)
Hnnn(ilat,jlat,ispin,jspin,iorb,jorb) = Hlso(is,js)
enddo
enddo
enddo
enddo
enddo
enddo
end function d_nlso2nnn
!
function c_nlso2nnn(Hlso,Nlat,Nspin,Norb) result(Hnnn)
integer :: Nlat,Nspin,Norb
complex(8),dimension(Nlat*Nspin*Norb,Nlat*Nspin*Norb) :: Hlso
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb) :: Hnnn
integer :: ilat,jlat
integer :: iorb,jorb
integer :: ispin,jspin
integer :: is,js
Hnnn=zero
do ilat=1,Nlat
do jlat=1,Nlat
do ispin=1,Nspin
do jspin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
is = index_stride_lso(ilat,ispin,iorb)
js = index_stride_lso(jlat,jspin,jorb)
Hnnn(ilat,jlat,ispin,jspin,iorb,jorb) = Hlso(is,js)
enddo
enddo
enddo
enddo
enddo
enddo
end function c_nlso2nnn
!
function d_nso2nn(Hso,Nspin,Norb) result(Hnn)
integer :: Nspin,Norb
real(8),dimension(Nspin*Norb,Nspin*Norb) :: Hso
real(8),dimension(Nspin,Nspin,Norb,Norb) :: Hnn
integer :: iorb,ispin,is
integer :: jorb,jspin,js
Hnn=zero
do ispin=1,Nspin
do jspin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
is = iorb + (ispin-1)*Norb !spin-orbit stride
js = jorb + (jspin-1)*Norb !spin-orbit stride
Hnn(ispin,jspin,iorb,jorb) = Hso(is,js)
enddo
enddo
enddo
enddo
end function d_nso2nn
!
function c_nso2nn(Hso,Nspin,Norb) result(Hnn)
integer :: Nspin,Norb
complex(8),dimension(Nspin*Norb,Nspin*Norb) :: Hso
complex(8),dimension(Nspin,Nspin,Norb,Norb) :: Hnn
integer :: iorb,ispin,is
integer :: jorb,jspin,js
Hnn=zero
do ispin=1,Nspin
do jspin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
is = iorb + (ispin-1)*Norb !spin-orbit stride
js = jorb + (jspin-1)*Norb !spin-orbit stride
Hnn(ispin,jspin,iorb,jorb) = Hso(is,js)
enddo
enddo
enddo
enddo
end function c_nso2nn
!+-----------------------------------------------------------------------------+!
!PURPOSE:
! reshape a matrix from the [Nlat][Nspin][Nspin][Norb][Norb] shape
! from/to the [Nlso][Nlso] shape.
! _nnn2nlso : from [Nlat][Nspin][Nspin][Norb][Norb] to [Nlso][Nlso]
! _nn2nso : from [Nspin][Nspin][Norb][Norb] to [Nso][Nso]
!+-----------------------------------------------------------------------------+!
function d_nnn2nlso(Hnnn,Nlat,Nspin,Norb) result(Hlso)
integer :: Nlat,Nspin,Norb
real(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb) :: Hnnn
real(8),dimension(Nlat*Nspin*Norb,Nlat*Nspin*Norb) :: Hlso
integer :: ilat,jlat
integer :: iorb,jorb
integer :: ispin,jspin
integer :: is,js
Hlso=zero
do ilat=1,Nlat
do jlat=1,Nlat
do ispin=1,Nspin
do jspin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
is = index_stride_lso(ilat,ispin,iorb)
js = index_stride_lso(jlat,jspin,jorb)
Hlso(is,js) = Hnnn(ilat,jlat,ispin,jspin,iorb,jorb)
enddo
enddo
enddo
enddo
enddo
enddo
end function d_nnn2nlso
!
function c_nnn2nlso(Hnnn,Nlat,Nspin,Norb) result(Hlso)
integer :: Nlat,Nspin,Norb
complex(8),dimension(Nlat,Nlat,Nspin,Nspin,Norb,Norb) :: Hnnn
complex(8),dimension(Nlat*Nspin*Norb,Nlat*Nspin*Norb) :: Hlso
integer :: ilat,jlat
integer :: iorb,jorb
integer :: ispin,jspin
integer :: is,js
Hlso=zero
do ilat=1,Nlat
do jlat=1,Nlat
do ispin=1,Nspin
do jspin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
is = index_stride_lso(ilat,ispin,iorb)
js = index_stride_lso(jlat,jspin,jorb)
Hlso(is,js) = Hnnn(ilat,jlat,ispin,jspin,iorb,jorb)
enddo
enddo
enddo
enddo
enddo
enddo
end function c_nnn2nlso
!
function d_nn2nso(Hnn,Nspin,Norb) result(Hso)
integer :: Nspin,Norb
real(8),dimension(Nspin,Nspin,Norb,Norb) :: Hnn
real(8),dimension(Nspin*Norb,Nspin*Norb) :: Hso
integer :: iorb,ispin,is
integer :: jorb,jspin,js
Hso=zero
do ispin=1,Nspin
do jspin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
is = iorb + (ispin-1)*Norb !spin-orbit stride
js = jorb + (jspin-1)*Norb !spin-orbit stride
Hso(is,js) = Hnn(ispin,jspin,iorb,jorb)
enddo
enddo
enddo
enddo
end function d_nn2nso
function c_nn2nso(Hnn,Nspin,Norb) result(Hso)
integer :: Nspin,Norb
complex(8),dimension(Nspin,Nspin,Norb,Norb) :: Hnn
complex(8),dimension(Nspin*Norb,Nspin*Norb) :: Hso
integer :: iorb,ispin,is
integer :: jorb,jspin,js
Hso=zero
do ispin=1,Nspin
do jspin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
is = iorb + (ispin-1)*Norb !spin-orbit stride
js = jorb + (jspin-1)*Norb !spin-orbit stride
Hso(is,js) = Hnn(ispin,jspin,iorb,jorb)
enddo
enddo
enddo
enddo
end function c_nn2nso
#if __GNUC__ > 6
!##################################################################
!##################################################################
! ROUTINES TO READ AND WRITE CLUSTER GREEN FUNCTION
! unformatted and formatted I/O
!##################################################################
!##################################################################
!+-------------------------------------------------------------------+
!PURPOSE : write overload for GFmatrix type (formatted)
!+-------------------------------------------------------------------+
subroutine write_formatted(dtv, unit, iotype, v_list, iostat, iomsg)
class(GFmatrix), intent(in) :: dtv
integer, intent(in) :: unit
integer, intent(out) :: iostat
character(*), intent(in) :: iotype
integer, intent(in) :: v_list(:)
integer :: Nexc,iexc,Ichan,ilat,jlat,iorb,ispin,istate
integer :: Nchan,Nstates
character(*), intent(inout) :: iomsg
!
!
Nstates = size(dtv%state)
write (unit, *,IOSTAT=iostat, IOMSG=iomsg) Nstates
do istate=1,Nstates
Nchan = size(dtv%state(istate)%channel)
write (unit, *,IOSTAT=iostat, IOMSG=iomsg) Nchan
do ichan=1,Nchan
write (unit, *,IOSTAT=iostat, IOMSG=iomsg) size(dtv%state(istate)%channel(ichan)%poles)
write (unit, *,IOSTAT=iostat, IOMSG=iomsg) dtv%state(istate)%channel(ichan)%poles
write (unit, *,IOSTAT=iostat, IOMSG=iomsg) dtv%state(istate)%channel(ichan)%weight
enddo
write (unit, *,IOSTAT=iostat, IOMSG=iomsg) "\n"
enddo
!
end subroutine write_formatted
!+-------------------------------------------------------------------+
!PURPOSE : read overload for GFmatrix type (formatted)
!+-------------------------------------------------------------------+
subroutine read_formatted(dtv, unit,iotype, v_list, iostat, iomsg)
class(GFmatrix), intent(inout) :: dtv
integer, intent(in) :: unit
integer, intent(out) :: iostat
character(*), intent(in) :: iotype
integer, intent(in) :: v_list(:)
character(*), intent(inout) :: iomsg
logical :: alloc
integer :: ichan,Nchan,Nlanc,istate,Nstates
!
read (unit,*,IOSTAT=iostat, IOMSG=iomsg) Nstates
call GFmatrix_allocate(dtv,Nstate=Nstates)
do istate=1,Nstates
read (unit,*,IOSTAT=iostat, IOMSG=iomsg) Nchan
call GFmatrix_allocate(dtv,istate=istate,Nchan=Nchan)
do ichan=1,Nchan
read (unit,*, IOSTAT=iostat, IOMSG=iomsg) Nlanc
call GFmatrix_allocate(dtv,istate=istate,ichan=ichan,Nexc=Nlanc)
read (unit, *, IOSTAT=iostat, IOMSG=iomsg) dtv%state(istate)%channel(ichan)%poles
read (unit, *, IOSTAT=iostat, IOMSG=iomsg) dtv%state(istate)%channel(ichan)%weight
enddo
enddo
!
end subroutine read_formatted
!+-------------------------------------------------------------------+
!PURPOSE : write overload for GFmatrix type (unformatted)
!+-------------------------------------------------------------------+
subroutine write_unformatted(dtv, unit, iostat, iomsg)
class(GFmatrix), intent(in) :: dtv
integer, intent(in) :: unit
integer, intent(out) :: iostat
integer :: Nexc,iexc,Ichan,ilat,jlat,iorb,ispin,istate
integer :: Nchan,Nstates
character(*), intent(inout) :: iomsg
!
!
Nstates = size(dtv%state)
write (unit, IOSTAT=iostat, IOMSG=iomsg) Nstates
do istate=1,Nstates
Nchan = size(dtv%state(istate)%channel)
write (unit, IOSTAT=iostat, IOMSG=iomsg) Nchan
do ichan=1,Nchan
write (unit, IOSTAT=iostat, IOMSG=iomsg) size(dtv%state(istate)%channel(ichan)%poles), dtv%state(istate)%channel(ichan)%poles, dtv%state(istate)%channel(ichan)%weight
enddo
enddo
!
end subroutine write_unformatted
!+-------------------------------------------------------------------+
!PURPOSE : read overload for GFmatrix type (unformatted)
!+-------------------------------------------------------------------+
subroutine read_unformatted(dtv, unit, iostat, iomsg)
class(GFmatrix), intent(inout) :: dtv
integer, intent(in) :: unit
integer, intent(out) :: iostat
character(*), intent(inout) :: iomsg
logical :: alloc
integer :: ichan,Nchan,Nlanc,istate,Nstates
!
read (unit, IOSTAT=iostat, IOMSG=iomsg) Nstates
call GFmatrix_allocate(dtv,Nstate=Nstates)
do istate=1,Nstates
read (unit, IOSTAT=iostat, IOMSG=iomsg) Nchan
call GFmatrix_allocate(dtv,istate=istate,Nchan=Nchan)
do ichan=1,Nchan
read (unit, IOSTAT=iostat, IOMSG=iomsg) Nlanc
call GFmatrix_allocate(dtv,istate=istate,ichan=ichan,Nexc=Nlanc)
read (unit, IOSTAT=iostat, IOMSG=iomsg) dtv%state(istate)%channel(ichan)%poles
read (unit, IOSTAT=iostat, IOMSG=iomsg) dtv%state(istate)%channel(ichan)%weight
enddo
enddo
!
end subroutine read_unformatted
#endif
!+-------------------------------------------------------------------+
!PURPOSE : Save cluster GF to file
!+-------------------------------------------------------------------+
subroutine save_gfprime(file,used,use_formatted)
character(len=*),optional :: file
character(len=256) :: file_
logical,optional :: used
logical :: used_
logical,optional :: use_formatted
logical :: use_formatted_
character(len=16) :: extension
integer :: unit_,Nchannel,Nexc,ichan,iexc,ilat,jlat,ispin,iorb,jorb
!
#if __GNUC__ > 6
if(.not.allocated(impGmatrix))stop "ed_gf_cluster ERROR: impGmatrix not allocated!"
used_=.false.;if(present(used))used_=used
use_formatted_=.false.;if(present(use_formatted))use_formatted_=use_formatted
extension=".restart";if(used_)extension=".used"
file_=str(str(file)//str(ed_file_suffix)//str(extension))
unit_=free_unit()
!
if(use_formatted_)then
open(unit_,file=str(file_),access='sequential')
else
open(unit_,file=str(file_),form='unformatted',access='sequential')
endif
!
do ilat=1,Nlat
do jlat=1,Nlat
do ispin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
if(use_formatted_)then
write(unit_,*)impGmatrix(ilat,jlat,ispin,ispin,iorb,jorb)
else
write(unit_)impGmatrix(ilat,jlat,ispin,ispin,iorb,jorb)
endif
enddo
enddo
enddo
enddo
enddo
close(unit_)
#else
print*,"Rear/write overloading requires Gfortran 6+"
#endif
end subroutine save_gfprime
!+-------------------------------------------------------------------+
!PURPOSE : Read cluster GF from file
!+-------------------------------------------------------------------+
subroutine read_gfprime(file,used,use_formatted)
character(len=*),optional :: file
character(len=256) :: file_
logical,optional :: used
logical :: used_
logical,optional :: use_formatted
logical :: use_formatted_
character(len=16) :: extension
integer :: unit_,Nchannel,Nexc,ichan,iexc,ilat,jlat,ispin,iorb,jorb
!
#if __GNUC__ > 6
if(allocated(impGmatrix))deallocate(impGmatrix)
allocate(impGmatrix(Nlat,Nlat,Nspin,Nspin,Norb,Norb))
used_=.false.;if(present(used))used_=used
use_formatted_=.false.;if(present(use_formatted))use_formatted_=use_formatted
extension=".restart";if(used_)extension=".used"
file_=str(str(file)//str(ed_file_suffix)//str(extension))
unit_=free_unit()
!
if(use_formatted_)then
open(unit_,file=str(file_),access='sequential')
else
open(unit_,file=str(file_),form='unformatted',access='sequential')
endif
!
rewind(unit_)
!
do ilat=1,Nlat
do jlat=1,Nlat
do ispin=1,Nspin
do iorb=1,Norb
do jorb=1,Norb
if(use_formatted_)then
read(unit_,*)impGmatrix(ilat,jlat,ispin,ispin,iorb,jorb)
else
read(unit_)impGmatrix(ilat,jlat,ispin,ispin,iorb,jorb)
endif
enddo
enddo
enddo
enddo
enddo
close(unit_)
#else
print*,"Rear/write overloading requires Gfortran 6+"
#endif
end subroutine read_gfprime
!##################################################################
!##################################################################
! ROUTINES TO SEARCH CHEMICAL POTENTIAL UP TO SOME ACCURACY
! can be used to fix any other *var so that *ntmp == nread
!##################################################################
!##################################################################
!+------------------------------------------------------------------+
!PURPOSE :
!+------------------------------------------------------------------+
subroutine ed_search_variable(var,ntmp,converged)
real(8),intent(inout) :: var
real(8),intent(in) :: ntmp
logical,intent(inout) :: converged
logical :: bool
real(8),save :: chich
real(8),save :: nold
real(8),save :: var_new
real(8),save :: var_old
real(8) :: var_sign
!
real(8) :: ndiff
integer,save :: count=0,totcount=0,i
integer :: unit
!
!check actual value of the density *ntmp* with respect to goal value *nread*
count=count+1
totcount=totcount+1
!
if(count==1)then
chich = ndelta !~0.2
inquire(file="var_compressibility.restart",EXIST=bool)
if(bool)then
open(free_unit(unit),file="var_compressibility.restart")
read(unit,*)chich
close(unit)
endif
var_old = var
endif
!
ndiff=ntmp-nread
!
!Get 'charge compressibility"
if(count>1)chich = (ntmp-nold)/(var-var_old)
!
!Add here controls on chich: not to be too small....
if(chich<0.1d0)chich=0.1d0*chich/abs(chich)
!
!update chemical potential
var_new = var - ndiff/chich
!
!
!re-define variables:
nold = ntmp
var_old = var
var = var_new
!
!Print information
write(LOGfile,"(A9,F16.9,A,F15.9)") "n = ",ntmp,"| instead of",nread
write(LOGfile,"(A9,ES16.9,A,ES16.9)")"dn = ",ndiff,"/",nerr
var_sign = (var-var_old)/abs(var-var_old)
if(var_sign>0d0)then
write(LOGfile,"(A9,ES16.9,A4)")"shift = ",ndiff/chich," ==>"
else
write(LOGfile,"(A9,ES16.9,A4)")"shift = ",ndiff/chich," <=="
end if
write(LOGfile,"(A9,F16.9)")"var = ",var
!
!Save info about search variable iteration:
open(free_unit(unit),file="search_variable_iteration_info"//reg(ed_file_suffix)//".ed",position="append")
if(count==1)write(unit,*)"#var,ntmp,ndiff"
write(unit,*)var,ntmp,ndiff
close(unit)
!
!If density is not converged set convergence to .false.
if(abs(ndiff)>nerr)converged=.false.
!
write(LOGfile,"(A18,I5)")"Search var count= ",count
write(LOGfile,"(A19,L2)")"Converged = ",converged
print*,""
!
open(free_unit(unit),file="var_compressibility.used")
write(unit,*)chich
close(unit)
!
end subroutine ed_search_variable
END MODULE ED_AUX_FUNX