TABLE OF CONTENTS


ABINIT/m_plowannier [ Modules ]

[ Top ] [ Modules ]

NAME

  m_plowannier

FUNCTION

COPYRIGHT

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

PARENTS

NOTES

SOURCE

20 #if defined HAVE_CONFIG_H
21 #include "config.h"
22 #endif
23 
24 #include "abi_common.h"
25 
26 MODULE m_plowannier
27 
28  use defs_basis
29  use m_errors
30  use m_abicore
31 
32  use m_io_tools,  only : open_file
33  use m_mpinfo,    only : proc_distrb_cycle
34 
35  implicit none
36 
37  private
38 
39  public :: init_plowannier
40  public :: copy_orbital
41  public :: compute_coeff_plowannier
42  public :: destroy_plowannier
43  public :: initialize_operwan
44  public :: destroy_operwan
45  public :: compute_oper_ks2wan
46  public :: normalization_plowannier
47  public :: print_operwan

m_plowannier/allocate_orbital [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  allocate_orbital

FUNCTION

  allocate an array of orbital_type

INPUTS

  lorbital1

OUTPUT

  lorbital2

PARENTS

      m_plowannier

CHILDREN

SOURCE

551 subroutine allocate_orbital(orbital1,orbital2,n1,n2,n3)
552 
553 
554 !This section has been created automatically by the script Abilint (TD).
555 !Do not modify the following lines by hand.
556 #undef ABI_FUNC
557 #define ABI_FUNC 'allocate_orbital'
558 !End of the abilint section
559 
560  implicit none
561 
562  !Arguments----------------
563  integer,intent(in) :: n1,n2,n3
564  type(orbital_type), intent(in) :: orbital1(n1,n2,n3)
565  type(orbital_type),intent(inout) :: orbital2(n1,n2,n3)
566 
567  !Local variable-----------
568  integer :: n4,n5,n6,n7
569  integer :: i,j,k,l
570 
571  do i = 1,n1
572    do j = 1,n2
573      do k = 1,n3
574        n4 = size(orbital1(i,j,k)%atom,1)
575        ABI_DATATYPE_ALLOCATE(orbital2(i,j,k)%atom,(n4))
576        do l = 1,n4
577          n5 = size(orbital1(i,j,k)%atom(l)%matl,1)
578          n6 = size(orbital1(i,j,k)%atom(l)%matl,2)
579          n7 = size(orbital1(i,j,k)%atom(l)%matl,3)
580          ABI_ALLOCATE(orbital2(i,j,k)%atom(l)%matl,(n5,n6,n7))
581        end do
582      end do
583    end do
584  end do
585 
586 end subroutine allocate_orbital

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)

PARENTS

      m_plowannier

CHILDREN

SOURCE

2408  subroutine compute_oper_ks2wan(wan,operks,operwan,option)
2409 
2410 #ifdef FC_INTEL
2411 !DEC$ NOOPTIMIZE
2412 #endif
2413 
2414 !This section has been created automatically by the script Abilint (TD).
2415 !Do not modify the following lines by hand.
2416 #undef ABI_FUNC
2417 #define ABI_FUNC 'compute_oper_ks2wan'
2418 !End of the abilint section
2419 
2420    implicit none
2421 
2422    !Arguments--------------------------
2423    type(plowannier_type), intent(in) :: wan
2424    type(operwan_type), intent(inout) :: operwan(:,:,:)
2425    complex(dpc), intent(in) :: operks(:,:,:,:)
2426    integer, intent(in) :: option
2427 
2428    !Local variables--------------------
2429    integer :: iatom1, iatom2, il1, il2, isppol, ispinor1, ispinor2, iband1, iband2, im1, im2
2430 
2431    ! ----------------------------------
2432    !Transformation KS2WAN
2433    ! ----------------------------------
2434 
2435 
2436    !!operation on reciprocal space
2437    do iatom1 = 1,wan%natom_wan
2438      do iatom2 = 1,wan%natom_wan
2439        do isppol = 1,wan%nsppol
2440          do il1 = 1,wan%nbl_atom_wan(iatom1)
2441            do il2 = 1,wan%nbl_atom_wan(iatom2)
2442              do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
2443                do im2 = 1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
2444                  do ispinor1 = 1,wan%nspinor
2445                    do ispinor2 = 1,wan%nspinor
2446                      !!sum over the bands
2447                      do iband1 = 1,wan%bandf_wan-wan%bandi_wan+1
2448                        do iband2 = 1,wan%bandf_wan-wan%bandi_wan+1
2449                          operwan(option,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)=&
2450                            operwan(option,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)&
2451                           + conjg(wan%psichi(option,iband2,iatom2)%atom(il2)%matl(im2,isppol,ispinor2))&
2452                    *operks(option,iband1,iband2,isppol)*wan%psichi(option,iband1,iatom1)%atom(il1)%matl(im1,isppol,ispinor1)
2453                        end do
2454                      end do
2455                    end do
2456                  end do
2457                end do
2458              end do
2459            end do
2460          end do
2461        end do
2462      end do
2463    end do
2464 
2465 end subroutine compute_oper_ks2wan

m_plowannier/copy_orbital [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  copy_orbital

FUNCTION

  Copy an array of orbital_type

INPUTS

  lorbital1

OUTPUT

  lorbital2

PARENTS

      m_plowannier

CHILDREN

SOURCE

486 subroutine copy_orbital(orbital1,orbital2,n1,n2,n3)
487 
488 
489 !This section has been created automatically by the script Abilint (TD).
490 !Do not modify the following lines by hand.
491 #undef ABI_FUNC
492 #define ABI_FUNC 'copy_orbital'
493 !End of the abilint section
494 
495  implicit none
496 
497  !Arguments----------------
498  integer,intent(in) :: n1,n2,n3
499  type(orbital_type), intent(in) :: orbital1(n1,n2,n3)
500  type(orbital_type),intent(inout) :: orbital2(n1,n2,n3)
501 
502  !Local variable-----------
503  integer :: n4,n5,n6,n7
504  integer :: i,j,k,l,m,p,q
505 
506  do i = 1,n1
507    do j = 1,n2
508      do k = 1,n3
509        n4 = size(orbital1(i,j,k)%atom,1)
510        do l = 1,n4
511          n5 = size(orbital1(i,j,k)%atom(l)%matl,1)
512          n6 = size(orbital1(i,j,k)%atom(l)%matl,2)
513          n7 = size(orbital1(i,j,k)%atom(l)%matl,3)
514          do m = 1,n5
515            do p = 1,n6
516              do q = 1,n7
517                orbital2(i,j,k)%atom(l)%matl(m,p,q) = orbital1(i,j,k)%atom(l)%matl(m,p,q)
518              end do
519            end do
520          end do
521        end do
522      end do
523    end do
524  end do
525 
526 end subroutine copy_orbital

m_plowannier/destroy_operwan [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  destroy_operwan

FUNCTION

  destroy operwan

INPUTS

  wan

OUTPUT

  operwan

PARENTS

      m_plowannier

CHILDREN

SOURCE

2353  subroutine destroy_operwan(wan,operwan)
2354 
2355 
2356 !This section has been created automatically by the script Abilint (TD).
2357 !Do not modify the following lines by hand.
2358 #undef ABI_FUNC
2359 #define ABI_FUNC 'destroy_operwan'
2360 !End of the abilint section
2361 
2362    implicit none
2363 
2364    !Arguments----------------------------------
2365    type(plowannier_type), intent(in) :: wan
2366    type(operwan_type), intent(inout) :: operwan(wan%nkpt,wan%natom_wan,wan%natom_wan)
2367 
2368    !Local variables----------------------------
2369    integer :: ikpt,iatom1,iatom2,il1,il2
2370 
2371 
2372    do ikpt = 1,wan%nkpt
2373      do iatom1 = 1,wan%natom_wan
2374        do iatom2 = 1,wan%natom_wan
2375          do il1 = 1,wan%nbl_atom_wan(iatom1)
2376            do il2 = 1,wan%nbl_atom_wan(iatom2)
2377              ABI_DEALLOCATE(operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl)
2378            end do
2379          end do
2380          ABI_DATATYPE_DEALLOCATE(operwan(ikpt,iatom1,iatom2)%atom)
2381        end do
2382      end do
2383    end do
2384  end subroutine destroy_operwan

m_plowannier/destroy_orbital [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  destroy_orbital

FUNCTION

  destroy an array of orbital_type

INPUTS

  lorbital1

OUTPUT

  lorbital2

PARENTS

      m_plowannier

CHILDREN

SOURCE

612 subroutine destroy_orbital(orbital2,n1,n2,n3)
613 
614 
615 !This section has been created automatically by the script Abilint (TD).
616 !Do not modify the following lines by hand.
617 #undef ABI_FUNC
618 #define ABI_FUNC 'destroy_orbital'
619 !End of the abilint section
620 
621  implicit none
622 
623  !Arguments----------------
624  integer,intent(in) :: n1,n2,n3
625  type(orbital_type),intent(inout) :: orbital2(n1,n2,n3)
626 
627  !Local variable-----------
628  integer :: n4
629  integer :: i,j,k,l
630 
631  do i = 1,n1
632    do j = 1,n2
633      do k = 1,n3
634        n4 = size(orbital2(i,j,k)%atom,1)
635        do l = 1,n4
636          ABI_DEALLOCATE(orbital2(i,j,k)%atom(l)%matl)
637        end do
638        ABI_DATATYPE_DEALLOCATE(orbital2(i,j,k)%atom)
639      end do
640    end do
641  end do
642 
643 end subroutine destroy_orbital

m_plowannier/destroy_plowannier [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  destroy_plowannier

FUNCTION

  deallocate variables

INPUTS

  wan

OUTPUT

PARENTS

      outscfcv

CHILDREN

SOURCE

2196  subroutine destroy_plowannier(wan)
2197 
2198 
2199 !This section has been created automatically by the script Abilint (TD).
2200 !Do not modify the following lines by hand.
2201 #undef ABI_FUNC
2202 #define ABI_FUNC 'destroy_plowannier'
2203 !End of the abilint section
2204 
2205  implicit none
2206 
2207 !Arguments-------------------------------------
2208  type(plowannier_type), intent(inout) :: wan
2209 !Local variables-------------------------------
2210  integer :: iatom,ikpt,iband,il
2211 
2212  do iatom=1,wan%natom_wan
2213    ABI_DEALLOCATE(wan%latom_wan(iatom)%lcalc)
2214    ABI_DEALLOCATE(wan%projector_wan(iatom)%lproj)
2215    ABI_DEALLOCATE(wan%nposition(iatom)%pos)
2216  enddo
2217  do iatom = 1,wan%natom_wan
2218    do il = 1,wan%nbl_atom_wan(iatom)
2219      ABI_DEALLOCATE(wan%psichi(1,1,iatom)%atom(il)%ph0phiint)
2220    end do
2221  end do
2222  do ikpt = 1,wan%nkpt
2223    do iband = wan%bandi_wan,wan%bandf_wan
2224      do iatom = 1,wan%natom_wan
2225        do il = 1,wan%nbl_atom_wan(iatom)
2226         ABI_DEALLOCATE(wan%psichi(ikpt,iband-wan%bandi_wan+1,iatom)%atom(il)%matl)
2227        end do
2228        ABI_DATATYPE_DEALLOCATE(wan%psichi(ikpt,iband-wan%bandi_wan+1,iatom)%atom)
2229      end do
2230    end do
2231  end do
2232 
2233  if (allocated(wan%kpt)) then
2234    ABI_DEALLOCATE(wan%kpt)
2235  end if
2236  if (allocated(wan%iatom_wan)) then
2237    ABI_DEALLOCATE(wan%iatom_wan)
2238  end if
2239  if (allocated(wan%nbl_atom_wan)) then
2240    ABI_DEALLOCATE(wan%nbl_atom_wan)
2241  end if
2242  if (allocated(wan%latom_wan)) then
2243    ABI_DATATYPE_DEALLOCATE(wan%latom_wan)
2244  end if
2245  if (allocated(wan%nbproj_atom_wan)) then
2246    ABI_DEALLOCATE(wan%nbproj_atom_wan)
2247  end if
2248  if (allocated(wan%projector_wan)) then
2249    ABI_DATATYPE_DEALLOCATE(wan%projector_wan)
2250  end if
2251  if (allocated(wan%position)) then
2252    ABI_DEALLOCATE(wan%position)
2253  end if
2254  if (allocated(wan%wtk)) then
2255    ABI_DEALLOCATE(wan%wtk)
2256  end if
2257  if (allocated(wan%acell)) then
2258    ABI_DEALLOCATE(wan%acell)
2259  end if
2260 
2261 
2262  if (allocated(wan%nposition)) then
2263    ABI_DATATYPE_DEALLOCATE(wan%nposition)
2264  end if
2265  if (allocated(wan%psichi)) then
2266    ABI_DATATYPE_DEALLOCATE(wan%psichi)
2267  end if
2268 
2269 
2270  end subroutine destroy_plowannier

m_plowannier/initialize_operwan [ Functions ]

[ Top ] [ m_plowannier ] [ Functions ]

NAME

  initialize_operwan

FUNCTION

  initialize operwan

INPUTS

  wan

OUTPUT

  operwan

PARENTS

      m_plowannier

CHILDREN

SOURCE

2295  subroutine initialize_operwan(wan,operwan)
2296 
2297 
2298 !This section has been created automatically by the script Abilint (TD).
2299 !Do not modify the following lines by hand.
2300 #undef ABI_FUNC
2301 #define ABI_FUNC 'initialize_operwan'
2302 !End of the abilint section
2303 
2304    implicit none
2305 
2306    !Arguments----------------------------------
2307    type(plowannier_type), intent(in) :: wan
2308    type(operwan_type), intent(inout) :: operwan(wan%nkpt,wan%natom_wan,wan%natom_wan)
2309 
2310    !Local variables----------------------------
2311    integer :: ikpt,iatom1,iatom2,il1,il2,n1,n2
2312 
2313    do ikpt = 1,wan%nkpt
2314      do iatom1 = 1,wan%natom_wan
2315        do iatom2 = 1,wan%natom_wan
2316          ABI_DATATYPE_ALLOCATE(operwan(ikpt,iatom1,iatom2)%atom,(wan%nbl_atom_wan(iatom1),wan%nbl_atom_wan(iatom2)))
2317          do il1 = 1,wan%nbl_atom_wan(iatom1)
2318            do il2 = 1,wan%nbl_atom_wan(iatom2)
2319              n1=2*wan%latom_wan(iatom1)%lcalc(il1)+1
2320              n2=2*wan%latom_wan(iatom2)%lcalc(il2)+1
2321    ABI_ALLOCATE(operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl,(n1,n2,wan%nsppol,wan%nspinor,wan%nspinor))
2322              operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl = zero
2323            end do
2324          end do
2325        end do
2326      end do
2327    end do
2328 
2329  end subroutine initialize_operwan

m_plowannier/latom_wan_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  latom_wan_type

FUNCTION

SOURCE

60 type, public :: latom_wan_type
61 
62   integer, allocatable :: lcalc(:)
63   ! array of the l we want to compute the psichi with
64 
65 end type latom_wan_type

m_plowannier/lorbital2_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  lorbital2_type

FUNCTION

SOURCE

153 type, public :: lorbital2_type
154 
155   complex(dpc), allocatable :: matl(:,:,:,:,:)
156   ! size (2l1+1,2l2+1,nspppol,nspinor,nspinor)
157 
158   real(dp), allocatable :: ph0phiint(:)
159   ! size (nproj), stocks the value of ph0phiint we may want
160 
161 
162 end type lorbital2_type

m_plowannier/lorbital_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  lorbital_type

FUNCTION

SOURCE

113 type, public :: lorbital_type
114 
115   complex(dpc), allocatable :: matl(:,:,:)
116   !details for different m
117 
118   real(dp), allocatable :: ph0phiint(:)
119   ! stocks the values for each projector of the l considered
120   ! size(total number of projectors for this l)
121 
122 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, option (=0 if we normalize over all k)

OUTPUT

  wan itself is modified

PARENTS

      m_plowannier

CHILDREN

SOURCE

2491 subroutine normalization_plowannier(wan)
2492 
2493 
2494   use m_matrix, only : invsqrt_matrix
2495 
2496 !This section has been created automatically by the script Abilint (TD).
2497 !Do not modify the following lines by hand.
2498 #undef ABI_FUNC
2499 #define ABI_FUNC 'normalization_plowannier'
2500 !End of the abilint section
2501 
2502   implicit none
2503 !Arguments------------------
2504 
2505 !This section has been created automatically by the script Abilint (TD).
2506 !Do not modify the following lines by hand.
2507 #undef ABI_FUNC
2508 #define ABI_FUNC 'normalization_plowannier'
2509 !End of the abilint section
2510 
2511   type(plowannier_type), intent(inout) :: wan
2512 
2513 !Local----------------------
2514   complex(dpc), allocatable :: operks(:,:,:,:)
2515   type(operwan_type), allocatable :: operwan(:,:,:)
2516   complex(dpc), allocatable :: operwansquare(:,:,:,:)
2517   complex(dpc), allocatable :: tmp_operwansquare(:,:)
2518   integer :: ikpt, iband, iband1, iband2, isppol,  ispinor1, ispinor2, iatom1, &
2519 &  iatom2, il1, il2, im1, im2, index_c, index_l, n1,n2,n3
2520   type(orbital_type), allocatable :: psichinormalized(:,:,:)
2521   !character(len = 50) :: mat_writing2
2522   !character(len = 5000) :: mat_writing
2523   character(len = 500) :: message
2524 
2525 
2526   !First, creation of the ks identity operator
2527   ABI_ALLOCATE(operks,(wan%nkpt,wan%bandf_wan-wan%bandi_wan+1,wan%bandf_wan-wan%bandi_wan+1,wan%nsppol))
2528   operks = czero
2529   do iband1 = 1,wan%bandf_wan-wan%bandi_wan+1
2530     do iband2 = 1,wan%bandf_wan-wan%bandi_wan+1
2531       if (iband1.eq.iband2) then
2532         do ikpt = 1,wan%nkpt
2533           do isppol= 1,wan%nsppol
2534             operks(ikpt,iband1,iband2,isppol) = cone
2535           end do
2536         end do
2537       end if
2538     end do
2539   end do
2540 
2541 
2542   !Allocation of operwan
2543   ABI_DATATYPE_ALLOCATE(operwan,(wan%nkpt,wan%natom_wan,wan%natom_wan))
2544   call initialize_operwan(wan,operwan)
2545 
2546 
2547 
2548   !Computation of the overlap
2549   do ikpt = 1,wan%nkpt
2550     call compute_oper_ks2wan(wan,operks,operwan,ikpt)
2551   end do
2552 
2553 
2554 
2555   !transform the operwan into an inversible matrix
2556   !!operwansquare is the overlap square matrix (wan%size_wan * wan%size_wan)
2557   ABI_ALLOCATE(operwansquare,(wan%nkpt,wan%nsppol,wan%nspinor*wan%size_wan,wan%nspinor*wan%size_wan))
2558 
2559   operwansquare = czero
2560 
2561   n1=size(wan%psichi,1)
2562   n2=size(wan%psichi,2)
2563   n3=size(wan%psichi,3)
2564   ABI_DATATYPE_ALLOCATE(psichinormalized,(n1,n2,n3))
2565   call allocate_orbital(wan%psichi,psichinormalized,n1,n2,n3)
2566   call copy_orbital(wan%psichi,psichinormalized,n1,n2,n3)
2567 
2568   do isppol = 1,wan%nsppol
2569     do ispinor1 = 1,wan%nspinor
2570       do ispinor2 = 1,wan%nspinor
2571         index_l = 0 ! index_l is set to 0 at the beginning
2572         do iatom1 = 1,wan%natom_wan
2573           do il1 = 1,wan%nbl_atom_wan(iatom1)
2574             do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
2575               index_l = index_l + 1 ! the line changes
2576               index_c = 1 ! counter_c is set to one each time the line changes
2577               do iatom2 = 1,wan%natom_wan
2578                 do il2 = 1,wan%nbl_atom_wan(iatom2)
2579                   do im2 = 1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
2580                     do ikpt = 1,wan%nkpt
2581          operwansquare(ikpt,isppol,index_l+wan%size_wan*(ispinor1-1),index_c+wan%size_wan*(ispinor2-1))  &
2582 &                 = operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)
2583                     end do
2584                     index_c = index_c + 1
2585                   end do !im2
2586                 end do !il2
2587               end do ! iatom2 (the line changes)
2588             end do ! im1
2589           end do ! il1
2590         end do !iatom1
2591       end do
2592     end do
2593   end do
2594 
2595 
2596 
2597 
2598   !!Write the overlap matrix for ikpt = 1 in a nice shape
2599 !    do isppol = 1,wan%nsppol
2600 !      do il1 = 1,size(operwansquare,3) !! dummy variable without any meaning
2601 !        write(mat_writing,'(a,i0,i0)') 'Overlap matrix before orthonormalization 1 ',isppol,il1
2602 !        do il2 = 1,size(operwansquare,4)
2603 !         write(mat_writing2,'(F10.6)') real(operwansquare(1,isppol,il1,il2))
2604 !          mat_writing = trim(mat_writing)//trim(mat_writing2)
2605 !        end do
2606 !        print*,trim(mat_writing)
2607 !      end do
2608 !    end do
2609 
2610 
2611 
2612   !take the square root inverse of operwansquare for normalization purposes
2613   ABI_ALLOCATE(tmp_operwansquare,(wan%nspinor*wan%size_wan,wan%nspinor*wan%size_wan))
2614   do isppol = 1,wan%nsppol
2615     do ikpt = 1,wan%nkpt
2616       tmp_operwansquare(:,:)=operwansquare(ikpt,isppol,:,:)
2617       call invsqrt_matrix(tmp_operwansquare,wan%nspinor*wan%size_wan)
2618       operwansquare(ikpt,isppol,:,:)=tmp_operwansquare(:,:)
2619     end do
2620   end do
2621   ABI_DEALLOCATE(tmp_operwansquare)
2622 
2623   do ikpt = 1,wan%nkpt
2624     do iband = 1,wan%bandf_wan-wan%bandi_wan+1
2625       do iatom1 = 1,wan%natom_wan
2626         do il1 = 1,wan%nbl_atom_wan(iatom1)
2627           psichinormalized(ikpt,iband,iatom1)%atom(il1)%matl = czero
2628         end do
2629       end do
2630     end do
2631   end do
2632 
2633 
2634 
2635   ! compute the new psichi normalized
2636   do isppol = 1,wan%nsppol
2637     do ispinor1 = 1,wan%nspinor
2638       do iband = 1,wan%bandf_wan-wan%bandi_wan+1
2639         index_l = 0
2640         do iatom1 = 1,wan%natom_wan
2641           do il1 = 1,wan%nbl_atom_wan(iatom1)
2642             do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
2643               ! sum
2644               do ispinor2 = 1,wan%nspinor
2645                 index_l = index_l + 1 ! the line changes
2646                 index_c = 1 ! when the line changes, index_c is set to 1
2647                 do iatom2 = 1,wan%natom_wan
2648                   do il2 = 1,wan%nbl_atom_wan(iatom2)
2649                     do im2 = 1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
2650                       do ikpt = 1,wan%nkpt
2651                         psichinormalized(ikpt,iband,iatom1)%atom(il1)%matl(im1,isppol,ispinor1) =&
2652                           psichinormalized(ikpt,iband,iatom1)%atom(il1)%matl(im1,isppol,ispinor1) +&
2653                           wan%psichi(ikpt,iband,iatom2)%atom(il2)%matl(im2,isppol,ispinor2)*&
2654                           operwansquare(ikpt,isppol,index_l+wan%size_wan*(ispinor1-1),index_c+&
2655                           wan%size_wan*(ispinor2-1))
2656                       end do
2657                       index_c = index_c + 1
2658                     end do !im2
2659                   end do !il2
2660                 end do !iatom2
2661               end do ! ispinor2
2662             end do !im1
2663           end do ! il1
2664         end do ! iatom1
2665       end do ! iband
2666     end do ! ispinor1
2667   end do !isppol
2668   ! copy the new psichi normalized
2669   call copy_orbital(psichinormalized,wan%psichi,n1,n2,n3)
2670   call destroy_orbital(psichinormalized,n1,n2,n3)
2671   ABI_DATATYPE_DEALLOCATE(psichinormalized)
2672 
2673   call destroy_operwan(wan,operwan)
2674   ABI_DATATYPE_DEALLOCATE(operwan)
2675 
2676 
2677 
2678 
2679 
2680 
2681 
2682 
2683 !!
2684 !!  !-------------------------------------------------------------
2685 !!  !check if the new norm is one
2686 
2687   ABI_DATATYPE_ALLOCATE(operwan,(wan%nkpt,wan%natom_wan,wan%natom_wan))
2688   call initialize_operwan(wan,operwan)
2689   do ikpt = 1,wan%nkpt
2690     call compute_oper_ks2wan(wan,operks,operwan,ikpt)
2691   end do
2692 
2693 
2694   do isppol = 1,wan%nsppol
2695     do ikpt = 1,wan%nkpt
2696       do iatom1 = 1,wan%natom_wan
2697         do iatom2 = 1,wan%natom_wan
2698           do il1 = 1,wan%nbl_atom_wan(iatom1)
2699             do il2 = 1,wan%nbl_atom_wan(iatom2)
2700               do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
2701                 do im2 = 1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
2702                   do ispinor1 = 1,wan%nspinor
2703                     do ispinor2 = 1,wan%nspinor
2704                       if (iatom1.eq.iatom2 .and. il1.eq.il2 .and. im1.eq.im2 .and. ispinor1.eq.ispinor2) then
2705                         if (abs(cmplx(1.0,0.0,dpc)-&
2706              &             operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)) > 1d-8) then
2707                           write(message,'(a,i0,a,F18.11)') 'Normalization error for ikpt =',ikpt,&
2708          &                         ' on diag, value = ',&
2709           &                        abs(operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2))
2710                           MSG_ERROR(message)
2711                         end if
2712                       else
2713                         if (abs(operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)) > 1d-8) then
2714                           write(message,'(a,i0,a,F10.3)') 'Normalization error for ikpt =',ikpt,&
2715          &                         ' not on diag, value = ',&
2716          &                           abs(operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2))
2717                          MSG_ERROR(message)
2718                         end if
2719                       end if
2720                     end do
2721                   end do
2722                 end do
2723               end do
2724             end do
2725           end do
2726         end do
2727       end do
2728     end do
2729   end do
2730 
2731 
2732   !!Uncomment to print the overlap matrix in the log file (for ikpt = 1)
2733 !  do isppol = 1,wan%nsppol
2734 !    do ispinor1 = 1,wan%nspinor
2735 !      do ispinor2 = 1,wan%nspinor
2736 !        index_l = 0 ! index_l is set to 0 at the beginning
2737 !        do iatom1 = 1,wan%natom_wan
2738 !          do il1 = 1,wan%nbl_atom_wan(iatom1)
2739 !            do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
2740 !              index_l = index_l + 1 ! the line changes
2741 !              index_c = 1 ! counter_c is set to one each time the line changes
2742 !              do iatom2 = 1,wan%natom_wan
2743 !                do il2 = 1,wan%nbl_atom_wan(iatom2)
2744 !                  do im2 = 1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
2745 !                    do ikpt = 1,wan%nkpt
2746 !                      operwansquare(ikpt,isppol,index_l+wan%size_wan*(ispinor1-1),&
2747 !            &index_c+wan%size_wan*(ispinor2-1)) = operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,ispinor1,ispinor2)
2748 !                    end do
2749 !                    index_c = index_c + 1
2750 !                  end do !im2
2751 !                end do !il2
2752 !              end do ! iatom2 (the line changes)
2753 !            end do ! im1
2754 !          end do ! il1
2755 !        end do !iatom1
2756 !      end do
2757 !    end do
2758 !  end do
2759 
2760   !!Print the overlap matrix in a nice shape
2761 !  do isppol = 1,wan%nsppol
2762 !    do il1 = 1,size(operwansquare,3) !! dummy variable without any meaning
2763 !      write(mat_writing,'(a,i0,i0)') 'Overlap matrix after orthonormalization ',isppol,il1
2764 !      do il2 = 1,size(operwansquare,4)
2765 !        write(mat_writing2,'(F10.6)') real(operwansquare(9,isppol,il1,il2))
2766 !        mat_writing = trim(mat_writing)//trim(mat_writing2)
2767 !      end do
2768 !      print*,trim(mat_writing)
2769 !    end do
2770 !  end do
2771 
2772 
2773 
2774 
2775 !!  !----------------------------------------------------------------
2776   ABI_DEALLOCATE(operwansquare)
2777   ABI_DEALLOCATE(operks)
2778   call destroy_operwan(wan,operwan)
2779   ABI_DATATYPE_DEALLOCATE(operwan)
2780 
2781 
2782 end subroutine normalization_plowannier

m_plowannier/operwan_realspace_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  operwan_realspace_type

FUNCTION

SOURCE

191 type, public :: operwan_realspace_type
192 
193   type(operwan_type), allocatable :: position(:,:)
194   ! size (number of positions chosen for atom1, number of positions chosen for atom2)
195 
196 end type operwan_realspace_type

m_plowannier/operwan_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  operwan_type

FUNCTION

SOURCE

174 type, public :: operwan_type
175 
176   type(lorbital2_type), allocatable :: atom(:,:)
177   ! l chosen for each on of both atoms
178 
179 end type operwan_type

m_plowannier/orbital_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  orbital_type

FUNCTION

SOURCE

134 type, public :: orbital_type
135 
136   type(lorbital_type), allocatable :: atom(:)
137   ! details of the psichi coefficients for each atom
138   ! size of number of l chosen
139 
140 end type orbital_type

m_plowannier/plowannier_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  plowannier_type

FUNCTION

SOURCE

208 type, public :: plowannier_type
209 
210   integer :: nkpt
211   ! number of k points in Brillouin zone
212 
213   integer :: bandi_wan
214   ! energy band minimum considered
215 
216   integer :: bandf_wan
217   ! energy band maximum considered
218 
219   integer :: natom_wan
220   ! number of atoms in a cell
221 
222   integer :: size_wan
223   ! sum of all the m possible for every atom considered
224 
225   integer, allocatable :: iatom_wan(:)
226   ! array of each atom's type
227 
228   integer, allocatable :: nbl_atom_wan(:)
229   ! array of the number of l considered for each atom
230 
231   type(latom_wan_type), allocatable :: latom_wan(:)
232   ! for each atom, it contains an array of the l we are interested in
233 
234   integer, allocatable :: nbproj_atom_wan(:)
235   ! array of the number of projectors considered for each atom
236 
237   type(projector_wan_type), allocatable :: projector_wan(:)
238   ! for each atom, it contains an array of the projectors we are interested in
239 
240   type(position_wan_type), allocatable :: nposition(:)
241   ! array of the number of position considered for each atom
242 
243   integer :: nsppol
244   ! number of polarization
245 
246   integer :: nspinor
247   ! number of spinorial components
248 
249   type(orbital_type), allocatable :: psichi(:,:,:)
250   ! arrays of psichi
251 
252   integer, allocatable :: position(:,:)
253   ! size natom,3, gives the position of the cell for this atom (rprim coordinates)
254 
255   real,allocatable :: kpt(:,:)
256   ! gives the coordinates in the BZ of the kpoint
257   ! size (3,nkpt)
258 
259   real,allocatable :: wtk(:)
260   !weight of each kpoint
261 
262   real,allocatable :: acell(:)
263   !size of the cell
264 
265 end type plowannier_type

m_plowannier/position_wan_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  position_wan_type

FUNCTION

SOURCE

 96 type, public :: position_wan_type
 97 
 98   integer, allocatable :: pos(:,:)
 99   ! size (number of position,3)
100 
101 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

PARENTS

      m_plowannier

CHILDREN

SOURCE

2809 subroutine print_operwan(wan,operwan,name,convert)
2810 
2811 !Arguments----------------------------------
2812 
2813 !This section has been created automatically by the script Abilint (TD).
2814 !Do not modify the following lines by hand.
2815 #undef ABI_FUNC
2816 #define ABI_FUNC 'print_operwan'
2817 !End of the abilint section
2818 
2819   type(operwan_type),intent(in) :: operwan(:,:,:)
2820   type(plowannier_type), intent(in) :: wan
2821   character(len=*), intent(in) :: name
2822   real(dp), intent(in) :: convert
2823 
2824 !Local variables----------------------------
2825   integer :: iatom1,iatom2,pos1,pos2,il1,il2,im1,im2,isppol,ikpt,unt
2826   real(dp) :: sum
2827   character(len = 500) :: str1,str2,msg
2828 
2829   if (open_file(name, msg, newunit=unt) /= 0) then
2830     MSG_ERROR(msg)
2831   end if
2832 
2833   write(unt,'(a)') '\documentclass[11pt,a4paper,landscape]{article}'
2834  write(unt,'(a)') '\usepackage[T1]{fontenc}'
2835   write(unt,'(a)') '\usepackage{geometry,tabularx,graphicx}'
2836   write(unt,'(a)') '\geometry{left=0.5cm,right=0.5cm}'
2837 
2838   write(unt,'(a)') '\begin{document}'
2839   write(unt,'(a)') '\noindent'
2840 
2841 !  write(unt,'(a,i0,a,F7.3,a)') "% ",wan%natom_wan," atom in a ",wan%acell(1)," cell"
2842 !  write(unt,'(a)') "% atom isppol proj"
2843 
2844 
2845 do isppol = 1,wan%nsppol
2846   write(unt,'(a)') '\begin{figure}'
2847   write(unt,'(a)') '\resizebox{\linewidth}{!}{%'
2848   write(unt,'(a)') '$ \left('
2849 
2850 
2851   write(str1,'(a)') '\begin{array}{'
2852   do iatom1 = 1,wan%natom_wan
2853     do pos1 = 1,size(wan%nposition(iatom1)%pos,1)
2854       do il1 = 1,wan%nbl_atom_wan(iatom1)
2855         do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
2856           str1 = trim(str1)//"c"
2857         end do
2858         if (iatom1 .ne. wan%natom_wan .or. il1 .ne. wan%nbl_atom_wan(iatom1) .or. pos1 .ne. size(wan%nposition(iatom1)%pos,1)) then
2859           str1 = trim(str1)//'|'
2860         end if
2861       end do
2862     end do
2863   end do
2864   str1 = trim(str1)//'}'
2865   write(unt,'(a)') str1
2866 
2867 
2868   do iatom1 = 1,wan%natom_wan
2869     do pos1 = 1,size(wan%nposition(iatom1)%pos,1)
2870       do il1 = 1,wan%nbl_atom_wan(iatom1)
2871         do im1 = 1,2*wan%latom_wan(iatom1)%lcalc(il1)+1
2872           write(str1,'(a)') ""
2873           do iatom2 = 1,wan%natom_wan
2874             do pos2 = 1,size(wan%nposition(iatom2)%pos,1)
2875               do il2 = 1,wan%nbl_atom_wan(iatom2)
2876                 do im2 = 1,2*wan%latom_wan(iatom2)%lcalc(il2)+1
2877                   sum = 0
2878                   do ikpt = 1,wan%nkpt
2879                     sum = sum + convert*real(operwan(ikpt,iatom1,iatom2)%atom(il1,il2)%matl(im1,im2,isppol,1,1)*wan%wtk(ikpt)*&
2880         &      exp(cmplx(0.0,1.0)*two_pi*(wan%kpt(1,ikpt)*(  &
2881         &           wan%nposition(iatom1)%pos(pos1,1)-wan%nposition(iatom2)%pos(pos2,1))+ &
2882         &           wan%kpt(2,ikpt)*(wan%nposition(iatom1)%pos(pos1,2)-wan%nposition(iatom2)%pos(pos2,2))+ &
2883         &           wan%kpt(3,ikpt)*(wan%nposition(iatom1)%pos(pos1,3)-wan%nposition(iatom2)%pos(pos2,3)))))
2884                   end do
2885                   write(str2,'(F10.6)') real(sum)
2886                   if ( len_trim(str1) .ge. 2) then
2887                     str1 = trim(str1)//"&"//trim(str2)
2888                   else
2889                     str1 = trim(str2)
2890                   end if
2891                   if (iatom2 .eq. wan%natom_wan .and. il2 .eq. wan%nbl_atom_wan(iatom2) .and. im2 &
2892    &                       .eq. 2*wan%latom_wan(iatom2)%lcalc(il2)+1 .and. pos2 .eq. size(wan%nposition(iatom2)%pos,1)) then
2893                     str1 = trim(str1)//'\\'
2894                   end if
2895                 end do
2896               end do
2897             end do
2898           end do
2899           write(unt,'(a)') trim(str1)
2900         end do
2901         if (iatom1 .ne. wan%natom_wan .or. il1 .ne. wan%nbl_atom_wan(iatom1) .or. pos1 .ne. size(wan%nposition(iatom1)%pos,1)) then
2902           write(unt,'(a)') '\hline'
2903         end if
2904       end do
2905     end do
2906   end do
2907 
2908 
2909 
2910   write(unt,'(a)') '\end{array} \right) $ }'
2911 
2912   if (name(len_trim(name)-2:len(trim(name))) == 'gen') then
2913     write(unt,'(a,i0,a,F7.3,a)')     "\caption{Energy matrix in real space for isppol = ", &
2914 &    isppol," in a ",wan%acell(1), " a.u. cell}"
2915   end if
2916 
2917   if (name(len_trim(name)-2:len(trim(name))) == 'occ') then
2918     write(unt,'(a,i0,a,F7.3,a)')     "\caption{Occupation matrix in real space for isppol = ",&
2919 &    isppol," in a ",wan%acell(1), " a.u. cell}"
2920   end if
2921 
2922 
2923 
2924   write(unt,'(a)') '\end{figure}'
2925   write(unt,'(a)') '\end{document}'
2926 
2927 end do
2928   close(unt)
2929 
2930 end subroutine print_operwan
2931 
2932 END MODULE m_plowannier

m_plowannier/projector_wan_type [ Types ]

[ Top ] [ m_plowannier ] [ Types ]

NAME

  projector_wan_type

FUNCTION

SOURCE

78 type, public :: projector_wan_type
79 
80   integer, allocatable :: lproj(:)
81   ! gives the list of the projector chosen
82 
83 end type projector_wan_type