TABLE OF CONTENTS


ABINIT/getdim_nloc [ Functions ]

[ Top ] [ Functions ]

NAME

 getdim_nloc

FUNCTION

 Determine the dimensions of arrays that contain
 the definition of non-local projectors : ekb, ffspl, indlmn

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (XG)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

  mixalch(npspalch,ntypalch,nimage)=alchemical mixing coefficients
  nimage=number of images
  npsp=number of pseudopotentials
  npspalch=number of pseudopotentials for alchemical purposes
  ntypat=number of types of pseudo atoms
  ntypalch=number of types of alchemical pseudo atoms
  pspheads(npsp)=<type pspheader_type>all the important information from the
   pseudopotential file headers, as well as the psp file names

OUTPUT

  lmnmax=maximum number of l,m,n projectors, not taking into account the spin-orbit
  lmnmaxso=maximum number of l,m,n projectors, taking into account the spin-orbit
  lnmax=maximum number of l,n projectors, not taking into account the spin-orbit
  lnmaxso=maximum number of l,n projectors, taking into account the spin-orbit

PARENTS

      m_psps,memory_eval

CHILDREN

      wrtout

SOURCE

 41 #if defined HAVE_CONFIG_H
 42 #include "config.h"
 43 #endif
 44 
 45 #include "abi_common.h"
 46 
 47 subroutine getdim_nloc(lmnmax,lmnmaxso,lnmax,lnmaxso,mixalch,nimage,npsp,npspalch,&
 48 & ntypat,ntypalch,pspheads)
 49 
 50  use defs_basis
 51  use defs_datatypes
 52  use m_errors
 53  use m_profiling_abi
 54 
 55 !This section has been created automatically by the script Abilint (TD).
 56 !Do not modify the following lines by hand.
 57 #undef ABI_FUNC
 58 #define ABI_FUNC 'getdim_nloc'
 59  use interfaces_14_hidewrite
 60 !End of the abilint section
 61 
 62  implicit none
 63 
 64 !Arguments ------------------------------------
 65 !scalars
 66  integer,intent(in) :: nimage,npsp,npspalch,ntypalch,ntypat
 67  integer,intent(out) :: lmnmax,lmnmaxso,lnmax,lnmaxso
 68 !arrays
 69  real(dp),intent(in) :: mixalch(npspalch,ntypalch,nimage)
 70  type(pspheader_type),intent(in) :: pspheads(npsp)
 71 
 72 !Local variables-------------------------------
 73 !scalars
 74  integer :: ilang,ipsp,ipspalch,itypalch,itypat,ntyppure
 75 !integer :: llmax
 76  character(len=500) :: message
 77 !arrays
 78  integer,allocatable :: lmnproj_typat(:),lmnprojso_typat(:),lnproj_typat(:)
 79  integer,allocatable :: lnprojso_typat(:),nproj_typat(:,:),nprojso_typat(:,:)
 80 
 81 ! *************************************************************************
 82 
 83 !write(std_out,*)' getdim_nloc: 'pspheads(1)%nproj(0:3)=',pspheads(1)%nproj(0:3)
 84 
 85  ABI_ALLOCATE(lmnproj_typat,(ntypat))
 86  ABI_ALLOCATE(lmnprojso_typat,(ntypat))
 87  ABI_ALLOCATE(lnproj_typat,(ntypat))
 88  ABI_ALLOCATE(lnprojso_typat,(ntypat))
 89  ABI_ALLOCATE(nproj_typat,(0:3,ntypat))
 90  ABI_ALLOCATE(nprojso_typat,(3,ntypat))
 91  lmnproj_typat(:)=0 ; lmnprojso_typat(:)=0
 92  lnproj_typat(:)=0 ; lnprojso_typat(:)=0
 93  nproj_typat(:,:)=0 ; nprojso_typat(:,:)=0
 94 
 95  ntyppure=ntypat-ntypalch
 96 
 97 !For each type of pseudo atom, compute the number of projectors
 98 !First, pure pseudo atoms
 99  if(ntyppure>0)then
100    do itypat=1,ntyppure
101      nproj_typat(0:3,itypat)=pspheads(itypat)%nproj(0:3)
102      nprojso_typat(:,itypat)=pspheads(itypat)%nprojso(:)
103    end do
104  end if
105 
106 !Then, alchemical pseudo atoms
107  if(ntypalch>0)then
108    do itypat=ntyppure+1,ntypat
109      itypalch=itypat-ntyppure
110      do ipsp=ntyppure+1,npsp
111        ipspalch=ipsp-ntyppure
112 !      If there is some mixing, must accumulate the projectors
113        if(sum(abs(mixalch(ipspalch,itypalch,:)))>tol10)then
114          nproj_typat(0:3,itypat)=nproj_typat(0:3,itypat)+pspheads(ipsp)%nproj(0:3)
115          nprojso_typat(:,itypat)=nprojso_typat(:,itypat)+pspheads(ipsp)%nprojso(:)
116        end if
117      end do
118    end do
119  end if
120 
121 !Now that the number of projectors is known, accumulate the dimensions
122  do itypat=1,ntypat
123    do ilang=0,3
124      lnproj_typat(itypat)=lnproj_typat(itypat)+nproj_typat(ilang,itypat)
125      lmnproj_typat(itypat)=lmnproj_typat(itypat)+nproj_typat(ilang,itypat)*(2*ilang+1)
126    end do
127    lnprojso_typat(itypat)=lnproj_typat(itypat)
128    lmnprojso_typat(itypat)=lmnproj_typat(itypat)
129    do ilang=1,3
130      lnprojso_typat(itypat)=lnprojso_typat(itypat)+nprojso_typat(ilang,itypat)
131      lmnprojso_typat(itypat)=lmnprojso_typat(itypat)+nprojso_typat(ilang,itypat)*(2*ilang+1)
132    end do
133  end do
134 
135 !Compute the maximal bounds, at least equal to 1, even for local psps
136  lmnmax=1;lmnmaxso=1;lnmax=1;lnmaxso=1
137  do itypat=1,ntypat
138    lmnmax  =max(lmnmax  ,lmnproj_typat  (itypat))
139    lmnmaxso=max(lmnmaxso,lmnprojso_typat(itypat))
140    lnmax   =max(lnmax   ,lnproj_typat   (itypat))
141    lnmaxso =max(lnmaxso ,lnprojso_typat (itypat))
142  end do
143 !The initial coding (below) was not totally portable (MT 110215)
144 !lmnmax=max(maxval(lmnproj_typat(1:ntypat)),1)
145 !lmnmaxso=max(maxval(lmnprojso_typat(1:ntypat)),1)
146 !lnmax=max(maxval(lnproj_typat(1:ntypat)),1)
147 !lnmaxso=max(maxval(lnprojso_typat(1:ntypat)),1)
148 
149  if(maxval(lmnproj_typat(1:ntypat))==0)then
150    write(message, '(3a)' )&
151 &   'Despite there is only a local part to pseudopotential(s),',ch10,&
152 &   'lmnmax and lnmax are set to 1.'
153    MSG_COMMENT(message)
154  end if
155 
156 !XG040806 : These lines make modifications of lnmax and lmnmax
157 !that are unjustified in many cases, according to the many tests cases
158 !where they produce a changes, while the test case was working properly.
159 !One should understand better the needs, and code more appropriate changes ...
160 !lnmax/lmnmax has to be bigger than 1+lmax (for compatibility reasons)
161 !llmax=maxval(pspheads(1:ntypat)%lmax)+1 ! And this line might have trouble with HP compiler
162 !if (lnmax   <llmax) lnmax=llmax
163 !if (lnmaxso <llmax) lnmaxso=llmax
164 !if (lmnmax  <llmax) lmnmax=llmax
165 !if (lmnmaxso<llmax) lmnmaxso=llmax
166 
167  write(message, '(a,a,i4,a,i4,3a,i4,a,i4,a)' ) ch10,&
168 & ' getdim_nloc : deduce lmnmax  =',lmnmax,', lnmax  =',lnmax,',',ch10,&
169 & '                      lmnmaxso=',lmnmaxso,', lnmaxso=',lnmaxso,'.'
170  call wrtout(std_out,message,'COLL')
171 
172  ABI_DEALLOCATE(lmnproj_typat)
173  ABI_DEALLOCATE(lmnprojso_typat)
174  ABI_DEALLOCATE(lnproj_typat)
175  ABI_DEALLOCATE(lnprojso_typat)
176  ABI_DEALLOCATE(nproj_typat)
177  ABI_DEALLOCATE(nprojso_typat)
178 
179 end subroutine getdim_nloc