TABLE OF CONTENTS


ABINIT/wvl_setBoxGeometry [ Functions ]

[ Top ] [ Functions ]

NAME

 wvl_setBoxGeometry

FUNCTION

 When wavelets are used, the box definition needs to be changed.
 The box size is recomputed knowing some psp informations such as
 the radius for coarse and fine grid. Then, the atoms are translated
 to be included in the new box. Finally the FFT grid is computed using
 the fine wavelet mesh and a buffer characteristic of used wavelets plus
 a buffer used to be multiple of 2, 3 or 5.

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DC)
 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

  psps <type(pseudopotential_type)>=variables related to pseudopotentials
  radii= the radii for each type of atoms, giving the fine and the coarse
         grid.

OUTPUT

  rprimd(3,3)=dimensional primitive translations in real space (bohr)

SIDE EFFECTS

  wvl <type(wvl_internal_type)>=internal variables used by wavelets, describing
                             the box are set.
  xred(3,natom)=reduced dimensionless atomic coordinates

PARENTS

      gstate,wvl_memory,wvl_wfsinp_reformat

CHILDREN

      mkrdim,nullify_locreg_descriptors,system_size,wrtout,xcart2xred
      xred2xcart

SOURCE

 43 #if defined HAVE_CONFIG_H
 44 #include "config.h"
 45 #endif
 46 
 47 #include "abi_common.h"
 48 
 49 subroutine wvl_setBoxGeometry(prtvol, radii, rprimd, xred, wvl, wvl_crmult, wvl_frmult)
 50 
 51  use defs_basis
 52  use defs_wvltypes
 53  use m_profiling_abi
 54  use m_errors
 55 
 56  use m_geometry,      only : mkrdim, xcart2xred, xred2xcart
 57 #if defined HAVE_BIGDFT
 58  use BigDFT_API, only: system_size,nullify_locreg_descriptors
 59 #endif
 60 
 61 !This section has been created automatically by the script Abilint (TD).
 62 !Do not modify the following lines by hand.
 63 #undef ABI_FUNC
 64 #define ABI_FUNC 'wvl_setBoxGeometry'
 65  use interfaces_14_hidewrite
 66 !End of the abilint section
 67 
 68  implicit none
 69 
 70 !Arguments ------------------------------------
 71 !scalars
 72  integer,intent(in) :: prtvol
 73  real(dp), intent(in) :: wvl_crmult, wvl_frmult
 74  type(wvl_internal_type), intent(inout) :: wvl
 75 !arrays
 76  real(dp),intent(in) :: radii(:,:)
 77  real(dp),intent(inout) :: rprimd(3,3),xred(:,:)
 78 
 79 !Local variables-------------------------------
 80 #if defined HAVE_BIGDFT
 81 !scalars
 82  integer :: ii
 83  logical,parameter :: OCLconv=.false.
 84  character(len=500) :: message
 85 !arrays
 86  real(dp) :: rprim(3,3),acell(3)
 87  real(dp),allocatable :: xcart(:,:)
 88 #endif
 89 
 90 ! *********************************************************************
 91 
 92 #if defined HAVE_BIGDFT
 93  if (prtvol == 0) then
 94    write(message, '(a,a,a,a)' ) ch10,&
 95 &   ' wvl_setBoxGeometry : Changing the box for wavelets computation.'
 96    call wrtout(std_out,message,'COLL')
 97  end if
 98 
 99 !Store xcart for each atom
100  ABI_ALLOCATE(xcart,(3, wvl%atoms%astruct%nat))
101  call xred2xcart(wvl%atoms%astruct%nat, rprimd, xcart, xred)
102 
103  call nullify_locreg_descriptors(wvl%Glr)
104  call system_size(wvl%atoms, xcart, radii, wvl_crmult, &
105 & wvl_frmult, wvl%h(1), wvl%h(2), wvl%h(3), OCLconv, wvl%Glr, wvl%shift)
106 
107  acell(:) = wvl%atoms%astruct%cell_dim(:)
108 
109  if (prtvol == 0) then
110    write(message, '(a,3F12.6)' ) &
111 &   '  | acell is now:         ', acell
112    call wrtout(std_out,message,'COLL')
113    write(message, '(a,2I5,a,a,2I5,a,a,2I5)' ) &
114 &   '  | nfl1, nfu1:           ', wvl%Glr%d%nfl1, wvl%Glr%d%nfu1, ch10, &
115 &   '  | nfl2, nfu2:           ', wvl%Glr%d%nfl2, wvl%Glr%d%nfu2, ch10, &
116 &   '  | nfl3, nfu3:           ', wvl%Glr%d%nfl3, wvl%Glr%d%nfu3
117    call wrtout(std_out,message,'COLL')
118  end if
119 
120 !Change the metric to orthogonal one
121  rprim(:, :) = real(0., dp)
122  do ii = 1, 3, 1
123    rprim(ii,ii) = real(1., dp)
124  end do
125  call mkrdim(acell, rprim, rprimd)
126 
127 !Save shifted atom positions into xred
128  call xcart2xred(wvl%atoms%astruct%nat, rprimd, xcart, xred)
129  ABI_DEALLOCATE(xcart)
130 
131  if (prtvol == 0) then
132    write(message, '(a,3I12)' ) &
133 &   '  | box size for datas:   ', wvl%Glr%d%n1i, wvl%Glr%d%n2i, wvl%Glr%d%n3i
134    call wrtout(std_out,message,'COLL')
135    write(message, '(a,3I12)' ) &
136 &   '  | box size for wavelets:', wvl%Glr%d%n1, wvl%Glr%d%n2, wvl%Glr%d%n3
137    call wrtout(std_out,message,'COLL')
138  end if
139 
140 #else
141  BIGDFT_NOTENABLED_ERROR()
142  if (.false.) write(std_out,*)  prtvol,wvl_crmult,wvl_frmult,wvl%h(1),&
143 & radii(1,1),rprimd(1,1),xred(1,1)
144 #endif
145 
146 end subroutine wvl_setBoxGeometry