TABLE OF CONTENTS


ABINIT/m_plowannier [ Modules ]

[ Top ] [ Modules ]

NAME

  m_plowannier

FUNCTION

COPYRIGHT

 Copyright (C) 2006-2024 ABINIT group (BAmadon,AGerossier,ROuterovitch)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .

NOTES

SOURCE

17 #if defined HAVE_CONFIG_H
18 #include "config.h"
19 #endif
20 
21 #include "abi_common.h"
22 
23 MODULE m_plowannier
24 
25 
26 #ifndef HAVE_CRPA_OPTIM
27 #ifdef FC_INTEL
28 #warning "optimization of m_plowannier is deactivated on intel fortran"
29 !DEC$ NOOPTIMIZE
30 #endif
31 #endif
32 
33  use defs_basis
34  use m_errors
35  use m_abicore
36  use m_dtset
37  use m_dtfil
38  use defs_wvltypes
39  use m_xmpi
40 
41  use defs_datatypes, only : pseudopotential_type
42  use defs_abitypes, only : MPI_type
43  use m_io_tools,  only : open_file
44  use m_mpinfo,    only : proc_distrb_cycle
45  use m_crystal, only : crystal_t
46  use m_pawtab, only : pawtab_type
47  use m_pawcprj, only : pawcprj_type,pawcprj_alloc,pawcprj_get,pawcprj_free
48  use m_pawrad, only : pawrad_type, simp_gen
49 
50 
51  implicit none
52 
53  private
54 
55  public :: init_plowannier
56  public :: copy_orbital
57  public :: compute_coeff_plowannier
58  public :: destroy_plowannier
59  public :: print_plowannier
60  public :: get_plowannier
61  public :: fullbz_plowannier
62  public :: initialize_operwan
63  public :: destroy_operwan
64  public :: zero_operwan
65  public :: compute_oper_ks2wan
66  public :: normalization_plowannier
67  public :: print_operwan
68  public :: init_operwan_realspace
69  public :: reduce_operwan_realspace
70  public :: destroy_operwan_realspace
71  public :: zero_operwan_realspace
72  public :: compute_oper_wank2realspace

m_plowannier/allocate_orbital [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  allocate_orbital

FUNCTION

  allocate an array of orbital_type

INPUTS

  lorbital1

OUTPUT

  lorbital2

SOURCE

581 subroutine allocate_orbital(orbital1,orbital2,n1,n2,n3)
582 
583  !Arguments----------------
584  integer,intent(in) :: n1,n2,n3
585  type(orbital_type), intent(in) :: orbital1(n1,n2,n3)
586  type(orbital_type),intent(inout) :: orbital2(n1,n2,n3)
587 
588  !Local variable-----------
589  integer :: n4,n5,n6,n7
590  integer :: i,j,k,l
591 
592  do i = 1,n1
593    do j = 1,n2
594      do k = 1,n3
595        n4 = size(orbital1(i,j,k)%atom,1)
596        ABI_MALLOC(orbital2(i,j,k)%atom,(n4))
597        do l = 1,n4
598          n5 = size(orbital1(i,j,k)%atom(l)%matl,1)
599          n6 = size(orbital1(i,j,k)%atom(l)%matl,2)
600          n7 = size(orbital1(i,j,k)%atom(l)%matl,3)
601          ABI_MALLOC(orbital2(i,j,k)%atom(l)%matl,(n5,n6,n7))
602        end do
603      end do
604    end do
605  end do
606 
607 end subroutine allocate_orbital

m_plowannier/atom_index_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  atom_index_type_type

FUNCTION

SOURCE

216 type, public :: atom_index_type
217 
218   type(operwan_type), allocatable :: position(:,:)
219   ! size (number of positions chosen for atom1, number of positions chosen for atom2)
220 
221 end type atom_index_type

m_plowannier/compute_oper_ks2wan [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  compute_oper_ks2wan

FUNCTION

  transform ks operator into wan one

INPUTS

  wan,operks,option

OUTPUT

  if option = ikpt, gives the wan operator in reciprocal space (for each k)

SOURCE

2833  subroutine compute_oper_ks2wan(wan,operks,operwan,option)
2834 
2835    !Arguments--------------------------
2836    type(plowannier_type), intent(in) :: wan
2837    type(operwan_type), intent(inout) :: operwan(:,:,:)
2838    complex(dpc), intent(in) :: operks(:,:,:,:)
2839    integer, intent(in) :: option
2840 
2841    !Local variables--------------------
2842    integer :: iatom1, iatom2, il1, il2, isppol, ispinor1, ispinor2, iband1, iband2, im1, im2
2843 
2844    ! ----------------------------------
2845    !Transformation KS2WAN
2846    ! ----------------------------------
2847 
2848 
2849    !!operation on reciprocal space
2850    do iatom1 = 1,wan%natom_wan
2851      do iatom2 = 1,wan%natom_wan
2852        do isppol = 1,wan%nsppol
2853          do il1 = 1,wan%nbl_atom_wan(iatom1)
2854            do il2 = 1,wan%nbl_atom_wan(iatom2)
2855              do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
2856                do im2 = 1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
2857                  do ispinor1 = 1,wan%nspinor
2858                    do ispinor2 = 1,wan%nspinor
2859                      !!sum over the bands
2860                      do iband1 = 1,wan%bandf_wan-wan%bandi_wan+1
2861                        do iband2 = 1,wan%bandf_wan-wan%bandi_wan+1
2862                          operwan(option,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)=&
2863                            operwan(option,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)&
2864                           + conjg(wan%psichi(option,iband2,iatom2)%atom(il2)%matl(im2,isppol,ispinor2))&
2865                    *operks(option,iband1,iband2,isppol)*wan%psichi(option,iband1,iatom1)%atom(il1)%matl(im1,isppol,ispinor1)
2866                        end do
2867                      end do
2868                     ! write(6,*) "operwan",im1,im2,operwan(option,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)
2869                    end do
2870                  end do
2871                end do
2872              end do
2873            end do
2874          end do
2875        end do
2876      end do
2877    end do
2878 
2879 end subroutine compute_oper_ks2wan

m_plowannier/compute_oper_wank2realspace [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  compute_operwan_wanK2realspace

FUNCTION

 Compute an operator from WannierK space to real space

INPUTS

  wan,operwan, operwan_realspace

OUTPUT

 operwan_realspace

SOURCE

3641 subroutine compute_oper_wank2realspace(wan,operwan,operwan_realspace)
3642 
3643 !Arguments----------------------------------
3644   type(operwan_realspace_type),intent(inout) :: operwan_realspace
3645   type(operwan_type),intent(in) :: operwan(:,:,:)
3646   type(plowannier_type), intent(in) :: wan
3647 
3648 
3649 !Local variables----------------------------
3650   integer :: isppol,iatom1,pos1,il1,im1,iatom2,pos2,il2,im2,ikpt,ispinor1,ispinor2
3651 
3652 
3653   do isppol = 1,wan%nsppol
3654     do ispinor1=1,wan%nspinor
3655       do ispinor2=1,wan%nspinor
3656         do iatom1 = 1,wan%natom_wan
3657           do pos1 = 1,size(wan%nposition(iatom1)%pos,1)
3658             do il1 = 1,wan%nbl_atom_wan(iatom1)
3659               do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
3660                 do iatom2 = 1,wan%natom_wan
3661                   do pos2 = 1,size(wan%nposition(iatom2)%pos,1)
3662                     do il2 = 1,wan%nbl_atom_wan(iatom2)
3663                       do im2 = 1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
3664                        !sum over ikpt
3665                         do ikpt = 1,wan%nkpt
3666                           operwan_realspace%atom_index(iatom1,iatom2)%position(pos1,pos2)%&
3667                             &atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2) =&
3668                             operwan_realspace%atom_index(iatom1,iatom2)%position(pos1,pos2)%&
3669                             &atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)&
3670                             + real(operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)&
3671                             * wan%wtk(ikpt) * exp( cmplx(0.0,1.0) * two_pi * ( &
3672                             wan%kpt(1,ikpt) * ( wan%nposition(iatom1)%pos(pos1,1) - wan%nposition(iatom2)%pos(pos2,1) )+&
3673                             wan%kpt(2,ikpt) * ( wan%nposition(iatom1)%pos(pos1,2) - wan%nposition(iatom2)%pos(pos2,2) )+&
3674                             wan%kpt(3,ikpt) * ( wan%nposition(iatom1)%pos(pos1,3) - wan%nposition(iatom2)%pos(pos2,3)))))
3675                         end do
3676                        !end of the sum
3677                       end do
3678                     enddo
3679                   enddo
3680                 end do
3681               end do
3682             end do
3683           end do
3684         end do
3685       end do
3686     end do
3687   end do
3688 end subroutine compute_oper_wank2realspace

m_plowannier/copy_orbital [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  copy_orbital

FUNCTION

  Copy an array of orbital_type

INPUTS

  lorbital1

OUTPUT

  lorbital2

SOURCE

530 subroutine copy_orbital(orbital1,orbital2,n1,n2,n3)
531 
532  !Arguments----------------
533  integer,intent(in) :: n1,n2,n3
534  type(orbital_type), intent(in) :: orbital1(n1,n2,n3)
535  type(orbital_type),intent(inout) :: orbital2(n1,n2,n3)
536 
537  !Local variable-----------
538  integer :: n4,n5,n6,n7
539  integer :: i,j,k,l,m,p,q
540 
541  do i = 1,n1
542    do j = 1,n2
543      do k = 1,n3
544        n4 = size(orbital1(i,j,k)%atom,1)
545        do l = 1,n4
546          n5 = size(orbital1(i,j,k)%atom(l)%matl,1)
547          n6 = size(orbital1(i,j,k)%atom(l)%matl,2)
548          n7 = size(orbital1(i,j,k)%atom(l)%matl,3)
549          do m = 1,n5
550            do p = 1,n6
551              do q = 1,n7
552                orbital2(i,j,k)%atom(l)%matl(m,p,q) = orbital1(i,j,k)%atom(l)%matl(m,p,q)
553              end do
554            end do
555          end do
556        end do
557      end do
558    end do
559  end do
560 
561 end subroutine copy_orbital

m_plowannier/destroy_operwan [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  destroy_operwan

FUNCTION

  destroy operwan

INPUTS

  wan

OUTPUT

  operwan

SOURCE

2741  subroutine destroy_operwan(wan,operwan)
2742 
2743    !Arguments----------------------------------
2744    type(plowannier_type), intent(in) :: wan
2745    type(operwan_type), intent(inout) :: operwan(wan%nkpt,wan%natom_wan,wan%natom_wan)
2746 
2747    !Local variables----------------------------
2748    integer :: ikpt,iatom1,iatom2,il1,il2
2749 
2750 
2751    do ikpt = 1,wan%nkpt
2752      do iatom1 = 1,wan%natom_wan
2753        do iatom2 = 1,wan%natom_wan
2754          do il1 = 1,wan%nbl_atom_wan(iatom1)
2755            do il2 = 1,wan%nbl_atom_wan(iatom2)
2756              ABI_FREE(operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl)
2757            end do
2758          end do
2759          ABI_FREE(operwan(ikpt,iatom1,iatom2)%atom)
2760        end do
2761      end do
2762    end do
2763  end subroutine destroy_operwan

m_plowannier/destroy_operwan_realspace [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  destroy_operwan_realspace

FUNCTION

  Destroy an operwan_realspace type variable

INPUTS

  wan, operwan_realspace

OUTPUT

 operwan_realspace

SOURCE

3542 subroutine destroy_operwan_realspace(wan,operwan_realspace)
3543 
3544 !Arguments----------------------------------
3545   type(operwan_realspace_type),intent(inout) :: operwan_realspace
3546   type(plowannier_type), intent(in) :: wan
3547 
3548 !Local variables----------------------------
3549   integer :: iatom1,iatom2,pos1,pos2,il1,il2
3550 
3551 
3552  do iatom1 = 1,wan%natom_wan
3553    do iatom2 = 1,wan%natom_wan
3554      do pos1 = 1,size(wan%nposition(iatom1)%pos,1)
3555        do pos2 = 1,size(wan%nposition(iatom2)%pos,1)
3556          do il1 = 1,wan%nbl_atom_wan(iatom1)
3557            do il2 = 1,wan%nbl_atom_wan(iatom2)
3558             ABI_FREE(operwan_realspace%atom_index(iatom1,iatom2)%position(pos1,pos2)%atom(il1,il2)%matl)
3559            end do
3560          end do
3561          ABI_FREE(operwan_realspace%atom_index(iatom1,iatom2)%position(pos1,pos2)%atom)
3562        end do
3563      end do
3564      ABI_FREE(operwan_realspace%atom_index(iatom1,iatom2)%position)
3565    end do
3566  end do
3567  ABI_FREE(operwan_realspace%atom_index)
3568 
3569 end subroutine destroy_operwan_realspace

m_plowannier/destroy_orbital [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  destroy_orbital

FUNCTION

  destroy an array of orbital_type

INPUTS

  lorbital1

OUTPUT

  lorbital2

SOURCE

628 subroutine destroy_orbital(orbital2,n1,n2,n3)
629 
630  !Arguments----------------
631  integer,intent(in) :: n1,n2,n3
632  type(orbital_type),intent(inout) :: orbital2(n1,n2,n3)
633 
634  !Local variable-----------
635  integer :: n4
636  integer :: i,j,k,l
637 
638  do i = 1,n1
639    do j = 1,n2
640      do k = 1,n3
641        n4 = size(orbital2(i,j,k)%atom,1)
642        do l = 1,n4
643          ABI_FREE(orbital2(i,j,k)%atom(l)%matl)
644        end do
645        ABI_FREE(orbital2(i,j,k)%atom)
646      end do
647    end do
648  end do
649 
650 end subroutine destroy_orbital

m_plowannier/destroy_plowannier [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  destroy_plowannier

FUNCTION

  deallocate variables

INPUTS

  wan

OUTPUT

SOURCE

2612  subroutine destroy_plowannier(wan)
2613 
2614 !Arguments-------------------------------------
2615  type(plowannier_type), intent(inout) :: wan
2616 !Local variables-------------------------------
2617  integer :: iatom,ikpt,iband,il
2618 
2619  do iatom=1,wan%natom_wan
2620    ABI_FREE(wan%latom_wan(iatom)%lcalc)
2621    ABI_FREE(wan%projector_wan(iatom)%lproj)
2622    ABI_FREE(wan%nposition(iatom)%pos)
2623  enddo
2624  do iatom = 1,wan%natom_wan
2625    do il = 1,wan%nbl_atom_wan(iatom)
2626      ABI_FREE(wan%psichi(1,1,iatom)%atom(il)%ph0phiint)
2627    end do
2628  end do
2629  do ikpt = 1,wan%nkpt
2630    do iband = wan%bandi_wan,wan%bandf_wan
2631      do iatom = 1,wan%natom_wan
2632        do il = 1,wan%nbl_atom_wan(iatom)
2633         ABI_FREE(wan%psichi(ikpt,iband-wan%bandi_wan+1,iatom)%atom(il)%matl)
2634        end do
2635        ABI_FREE(wan%psichi(ikpt,iband-wan%bandi_wan+1,iatom)%atom)
2636      end do
2637    end do
2638  end do
2639 
2640  if (allocated(wan%kpt)) then
2641    ABI_FREE(wan%kpt)
2642  end if
2643  if (allocated(wan%iatom_wan)) then
2644    ABI_FREE(wan%iatom_wan)
2645  end if
2646  if (allocated(wan%nbl_atom_wan)) then
2647    ABI_FREE(wan%nbl_atom_wan)
2648  end if
2649  if (allocated(wan%latom_wan)) then
2650    ABI_FREE(wan%latom_wan)
2651  end if
2652  if (allocated(wan%nbproj_atom_wan)) then
2653    ABI_FREE(wan%nbproj_atom_wan)
2654  end if
2655  if (allocated(wan%projector_wan)) then
2656    ABI_FREE(wan%projector_wan)
2657  end if
2658  if (allocated(wan%position)) then
2659    ABI_FREE(wan%position)
2660  end if
2661  if (allocated(wan%wtk)) then
2662    ABI_FREE(wan%wtk)
2663  end if
2664  if (allocated(wan%acell)) then
2665    ABI_FREE(wan%acell)
2666  end if
2667 
2668 
2669  if (allocated(wan%nposition)) then
2670    ABI_FREE(wan%nposition)
2671  end if
2672  if (allocated(wan%psichi)) then
2673    ABI_FREE(wan%psichi)
2674  end if
2675 
2676 
2677  end subroutine destroy_plowannier

m_plowannier/fullbz_plowannier [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  fullbz_plowannier

FUNCTION

  Reconstruct the pischis on the full BZ

INPUTS

 dtset,kmesh,cryst,wanibz

OUTPUT

 wanbz

SOURCE

2504  subroutine fullbz_plowannier(dtset,kmesh,cryst,pawang,wanibz,wanbz)
2505 
2506    use m_abicore
2507    use m_specialmsg, only : wrtout
2508    use defs_abitypes
2509    use m_bz_mesh, only : kmesh_t
2510    use m_crystal, only : crystal_t
2511    use m_pawang, only  : pawang_type
2512    implicit none
2513 
2514 !Arguments-------------------------
2515    type(plowannier_type),intent(inout) :: wanibz
2516    type(plowannier_type),intent(out) :: wanbz
2517    type(dataset_type),intent(in) :: dtset
2518    type(kmesh_t),intent(in) :: kmesh
2519    type(crystal_t),intent(in) :: cryst
2520    type(pawang_type),intent(in) :: pawang
2521 !Local variables----------------------
2522    character(len=500) :: msg
2523    integer :: sym,iatom,spin,ik_bz,iband,ibandc,il,ispinor,im,ik_ibz,isym,itim
2524    integer ::  at_indx,indx, iat,m1,m2,l
2525    real(dp) :: kbz(3),wtk(kmesh%nbz)
2526 !*****************************************************************************************
2527 
2528    wtk=one
2529    sym=kmesh%nbz/kmesh%nibz
2530    call init_plowannier(wanibz%bandf_wan,wanibz%bandi_wan,dtset%plowan_compute,&
2531      &dtset%plowan_iatom,dtset%plowan_it,dtset%plowan_lcalc,dtset%plowan_natom,&
2532      &dtset%plowan_nbl,dtset%plowan_nt,dtset%plowan_projcalc,dtset%acell_orig,&
2533      &kmesh%bz,dtset%nimage,kmesh%nbz,dtset%nspinor,dtset%nsppol,wtk,dtset%dmft_t2g,wanbz)
2534 
2535    write(msg,'(a)')" Reconstruction of the full Brillouin Zone using data.plowann in the IBZ"
2536    call wrtout(std_out,msg,'COLL');call wrtout(ab_out,msg,'COLL')
2537    if (cryst%nsym==1) then
2538      do ik_bz=1,kmesh%nbz
2539        do iband=wanbz%bandi_wan,wanbz%bandf_wan
2540          ibandc=iband-wanbz%bandi_wan+1
2541          do iatom=1,wanbz%natom_wan
2542            do il=1,wanbz%nbl_atom_wan(iatom)
2543              do im=1,2*wanbz%latom_wan(iatom)%lcalc(il)+1
2544                do spin=1,wanbz%nsppol
2545                  do ispinor=1,wanbz%nspinor
2546                    if (kmesh%tabi(ik_bz)==1) then
2547                      wanbz%psichi(ik_bz,ibandc,iatom)%atom(il)%matl(im,spin,ispinor)=&
2548                        &wanibz%psichi(kmesh%tab(ik_bz),ibandc,iatom)%atom(il)%matl(im,spin,ispinor)
2549                    else if (kmesh%tabi(ik_bz)==-1) then
2550                      wanbz%psichi(ik_bz,ibandc,iatom)%atom(il)%matl(im,spin,ispinor)=&
2551                        &conjg(wanibz%psichi(kmesh%tab(ik_bz),ibandc,iatom)%atom(il)%matl(im,spin,ispinor))
2552                    endif
2553                  enddo
2554                enddo
2555              enddo
2556            enddo
2557          enddo
2558        enddo
2559      enddo
2560    else if (cryst%nsym>1) then
2561     do ik_bz=1,kmesh%nbz
2562       call kmesh%get_BZ_item(ik_bz,kbz,ik_ibz,isym,itim)
2563       do iatom=1,wanbz%natom_wan
2564         indx=cryst%indsym(4,isym,wanibz%iatom_wan(iatom))
2565 !Link beetween full list and wan list of atom
2566         do iat=1,wanbz%natom_wan
2567           if (indx==wanbz%iatom_wan(iat))then
2568             at_indx=iat
2569           end if
2570         end do
2571 !
2572         do spin=1,wanibz%nsppol
2573           do ispinor=1,wanibz%nspinor
2574             do il=1,wanibz%nbl_atom_wan(iatom)
2575               l=wanibz%latom_wan(iatom)%lcalc(il)
2576               do m1=1,2*l+1
2577                 do m2=1,2*l+1
2578                   do iband=wanibz%bandi_wan,wanibz%bandf_wan
2579                     ibandc=iband-wanibz%bandi_wan+1
2580                     wanbz%psichi(ik_bz,ibandc,iatom)%atom(il)%matl(m1,spin,ispinor)=&
2581                       &wanbz%psichi(ik_bz,ibandc,iatom)%atom(il)%matl(m1,spin,ispinor)+&
2582                       &wanibz%psichi(ik_ibz,ibandc,at_indx)%atom(il)%matl(m2,spin,ispinor)&
2583                       &*pawang%zarot(m2,m1,l+1,isym)
2584                   end do
2585                 enddo
2586               enddo
2587             enddo
2588           enddo
2589         enddo
2590       enddo
2591     enddo
2592   end if
2593   call destroy_plowannier(wanibz)
2594 end subroutine fullbz_plowannier

m_plowannier/get_plowannier [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  get_plowannier

FUNCTION

  get the psichies (Wannier weights) from a data.plowann file

INPUTS

 wan

OUTPUT

 wan

SOURCE

2408  subroutine get_plowannier(wan_in,wan_out,dtset)
2409 
2410  use m_abicore
2411  use defs_abitypes
2412  use m_io_tools,  only : open_file
2413  use m_specialmsg, only : wrtout
2414  implicit none
2415 
2416  !Arguments-------------------------
2417  type(plowannier_type),intent(inout) :: wan_in
2418  type(plowannier_type),intent(inout) :: wan_out
2419  type(dataset_type),intent(in) :: dtset
2420  !Local variables-------------------
2421  character(len=500) :: msg
2422  integer :: unt,iatom,spin,ikpt,iband,ibandc,il,ispinor,im,dummy,natom,bandi,bandf,nbl,nspin,nkpt
2423   integer :: t2g
2424  real(dp) ::xx,yy
2425 
2426 t2g=dtset%dmft_t2g
2427 
2428 
2429  !Opening of the data.plowann file
2430  if (open_file('data.plowann',msg,newunit=unt,form='formatted',status='old') /= 0) then
2431   ABI_ERROR(msg)
2432  end if
2433  rewind(unt)
2434 
2435 
2436  !Reading of the header of data.plowann
2437  read(unt,'(a22,i2)') msg, natom
2438  read(unt,*)
2439  read(unt,'(a7,2i4)') msg, bandi,bandf
2440  read(unt,'(a26,i2)') msg, nbl
2441  do iatom=1,wan_in%natom_wan
2442    read(unt,*)
2443  enddo
2444  read(unt,'(a16,i2)') msg, nspin
2445  read(unt,'(a19,i4)') msg, nkpt
2446 
2447  !Testing the header
2448  if (natom /= wan_in%natom_wan .OR.&
2449 & nbl/= sum(wan_in%nbl_atom_wan(:)) .OR. nspin /= wan_in%nsppol .OR. nkpt/=wan_in%nkpt ) then
2450    write(msg,'(a,3i3)')"Not the same atoms or bands in both datasets",natom,bandi,bandf
2451    ABI_ERROR(msg)
2452  endif
2453 
2454  call init_plowannier(bandf,bandi,dtset%plowan_compute,&
2455      &dtset%plowan_iatom,dtset%plowan_it,dtset%plowan_lcalc,dtset%plowan_natom,&
2456      &dtset%plowan_nbl,dtset%plowan_nt,dtset%plowan_projcalc,dtset%acell_orig,&
2457      &dtset%kptns,dtset%nimage,dtset%nkpt,dtset%nspinor,dtset%nsppol,dtset%wtk,dtset%dmft_t2g,wan_out)
2458 
2459  call destroy_plowannier(wan_in)
2460 
2461  write(msg,'(a)')"Reading of the Wannier weights from data.plowann"
2462  call wrtout(std_out,msg,'COLL')
2463  call wrtout(ab_out,msg,'COLL')
2464  !Reading of the psichis
2465  do ikpt=1,wan_out%nkpt
2466    read(unt,*)
2467     do spin=1,wan_out%nsppol
2468       do ispinor=1,wan_out%nspinor
2469         do iband=wan_out%bandi_wan,wan_out%bandf_wan
2470           ibandc=iband-wan_out%bandi_wan+1
2471           read(unt,*)
2472           do iatom=1,wan_out%natom_wan
2473             do il=1,wan_out%nbl_atom_wan(iatom)
2474               do im=1,2*wan_out%latom_wan(iatom)%lcalc(il)+1
2475                 read(unt,'(8x,3i3,2x,2f23.15)')dummy,dummy,dummy,xx,yy
2476                 wan_out%psichi(ikpt,ibandc,iatom)%atom(il)%matl(im,spin,ispinor)=cmplx(xx,yy)
2477               enddo!m
2478             enddo!l
2479           enddo!atom
2480         enddo!band
2481       enddo!spinor
2482     enddo!spin
2483   enddo!k-point
2484  close(unt)
2485 end subroutine get_plowannier

m_plowannier/init_operwan_realspace [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  init_operwan_realspace

FUNCTION

  Initialize an operwan_realspace type variable

INPUTS

  wan, operwan_realspace

OUTPUT

 operwan_realspace

SOURCE

3360 subroutine init_operwan_realspace(wan,oprs)
3361 
3362 !Arguments----------------------------------
3363   type(operwan_realspace_type),intent(inout) :: oprs
3364   type(plowannier_type), intent(in) :: wan
3365 
3366 !Local variables----------------------------
3367   integer :: i1,i2,n1,n2,p1,p2,l1,l2,sp,pi
3368 
3369  !variable names is shorten to achieve not too long line lenght
3370   sp=wan%nsppol
3371   pi=wan%nspinor
3372   ABI_MALLOC(oprs%atom_index,(wan%natom_wan,wan%natom_wan))
3373   do i1 = 1,wan%natom_wan
3374     do i2 = 1,wan%natom_wan
3375       n1=size(wan%nposition(i1)%pos,1)
3376       n2=size(wan%nposition(i2)%pos,1)
3377       ABI_MALLOC(oprs%atom_index(i1,i2)%position,(n1,n2))
3378       do p1 = 1,size(wan%nposition(i1)%pos,1)
3379         do p2 = 1,size(wan%nposition(i2)%pos,1)
3380           n1=wan%nbl_atom_wan(i1)
3381           n2=wan%nbl_atom_wan(i2)
3382           ABI_MALLOC(oprs%atom_index(i1,i2)%position(p1,p2)%atom,(n1,n2))
3383           do l1 = 1,wan%nbl_atom_wan(i1)
3384             do l2 = 1,wan%nbl_atom_wan(i2)
3385               n1=2*wan%latom_wan(i1)%lcalc(l1)+1
3386               n2=2*wan%latom_wan(i2)%lcalc(l2)+1
3387               ABI_MALLOC(oprs%atom_index(i1,i2)%position(p1,p2)%atom(l1,l2)%matl,(n1,n2,sp,pi,pi))
3388               oprs%atom_index(i1,i2)%position(p1,p2)%atom(l1,l2)%matl = czero
3389             end do
3390          end do
3391        end do
3392      end do
3393    end do
3394  end do
3395 
3396 end subroutine init_operwan_realspace

m_plowannier/initialize_operwan [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  initialize_operwan

FUNCTION

  initialize operwan

INPUTS

  wan

OUTPUT

  operwan

SOURCE

2697  subroutine initialize_operwan(wan,operwan)
2698 
2699    !Arguments----------------------------------
2700    type(plowannier_type), intent(in) :: wan
2701    type(operwan_type), intent(inout) :: operwan(wan%nkpt,wan%natom_wan,wan%natom_wan)
2702 
2703    !Local variables----------------------------
2704    integer :: ikpt,iatom1,iatom2,il1,il2,n1,n2
2705 
2706    do ikpt = 1,wan%nkpt
2707      do iatom1 = 1,wan%natom_wan
2708        do iatom2 = 1,wan%natom_wan
2709          ABI_MALLOC(operwan(ikpt,iatom1,iatom2)%atom,(wan%nbl_atom_wan(iatom1),wan%nbl_atom_wan(iatom2)))
2710          do il1 = 1,wan%nbl_atom_wan(iatom1)
2711            do il2 = 1,wan%nbl_atom_wan(iatom2)
2712              n1=2*wan%latom_wan(iatom1)%lcalc(il1)+1
2713              n2=2*wan%latom_wan(iatom2)%lcalc(il2)+1
2714    ABI_MALLOC(operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl,(n1,n2,wan%nsppol,wan%nspinor,wan%nspinor))
2715              operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl = zero
2716            end do
2717          end do
2718        end do
2719      end do
2720    end do
2721 
2722  end subroutine initialize_operwan

m_plowannier/latom_wan_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  latom_wan_type

FUNCTION

SOURCE

85 type, public :: latom_wan_type
86 
87   integer, allocatable :: lcalc(:)
88   ! array of the l we want to compute the psichi with
89 
90 end type latom_wan_type

m_plowannier/lorbital2_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  lorbital2_type

FUNCTION

SOURCE

178 type, public :: lorbital2_type
179 
180   complex(dpc), allocatable :: matl(:,:,:,:,:)
181   ! size (2l1+1,2l2+1,nspppol,nspinor,nspinor)
182 
183   real(dp), allocatable :: ph0phiint(:)
184   ! size (nproj), stocks the value of ph0phiint we may want
185 
186 
187 end type lorbital2_type

m_plowannier/lorbital_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  lorbital_type

FUNCTION

SOURCE

138 type, public :: lorbital_type
139 
140   complex(dpc), allocatable :: matl(:,:,:)
141   !details for different m
142 
143   real(dp), allocatable :: ph0phiint(:)
144   ! stocks the values for each projector of the l considered
145   ! size(total number of projectors for this l)
146 
147 end type lorbital_type

m_plowannier/normalization_plowannier [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  normalization_plowannier

FUNCTION

  Use compute_oper_ks2wan to calculate overlap and do the normalization for the wan%psichi coefficients

INPUTS

  wan, opt (=0 normalize k-point by k-point; =1 normalize the sum over k)

OUTPUT

  wan itself is modified

SOURCE

2900 subroutine normalization_plowannier(wan,opt)
2901 
2902   use m_matrix, only : invsqrt_matrix
2903 
2904 !Arguments------------------
2905   type(plowannier_type), intent(inout) :: wan
2906   integer, intent(in) :: opt
2907 !Local----------------------
2908   complex(dpc), allocatable :: operks(:,:,:,:)
2909   type(operwan_type), allocatable :: operwan(:,:,:)
2910   complex(dpc), allocatable :: operwansquare(:,:,:,:)
2911   complex(dpc), allocatable :: tmp_operwansquare(:,:)
2912   integer :: ikpt, iband, iband1, iband2, isppol,  ispinor1, ispinor2, iatom1,nb_zeros_tot
2913   integer :: iatom2, il1, il2, im1, im2, index_c, index_l, n1,n2,n3, nkpt,nb_of_zeros
2914   type(orbital_type), allocatable :: psichinormalized(:,:,:)
2915   !character(len = 50) :: mat_writing2
2916   !character(len = 5000) :: mat_writing
2917   character(len = 500) :: message
2918 
2919   !Initialize nkpt (wan%nkpt if opt=0, 1 if opt=1)
2920   if (opt==1) then
2921     nkpt=1
2922   else
2923     nkpt=wan%nkpt
2924   end if
2925 
2926   !First, creation of the ks identity operator
2927   ABI_MALLOC(operks,(wan%nkpt,wan%bandf_wan-wan%bandi_wan+1,wan%bandf_wan-wan%bandi_wan+1,wan%nsppol))
2928   operks = czero
2929   do iband1 = 1,wan%bandf_wan-wan%bandi_wan+1
2930     do iband2 = 1,wan%bandf_wan-wan%bandi_wan+1
2931       if (iband1.eq.iband2) then
2932         do ikpt = 1,wan%nkpt
2933           do isppol= 1,wan%nsppol
2934             operks(ikpt,iband1,iband2,isppol) = cone
2935           end do
2936         end do
2937       end if
2938     end do
2939   end do
2940 
2941 
2942   !Allocation of operwan
2943   ABI_MALLOC(operwan,(wan%nkpt,wan%natom_wan,wan%natom_wan))
2944   call initialize_operwan(wan,operwan)
2945 
2946 
2947 
2948   !Computation of the overlap
2949   do ikpt = 1,wan%nkpt
2950     call compute_oper_ks2wan(wan,operks,operwan,ikpt)
2951   end do
2952 
2953 
2954 
2955   !transform the operwan into an inversible matrix
2956   !!operwansquare is the overlap square matrix (wan%size_wan * wan%size_wan)
2957   ABI_MALLOC(operwansquare,(wan%nkpt,wan%nsppol,wan%nspinor*wan%size_wan,wan%nspinor*wan%size_wan))
2958 
2959   operwansquare = czero
2960 
2961   n1=size(wan%psichi,1)
2962   n2=size(wan%psichi,2)
2963   n3=size(wan%psichi,3)
2964   ABI_MALLOC(psichinormalized,(n1,n2,n3))
2965   call allocate_orbital(wan%psichi,psichinormalized,n1,n2,n3)
2966   call copy_orbital(wan%psichi,psichinormalized,n1,n2,n3)
2967 
2968   do isppol = 1,wan%nsppol
2969     do ispinor1 = 1,wan%nspinor
2970       do ispinor2 = 1,wan%nspinor
2971         index_l = 0 ! index_l is set to 0 at the beginning
2972         do iatom1 = 1,wan%natom_wan
2973           do il1 = 1,wan%nbl_atom_wan(iatom1)
2974             do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
2975               index_l = index_l + 1 ! the line changes
2976               index_c = 1 ! counter_c is set to one each time the line changes
2977               do iatom2 = 1,wan%natom_wan
2978                 do il2 = 1,wan%nbl_atom_wan(iatom2)
2979                   do im2 = 1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
2980                     do ikpt = 1,wan%nkpt
2981                       if (opt==1) then ! operwansquare is a sum of operwan over k points
2982                         operwansquare(1,isppol,index_l+wan%size_wan*(ispinor1-1),index_c+wan%size_wan*(ispinor2-1))  &
2983 &                        =operwansquare(1,isppol,index_l+wan%size_wan*(ispinor1-1),index_c+wan%size_wan*(ispinor2-1))+ &
2984 &                        (operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)*wan%wtk(ikpt))
2985                       else !operwansquare is a matrix with a k-point dimension
2986                         operwansquare(ikpt,isppol,index_l+wan%size_wan*(ispinor1-1),index_c+wan%size_wan*(ispinor2-1))  &
2987                           &= operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)
2988                       end if
2989                     end do
2990                     index_c = index_c + 1
2991                   end do !im2
2992                 end do !il2
2993               end do ! iatom2 (the line changes)
2994             end do ! im1
2995           end do ! il1
2996         end do !iatom1
2997       end do
2998     end do
2999   end do
3000 
3001 
3002 
3003 
3004   !!Write the overlap matrix for ikpt = 1 in a nice shape
3005 !    do isppol = 1,wan%nsppol
3006 !      do il1 = 1,size(operwansquare,3) !! dummy variable without any meaning
3007 !        write(mat_writing,'(a,i0,i0)') 'Overlap matrix before orthonormalization 1 ',isppol,il1
3008 !        do il2 = 1,size(operwansquare,4)
3009 !         write(mat_writing2,'(F10.6)') real(operwansquare(1,isppol,il1,il2))
3010 !          mat_writing = trim(mat_writing)//trim(mat_writing2)
3011 !        end do
3012 !        print*,trim(mat_writing)
3013 !      end do
3014 !    end do
3015 
3016 
3017 
3018   !take the square root inverse of operwansquare for normalization purposes
3019   nb_zeros_tot=0
3020   ABI_MALLOC(tmp_operwansquare,(wan%nspinor*wan%size_wan,wan%nspinor*wan%size_wan))
3021   do isppol = 1,wan%nsppol
3022     do ikpt = 1,nkpt
3023       write(std_out,*)"ikpt = ", ikpt
3024       tmp_operwansquare(:,:)=operwansquare(ikpt,isppol,:,:)
3025       call invsqrt_matrix(tmp_operwansquare,wan%nspinor*wan%size_wan,nb_of_zeros)
3026       operwansquare(ikpt,isppol,:,:)=tmp_operwansquare(:,:)
3027       nb_zeros_tot=nb_zeros_tot+nb_of_zeros
3028     end do
3029   end do
3030   ABI_FREE(tmp_operwansquare)
3031 
3032   do ikpt = 1,wan%nkpt
3033     do iband = 1,wan%bandf_wan-wan%bandi_wan+1
3034       do iatom1 = 1,wan%natom_wan
3035         do il1 = 1,wan%nbl_atom_wan(iatom1)
3036           psichinormalized(ikpt,iband,iatom1)%atom(il1)%matl = czero
3037         end do
3038       end do
3039     end do
3040   end do
3041 
3042 
3043 
3044   ! compute the new psichi normalized
3045   do isppol = 1,wan%nsppol
3046     do ispinor1 = 1,wan%nspinor
3047       do iband = 1,wan%bandf_wan-wan%bandi_wan+1
3048         index_l = 0
3049         do iatom1 = 1,wan%natom_wan
3050           do il1 = 1,wan%nbl_atom_wan(iatom1)
3051             do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
3052               ! sum
3053               do ispinor2 = 1,wan%nspinor
3054                 index_l = index_l + 1 ! the line changes
3055                 index_c = 1 ! when the line changes, index_c is set to 1
3056                 do iatom2 = 1,wan%natom_wan
3057                   do il2 = 1,wan%nbl_atom_wan(iatom2)
3058                     do im2 = 1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
3059                       do ikpt = 1,wan%nkpt
3060                         if (opt==1)then ! all the psichi are normalized with the same operwansquare
3061                           psichinormalized(ikpt,iband,iatom1)%atom(il1)%matl(im1,isppol,ispinor1) =&
3062 &                           psichinormalized(ikpt,iband,iatom1)%atom(il1)%matl(im1,isppol,ispinor1) +&
3063 &                           wan%psichi(ikpt,iband,iatom2)%atom(il2)%matl(im2,isppol,ispinor2)*&
3064 &                           operwansquare(1,isppol,index_l+wan%size_wan*(ispinor1-1),index_c+&
3065 &                           wan%size_wan*(ispinor2-1))
3066                         else ! each psichi is normalized with his own operwansquare
3067                           psichinormalized(ikpt,iband,iatom1)%atom(il1)%matl(im1,isppol,ispinor1) =&
3068 &                           psichinormalized(ikpt,iband,iatom1)%atom(il1)%matl(im1,isppol,ispinor1) +&
3069 &                           wan%psichi(ikpt,iband,iatom2)%atom(il2)%matl(im2,isppol,ispinor2)*&
3070 &                           operwansquare(ikpt,isppol,index_l+wan%size_wan*(ispinor1-1),index_c+&
3071 &                           wan%size_wan*(ispinor2-1))
3072                         end if
3073                       end do
3074                       index_c = index_c + 1
3075                     end do !im2
3076                   end do !il2
3077                 end do !iatom2
3078               end do ! ispinor2
3079             end do !im1
3080           end do ! il1
3081         end do ! iatom1
3082       end do ! iband
3083     end do ! ispinor1
3084   end do !isppol
3085   ! copy the new psichi normalized
3086   call copy_orbital(psichinormalized,wan%psichi,n1,n2,n3)
3087   call destroy_orbital(psichinormalized,n1,n2,n3)
3088   ABI_FREE(psichinormalized)
3089 
3090   call destroy_operwan(wan,operwan)
3091   ABI_FREE(operwan)
3092 
3093 
3094 
3095 
3096 
3097 
3098 
3099 
3100 !!
3101 !!  !-------------------------------------------------------------
3102 !!  !check if the new norm is one
3103 
3104   ABI_MALLOC(operwan,(wan%nkpt,wan%natom_wan,wan%natom_wan))
3105   call initialize_operwan(wan,operwan)
3106   do ikpt = 1,wan%nkpt
3107     call compute_oper_ks2wan(wan,operks,operwan,ikpt)
3108   end do
3109 
3110 
3111   do isppol = 1,wan%nsppol
3112     do ikpt = 1,wan%nkpt
3113       do iatom1 = 1,wan%natom_wan
3114         do iatom2 = 1,wan%natom_wan
3115           do il1 = 1,wan%nbl_atom_wan(iatom1)
3116             do il2 = 1,wan%nbl_atom_wan(iatom2)
3117               do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
3118                 do im2 = 1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
3119                   do ispinor1 = 1,wan%nspinor
3120                     do ispinor2 = 1,wan%nspinor
3121                       if (opt==0 .and. nb_zeros_tot==0) then
3122                         if (iatom1.eq.iatom2 .and. il1.eq.il2 .and. im1.eq.im2 .and. ispinor1.eq.ispinor2) then
3123                           if (abs(cmplx(1.0,0.0,dpc)-&
3124                             &operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%&
3125                             &matl(im1,im2,isppol,ispinor1,ispinor2)) > 1d-8) then
3126                             write(message,'(a,i0,a,F18.11)') 'Normalization error for ikpt =',ikpt,&
3127                               &' on diag, value = ',&
3128                               &abs(operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2))
3129                             ABI_ERROR(message)
3130                           end if
3131                         else
3132                           if (abs(operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)) > 1d-8) then
3133                             write(message,'(a,i0,a,F10.3)') 'Normalization error for ikpt =',ikpt,&
3134                               &' not on diag, value = ',&
3135                               &abs(operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2))
3136                             ABI_ERROR(message)
3137                           end if
3138                         end if
3139                       end if
3140                     end do
3141                   end do
3142                 end do
3143               end do
3144             end do
3145           end do
3146         end do
3147       end do
3148     end do
3149   end do
3150   if (opt==0 .and. nb_zeros_tot/=0) then
3151     write(message,'(a,i2,a)')"The matrix inversion detects ",nb_zeros_tot,&
3152     " zero(s) on the diagonals. Take results with caution or modify nkpt and/or bands for plowan"
3153     ABI_COMMENT(message)
3154   end if
3155 
3156   !!Uncomment to print the overlap matrix in the log file (for ikpt = 1)
3157 !  do isppol = 1,wan%nsppol
3158 !    do ispinor1 = 1,wan%nspinor
3159 !      do ispinor2 = 1,wan%nspinor
3160 !        index_l = 0 ! index_l is set to 0 at the beginning
3161 !        do iatom1 = 1,wan%natom_wan
3162 !          do il1 = 1,wan%nbl_atom_wan(iatom1)
3163 !            do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
3164 !              index_l = index_l + 1 ! the line changes
3165 !              index_c = 1 ! counter_c is set to one each time the line changes
3166 !              do iatom2 = 1,wan%natom_wan
3167 !                do il2 = 1,wan%nbl_atom_wan(iatom2)
3168 !                  do im2 = 1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
3169 !                    do ikpt = 1,wan%nkpt
3170 !                      operwansquare(ikpt,isppol,index_l+wan%size_wan*(ispinor1-1),&
3171 !            &index_c+wan%size_wan*(ispinor2-1)) = operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)
3172 !                    end do
3173 !                    index_c = index_c + 1
3174 !                  end do !im2
3175 !                end do !il2
3176 !              end do ! iatom2 (the line changes)
3177 !            end do ! im1
3178 !          end do ! il1
3179 !        end do !iatom1
3180 !      end do
3181 !    end do
3182 !  end do
3183 
3184   !!Print the overlap matrix in a nice shape
3185 !  do isppol = 1,wan%nsppol
3186 !    do il1 = 1,size(operwansquare,3) !! dummy variable without any meaning
3187 !      write(mat_writing,'(a,i0,i0)') 'Overlap matrix after orthonormalization ',isppol,il1
3188 !      do il2 = 1,size(operwansquare,4)
3189 !        write(mat_writing2,'(F10.6)') real(operwansquare(9,isppol,il1,il2))
3190 !        mat_writing = trim(mat_writing)//trim(mat_writing2)
3191 !      end do
3192 !      print*,trim(mat_writing)
3193 !    end do
3194 !  end do
3195 
3196 
3197 
3198 
3199 !!  !----------------------------------------------------------------
3200   ABI_FREE(operwansquare)
3201   ABI_FREE(operks)
3202   call destroy_operwan(wan,operwan)
3203   ABI_FREE(operwan)
3204 
3205 
3206 end subroutine normalization_plowannier

m_plowannier/operwan_realspace_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  operwan_realspace_type

FUNCTION

SOURCE

233 type, public :: operwan_realspace_type
234 
235   type(atom_index_type), allocatable :: atom_index(:,:)
236   ! size (number of atom, number of atom)
237 
238 end type operwan_realspace_type

m_plowannier/operwan_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  operwan_type

FUNCTION

SOURCE

199 type, public :: operwan_type
200 
201   type(lorbital2_type), allocatable :: atom(:,:)
202   ! l chosen for each on of both atoms
203 
204 end type operwan_type

m_plowannier/orbital_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  orbital_type

FUNCTION

SOURCE

159 type, public :: orbital_type
160 
161   type(lorbital_type), allocatable :: atom(:)
162   ! details of the psichi coefficients for each atom
163   ! size of number of l chosen
164 
165 end type orbital_type

m_plowannier/plowannier_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  plowannier_type

FUNCTION

SOURCE

250 type, public :: plowannier_type
251 
252   integer :: nkpt
253   ! number of k points in Brillouin zone
254 
255   integer :: bandi_wan
256   ! energy band minimum considered
257 
258   integer :: bandf_wan
259   ! energy band maximum considered
260 
261   integer :: natom_wan
262   ! number of atoms (used to compute Wannier functions)
263 
264   integer :: size_wan
265   ! sum of all the m possible for every atom considered
266 
267   integer, allocatable :: iatom_wan(:)
268   ! array of each atom (used to compute Wannier functions)
269 
270   integer, allocatable :: nbl_atom_wan(:)
271   ! array of the number of l considered for each atom
272 
273   type(latom_wan_type), allocatable :: latom_wan(:)
274   ! for each atom, it contains an array of the l we are interested in
275 
276   integer, allocatable :: nbproj_atom_wan(:)
277   ! array of the number of projectors considered for each atom
278 
279   type(projector_wan_type), allocatable :: projector_wan(:)
280   ! for each atom, it contains an array of the projectors we are interested in
281 
282   type(position_wan_type), allocatable :: nposition(:)
283   ! array of the number of position considered for each atom
284 
285   integer :: nsppol
286   ! number of polarization
287 
288   integer :: nspinor
289   ! number of spinorial components
290 
291   type(orbital_type), allocatable :: psichi(:,:,:)
292   ! arrays of psichi
293 
294   integer, allocatable :: position(:,:)
295   ! size natom,3, gives the position of the cell for this atom (rprim coordinates)
296 
297   real(dp),allocatable :: kpt(:,:)
298   ! gives the coordinates in the BZ of the kpoint
299   ! size (3,nkpt)
300 
301   real(dp),allocatable :: wtk(:)
302   !weight of each kpoint
303 
304   real(dp),allocatable :: acell(:)
305   !size of the cell
306 
307 end type plowannier_type

m_plowannier/position_wan_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  position_wan_type

FUNCTION

SOURCE

121 type, public :: position_wan_type
122 
123   integer, allocatable :: pos(:,:)
124   ! size (number of position,3)
125 
126 end type position_wan_type

m_plowannier/print_operwan [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  print_operwan

FUNCTION

  Print the Wannier operator (real space) in a latex file

INPUTS

  wan, operwan, name

OUTPUT

SOURCE

3228 subroutine print_operwan(wan,operwan,name,convert)
3229 
3230 !Arguments----------------------------------
3231   type(operwan_type),intent(in) :: operwan(:,:,:)
3232   type(plowannier_type), intent(in) :: wan
3233   character(len=*), intent(in) :: name
3234   real(dp), intent(in) :: convert
3235 
3236 !Local variables----------------------------
3237   integer :: iatom1,iatom2,pos1,pos2,il1,il2,im1,im2,isppol,ikpt,unt
3238   real(dp) :: sum
3239   character(len = 500) :: str1,str2,msg
3240 
3241   if (open_file(name, msg, newunit=unt) /= 0) then
3242     ABI_ERROR(msg)
3243   end if
3244 
3245   write(unt,'(a)') '\documentclass[11pt,a4paper,landscape]{article}'
3246  write(unt,'(a)') '\usepackage[T1]{fontenc}'
3247   write(unt,'(a)') '\usepackage{geometry,tabularx,graphicx}'
3248   write(unt,'(a)') '\geometry{left=0.5cm,right=0.5cm}'
3249 
3250   write(unt,'(a)') '\begin{document}'
3251   write(unt,'(a)') '\noindent'
3252 
3253 !  write(unt,'(a,i0,a,F7.3,a)') "% ",wan%natom_wan," atom in a ",wan%acell(1)," cell"
3254 !  write(unt,'(a)') "% atom isppol proj"
3255 
3256 
3257 do isppol = 1,wan%nsppol
3258   write(unt,'(a)') '\begin{figure}'
3259   write(unt,'(a)') '\resizebox{\linewidth}{!}{%'
3260   write(unt,'(a)') '$ \left('
3261 
3262 
3263   write(str1,'(a)') '\begin{array}{'
3264   do iatom1 = 1,wan%natom_wan
3265     do pos1 = 1,size(wan%nposition(iatom1)%pos,1)
3266       do il1 = 1,wan%nbl_atom_wan(iatom1)
3267         do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
3268           str1 = trim(str1)//"c"
3269         end do
3270         if (iatom1 .ne. wan%natom_wan .or. il1 .ne. wan%nbl_atom_wan(iatom1) .or. pos1 .ne. size(wan%nposition(iatom1)%pos,1)) then
3271           str1 = trim(str1)//'|'
3272         end if
3273       end do
3274     end do
3275   end do
3276   str1 = trim(str1)//'}'
3277   write(unt,'(a)') str1
3278 
3279 
3280   do iatom1 = 1,wan%natom_wan
3281     do pos1 = 1,size(wan%nposition(iatom1)%pos,1)
3282       do il1 = 1,wan%nbl_atom_wan(iatom1)
3283         do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
3284           write(str1,'(a)') ""
3285           do iatom2 = 1,wan%natom_wan
3286             do pos2 = 1,size(wan%nposition(iatom2)%pos,1)
3287               do il2 = 1,wan%nbl_atom_wan(iatom2)
3288                 do im2 = 1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
3289                   sum = 0
3290                   do ikpt = 1,wan%nkpt
3291                     sum = sum + convert*real(operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,1,1)*wan%wtk(ikpt)*&
3292         &      exp(cmplx(0.0,1.0)*two_pi*(wan%kpt(1,ikpt)*(  &
3293         &           wan%nposition(iatom1)%pos(pos1,1)-wan%nposition(iatom2)%pos(pos2,1))+ &
3294         &           wan%kpt(2,ikpt)*(wan%nposition(iatom1)%pos(pos1,2)-wan%nposition(iatom2)%pos(pos2,2))+ &
3295         &           wan%kpt(3,ikpt)*(wan%nposition(iatom1)%pos(pos1,3)-wan%nposition(iatom2)%pos(pos2,3)))))
3296                   end do
3297                   write(str2,'(F10.6)') real(sum)
3298                   if ( len_trim(str1) .ge. 2) then
3299                     str1 = trim(str1)//"&"//trim(str2)
3300                   else
3301                     str1 = trim(str2)
3302                   end if
3303                   if (iatom2 .eq. wan%natom_wan .and. il2 .eq. wan%nbl_atom_wan(iatom2) .and. im2 &
3304    &                       .eq. 2*wan%latom_wan(iatom2)%lcalc(il2)+1 .and. pos2 .eq. size(wan%nposition(iatom2)%pos,1)) then
3305                     str1 = trim(str1)//'\\'
3306                   end if
3307                 end do
3308               end do
3309             end do
3310           end do
3311           write(unt,'(a)') trim(str1)
3312         end do
3313         if (iatom1 .ne. wan%natom_wan .or. il1 .ne. wan%nbl_atom_wan(iatom1) .or. pos1 .ne. size(wan%nposition(iatom1)%pos,1)) then
3314           write(unt,'(a)') '\hline'
3315         end if
3316       end do
3317     end do
3318   end do
3319 
3320 
3321 
3322   write(unt,'(a)') '\end{array} \right) $ }'
3323 
3324   if (name(len_trim(name)-2:len(trim(name))) == 'gen') then
3325     write(unt,'(a,i0,a,F7.3,a)')     "\caption{Energy matrix in real space for isppol = ", &
3326 &    isppol," in a ",wan%acell(1), " a.u. cell}"
3327   end if
3328 
3329   if (name(len_trim(name)-2:len(trim(name))) == 'occ') then
3330     write(unt,'(a,i0,a,F7.3,a)')     "\caption{Occupation matrix in real space for isppol = ",&
3331 &    isppol," in a ",wan%acell(1), " a.u. cell}"
3332   end if
3333 
3334 
3335 
3336   write(unt,'(a)') '\end{figure}'
3337   write(unt,'(a)') '\end{document}'
3338 
3339 end do
3340   close(unt)
3341 
3342 end subroutine print_operwan

m_plowannier/print_plowannier [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  print_plowannier

FUNCTION

  print the wannier weight (psichi) on a forlb.ovlp file

INPUTS

 dtset%typat,wan

OUTPUT

SOURCE

2336  subroutine print_plowannier(wan)
2337 
2338  use m_abicore
2339  use m_io_tools,  only : open_file
2340  use m_specialmsg, only : wrtout
2341  implicit none
2342 
2343  !Arguments-------------------------
2344  type(plowannier_type),intent(in) :: wan
2345  !Local variables-------------------
2346  character(len=500) :: msg
2347  integer :: unt,iatom,spin,ikpt,iband,ibandc,il,ispinor,im
2348 
2349  !Creation of the data.plowann file
2350  if (open_file('data.plowann',msg,newunit=unt,form='formatted',status='replace') /= 0) then
2351   ABI_ERROR(msg)
2352  end if
2353  rewind(unt)
2354 
2355  write(msg,'(2a)') ch10,' Print the psichi coefficients in data.plowann'
2356  call wrtout(std_out,msg,'COLL') ; call wrtout(ab_out,msg,'COLL')
2357 
2358  !Header of the file data.plowann
2359  write(unt,'(a22,i2)')"Total number of atom =", wan%natom_wan
2360  write(unt,*)"List of atoms", wan%iatom_wan(:)
2361  write(unt,'(a7,2i4)')"Bands =",wan%bandi_wan,wan%bandf_wan
2362  write(unt,'(a26,i2)')"Total number of orbitals =",sum(wan%nbl_atom_wan(:))
2363  do iatom=1,wan%natom_wan
2364    write(unt,'(a17,i2,a3,4i2)')"Orbitals for atom",wan%iatom_wan(iatom)," = ",wan%latom_wan(iatom)%lcalc(:)
2365  enddo
2366  write(unt,'(a16,i2)')"Number of spin =",wan%nsppol
2367  write(unt,'(a19,i4)')"Number of k-points=",wan%nkpt
2368  do ikpt=1,wan%nkpt
2369    write(unt,'(a,2x,i4)')"ikpt =",ikpt
2370     do spin=1,wan%nsppol
2371       do ispinor=1,wan%nspinor
2372         do iband=wan%bandi_wan,wan%bandf_wan
2373           ibandc=iband-wan%bandi_wan+1
2374           write(unt,'(2x,a,2x,i2,2x,i2)')"iband =",iband
2375           do iatom=1,wan%natom_wan
2376             do il=1,wan%nbl_atom_wan(iatom)
2377               do im=1,2*wan%latom_wan(iatom)%lcalc(il)+1
2378                 write(unt,'(8x,3i3,2x,2f23.15)')iatom,wan%latom_wan(iatom)%lcalc(il),im,&
2379                 &real(wan%psichi(ikpt,ibandc,iatom)%atom(il)%matl(im,spin,ispinor))&
2380                 &,aimag(wan%psichi(ikpt,ibandc,iatom)%atom(il)%matl(im,spin,ispinor))
2381               enddo!m
2382             enddo!l
2383           enddo!atom
2384         enddo!band
2385       enddo!spinor
2386     enddo!spin
2387   enddo!k-point
2388  close(unt)
2389  end subroutine print_plowannier

m_plowannier/projector_wan_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  projector_wan_type

FUNCTION

SOURCE

103 type, public :: projector_wan_type
104 
105   integer, allocatable :: lproj(:)
106   ! gives the list of the projector chosen
107 
108 end type projector_wan_type

m_plowannier/reduce_operwan_realspace [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  reduce_operwan_realspace

FUNCTION

  reduce a table of operwan_realspace type

INPUTS

  wan,rhot1,npwx,nibz,comm,nbz,nsppol

OUTPUT

 rhot1

SOURCE

3414 subroutine reduce_operwan_realspace(wan,rhot1,npwx,nibz,comm,nbz,nsppol)
3415 
3416 
3417   use m_xmpi, only : xmpi_barrier,xmpi_sum
3418 !Arguments---------------------------------------
3419   type(plowannier_type),intent(in) :: wan
3420   integer, intent(in) :: npwx,nibz,comm,nbz,nsppol
3421   type(operwan_realspace_type),target,intent(inout) :: rhot1(npwx,nibz)
3422 !Local variables----------------------------------
3423   complex(dpc),allocatable ::  buffer(:)
3424   integer :: dim,pwx,ibz, spin, ispinor1, ispinor2, iatom1, iatom2, pos1, pos2
3425   integer :: il1, il2, im1, im2, nnn, ierr
3426   complex(dpc),pointer :: oper_ptr(:,:,:,:,:)
3427 
3428 
3429    dim=0
3430      do pwx=1,npwx
3431      do ibz=1,nibz
3432        do spin=1,wan%nsppol
3433        do ispinor1=1,wan%nspinor
3434        do ispinor2=1,wan%nspinor
3435          do iatom1=1,wan%natom_wan
3436          do iatom2=1,wan%natom_wan
3437            do pos1=1,size(wan%nposition(iatom1)%pos,1)
3438            do pos2=1,size(wan%nposition(iatom2)%pos,1)
3439              do il1=1,wan%nbl_atom_wan(iatom1)
3440              do il2=1,wan%nbl_atom_wan(iatom2)
3441                do im1=1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
3442                do im2=1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
3443      dim=dim+1
3444                enddo!im2
3445                enddo!im1
3446              enddo!il2
3447              enddo!il1
3448            enddo!pos2
3449            enddo!pos1
3450          enddo!iatom2
3451          enddo!iatom1
3452        enddo!ispinor2
3453        enddo!ispinor1
3454        enddo!spin
3455      enddo!ibz
3456      enddo!pwx
3457      ABI_MALLOC(buffer,(dim))
3458      nnn=0
3459      do pwx=1,npwx
3460      do ibz=1,nibz
3461        do iatom1=1,wan%natom_wan
3462        do iatom2=1,wan%natom_wan
3463          do pos1=1,size(wan%nposition(iatom1)%pos,1)
3464          do pos2=1,size(wan%nposition(iatom2)%pos,1)
3465            do il1=1,wan%nbl_atom_wan(iatom1)
3466            do il2=1,wan%nbl_atom_wan(iatom2)
3467              oper_ptr=>rhot1(pwx,ibz)%atom_index(iatom1,iatom2)%position(pos1,pos2)%atom(il1,il2)%matl
3468              do im1=1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
3469              do im2=1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
3470                do spin=1,wan%nsppol
3471                  do ispinor1=1,wan%nspinor
3472                  do ispinor2=1,wan%nspinor
3473      nnn=nnn+1
3474      buffer(nnn)=oper_ptr(im1,im2,spin,ispinor1,ispinor2)
3475                  enddo!ispinor2
3476                  enddo!ispinor1
3477                enddo!spin
3478              enddo!im2
3479              enddo!im1
3480            enddo!il2
3481            enddo!il1
3482          enddo!pos2
3483          enddo!pos1
3484        enddo!iatom2
3485        enddo!iatom1
3486      enddo!ibz
3487      enddo!pwx
3488      call xmpi_barrier(comm)
3489      call xmpi_sum(buffer,comm,ierr)
3490      call xmpi_barrier(comm)
3491      buffer=buffer/nbz/nsppol
3492      nnn=0
3493      do pwx=1,npwx
3494      do ibz=1,nibz
3495          do iatom1=1,wan%natom_wan
3496          do iatom2=1,wan%natom_wan
3497            do pos1=1,size(wan%nposition(iatom1)%pos,1)
3498            do pos2=1,size(wan%nposition(iatom2)%pos,1)
3499              do il1=1,wan%nbl_atom_wan(iatom1)
3500              do il2=1,wan%nbl_atom_wan(iatom2)
3501                oper_ptr=>rhot1(pwx,ibz)%atom_index(iatom1,iatom2)%position(pos1,pos2)%atom(il1,il2)%matl
3502                do im1=1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
3503                do im2=1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
3504                  do spin=1,wan%nsppol
3505                    do ispinor1=1,wan%nspinor
3506                    do ispinor2=1,wan%nspinor
3507       nnn=nnn+1
3508       oper_ptr(im1,im2,spin,ispinor1,ispinor2)=buffer(nnn)
3509                enddo!im2
3510                enddo!im1
3511              enddo!il2
3512              enddo!il1
3513            enddo!pos2
3514            enddo!pos1
3515          enddo!iatom2
3516          enddo!iatom1
3517        enddo!ispinor2
3518        enddo!ispinor1
3519        enddo!spin
3520      enddo!ibz
3521      enddo!pwx
3522      ABI_FREE(buffer)
3523 
3524 
3525 end subroutine reduce_operwan_realspace

m_plowannier/zero_operwan [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  zero_operwan

FUNCTION

  zero operwan

INPUTS

  wan

OUTPUT

  operwan

SOURCE

2781  subroutine zero_operwan(wan,operwan)
2782 
2783    implicit none
2784 
2785    !Arguments----------------------------------
2786    type(plowannier_type), intent(in) :: wan
2787    type(operwan_type), intent(inout) :: operwan(wan%nkpt,wan%natom_wan,wan%natom_wan)
2788 
2789    !Local variables----------------------------
2790    integer :: ikpt,  iatom1, iatom2, il1, il2, isppol, ispinor1, ispinor2,  im1, im2
2791 
2792 
2793    do ikpt = 1,wan%nkpt
2794      do isppol = 1,wan%nsppol
2795        do iatom1 = 1,wan%natom_wan
2796          do iatom2 = 1,wan%natom_wan
2797            do il1 = 1,wan%nbl_atom_wan(iatom1)
2798              do il2 = 1,wan%nbl_atom_wan(iatom2)
2799                do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
2800                  do im2 = 1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
2801                    do ispinor1 = 1,wan%nspinor
2802                      do ispinor2 = 1,wan%nspinor
2803                        operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)=czero
2804                      end do
2805                    end do
2806                  end do
2807                end do
2808              end do
2809            end do
2810          end do
2811        end do
2812      end do
2813    end do
2814 
2815  end subroutine zero_operwan

m_plowannier/zero_operwan_realspace [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  zero_operwan_realspace

FUNCTION

  Set an operwan_realspace to zero

INPUTS

  wan, operwan_realspace

OUTPUT

 operwan_realspace

SOURCE

3587 subroutine zero_operwan_realspace(wan,operwan_realspace)
3588 
3589 !Arguments----------------------------------
3590   type(operwan_realspace_type),intent(inout) :: operwan_realspace
3591   type(plowannier_type), intent(in) :: wan
3592 
3593 !Local variables----------------------------
3594   integer :: isppol,iatom1,iatom2,pos1,pos2,il1,il2,im1,im2,ispinor1,ispinor2
3595 
3596 
3597   do isppol = 1,wan%nsppol
3598     do ispinor1=1,wan%nspinor
3599       do ispinor2=1,wan%nspinor
3600         do iatom1 = 1,wan%natom_wan
3601           do pos1 = 1,size(wan%nposition(iatom1)%pos,1)
3602             do il1 = 1,wan%nbl_atom_wan(iatom1)
3603               do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
3604                 do iatom2 = 1,wan%natom_wan
3605                   do pos2 = 1,size(wan%nposition(iatom2)%pos,1)
3606                     do il2 = 1,wan%nbl_atom_wan(iatom2)
3607                       do im2 = 1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
3608                         operwan_realspace%atom_index(iatom1,iatom2)%position(pos1,pos2)%&
3609                           &atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)=czero
3610                       enddo
3611                     enddo
3612                   enddo
3613                 enddo
3614               enddo
3615             enddo
3616           enddo
3617         enddo
3618       enddo
3619     enddo
3620   enddo
3621 end subroutine zero_operwan_realspace