TABLE OF CONTENTS


ABINIT/m_pawang [ Modules ]

[ Top ] [ Modules ]

NAME

  m_pawang

FUNCTION

  This module contains the definition of the pawang_type structured datatype,
  as well as related functions and methods.
  pawang_type variables define ANGular mesh discretization of PAW augmentation
  regions and related data.

COPYRIGHT

 Copyright (C) 2013-2024 ABINIT group (MT, FJ, BA)
 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

  FOR DEVELOPPERS: in order to preserve the portability of libPAW library,
  please consult ~abinit/src/??_libpaw/libpaw-coding-rules.txt

SOURCE

23 #include "libpaw.h"
24 
25 MODULE m_pawang
26 
27  USE_DEFS
28  USE_MSG_HANDLING
29  USE_MPI_WRAPPERS
30  USE_MEMORY_PROFILING
31 
32  use m_paw_sphharm, only : initylmr, ylm_angular_mesh, mat_mlms2jmj, mat_slm2ylm, &
33 &                          lsylm, realgaunt, nablarealgaunt
34 
35  implicit none
36 
37  private
38 
39 !Public procedures.
40  public :: pawang_init       ! Constructor
41  public :: pawang_free       ! Free memory

m_pawang/pawang_free [ Functions ]

[ Top ] [ m_pawang ] [ Functions ]

NAME

 pawang_free

FUNCTION

  Free all dynamic memory and reset all flags stored in a pawang datastructure

SIDE EFFECTS

  Pawang <type(pawang_type)>=ANGular mesh discretization and related data

SOURCE

311 subroutine pawang_free(Pawang)
312 
313 !Arguments ------------------------------------
314 !scalars
315  type(Pawang_type),intent(inout) :: Pawang
316 
317 ! *************************************************************************
318 
319  !@Pawang_type
320  if (allocated(pawang%angwgth))    then
321    LIBPAW_DEALLOCATE(pawang%angwgth)
322  end if
323  if (allocated(pawang%anginit))    then
324    LIBPAW_DEALLOCATE(pawang%anginit)
325  end if
326  if (allocated(pawang%zarot))      then
327    LIBPAW_DEALLOCATE(pawang%zarot)
328  end if
329  if (allocated(pawang%gntselect))  then
330    LIBPAW_DEALLOCATE(pawang%gntselect)
331  end if
332   if (allocated(pawang%nablagntselect))  then
333    LIBPAW_DEALLOCATE(pawang%nablagntselect)
334  end if
335  if (allocated(pawang%realgnt))    then
336    LIBPAW_DEALLOCATE(pawang%realgnt)
337  end if
338  if (allocated(pawang%nablarealgnt))    then
339    LIBPAW_DEALLOCATE(pawang%nablarealgnt)
340  end if
341  if (allocated(pawang%ylmr))       then
342    LIBPAW_DEALLOCATE(pawang%ylmr)
343  end if
344  if (allocated(pawang%ylmrgr))     then
345    LIBPAW_DEALLOCATE(pawang%ylmrgr)
346  end if
347  if (allocated(pawang%ls_ylm))     then
348    LIBPAW_DEALLOCATE(pawang%ls_ylm)
349  end if
350 
351  pawang%angl_size =0
352  pawang%ylm_size  =0
353  pawang%use_ls_ylm=0
354  pawang%l_max=-1
355  pawang%l_size_max=-1
356  pawang%gnt_option=-1
357  pawang%nabgnt_option=-1
358  pawang%ngnt=0
359 
360 end subroutine pawang_free

m_pawang/pawang_init [ Functions ]

[ Top ] [ m_pawang ] [ Functions ]

NAME

 pawang_init

FUNCTION

  Initialize a pawang datastructure

INPUTS

  gnt_option=flag activated if pawang%gntselect and pawang%realgnt have to be allocated
             also determine the size of these pointers
  nabgnt_option=flag activated if pawang%nablagntselect and pawang%nablarealgnt have to be allocated
  lmax=maximum value of angular momentum l
  nphi,ntheta=dimensions of paw angular mesh
  nsym=number of symetries
  ngrad2_ylm=order of spherical harmonics gradients to be stored (0, 1 or 2)
  use_angular_grid=flag activated if angular grid data have to be allocated
                   (pawang%angwgth, pawang%anginit)
  use_ylm=flag activated if spherical harmonics have to be allocated and computed
          (pawang%ylmr, pawang%ylmrgr)
  use_ls_ylm=flag activated if LS operator has to be allocated and computed
          (pawang%ls_ylm)

OUTPUT

  Pawang <type(pawang_type)>=ANGular mesh discretization and related data

SOURCE

190 subroutine pawang_init(Pawang,gnt_option,nabgnt_option,lmax,nphi,ntheta,nsym,ngrad2_ylm,&
191 &                      use_angular_grid,use_ylm,use_ls_ylm)
192 
193 !Arguments ------------------------------------
194 !scalars
195  integer,intent(in) :: gnt_option,nabgnt_option,lmax,nphi,nsym,ntheta,ngrad2_ylm
196  integer,intent(in) :: use_angular_grid,use_ylm,use_ls_ylm
197  type(Pawang_type),intent(inout) :: Pawang
198 
199 !Local variables-------------------------------
200 !scalars
201  integer :: ll,sz1,sz2
202 !arrays
203  real(dp),allocatable :: rgnt_tmp(:)
204  real(dp),allocatable :: nablargnt_tmp(:)
205 
206 ! *************************************************************************
207 
208  !@Pawang_type
209 
210  Pawang%l_max=lmax+1
211  Pawang%l_size_max=2*Pawang%l_max-1
212  Pawang%nsym=nsym
213 
214  if (use_angular_grid==1) then
215    Pawang%nphi=nphi
216    Pawang%ntheta=ntheta
217    call ylm_angular_mesh(Pawang%ntheta,Pawang%nphi,Pawang%angl_size,Pawang%anginit,Pawang%angwgth)
218  else
219    Pawang%nphi=0
220    Pawang%ntheta=0
221    Pawang%angl_size=0
222  end if
223 
224  if (use_ylm>0.and.Pawang%angl_size>0) then
225    ll=Pawang%l_size_max+1
226    Pawang%ylm_size=ll**2
227    LIBPAW_ALLOCATE(Pawang%ylmr,(Pawang%ylm_size,Pawang%angl_size))
228    if (ngrad2_ylm==2) then
229      LIBPAW_ALLOCATE(pawang%ylmrgr,(9,Pawang%ylm_size,Pawang%angl_size))
230      call initylmr(ll,0,pawang%angl_size,pawang%angwgth,3,pawang%anginit,pawang%ylmr,&
231 &                  ylmr_gr=pawang%ylmrgr)
232    else if (ngrad2_ylm==1) then
233        LIBPAW_ALLOCATE(Pawang%ylmrgr,(3,Pawang%ylm_size,Pawang%angl_size))
234        call initylmr(ll,0,pawang%angl_size,pawang%angwgth,2,pawang%anginit,pawang%ylmr,&
235 &                    ylmr_gr=pawang%ylmrgr)
236    else
237      call initylmr(ll,0,pawang%angl_size,pawang%angwgth,1,pawang%anginit,pawang%ylmr)
238    end if
239  else
240    Pawang%ylm_size=0
241  end if
242 
243  Pawang%gnt_option=gnt_option
244  if (Pawang%gnt_option==1.or.Pawang%gnt_option==2) then
245    if (Pawang%gnt_option==1) then
246      sz1=(Pawang%l_size_max)**2
247      sz2=((Pawang%l_max**2)*(Pawang%l_max**2+1))/2
248      LIBPAW_ALLOCATE(rgnt_tmp,(sz1*sz2))
249      LIBPAW_ALLOCATE(pawang%gntselect,(sz1,sz2))
250      call realgaunt(Pawang%l_max,Pawang%ngnt,Pawang%gntselect,rgnt_tmp)
251    else if (Pawang%gnt_option==2) then
252      sz1=(2*Pawang%l_size_max-1)**2
253      sz2=((Pawang%l_size_max)**2*(Pawang%l_size_max**2+1))/2
254      LIBPAW_ALLOCATE(rgnt_tmp,(sz1*sz2))
255      LIBPAW_ALLOCATE(pawang%gntselect,(sz1,sz2))
256      call realgaunt(Pawang%l_size_max,Pawang%ngnt,Pawang%gntselect,rgnt_tmp)
257    end if
258    if (allocated(pawang%realgnt))  then
259      LIBPAW_DEALLOCATE(pawang%realgnt)
260    end if
261    LIBPAW_ALLOCATE(Pawang%realgnt,(Pawang%ngnt))
262    Pawang%realgnt(1:Pawang%ngnt)=rgnt_tmp(1:Pawang%ngnt)
263    LIBPAW_DEALLOCATE(rgnt_tmp)
264  end if
265 
266  Pawang%nabgnt_option=nabgnt_option
267  if (Pawang%nabgnt_option==1) then
268    sz1=(Pawang%l_size_max)**2
269    sz2=(Pawang%l_max)**2
270    LIBPAW_ALLOCATE(nablargnt_tmp,(sz1*sz2*sz2))
271    LIBPAW_ALLOCATE(pawang%nablagntselect,(sz1,sz2,sz2))
272    call nablarealgaunt(pawang%l_size_max,pawang%l_max, &
273 &                      pawang%nnablagnt,pawang%nablagntselect,nablargnt_tmp)
274    if (allocated(pawang%nablarealgnt)) then
275      LIBPAW_DEALLOCATE(pawang%nablarealgnt)
276    end if
277    LIBPAW_ALLOCATE(pawang%nablarealgnt,(pawang%nnablagnt))
278    Pawang%nablarealgnt(1:Pawang%nnablagnt)=nablargnt_tmp(1:Pawang%nnablagnt)
279    LIBPAW_DEALLOCATE(nablargnt_tmp)
280  end if
281 
282  Pawang%use_ls_ylm=use_ls_ylm
283  if (use_ls_ylm>0) then
284    LIBPAW_ALLOCATE(pawang%ls_ylm,(2,Pawang%l_max**2*(Pawang%l_max**2+1)/2,2))
285    call lsylm(pawang%ls_ylm,lmax)
286  end if
287 
288  if (nsym>0) then
289    if (.not.allocated(Pawang%zarot)) then
290      LIBPAW_ALLOCATE(Pawang%zarot,(Pawang%l_size_max,Pawang%l_size_max,Pawang%l_max,nsym))
291    end if
292  end if
293 
294 end subroutine pawang_init

m_pawang/pawang_type [ Types ]

[ Top ] [ m_pawang ] [ Types ]

NAME

 pawang_type

FUNCTION

 For PAW: ANGular mesh discretization of PAW augmentation regions and related data.

SOURCE

 57  type,public :: pawang_type
 58 
 59 !Integer scalars
 60 
 61   integer :: angl_size=0
 62    ! Dimension of paw angular mesh
 63    ! angl_size=ntheta*nphi
 64 
 65   integer :: l_max=-1
 66    ! Maximum value of angular momentum l+1
 67 
 68   integer :: l_size_max=-1
 69    ! Maximum value of angular momentum +1
 70    ! leading to non-zero Gaunt coefficients
 71    ! l_size_max = 2*l_max-1
 72 
 73   integer :: ngnt=0
 74    ! Number of non zero Gaunt coefficients
 75 
 76   integer :: nnablagnt=0
 77    ! Number of non zero Gaunt coefficient derivatives
 78 
 79   integer :: ntheta, nphi
 80    ! Dimensions of paw angular mesh
 81 
 82   integer :: nsym
 83    ! Number of symmetry elements in space group
 84 
 85   integer :: nabgnt_option=-1
 86    ! Option for nablarealgaunt coefficients:
 87    ! nabgnt_option==0, nablarealgaunt coeffs are not computed (and not allocated)
 88    ! nabgnt_option==1, nablarealgaunt coeffs are computed up to l_max
 89 
 90   integer :: gnt_option=-1
 91    ! Option for Gaunt coefficients:
 92    !  gnt_option==0, Gaunt coeffs are not computed (and not allocated)
 93    !  gnt_option==1, Gaunt coeffs are computed up to 2*l_size_max-1
 94    !  gnt_option==2, Gaunt coeffs are computed up to l_size_max
 95 
 96   integer :: use_ls_ylm=0
 97    ! Flag: use_ls_ylm=1 if pawang%ls_ylm is allocated
 98 
 99   integer :: ylm_size=0
100    ! Size of ylmr/ylmrgr arrays
101 
102 !Integer arrays
103 
104   integer, allocatable :: gntselect(:,:)
105    ! gntselect(l_size_max**2,l_max**2*(l_max**2+1)/2)
106    ! Selection rules for Gaunt coefficients stored as (LM,ij) where ij is in packed form.
107    ! (if gntselect>0, Gaunt coeff. is non-zero)
108 
109   integer, allocatable :: nablagntselect(:,:,:)
110   ! nablagntselect(l_size_max**2,l_max**2,l_max**2)
111   ! Selection rules for nablaGaunt coefficients
112   ! (if nablagntselect>0, nablGaunt coeff. is non-zero)
113 
114 !Real (real(dp)) arrays
115 
116   real(dp), allocatable :: anginit(:,:)
117    ! anginit(3,angl_size)
118    ! For each point of the angular mesh, gives the coordinates
119    ! of the corresponding point on an unitary sphere
120    ! Not used in present version (5.3)
121 
122   real(dp), allocatable :: angwgth(:)
123    ! angwgth(angl_size)
124    ! For each point of the angular mesh, gives the weight
125    ! of the corresponding point on an unitary sphere
126 
127   real(dp), allocatable :: ls_ylm(:,:,:)
128    ! ls_ylm(2,l_max**2*(l_max**2+1)/2,2)
129    ! LS operator in the real spherical harmonics basis
130    ! ls_ylm(ilm1m2,ispin)= <sigma, y_lm1| LS |y_lm2, sigma_prime>
131 
132   real(dp), allocatable :: nablarealgnt(:)
133    ! nablarealgnt(2,nnablagnt)
134    ! Non zero real nablaGaunt coefficients
135 
136   real(dp), allocatable :: realgnt(:)
137    ! realgnt(ngnt)
138    ! Non zero real Gaunt coefficients
139 
140   real(dp), allocatable :: ylmr(:,:)
141    ! ylmr(ylm_size,angl_size)
142    ! Real Ylm calculated in real space
143 
144   real(dp), allocatable :: ylmrgr(:,:,:)
145    ! ylmrgr(1:3,ylm_size,angl_size)
146    ! First gradients of real Ylm calculated in real space (cart. coordinates)
147 
148   real(dp), allocatable :: zarot(:,:,:,:)
149    !  zarot(l_size_max,l_size_max,l_max,nsym)
150    !  Coeffs of the transformation of real spherical
151    !  harmonics under the symmetry operations symrec.
152 
153  end type pawang_type