TABLE OF CONTENTS


ABINIT/m_lmn_indices [ Modules ]

[ Top ] [ Modules ]

NAME

  m_lmn_indices

FUNCTION

  This module provides tools to calculate tables commonly used to iterate
  over the the (l,m,n) channels of the PAW partial waves.

COPYRIGHT

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

SOURCE

18 #if defined HAVE_CONFIG_H
19 #include "config.h"
20 #endif
21 
22 #include "abi_common.h"
23 
24 MODULE m_lmn_indices
25 
26  use defs_basis
27  use m_profiling_abi
28  use m_errors
29 
30  implicit none
31 
32  private
33 
34 ! Public procedures.
35  public :: ilm2lm          ! Returns the value of l and m in from the index ilm
36  public :: make_indlmn     ! Calculates the indlmn(6,lmn_size) table giving l,m,n,lm,ln,spin for i=lmn.
37  public :: make_indklmn    ! Calculates the indklmn(8,lmn2_size) table giving
38                            !   klm, kln, abs(il-jl), (il+jl), ilm and jlm, ilmn and jlmn for each symmetric klmn=(ilmn,jlmn)
39  public :: make_kln2ln     ! Calculates the kln2ln(6,ln2_size) table giving
40                            !   il, jl ,in, jn, iln, jln for each symmetric kln=(iln,jln)
41  public :: make_klm2lm     ! Calculates the klm2lm(6,lm2_size) table giving
42                            !   il, jl ,im, jm, ilm, jlm for each symmetric klm=(ilm,jlm)
43  public :: klmn2ijlmn      ! Calculates ilmn and jlmn from klmn.
44  public :: make_indln      ! Calculates indln(2,ln_size) giving l and n for i=ln
45  public :: uppert_index    ! The sequential index of an element in the upper triangle of a matrix 
46 
47 CONTAINS  !========================================================================================

m_lmn_indices/ilm2lm [ Functions ]

[ Top ] [ m_lmn_indices ] [ Functions ]

NAME

  ilm2lm

FUNCTION

  Returns the value of l and m in from the index ilm

INPUTS

  ilm=The contracted index for (l,m).

OUTPUT

  ll=The angular momentun defined in [0,1,2,3,4].
  mm=The magnetic number in the interval [-l,....+l].

PARENTS

CHILDREN

SOURCE

 70 subroutine ilm2lm(ilm,ll,mm)
 71 
 72 
 73 !This section has been created automatically by the script Abilint (TD).
 74 !Do not modify the following lines by hand.
 75 #undef ABI_FUNC
 76 #define ABI_FUNC 'ilm2lm'
 77 !End of the abilint section
 78 
 79  implicit none
 80 
 81 !Arguments ------------------------------------
 82 !scalars
 83  integer,intent(in) :: ilm
 84  integer,intent(out) :: ll,mm
 85 
 86 !Local variables-------------------------------
 87  integer :: ii
 88 ! *********************************************************************
 89 
 90  if (ilm<1) then 
 91    MSG_ERROR("Wrong ilm")
 92  end if
 93 
 94  ll = -1
 95  do ii=0,100
 96   if ( (ii+1)**2 >= ilm) then
 97    ll = ii
 98    EXIT
 99   end if
100  end do
101 
102  mm = ilm - ll**2 -ll-1
103 
104  if (ll==-1) then
105    MSG_ERROR("l>100 not programmed!")
106  end if
107 
108 end subroutine ilm2lm

m_lmn_indices/klmn2ijlmn [ Functions ]

[ Top ] [ m_lmn_indices ] [ Functions ]

NAME

  klmn2ijlmn

FUNCTION

  Find ilmn and jlmn from klmn and lmn_size.

INPUTS

  lmn_size=Number of (l,m,n) elements for the PAW basis set.
  klmn=The index corresponding to (ilmn,jlmn) in packed form.

OUTPUT

  jlmn, ilmn=The two symmetrix indices corresponding to klmn. NB: jlmn >= ilmn

PARENTS

      m_lmn_indices,m_paw_slater

CHILDREN

SOURCE

495 subroutine klmn2ijlmn(klmn,lmn_size,ilmn,jlmn)
496 
497 
498 !This section has been created automatically by the script Abilint (TD).
499 !Do not modify the following lines by hand.
500 #undef ABI_FUNC
501 #define ABI_FUNC 'klmn2ijlmn'
502 !End of the abilint section
503 
504  implicit none
505 
506 !Arguments ------------------------------------
507 !scalars
508  integer,intent(in) :: klmn,lmn_size
509  integer,intent(out) :: ilmn,jlmn
510 
511 !Local variables ------------------------------
512 !scalars
513  integer :: ii,jj,k0
514 !************************************************************************
515 
516  ilmn=-1; jlmn=-1
517 
518  do jj=1,lmn_size
519    k0=jj*(jj-1)/2
520    do ii=1,jj
521      if (klmn==ii+k0) then
522        ilmn = ii
523        jlmn = jj; RETURN
524      end if
525    end do
526  end do
527 
528  MSG_BUG("Not able to found ilmn and jlmn")
529 
530 end subroutine klmn2ijlmn

m_lmn_indices/make_indklmn [ Functions ]

[ Top ] [ m_lmn_indices ] [ Functions ]

NAME

  make_indklmn

FUNCTION

  Performs the setup of the indklmn table for PAW calculations.
  Compute the indklmn indexes giving klm, kln, abs(il-jl) and (il+jl), ilm and jlm, ilmn and jlmn
  for each klmn=(ilmn,jlmn) with jlmn >= ilmn

INPUTS

  lcutdens=Maximum l for densities/potentials moments computations
  lmn_size=Number of (l,m,n) elements for the PAW basis set.
  lmn2_size=Number of elements in the symmetric basis set: lmn2_size=lmn_size*(lmn_size+1)/2
  indlmn(6,lmn_size)=Array giving l,m,n,lm,ln,spin for i=lmn.

OUTPUT

  indklmn(6,lmn2_size)=Array giving klm, kln, abs(il-jl), (il+jl), ilm and jlm, ilmn and jlmn
    for each klmn=(ilmn,jlmn). Note: ilmn=(il,im,in) and ilmn<=jlmn
  klm_diag(lmn2_size)=1 il==jl and im==jm, 0 otherwise.

PARENTS

      m_atom

CHILDREN

SOURCE

213 subroutine make_indklmn(lcutdens,lmn_size,lmn2_size,indlmn,indklmn,klm_diag)
214 
215 
216 !This section has been created automatically by the script Abilint (TD).
217 !Do not modify the following lines by hand.
218 #undef ABI_FUNC
219 #define ABI_FUNC 'make_indklmn'
220 !End of the abilint section
221 
222  implicit none
223 
224 !Arguments ------------------------------------
225 !scalars
226  integer,intent(in) :: lmn2_size,lmn_size
227  integer,intent(in) ::  lcutdens
228 !scalars
229  integer,intent(in) :: indlmn(6,lmn_size)
230  integer,intent(out) :: indklmn(8,lmn2_size)
231  integer,intent(out) :: klm_diag(lmn2_size)
232 
233 !Local variables ------------------------------
234 !scalars
235  integer :: i0lm,i0ln,il,ilm,ilmn,iln
236  integer :: j0lm,j0lmn,j0ln,jl,jlm,jlmn,jln,klmn
237 
238 !************************************************************************
239 
240  klm_diag=0
241 
242  do jlmn=1,lmn_size
243   jl= indlmn(1,jlmn); jlm=indlmn(4,jlmn); jln=indlmn(5,jlmn)
244   j0lmn=jlmn*(jlmn-1)/2
245   j0lm =jlm *(jlm -1)/2
246   j0ln =jln *(jln -1)/2
247   do ilmn=1,jlmn
248    il=indlmn(1,ilmn); ilm=indlmn(4,ilmn); iln=indlmn(5,ilmn)
249    klmn=j0lmn+ilmn
250    if (ilm<=jlm) then
251     indklmn(1,klmn)=j0lm+ilm     !klm
252    else
253     i0lm=ilm*(ilm-1)/2
254     indklmn(1,klmn)=i0lm+jlm
255    end if
256    if (iln<=jln) then
257     indklmn(2,klmn)=j0ln+iln     !kln
258    else
259     i0ln=iln*(iln-1)/2
260     indklmn(2,klmn)=i0ln+jln
261    end if
262    !MG This is not safe, what happens if lcutdens < |il-jl|?
263    indklmn(3,klmn)=MIN(ABS(il-jl),lcutdens)  ! abs(li-lj) NB this is a l-value, not an index >=1.
264    indklmn(4,klmn)=MIN(il+jl,lcutdens)       ! abs(li+lj) NB this is a l-value, not an index >=1.
265 
266    indklmn(5,klmn)=ilm                       ! ilm
267    indklmn(6,klmn)=jlm                       ! jlm
268 
269    indklmn(7,klmn)=ilmn                      ! ilmn
270    indklmn(8,klmn)=jlmn                      ! jlmn
271 
272 
273    if (ilm==jlm) klm_diag(klmn)=1
274   end do
275  end do
276 
277 end subroutine make_indklmn

m_lmn_indices/make_indlmn [ Functions ]

[ Top ] [ m_lmn_indices ] [ Functions ]

NAME

  make_indlmn

FUNCTION

  Performs the setup of the indlmn table for PAW calculations (indices for (l,m,n) basis)

INPUTS

  ln_size= Total number of nl components
  lmn_size= Second dimension in indlmn. Total number of (l,m,n) components for this atom.
  orbitals(ln_size)=Give the value of l for each element of the augmented basis set.

OUTPUT

  indlmn(6,lmn_size)=array giving l,m,n,lm,ln,s for i=lmn

PARENTS

      m_atom

CHILDREN

SOURCE

135 subroutine make_indlmn(ln_size,lmn_size,orbitals,indlmn)
136 
137 
138 !This section has been created automatically by the script Abilint (TD).
139 !Do not modify the following lines by hand.
140 #undef ABI_FUNC
141 #define ABI_FUNC 'make_indlmn'
142 !End of the abilint section
143 
144  implicit none
145 
146 !Arguments ------------------------------------
147 !scalars
148  integer,intent(in) :: ln_size,lmn_size
149  integer,intent(in) :: orbitals(ln_size)
150 !scalars
151  integer,intent(out) :: indlmn(6,lmn_size)
152 
153 !Local variables ------------------------------
154 !scalars
155  integer :: ilmn,ib,il,iln,ilm
156 !arrays
157  integer,allocatable :: nprj(:)
158 
159 !************************************************************************
160 
161  ABI_ALLOCATE(nprj,(0:MAXVAL(orbitals)))
162 
163  ilmn=0; iln=0; nprj=0
164  do ib=1,ln_size
165    il=orbitals(ib)
166    nprj(il)=nprj(il)+1
167    iln=iln+1
168    do ilm=1,2*il+1
169      indlmn(1,ilmn+ilm)=il           ! l
170      indlmn(2,ilmn+ilm)=ilm-(il+1)   ! m
171      indlmn(3,ilmn+ilm)=nprj(il)     ! n
172      indlmn(4,ilmn+ilm)=il*il+ilm    ! lm index
173      indlmn(5,ilmn+ilm)=iln          ! ln index
174      indlmn(6,ilmn+ilm)=1            ! spin (not yet used!)
175    end do
176    ilmn=ilmn+2*il+1
177  end do
178 
179  ABI_DEALLOCATE(nprj)
180 
181 end subroutine make_indlmn

m_lmn_indices/make_indln [ Functions ]

[ Top ] [ m_lmn_indices ] [ Functions ]

NAME

  make_indln

FUNCTION

  Performs the setup of the indln table for PAW calculations.
  Compute the indln indexes giving ilmn=(l,m,n)

INPUTS

  lmn_size=Number of (l,m,n) elements for the PAW basis set.
  ln_size=Number of (l,n) elements
  indlmn(6,lmn_size)=Array giving l,m,n,lm,ln,spin for i=lmn.

OUTPUT

  indln(2,ln_size)=Array giving l and n for i=ln

PARENTS

      m_paw_slater

CHILDREN

SOURCE

558 subroutine make_indln(lmn_size,ln_size,indlmn,indln)
559 
560 
561 !This section has been created automatically by the script Abilint (TD).
562 !Do not modify the following lines by hand.
563 #undef ABI_FUNC
564 #define ABI_FUNC 'make_indln'
565 !End of the abilint section
566 
567  implicit none
568 
569 !Arguments ------------------------------------
570 !scalars
571  integer,intent(in) :: lmn_size,ln_size
572 !arrays
573  integer,intent(in) :: indlmn(6,lmn_size)
574  integer,intent(out) :: indln(2,ln_size)
575 
576 !Local variables ------------------------------
577 !scalars
578  integer :: ilmn,ll,nn,ll_old,nn_old,ii
579 !************************************************************************
580 
581  ii=0; ll_old=-1; nn_old=-1
582  do ilmn=1,lmn_size
583   ll=indlmn(1,ilmn)
584   nn=indlmn(3,ilmn)
585   if (ll/=ll_old.or.nn/=nn_old) then
586    ll_old=ll
587    nn_old=nn
588    ii=ii+1
589    indln(1,ii)=ll
590    indln(2,ii)=nn
591   end if
592  end do
593 
594  ABI_CHECK(ii==ln_size,"ii/=ln_size")
595 
596 end subroutine make_indln

m_lmn_indices/make_klm2lm [ Functions ]

[ Top ] [ m_lmn_indices ] [ Functions ]

NAME

  make_klm2lm

FUNCTION

  Performs the setup of the klm2lm table for PAW calculations.

INPUTS

  lmn_size=Number of (l,m,n) elements for the PAW basis set.
  lmn2_size=Number of elements in the symmetric basis set: lmn2_size=lmn_size*(lmn_size+1)/2
  lm2_size)=Number of (l.m) elements in the symmetric basis set.
  indlmn(6,lmn_size)=Array giving l,m,n,lm,ln,spin for i=lmn.
  indklmn(8,lmn2_size)=Array giving klm, kln, abs(il-jl), (il+jl), ilm and jlm
   for each klmn=(ilmn,jlmn). Note: ilmn=(il,im,in) and ilmn<=jlmn

OUTPUT

  klm2lm(6,lm2_size)=Table giving il, jl ,im, jm, ilm, jlm for each klm=(ilm,jlm)
  where ilm=(il,im) and ilm<=jlm. NB: klm2lm is an application and not a bijection.

NOTES

  klm2lm can be calculated easily if we assume that all (l,m) channels
  are ordered by increasing l and m. This is the standard convention
  used in most of the PAW datasets. This routines, howevever, works
  works also in the unlikely case in with (l,m) are not ordered.

PARENTS

      m_paw_slater

CHILDREN

SOURCE

400 subroutine make_klm2lm(lmn_size,lmn2_size,lm2_size,indlmn,indklmn,klm2lm)
401 
402 
403 !This section has been created automatically by the script Abilint (TD).
404 !Do not modify the following lines by hand.
405 #undef ABI_FUNC
406 #define ABI_FUNC 'make_klm2lm'
407 !End of the abilint section
408 
409  implicit none
410 
411 !Arguments ------------------------------------
412 !scalars
413  integer,intent(in) :: lmn2_size,lmn_size,lm2_size
414 !arrays
415  integer,intent(in) :: indlmn(6,lmn_size)
416  integer,intent(in) :: indklmn(8,lmn2_size)
417  integer,intent(out) :: klm2lm(6,lm2_size)
418 
419 !Local variables ------------------------------
420 !scalars
421  integer :: il,ilm,ilmn,jl,jlm,jlmn,im,jm
422  integer :: klmn,klm,klm_old !,iklm
423 
424 !************************************************************************
425 
426  klm2lm = 0
427  klm_old = -1
428  do klmn=1,lmn2_size
429    klm = indklmn(1,klmn)
430    if (klm /= klm_old) then
431      klm_old = klm
432 
433      ilm = indklmn(5,klmn)
434      jlm = indklmn(6,klmn)
435 
436      call klmn2ijlmn(klmn,lmn_size,ilmn,jlmn)
437 
438      il=indlmn(1,ilmn); im=indlmn(2,ilmn)
439      jl=indlmn(1,jlmn); jm=indlmn(2,jlmn)
440 
441      !shift to have the index instead of l or m.
442      klm2lm(1,klm) = il +1       ! il
443      klm2lm(2,klm) = jl +1       ! jl
444      klm2lm(3,klm) = im + il +1  ! im
445      klm2lm(4,klm) = jm + jl +1  ! jm
446      klm2lm(5,klm) = ilm         ! ilm
447      klm2lm(6,klm) = jlm         ! jlm
448    end if
449  end do !klmn
450 
451  if (ANY(klm2lm==0)) then
452    MSG_BUG("check klm2lm")
453  end if
454 
455 !DEBUG
456 #if 0
457  write(std_out,*)"Debugging make_klm2lm:"
458  do iklm=1,lm2_size
459   il  = klm2lm(1,iklm) -1      ! li
460   jl  = klm2lm(2,iklm) -1      ! lj
461   im  = klm2lm(3,iklm) -il -1  ! mi
462   jm  = klm2lm(4,iklm) -jl -1  ! mj
463   ilm = klm2lm(5,iklm)         ! ilm
464   jlm = klm2lm(6,iklm)         ! jlm
465   write(std_out,'(i3,2(a,2i3,a))')"iklm ",iklm," l m (",il,im,")"," l m (",jl,jm,")"
466  end do
467 #endif
468 
469 end subroutine make_klm2lm

m_lmn_indices/make_kln2ln [ Functions ]

[ Top ] [ m_lmn_indices ] [ Functions ]

NAME

  make_kln2ln

FUNCTION

  Performs the setup of the kln2ln table for PAW calculations.

INPUTS

  lmn_size=Number of (l,m,n) elements for the PAW basis set.
  lmn2_size=Number of elements in the symmetric basis set: lmn2_size=lmn_size*(lmn_size+1)/2
  ln2_size=Number of symmetric (l,n) channels i.e. ln_size*(ln_size+1)/2
  indlmn(6,lmn_size)=Array giving l,m,n,lm,ln,spin for i=lmn.
  indklmn(8,lmn2_size)=Array giving klm, kln, abs(il-jl), (il+jl), ilm and jlm, ilmn and jlmn
   for each klmn=(ilmn,jlmn). Note: ilmn=(il,im,in) and ilmn<=jlmn

OUTPUT

  kln2ln(6,ln2_size)=Table giving il, jl ,in, jn, iln, jln for each kln=(iln,jln)
  where iln=(il,in) and iln<=jln. NB: kln2ln is an application and not a bijection

PARENTS

      m_atom,m_paw_slater

CHILDREN

SOURCE

308 subroutine make_kln2ln(lmn_size,lmn2_size,ln2_size,indlmn,indklmn,kln2ln)
309 
310 
311 !This section has been created automatically by the script Abilint (TD).
312 !Do not modify the following lines by hand.
313 #undef ABI_FUNC
314 #define ABI_FUNC 'make_kln2ln'
315 !End of the abilint section
316 
317  implicit none
318 
319 !Arguments ------------------------------------
320 !scalars
321  integer,intent(in) :: lmn2_size,lmn_size,ln2_size
322 !arrays
323  integer,intent(in) :: indlmn(6,lmn_size)
324  integer,intent(in) :: indklmn(8,lmn2_size)
325  integer,intent(out) :: kln2ln(6,ln2_size)
326 
327 !Local variables ------------------------------
328 !scalars
329  integer :: il,in,ilmn,iln
330  integer :: jl,jn,jlmn,jln
331  integer :: klmn,kln,kln_old
332 
333 !************************************************************************
334 
335  kln2ln = 0
336 
337  kln_old = -1
338  do klmn=1,lmn2_size
339    kln = indklmn(2,klmn)
340    if (kln /= kln_old) then
341      kln_old = kln
342 
343      call klmn2ijlmn(klmn,lmn_size,ilmn,jlmn)
344 
345      il= indlmn(1,ilmn); iln=indlmn(5,ilmn)
346      jl= indlmn(1,jlmn); jln=indlmn(5,jlmn)
347      !$in = indklmn(9 ,klmn) !=indlmn(3,ilmn)
348      !$jn = indklmn(10,klmn) !=indlmn(3,jlmn)
349 
350      in = indlmn(3,ilmn)
351      jn = indlmn(3,jlmn)
352 
353      ! Shift il to have the index instead of the value of l.
354      kln2ln(1,kln) = il +1
355      kln2ln(2,kln) = jl +1
356      kln2ln(3,kln) = in
357      kln2ln(4,kln) = jn
358      kln2ln(5,kln) = iln
359      kln2ln(6,kln) = jln
360    end if
361  end do !klmn
362 
363 end subroutine make_kln2ln

m_lmn_indices/uppert_index [ Functions ]

[ Top ] [ m_lmn_indices ] [ Functions ]

NAME

   uppert_index

FUNCTION

  Helper function returning the sequential index of an element in the upper triangle of a matrix 
  given the row-index ii and the column-index jj. If ii>jj the index of the element a_{jj,ii} is returned.

INPUTS

  ii=Row index
  jj=column index

PARENTS

SOURCE

617 function uppert_index(ii,jj)
618 
619 
620 !This section has been created automatically by the script Abilint (TD).
621 !Do not modify the following lines by hand.
622 #undef ABI_FUNC
623 #define ABI_FUNC 'uppert_index'
624 !End of the abilint section
625 
626  implicit none
627 
628 !Arguments ------------------------------------
629 !scalars
630  integer,intent(in) :: ii,jj
631  integer :: uppert_index
632 !scalars
633 
634 !************************************************************************
635 
636  if (jj>=jj) then
637    uppert_index = ii + jj*(jj-1)/2
638  else
639    uppert_index = jj + ii*(ii-1)/2
640  end if
641 
642 end function uppert_index