File: bse_basic_structure.f90

package info (click to toggle)
espresso 6.7-4
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 311,068 kB
  • sloc: f90: 447,429; ansic: 52,566; sh: 40,631; xml: 37,561; tcl: 20,077; lisp: 5,923; makefile: 4,503; python: 4,379; perl: 1,219; cpp: 761; fortran: 618; java: 568; awk: 128
file content (920 lines) | stat: -rw-r--r-- 27,307 bytes parent folder | download | duplicates (3)
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
MODULE bse_basic_structures
!this module describes the basis structures 
!which are obtained from the DFT and GW code
  USE kinds, ONLY : DP

  REAL(kind=DP), ALLOCATABLE :: vg_q(:) ! contains the elements V(G) of the Coulomb potential obtained upon integration over q

  COMPLEX(kind=DP), ALLOCATABLE :: u_trans(:,:,:)!unitarian transformation from bloch wfcs to wannier'

  TYPE wannier_o
! this structures contains the overlap between the wannier square modules
     integer :: numb_v ! number of valence bands for the two spin channels 
     real(kind=dp),dimension(:,:), pointer ::o(:,:) ! overlap matrix (numb_v*numb_v)   
  END TYPE

  TYPE ii_mat
     integer :: numb_v! number of valence bands for the two spin channels 
     integer :: np_max ! maximum number of overlapping wannier orbitals for a given v
     integer, dimension (:,:), pointer :: iimat(:,:) ! (np_max,numb_v) the rows of this matrix contain for each iv, 
                                                     !the set of jv for which   o_mat(iv,jv)>=s_bse
 
  END TYPE

  TYPE vww_prod
!this type contains the v*wv*wv' products needed for the exchange part of the
!direct interaction term of the excitonic Hamiltonian
     integer :: numb_v! number of valence bands for the two spin channels 
     integer :: npw ! number of plane wave per processor
     integer :: np_max ! maximum number of overlapping wannier orbitals for a given v
     complex(kind=dp), dimension (:,:,:), pointer :: vww(:,:,:) ! v*wv*wv' product in G space (npw,np_max,numb_v)

  END TYPE
  
  TYPE bse_z
! this type contains the z terms to build up the Wc term of the excitonic Hamiltonian
!z_beta_v_v'=(v*phi_beta)*wv*Wv'
     integer :: numb_v! number of valence bands for the two spin channels 
     integer :: np_max ! maximum number of overlapping wannier orbitals for a given v
     integer :: numw_prod ! dimension of the polarizability basis
     real(kind=dp), dimension (:,:,:), pointer :: z(:,:,:) ! v*phi_beta*wv*wv' product (numw_prod,np_max,numb_v)
  END TYPE

  TYPE v_state
! this type contains the valence states wavefunctions and single particle energies
!
  integer :: nspin  ! number of spin channels
  integer :: numb_v(2) ! number valence state
  integer :: npw ! number of g-vectors per processor
  real(kind=dp), dimension (:,:),pointer  :: esp(:,:) ! single particle energies (numb_v,nspin) 
  complex(kind=dp), dimension(:,:,:), pointer :: wfn(:,:,:) ! wave function in G space (npw,numb_v,nspin)  
  integer ::gstart

  END TYPE

  TYPE v_state_r
! this type contains the valence states wfns in real space on the dual grid 
  integer :: nspin  ! number of spin channels
  integer :: numb_v(2) ! number of valence states
  integer :: nrxxt ! number of r points per processor
  real(kind=dp), dimension(:,:,:), pointer :: wfnrt(:,:,:) ! wave function in r-spce (dual grid) (nrxxt,numb_v,nspin)

  END TYPE

  TYPE c_state
! this type contains the valence states wavefunctions and single particle energies
!
  integer :: nspin  ! number of spin channels
  integer :: numb_c ! number valence state
  integer :: npw ! number of g-vectors per processor
  real(kind=dp), dimension (:,:),pointer  :: esp(:) ! single particle energies (numb_c) 
  complex(kind=dp), dimension(:,:), pointer :: wfn(:,:) ! wave function in G space (npw,numb_c)  
  integer ::gstart

  END TYPE

  TYPE c_state_r
! this type contains the valence states wfns in real space on the dual grid 
  integer :: nspin  ! number of spin channels
  integer :: numb_c ! number of valence states
  integer :: nrxxt ! number of r points per processor
  real(kind=dp), dimension(:,:), pointer :: wfnrt(:,:) ! wave function in r-spce (dual grid) (nrxxt,numb_c)

  END TYPE

  CONTAINS

      subroutine initialize_v_state_r(v_wfnr)
      implicit none
      type(v_state_r) :: v_wfnr
      nullify(v_wfnr%wfnrt)
      return
      end subroutine

      subroutine initialize_v_state(v_wfn)
      implicit none
      type(v_state) :: v_wfn
      nullify(v_wfn%wfn)
      nullify(v_wfn%esp)
      return
      end subroutine

      subroutine initialize_c_state_r(c_wfnr)
      implicit none
      type(c_state_r) :: c_wfnr
      nullify(c_wfnr%wfnrt)
      return
      end subroutine

      subroutine initialize_c_state(c_wfn)
      implicit none
      type(c_state) :: c_wfn
      nullify(c_wfn%wfn)
      nullify(c_wfn%esp)
      return
      end subroutine
  
      subroutine initialize_wannier_o(o)
      implicit none
      type(wannier_o) :: o
      nullify(o%o)
      return
      end subroutine
  
      subroutine initialize_imat(iimat)
      implicit none
      type(ii_mat) :: iimat
      nullify(iimat%iimat)
      return
      end subroutine

      subroutine initialize_vww_prod(vww)
      implicit none
      type(vww_prod) :: vww
      nullify(vww%vww)
      return
      end subroutine

      subroutine initialize_bse_z(z)
      implicit none
      type(bse_z) :: z
      nullify(z%z)
      return
      end subroutine

      subroutine free_v_state_r(v_wfnr)
      implicit none
      type(v_state_r) :: v_wfnr
      if(associated(v_wfnr%wfnrt)) deallocate (v_wfnr%wfnrt)
      nullify(v_wfnr%wfnrt)
      return
      end subroutine

      subroutine free_v_state(v_wfn)
      implicit none
      type(v_state) :: v_wfn
      if(associated(v_wfn%wfn)) deallocate (v_wfn%wfn)
      nullify(v_wfn%wfn)
      if(associated(v_wfn%esp)) deallocate (v_wfn%esp)
      nullify(v_wfn%esp)
      return
      end subroutine

      subroutine free_c_state_r(c_wfnr)
      implicit none
      type(c_state_r) :: c_wfnr
      if(associated(c_wfnr%wfnrt)) deallocate (c_wfnr%wfnrt)
      nullify(c_wfnr%wfnrt)
      return
      end subroutine

      subroutine free_c_state(c_wfn)
      implicit none
      type(c_state) :: c_wfn
      if(associated(c_wfn%wfn)) deallocate (c_wfn%wfn)
      nullify(c_wfn%wfn)
      if(associated(c_wfn%esp)) deallocate (c_wfn%esp)
      nullify(c_wfn%esp)
      return
      end subroutine
  
      subroutine free_wannier_o(o)
      implicit none
      type(wannier_o) :: o
      if(associated(o%o)) deallocate (o%o)
      nullify(o%o)
      return
      end subroutine
  
      subroutine free_imat(iimat)
      implicit none
      type(ii_mat) :: iimat
      if(associated(iimat%iimat)) deallocate (iimat%iimat)
      nullify(iimat%iimat)
      return
      end subroutine

      subroutine free_vww_prod(vww)
      implicit none
      type(vww_prod) :: vww
      if(associated(vww%vww)) deallocate (vww%vww)
      nullify(vww%vww)
      return
      end subroutine

      subroutine free_bse_z(z)
      implicit none
      type(bse_z) :: z
      if(associated(z%z)) deallocate (z%z)
      nullify(z%z)
      return
      end subroutine

      subroutine make_v_state(numb_v,v)
      use io_global, ONLY : stdout, ionode 
      USE gvect,                 ONLY : gstart
      USE lsda_mod,              ONLY : nspin
      use wavefunctions,  ONLY : evc
      use io_files,  ONLY : prefix, iunwfc, tmp_dir
      USE io_files, ONLY: nwordwfc
      USE wvfct,    ONLY : nbnd, npwx,npw,et
      use mp_world, ONLY : mpime
      USE mp,          ONLY :mp_barrier
      USE mp_world,             ONLY : world_comm

      implicit none

      type(v_state) :: v
      integer :: numb_v(2)
  
      integer :: is,ivmax,iv
      logical :: debug

      debug=.false.
     
      call start_clock('make_v_state')
     
      if(debug) then
         write(*,*) 'make_v_state: in, mpime=',mpime
         ! debug MARGHE
         write(*,*) 'nbnd=', nbnd
         write(*,*) 'numb_v(1)=', numb_v(1)
      endif


      v%nspin=nspin
      v%numb_v(:)=numb_v(:)
      v%npw=npw
      v%gstart=gstart


      allocate( evc( npwx, nbnd ) )
  
      if (nspin==1) then
         ivmax= v%numb_v(1)
      else 
         ivmax=max(v%numb_v(1),v%numb_v(2))
      endif
      


      allocate( v%wfn(v%npw,ivmax,v%nspin))
      allocate( v%esp(ivmax,v%nspin))


      do is=1,nspin
         call davcio(evc,2*nwordwfc,iunwfc,is,-1)
         do iv=1,v%numb_v(is)
            v%wfn(1:v%npw,1:v%numb_v(is),is)=evc(1:v%npw,1:v%numb_v(is))
         enddo  
            v%esp(1:v%numb_v(is),is)=et(1:v%numb_v(is),is)
      enddo

      deallocate(evc)

      if(debug) then
         write(*,*) 'make_v_state: out, mpime=',mpime
      endif

      call mp_barrier( world_comm )
      call stop_clock('make_v_state')

      return
      end subroutine

      subroutine make_c_state(numb_v,c)
      use io_global, ONLY : stdout, ionode 
      USE gvect,                 ONLY : gstart
      USE lsda_mod,              ONLY : nspin
      use wavefunctions,  ONLY : evc
      use io_files,  ONLY : prefix, iunwfc, tmp_dir
      USE io_files, ONLY: nwordwfc
      USE wvfct,    ONLY : nbnd, npwx,npw,et
      use mp_world, ONLY : mpime
      USE mp,          ONLY :mp_barrier
      USE mp_world,             ONLY : world_comm

      implicit none

      type(c_state) :: c
      integer :: numb_v(2)
  
      integer :: is,ic
      logical :: debug

      debug=.false.
     
      call start_clock('make_c_state')
     
      if(debug) then
         write(*,*) 'make_c_state: in, mpime=',mpime
         ! debug MARGHE
         write(*,*) 'nbnd=', nbnd
         write(*,*) 'numb_v(1)=', numb_v(1)
      endif


      c%nspin=nspin
      c%numb_c=nbnd-numb_v(1)
      c%npw=npw
      c%gstart=gstart


      allocate( evc( npwx, nbnd ) )
  
!      if (nspin==1) then
!         ivmax= v%numb_v(1)
!      else 
!         ivmax=max(v%numb_v(1),v%numb_v(2))
!      endif
      


      allocate( c%wfn(c%npw,c%numb_c))
      allocate( c%esp(c%numb_c))


      do is=1,nspin
         call davcio(evc,2*nwordwfc,iunwfc,is,-1)
         do ic=1,c%numb_c
            c%wfn(1:c%npw,1:c%numb_c)=evc(1:c%npw,numb_v(is)+1:nbnd)
         enddo  
            c%esp(1:c%numb_c)=et(numb_v(is)+1:nbnd,is)
      enddo

      deallocate(evc)

      if(debug) then
         write(*,*) 'make_c_state: out, mpime=',mpime
      endif

      call mp_barrier( world_comm )
      call stop_clock('make_c_state')

      return
      end subroutine

      subroutine c_times_cstate(v,cstate_in,cstate_out)
      ! this subroutine multiplies each line ic of the c_state vector by the ic real component of the v vector
      use kinds, only:DP
      !use bse_wannier, only: qpe_imin,qpe_imax

      implicit none

      type(c_state),intent(in) :: cstate_in
      type(c_state),intent(out) :: cstate_out

      integer :: ib
      real(kind=DP) :: v(cstate_in%numb_c)

      do ib=1,cstate_in%numb_c
         cstate_out%wfn(1:cstate_out%npw,ib)=cmplx(v(ib),0.d0)* cstate_out%wfn(1:cstate_in%npw,ib)
      enddo
     
      return
      end subroutine

 
      subroutine v_wfng_to_wfnr(vwfng,fc,vwfnr)
     !this subroutine FFT the valence wfns to real space in the dual grid

      USE kinds, ONLY : DP
      USE fft_custom_gwl
      USE bse_wannier, ONLY : dual_bse 
      USE wvfct,    ONLY : g2kin, npwx, npw, nbnd, nbndx
      USE io_global, ONLY : stdout, ionode, ionode_id
      USE mp_world, ONLY : mpime, nproc,world_comm
      USE mp_wave, ONLY : mergewf,splitwf
      USE mp,             ONLY : mp_sum
      USE gvect
      USE wavefunctions, ONLY :  psic

      implicit none

      type(v_state) vwfng
      type(v_state_r) vwfnr
      type(fft_cus) :: fc

      COMPLEX(kind=DP), allocatable :: vwfng_t(:,:,:)
      COMPLEX(kind=DP), ALLOCATABLE :: evc_g(:)

      integer :: ii,is
      integer ::ivmax  

      call start_clock('v_wfng_to_wfnr')
  
 
      if (vwfng%nspin==1) then
         ivmax= vwfng%numb_v(1)
      else 
         ivmax=max(vwfng%numb_v(1),vwfng%numb_v(2))
      endif


      allocate(vwfng_t(fc%npwt,ivmax,vwfng%nspin))

      vwfnr%nspin=vwfng%nspin
      vwfnr%nrxxt=fc%nrxxt       
      vwfnr%numb_v=vwfng%numb_v

      allocate(vwfnr%wfnrt(vwfnr%nrxxt,ivmax,vwfnr%nspin))
      
      allocate(evc_g(fc%ngmt_g ))

      if(fc%dual_t==4.d0) then
      do is=1,vwfng%nspin
         vwfng_t(1:fc%npwt,1:vwfng%numb_v(is),is)= vwfng%wfn(1:fc%npwt,1:vwfng%numb_v(is),is)
      enddo
      else
         do is=1,vwfng%nspin
           call reorderwfp_col(vwfng%numb_v(is),vwfng%npw,fc%npwt,vwfng%wfn(1,1,is),vwfng_t(1,1,is),vwfng%npw,&
                 & fc%npwt,ig_l2g,fc%ig_l2gt,fc%ngmt_g,mpime, nproc,world_comm )

           !do ii=1,vwfng%numb_v(is)
           !   call mergewf(vwfng%wfn(:,ii,is),evc_g,vwfng%npw,ig_l2g,mpime,nproc,ionode_id,intra_pool_comm)
           !   call splitwf(vwfng_t(:,ii,is),evc_g,fc%npwt,fc%ig_l2gt,mpime,nproc,ionode_id,intra_pool_comm)
           !enddo
         enddo
      endif

      do is=1,vwfng%nspin 
         do ii=1,vwfng%numb_v(is),2
            psic(1:fc%nrxxt)=(0.d0,0.d0)
            if (ii==vwfng%numb_v(is)) then
               psic(fc%nlt(1:fc%npwt))  = vwfng_t(1:fc%npwt,ii,is)
               psic(fc%nltm(1:fc%npwt)) = CONJG( vwfng_t(1:fc%npwt,ii,is) )
            else
               psic(fc%nlt(1:fc%npwt))=vwfng_t(1:fc%npwt,ii,is)+(0.d0,1.d0)*vwfng_t(1:fc%npwt,ii+1,is)
               psic(fc%nltm(1:fc%npwt)) =CONJG(vwfng_t(1:fc%npwt,ii,is))+(0.d0,1.d0)*CONJG(vwfng_t(1:fc%npwt,ii+1,is))
            endif
            CALL cft3t( fc, psic, fc%nr1t, fc%nr2t, fc%nr3t, fc%nrx1t, fc%nrx2t, fc%nrx3t, 2 )
            vwfnr%wfnrt(1:fc%nrxxt,ii,is)= DBLE(psic(1:fc%nrxxt))
            if(ii/=vwfng%numb_v(is)) vwfnr%wfnrt(1:fc%nrxxt,ii+1,is)= DIMAG(psic(1:fc%nrxxt))
         enddo
      enddo

      deallocate(evc_g)

      call stop_clock('v_wfng_to_wfnr')

      return
      end subroutine

      subroutine c_wfng_to_wfnr(cwfng,fc,cwfnr)
     !this subroutine FFT the conduction wfns to real space in the dual grid

      USE kinds, ONLY : DP
      USE fft_custom_gwl
      USE bse_wannier, ONLY : dual_bse, num_nbndv 
      USE wvfct,    ONLY : g2kin, npwx, npw, nbnd, nbndx
      USE io_global, ONLY : stdout, ionode, ionode_id
      USE mp_world, ONLY : mpime, nproc,world_comm
      USE mp_wave, ONLY : mergewf,splitwf
      USE mp,             ONLY : mp_sum
      USE gvect
      USE wavefunctions, ONLY :  psic

      implicit none

      type(c_state) cwfng
      type(c_state_r) cwfnr
      type(fft_cus) :: fc

      COMPLEX(kind=DP), allocatable :: cwfng_t(:,:)
      COMPLEX(kind=DP), ALLOCATABLE :: evc_g(:)

      integer :: ii,is
      integer ::icmax  

      call start_clock('c_wfng_to_wfnr')
  
 
!      if (vwfng%nspin==1) then
!         ivmax= vwfng%numb_v(1)
!      else 
!         ivmax=max(vwfng%numb_v(1),vwfng%numb_v(2))
!      endif
!       icmax=nbnd-num_nbndv(1)


      allocate(cwfng_t(fc%npwt,cwfng%numb_c))

      cwfnr%nrxxt=fc%nrxxt       
      cwfnr%numb_c=cwfng%numb_c

      allocate(cwfnr%wfnrt(cwfnr%nrxxt,cwfnr%numb_c))
      
      allocate(evc_g(fc%ngmt_g ))

      if(fc%dual_t==4.d0) then
         cwfng_t(1:fc%npwt,1:cwfng%numb_c)= cwfng%wfn(1:fc%npwt,1:cwfng%numb_c)
      else
         call reorderwfp_col(cwfng%numb_c,cwfng%npw,fc%npwt,cwfng%wfn(1,1),cwfng_t(1,1),cwfng%npw,&
                 & fc%npwt,ig_l2g,fc%ig_l2gt,fc%ngmt_g,mpime, nproc,world_comm )

      endif

      do ii=1,cwfng%numb_c,2
           psic(1:fc%nrxxt)=(0.d0,0.d0)
           if (ii==cwfng%numb_c) then
              psic(fc%nlt(1:fc%npwt))  = cwfng_t(1:fc%npwt,ii)
              psic(fc%nltm(1:fc%npwt)) = CONJG( cwfng_t(1:fc%npwt,ii) )
           else
              psic(fc%nlt(1:fc%npwt))=cwfng_t(1:fc%npwt,ii)+(0.d0,1.d0)*cwfng_t(1:fc%npwt,ii+1)
              psic(fc%nltm(1:fc%npwt)) =CONJG(cwfng_t(1:fc%npwt,ii))+(0.d0,1.d0)*CONJG(cwfng_t(1:fc%npwt,ii+1))
           endif
           CALL cft3t( fc, psic, fc%nr1t, fc%nr2t, fc%nr3t, fc%nrx1t, fc%nrx2t, fc%nrx3t, 2 )
           cwfnr%wfnrt(1:fc%nrxxt,ii)= DBLE(psic(1:fc%nrxxt))
           if(ii/=cwfng%numb_c) cwfnr%wfnrt(1:fc%nrxxt,ii+1)= DIMAG(psic(1:fc%nrxxt))
      enddo

      deallocate(evc_g)

      call stop_clock('c_wfng_to_wfnr')

      return
      end subroutine

      subroutine write_wfnr(wfnr)
      ! this subroutines writes on disk the type v_state_r for every processor
!      USE io_files,             ONLY : find_free_unit, prefix
      USE io_files,             ONLY :  tmp_dir,prefix
      USE mp_world,  ONLY : mpime
      implicit none

      INTEGER, EXTERNAL :: find_free_unit
      type(v_state_r) wfnr
      INTEGER :: iw, iunw,is
      CHARACTER(5) :: nproc

      iunw=find_free_unit()
      
      write(nproc,'(5i1)') &
           & mpime/10000,mod(mpime,10000)/1000,mod(mpime,1000)/100,mod(mpime,100)/10,mod(mpime,10)


      open( unit=iunw, file=trim(tmp_dir)//trim(prefix)//'.wfnr_t.'// nproc , status='unknown',form='unformatted')

      write(iunw) wfnr%numb_v  
      write(iunw) wfnr%nspin
      write(iunw) wfnr%nrxxt
      
      do is=1,wfnr%nspin
         do iw=1,wfnr%numb_v(is)
            write(iunw)  wfnr%wfnrt(1:wfnr%nrxxt,iw,is)
         enddo
      enddo 
      close(iunw)
      end subroutine

      subroutine read_wfnr(wfnr)
      ! this subroutines reads from disk the type v_state_r for every processor
!      USE io_files,             ONLY : find_free_unit, prefix
      USE io_files,             ONLY : tmp_dir, prefix
      USE mp_world,  ONLY : mpime
      implicit none

      INTEGER, EXTERNAL :: find_free_unit
      type(v_state_r) wfnr
      INTEGER :: iw, iunw,is
      CHARACTER(5) :: nproc

      iunw=find_free_unit()
      
      write(nproc,'(5i1)') &
           & mpime/10000,mod(mpime,10000)/1000,mod(mpime,1000)/100,mod(mpime,100)/10,mod(mpime,10)


      open( unit=iunw, file=trim(tmp_dir)//trim(prefix)//'.wfnr_t.'// nproc , status='old',form='unformatted')

      read(iunw) wfnr%numb_v  
      read(iunw) wfnr%nspin
      read(iunw) wfnr%nrxxt

      do is=1,wfnr%nspin
         do iw=1,wfnr%numb_v(is)
            read(iunw)  wfnr%wfnrt(1:wfnr%nrxxt,iw,is)
         enddo
      enddo
 
      close(iunw)
      end subroutine
	
      subroutine write_cwfnr(wfnr)
      ! this subroutines writes on disk the type v_state_r for every processor
!      USE io_files,             ONLY : find_free_unit, prefix
      USE io_files,             ONLY :  tmp_dir,prefix
      USE mp_world,  ONLY : mpime
      implicit none

      INTEGER, EXTERNAL :: find_free_unit
      type(c_state_r) wfnr
      INTEGER :: iw, iunw,is
      CHARACTER(5) :: nproc

      iunw=find_free_unit()
      
      write(nproc,'(5i1)') &
           & mpime/10000,mod(mpime,10000)/1000,mod(mpime,1000)/100,mod(mpime,100)/10,mod(mpime,10)


      open( unit=iunw, file=trim(tmp_dir)//trim(prefix)//'.cwfnr_t.'// nproc , status='unknown',form='unformatted')

      write(iunw) wfnr%numb_c  
      write(iunw) wfnr%nrxxt
      
      do iw=1,wfnr%numb_c
         write(iunw)  wfnr%wfnrt(1:wfnr%nrxxt,iw)
      enddo

      close(iunw)
      end subroutine

      subroutine read_cwfnr(wfnr)
      ! this subroutines reads from disk the type v_state_r for every processor
!      USE io_files,             ONLY : find_free_unit, prefix
      USE io_files,             ONLY : tmp_dir,prefix
      USE mp_world,  ONLY : mpime
      implicit none

      INTEGER, EXTERNAL :: find_free_unit
      type(c_state_r) wfnr
      INTEGER :: iw, iunw,is
      CHARACTER(5) :: nproc

      iunw=find_free_unit()
      
      write(nproc,'(5i1)') &
           & mpime/10000,mod(mpime,10000)/1000,mod(mpime,1000)/100,mod(mpime,100)/10,mod(mpime,10)


      open( unit=iunw, file=trim(tmp_dir)//trim(prefix)//'.cwfnr_t.'// nproc , status='old',form='unformatted')

      read(iunw) wfnr%numb_c  
      read(iunw) wfnr%nrxxt

      do iw=1,wfnr%numb_c
         read(iunw)  wfnr%wfnrt(1:wfnr%nrxxt,iw)
      enddo

      close(iunw)
      end subroutine
 

      subroutine read_omat(ispin,o)
      ! this subroutines reads the overlap matrix written by pw4gww
!      USE io_files,             ONLY : find_free_unit, prefix
      USE io_files,             ONLY : prefix,tmp_dir
      USE io_global,            ONLY : ionode, ionode_id
      USE mp,                   ONLY : mp_bcast
      USE kinds,                ONLY : DP
      USE mp_world,             ONLY : world_comm
      
      implicit none

      INTEGER, EXTERNAL :: find_free_unit

      type(wannier_o) :: o 
      integer ispin

      integer ii,iunu
      real(kind=DP) :: s_bse

      if(ionode) then
         iunu = find_free_unit()
         if (ispin==1) open(unit=iunu,file=trim(tmp_dir)//trim(prefix)//'.wbse1',status='old',form='unformatted')
         if (ispin==2) open(unit=iunu,file=trim(tmp_dir)//trim(prefix)//'.wbse2',status='old',form='unformatted')

         read(iunu) o%numb_v
         read(iunu) s_bse

         allocate(o%o(o%numb_v,o%numb_v))

         do ii=1,o%numb_v
            read(iunu) o%o(1:o%numb_v,ii)
         enddo
         close(iunu)
      endif
      
      CALL mp_bcast(o%numb_v, ionode_id , world_comm)
      if(.not.ionode) then     
        allocate(o%o(o%numb_v,o%numb_v))
      endif
      CALL mp_bcast(o%o, ionode_id, world_comm )

      return
      end subroutine

      subroutine read_iimat(iimat,ispin) 
      ! this subroutines reads the ii matrix written by pw4gww
!      USE io_files,             ONLY : find_free_unit, prefix
      USE io_files,             ONLY :  prefix, tmp_dir
      USE io_global,            ONLY : ionode, ionode_id
      USE mp,                   ONLY : mp_bcast
      USE mp_world,             ONLY : world_comm
      USE kinds,                ONLY : DP

      implicit none
      INTEGER, EXTERNAL :: find_free_unit
      type(ii_mat) :: iimat
      integer ispin

      real(kind=DP) :: s_bse
      integer       :: iv,iuni
      logical       :: debug

      debug=.false.

      if(ionode) then
        iuni = find_free_unit()
        if (ispin==1) open(unit=iuni,file=trim(tmp_dir)//trim(prefix)//'.iwwbse1',status='old',form='unformatted')
        if (ispin==2) open(unit=iuni,file=trim(tmp_dir)//trim(prefix)//'.iwwbse2',status='old',form='unformatted')
        read(iuni) iimat%numb_v
        read(iuni) s_bse
        read(iuni) iimat%np_max
        if(debug) then
          write(*,*) 'From read_iimat numb_v',iimat%numb_v
          write(*,*) 'From read_iimat s_bse', s_bse
          write(*,*) 'From read_iimat np_max',iimat%np_max
        endif
      endif

      CALL mp_bcast(iimat%numb_v, ionode_id, world_comm )
      CALL mp_bcast(iimat%np_max, ionode_id, world_comm )

      allocate(iimat%iimat(iimat%np_max,iimat%numb_v))

      if(ionode) then
        if(debug) then
           write(*,*) 'iimat matrix'  
        endif
        do iv=1, iimat%numb_v
           read(iuni) iimat%iimat(1:iimat%np_max,iv)
           if(debug) then
              write(*,*) 'iv=',iv, iimat%iimat(1:iimat%np_max,iv)
           endif
        enddo
        close(iuni)
      endif

      CALL mp_bcast(iimat%iimat, ionode_id, world_comm )

      return 
      end subroutine

      subroutine read_vww_prod(ispin,numb_v,npw,np_max,iimat,vww)
      !each processor reads the vww(G) written by pw4gww
      !be careful to check that the iimat that is passed to the subroutine is the related to the correct spin channel

!      USE io_files,             ONLY : find_free_unit, prefix,diropn
      USE io_files,             ONLY : prefix,diropn
      USE io_global, ONLY : stdout, ionode 

      implicit none
      INTEGER, EXTERNAL :: find_free_unit
      type(vww_prod) :: vww
      type(ii_mat)   :: iimat
      integer        :: numb_v,npw,np_max,ispin

      integer iv, ip, iungprod, ii,iundebug,i
      logical exst,debug

      debug=.false.
     
      if(debug) then
         iundebug = find_free_unit()
         open(iundebug,file='vww_bse.dat')
      endif      
 
      vww%numb_v=numb_v
      vww%npw=npw
      vww%np_max=np_max

      allocate(vww%vww(npw,np_max,numb_v))

      vww%vww(1:npw,1:np_max,1:numb_v)=dcmplx(0.d0,0.d0)
      
      iungprod = find_free_unit()
      if (ispin==1)  CALL diropn( iungprod, 'vww_bse1.',npw*2, exst)
      if (ispin==2)  CALL diropn( iungprod, 'vww_bse2.',npw*2, exst)

!      if(debug) then 
!         if(ionode) write(stdout,*) 'Read_vww_prod #1'
!      endif

      ii=0
      do iv=1,numb_v
         do ip=1, np_max 
            if(iimat%iimat(ip,iv)>0) then
!               if(debug) then 
!                  if(ionode) write(stdout,*) 'Read_vww_prod #', ii
!               endif
               ii=ii+1
               call davcio(vww%vww(:,ip,iv),npw*2,iungprod,ii,-1)
               if(debug) then
                  if(ionode) then
                     do i=1,npw
                        write(iundebug,*) vww%vww(i,ip,iv)   
                     enddo
                  endif
               endif
            endif
         enddo
      enddo

      close(iungprod)   
      if (debug) close(iundebug)   
      return
      end subroutine

      subroutine read_z(ispin,iimat,z)
      ! the ionode reads the z matrix and broadcast its value to the rest of the
      ! processors.
      !be careful to check that the iimat that is passed to the subroutine is the related to the correct spin channel


!      USE io_files,             ONLY : find_free_unit, prefix
      USE io_files,             ONLY :  prefix, tmp_dir
      USE io_global,            ONLY : ionode, ionode_id
      USE mp,                   ONLY : mp_bcast, mp_barrier
      USE mp_world, ONLY : world_comm
      USE kinds,                ONLY : DP
      USE io_global, ONLY : stdout,ionode


      implicit none
      INTEGER, EXTERNAL :: find_free_unit
      type(bse_z)    ::z 
      type(ii_mat)   :: iimat
!      integer        :: numw_prod
      integer        ::ispin

      real(kind=DP) :: s_bse
      integer       :: iv,iunz,ii

      logical debug

      debug=.false.

      if(ionode) then 
         iunz = find_free_unit()
         if(debug) then
           if(ionode) write(stdout,*) 'read_z ',trim(tmp_dir)//trim(prefix)//'.zbse1'
         endif


         if (ispin==1) open(unit=iunz,file=trim(tmp_dir)//trim(prefix)//'.zbse1',status='old',form='unformatted')
         if (ispin==2) open(unit=iunz,file=trim(tmp_dir)//trim(prefix)//'.zbse2',status='old',form='unformatted')
         read(iunz) z%numb_v
         read(iunz) s_bse 
         read(iunz) z%np_max
         read(iunz) z%numw_prod

         if(debug) then
           if(ionode) write(stdout,*) 'z%numb_v=', z%numb_v
           if(ionode) write(stdout,*) 's_bse=',s_bse
           if(ionode) write(stdout,*) 'z%np_max=',z%np_max
           if(ionode) write(stdout,*) 'z%numw_prod=', z%numw_prod
         endif

      endif 

      CALL mp_bcast(z%numb_v, ionode_id, world_comm )
      CALL mp_bcast(z%np_max, ionode_id, world_comm )
      CALL mp_bcast(z%numw_prod, ionode_id, world_comm )
      call mp_barrier(world_comm)
      
      allocate(z%z(z%numw_prod,z%np_max,z%numb_v))
     
      if(ionode) then
         do iv=1, z%numb_v
            do ii=1,z%np_max
               if(debug) then
                 if(ionode) write(stdout,*)'read_z, ii=',ii 
               endif
               if (iimat%iimat(ii,iv)>0) read(iunz) z%z(:,ii,iv)
            enddo
         enddo
      endif

      if(debug) then
        if(ionode) write(stdout,*) 'read_z #1'
      endif
!

      CALL mp_bcast(z%z, ionode_id, world_comm )
      call mp_barrier(world_comm)

      if(debug) then
        if(ionode) write(stdout,*) 'read_z #2'
      endif


      if(ionode) close(iunz)
      FLUSH( stdout )      
 
      return

      end subroutine

END MODULE