TABLE OF CONTENTS


defs_wvltypes/wvl_descr_psp_fill [ Functions ]

[ Top ] [ defs_wvltypes ] [ Functions ]

NAME

 wvl_descr_psp_fill

FUNCTION

 Read the radii from @pspunit, if any and fill values with default otherwise.

INPUTS

OUTPUT

PARENTS

      psp10in,psp2in,psp3in,pspatm

CHILDREN

      atomic_info,psp_from_data,wrtout

SOURCE

148 subroutine wvl_descr_psp_fill(gth_params, ipsp, ixc, nelpsp, nzatom, pspunit)
149 
150   use defs_datatypes
151   use m_errors
152 #if defined HAVE_BIGDFT
153   use BigDFT_API, only: atomic_info, UNINITIALIZED, psp_from_data
154 #endif
155 
156 !This section has been created automatically by the script Abilint (TD).
157 !Do not modify the following lines by hand.
158 #undef ABI_FUNC
159 #define ABI_FUNC 'wvl_descr_psp_fill'
160  use interfaces_14_hidewrite
161 !End of the abilint section
162 
163   implicit none
164 
165 !Arguments ------------------------------------
166   integer, intent(in) :: ipsp, pspunit, nzatom, nelpsp, ixc
167   type(pseudopotential_gth_type), intent(inout) :: gth_params
168 !Local variables-------------------------------
169 #if defined HAVE_BIGDFT
170   integer :: ios, ii, nzatom_, nelpsp_, npspcode_, ixc_
171   real(dp) :: ehomo, radfine
172   logical :: exists
173   character(len = 2) :: symbol
174   character(len=100) :: line
175   character(len=500) :: message
176 #endif
177 
178 ! ***************************************************************************
179 
180 #if defined HAVE_BIGDFT
181 
182   ! check if gth_params%psppar have been set
183  if (any(gth_params%psppar == UNINITIALIZED(1._dp))) then
184    call atomic_info(nzatom, nelpsp, symbol = symbol)
185    ixc_ = ixc
186    if (ixc_>0.and.ixc_< 10) ixc_=1
187    if (ixc_>0.and.ixc_>=10) ixc_=11
188    call psp_from_data(symbol, nzatom_, nelpsp_, npspcode_, ixc_, &
189 &   gth_params%psppar(:,:,ipsp), exists)
190    if(.not. exists) then
191      write(message,'(a,a,a,a)')ch10,&
192 &     "wvl_descr_psp_fill : bug, chemical element not found in BigDFT table",ch10,&
193 &     "Action: upgrade BigDFT table"
194      call wrtout(ab_out,message,'COLL')
195      MSG_BUG(message)
196    end if
197    gth_params%set(ipsp) = .true.
198  end if
199 
200   ! Try to read radii from pspunit
201  if (pspunit /= 0) then
202    read (pspunit, '(a100)', iostat = ios) line
203    if (ios /= 0) then
204      line=''
205    end if
206      !  We try to read the values from the pseudo.
207    read (line, *, iostat = ios) gth_params%radii_cf(ipsp, 1), gth_params%radii_cf(ipsp, 2), &
208 &   gth_params%radii_cf(ipsp, 3)
209    if (ios /= 0 .or. gth_params%radii_cf(ipsp, 3) < zero) then
210      read (line, *, iostat = ios) gth_params%radii_cf(ipsp, 1), gth_params%radii_cf(ipsp, 2)
211      gth_params%radii_cf(ipsp, 3) = UNINITIALIZED(gth_params%radii_cf(ipsp, 3))
212    end if
213    if (ios /= 0 .or. gth_params%radii_cf(ipsp, 1) < zero .or. gth_params%radii_cf(ipsp, 2) < zero) then
214      gth_params%radii_cf(ipsp, 1) = UNINITIALIZED(gth_params%radii_cf(ipsp, 1))
215      gth_params%radii_cf(ipsp, 2) = UNINITIALIZED(gth_params%radii_cf(ipsp, 2))
216      write(message, '(a,a,a,a,a,a,a)' ) '-', ch10,&
217 &     '- wvl_descr_psp_fill : COMMENT -',ch10,&
218 &     "-  the pseudo-potential does not include geometric informations,",ch10,&
219 &     '-  values have been computed.'
220      call wrtout(ab_out,message,'COLL')
221      call wrtout(std_out,  message,'COLL')
222    end if
223  end if
224 
225   ! Update radii.
226  if (gth_params%radii_cf(ipsp, 1) == UNINITIALIZED(gth_params%radii_cf(ipsp, 1))) then
227    call atomic_info(nzatom, nelpsp, ehomo = ehomo)
228      !assigning the radii by calculating physical parameters
229    gth_params%radii_cf(ipsp, 1)=1._dp/sqrt(abs(2._dp*ehomo))
230  end if
231  if (gth_params%radii_cf(ipsp, 2) == UNINITIALIZED(gth_params%radii_cf(ipsp, 2))) then
232    radfine = gth_params%psppar(0, 0, ipsp)
233    do ii = 1, 4, 1
234      if (gth_params%psppar(ii, 0, ipsp) /= zero) then
235        radfine = min(radfine, gth_params%psppar(ii, 0, ipsp))
236      end if
237    end do
238    gth_params%radii_cf(ipsp, 2)=radfine
239  end if
240  if (gth_params%radii_cf(ipsp, 3) == UNINITIALIZED(gth_params%radii_cf(ipsp, 3))) then
241    gth_params%radii_cf(ipsp, 3)=gth_params%radii_cf(ipsp, 2)
242  end if
243 
244   ! Set flag.
245  if (gth_params%radii_cf(ipsp, 1) >= 0.d0 .and. gth_params%radii_cf(ipsp, 2) >= 0.d0) then
246    gth_params%hasGeometry(ipsp) = .true.
247  else
248    gth_params%hasGeometry(ipsp) = .false.
249  end if
250 
251  write(message, '(a,f12.7,a,f12.7,a,f12.7)' )&
252 & '  radii_cf(1)=', gth_params%radii_cf(ipsp, 1),&
253 & '; radii_cf(2)=', gth_params%radii_cf(ipsp, 2),&
254 & '; rad_cov=', gth_params%radii_cf(ipsp, 3)
255  call wrtout(ab_out,message,'COLL')
256  call wrtout(std_out,  message,'COLL')
257 
258 ! Some consistency checks on radii, to be moved earlier, but need to update test refs.
259 ! gth_params%radii_cf(ipsp, 3) = min(gth_params%radii_cf(ipsp, 3), gth_params%radii_cf(ipsp, 2))
260 
261 ! This was one before
262 ! maxrad=zero
263 ! do ii=0,2,1
264 !   if (ii==1) maxrad=zero
265 !   if (gth_params%psppar(ii,0,ipsp)/=zero) maxrad=max(maxrad,gth_params%psppar(ii,0,ipsp))
266 ! end do
267 ! if (maxrad== zero) then
268 !   gth_params%radii_cf(ipsp,3)=zero
269 ! else
270 !   gth_params%radii_cf(ipsp,3)=max( &
271 !&    min(dtset%wvl_crmult*psps%gth_params%radii_cf(ipsp,1), &
272 !&    15._dp*maxrad)/dtset%wvl_frmult,gth_params%radii_cf(ipsp,2))
273 ! end if
274 
275 #else
276  BIGDFT_NOTENABLED_ERROR()
277  if (.false.) write(std_out,*) ipsp,pspunit,nzatom,nelpsp,ixc,gth_params%psppar
278 #endif
279 
280 end subroutine wvl_descr_psp_fill

defs_wvltypes/wvl_descr_psp_set [ Functions ]

[ Top ] [ defs_wvltypes ] [ Functions ]

NAME

 wvl_descr_psp_set

FUNCTION

 Defines the part of the wvl%atoms%-datastructure which
 depends on psps

INPUTS

 psps <type(pseudopotential_type)>=variables related to pseudopotentials
 dtset <type(dataset_type)>=all input variables for this dataset

OUTPUT

 nsppol=number of spin components
 spinat(3,natom)=spin polarization on each atom
 wvl <type(wvl_internal_type)> = wavelet type
                 | psppar   = The covalence radii for each pseudo 
                 | pspcod   = the format -or code- of psp generation
                 | iasctype = semicore code (see defs_datatype)
                 | nzatom   = charge of the nucleus
                 | nelpsp   = the ionic pseudo-charge
                 | natsc    = number of atoms with semicore 

PARENTS

      gstate

CHILDREN

      atomic_info,psp_from_data,wrtout

SOURCE

 32 #if defined HAVE_CONFIG_H
 33 #include "config.h"
 34 #endif
 35 
 36 #include "abi_common.h"
 37 
 38 subroutine wvl_descr_psp_set(filoccup, nsppol, psps, spinat, wvl)
 39 
 40  use m_profiling_abi
 41  use m_errors
 42 
 43  use defs_basis
 44  use defs_datatypes
 45  use defs_wvltypes
 46 #if defined HAVE_BIGDFT
 47  use BigDFT_API, only: aoig_set,UNINITIALIZED,dict_init,dict_free,dictionary, &
 48 &                      operator(//), bigdft_mpi, dict_set
 49  use BigDFT_API, only: psp_data_merge_to_dict, psp_dict_fill_all, atomic_info, &
 50 &                psp_dict_analyse, atomic_data_set_from_dict, merge_input_file_to_dict
 51 #endif
 52 
 53 !This section has been created automatically by the script Abilint (TD).
 54 !Do not modify the following lines by hand.
 55 #undef ABI_FUNC
 56 #define ABI_FUNC 'wvl_descr_psp_set'
 57 !End of the abilint section
 58 
 59   implicit none
 60 
 61 !Arguments ------------------------------------
 62 !scalars
 63   integer, intent(in)                    :: nsppol
 64   type(wvl_internal_type), intent(inout) :: wvl
 65   type(pseudopotential_type), intent(in) :: psps
 66   character(len = *), intent(in) :: filoccup
 67 !arrays
 68   real(dp),intent(in) :: spinat(:,:)
 69 
 70 !Local variables-------------------------------
 71 #if defined HAVE_BIGDFT
 72   integer :: ityp,pspcod
 73   logical :: exists
 74   real(dp) :: radii_cf(3)
 75   character(len=2) :: symbol
 76   character(len=27) :: filename
 77   type(dictionary), pointer :: dict
 78 #endif
 79 
 80 ! *********************************************************************
 81 
 82 #if defined HAVE_BIGDFT
 83 
 84 !We create the atoms_data structure, the part that is dependent from psp.
 85  do ityp=1,size(psps%pspcod)
 86    wvl%atoms%psppar(:,:,ityp)= psps%gth_params%psppar(:,:,ityp)
 87    wvl%atoms%npspcode(ityp)  = merge(psps%pspcod(ityp),7,psps%pspcod(ityp)/=17)
 88    wvl%atoms%ixcpsp(ityp)    = psps%pspxc(ityp)
 89    wvl%atoms%nzatom(ityp)    = int(psps%znucltypat(ityp))
 90    wvl%atoms%nelpsp(ityp)    = int(psps%ziontypat(ityp))
 91  end do
 92  wvl%atoms%donlcc = .false.
 93 
 94  call dict_init(dict)
 95  radii_cf = UNINITIALIZED(1._dp)
 96  do ityp = 1, wvl%atoms%astruct%ntypes, 1
 97    call atomic_info(wvl%atoms%nzatom(ityp), wvl%atoms%nelpsp(ityp), &
 98 &   symbol = symbol)
 99    write(wvl%atoms%astruct%atomnames(ityp), "(A)") symbol
100    filename = 'psppar.' // trim(symbol)
101    pspcod=psps%pspcod(ityp);if (pspcod==17) pspcod=7 ! PAW specificity
102    call psp_data_merge_to_dict(dict // filename, int(psps%znucltypat(ityp)), &
103 &   int(psps%ziontypat(ityp)), pspcod, psps%pspxc(ityp), &
104 &   psps%gth_params%psppar(0:4,0:6,ityp), radii_cf, &
105 &   UNINITIALIZED(1._dp), UNINITIALIZED(1._dp))
106    call psp_dict_fill_all(dict, trim(symbol), psps%pspxc(ityp))
107  end do
108 
109  call psp_dict_analyse(dict, wvl%atoms)
110  inquire(file = filoccup, exist = exists)
111  if (exists) then
112    call merge_input_file_to_dict(dict//"Atomic occupation",filoccup,bigdft_mpi)
113  end if
114  ! Need to add spinat in the dictionary, later
115  if (spinat(1,1) == zero .and. .false.) then
116    call dict_set(dict // "spinat", spinat(1,1))
117  end if
118  call atomic_data_set_from_dict(dict, "Atomic occupation", wvl%atoms, nsppol)
119  call dict_free(dict)
120 
121 #else
122  BIGDFT_NOTENABLED_ERROR()
123  if (.false.) write(std_out,*) nsppol,wvl%h(1),psps%npsp,filoccup,spinat(1,1)
124 #endif  
125 
126 end subroutine wvl_descr_psp_set