TABLE OF CONTENTS
- ABINIT/m_wvl_descr_psp
- defs_wvltypes/wvl_descr_atoms_set
- defs_wvltypes/wvl_descr_atoms_set_sym
- defs_wvltypes/wvl_descr_free
- defs_wvltypes/wvl_descr_psp_fill
- defs_wvltypes/wvl_descr_psp_set
ABINIT/m_wvl_descr_psp [ Modules ]
NAME
wvl_descr_psp
FUNCTION
COPYRIGHT
Copyright (C) 2008-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 .
PARENTS
CHILDREN
SOURCE
21 #if defined HAVE_CONFIG_H 22 #include "config.h" 23 #endif 24 25 #include "abi_common.h" 26 27 module m_wvl_descr_psp 28 29 use defs_basis 30 use m_errors 31 use m_abicore 32 33 implicit none 34 35 private
defs_wvltypes/wvl_descr_atoms_set [ Functions ]
[ Top ] [ defs_wvltypes ] [ Functions ]
NAME
wvl_descr_atoms_set
FUNCTION
Defines wvl%atoms% data structure
INPUTS
acell(3)=unit cell length scales (bohr) dtset <type(dataset_type)>=all input variables for this dataset
OUTPUT
wvl <type(wvl_internal_type)>= wavelet type | nat = number of atoms | ntypes = number of species | alat1 = acell(1) | alat2 = acell(2) | alat3 = acell(3) | iatype = types for atoms | lfrztyp = flag for the movement of atoms. | natpol = integer related to polarisation at the first step
PARENTS
gstate,wvl_memory
CHILDREN
allocate_atoms_nat,allocate_atoms_ntypes,astruct_set_n_atoms astruct_set_n_types,f_release_routine,f_routine
SOURCE
412 subroutine wvl_descr_atoms_set(acell, icoulomb, natom, ntypat, typat, wvl) 413 414 use defs_wvltypes 415 #if defined HAVE_BIGDFT 416 use BigDFT_API, only: atoms_data_null,f_routine,f_release_routine,& 417 & astruct_set_n_atoms,astruct_set_n_types,& 418 & allocate_atoms_nat,allocate_atoms_ntypes 419 #endif 420 421 !This section has been created automatically by the script Abilint (TD). 422 !Do not modify the following lines by hand. 423 #undef ABI_FUNC 424 #define ABI_FUNC 'wvl_descr_atoms_set' 425 !End of the abilint section 426 427 implicit none 428 429 !Arguments ------------------------------------ 430 !scalars 431 integer, intent(in) :: icoulomb, natom, ntypat 432 type(wvl_internal_type), intent(inout) :: wvl 433 !arrays 434 integer, intent(in) :: typat(natom) 435 real(dp), intent(in) :: acell(3) 436 437 !Local variables------------------------------- 438 !scalars 439 #if defined HAVE_BIGDFT 440 integer :: itype 441 #endif 442 443 ! ********************************************************************* 444 445 #if defined HAVE_BIGDFT 446 447 call f_routine(ABI_FUNC) 448 449 wvl%atoms=atoms_data_null() 450 451 !We create the atoms_data structure from this dataset 452 !to be used later in BigDFT routines. 453 if (icoulomb == 0) then 454 wvl%atoms%astruct%geocode = 'P' 455 else if (icoulomb == 1) then 456 wvl%atoms%astruct%geocode = 'F' 457 else if (icoulomb == 2) then 458 wvl%atoms%astruct%geocode = 'S' 459 end if 460 write(wvl%atoms%astruct%units, "(A)") "Bohr" 461 462 call astruct_set_n_atoms(wvl%atoms%astruct, natom) 463 call astruct_set_n_types(wvl%atoms%astruct, ntypat) 464 465 do itype = 1, ntypat, 1 466 write(wvl%atoms%astruct%atomnames(itype), "(A,I2)") "At. type", itype 467 end do 468 wvl%atoms%astruct%cell_dim(1) = acell(1) 469 wvl%atoms%astruct%cell_dim(2) = acell(2) 470 wvl%atoms%astruct%cell_dim(3) = acell(3) 471 wvl%atoms%astruct%iatype = typat 472 473 wvl%atoms%astruct%sym%symObj = 0 474 475 call allocate_atoms_nat(wvl%atoms) 476 call allocate_atoms_ntypes(wvl%atoms) 477 478 call f_release_routine() 479 480 #else 481 BIGDFT_NOTENABLED_ERROR() 482 if (.false.) write(std_out,*) icoulomb,natom,ntypat,wvl%h(1),typat(1),acell(1) 483 #endif 484 485 end subroutine wvl_descr_atoms_set
defs_wvltypes/wvl_descr_atoms_set_sym [ Functions ]
[ Top ] [ defs_wvltypes ] [ Functions ]
NAME
wvl_descr_atoms_set_sym
FUNCTION
Add symmetry information to wvl%atoms data structure.
INPUTS
acell(3)=unit cell length scales (bohr) dtset <type(dataset_type)>=all input variables for this dataset
OUTPUT
wvl <type(wvl_internal_type)>= wavelet type | nat = number of atoms | ntypes = number of species | alat1 = acell(1) | alat2 = acell(2) | alat3 = acell(3) | iatype = types for atoms | lfrztyp = flag for the movement of atoms. | natpol = integer related to polarisation at the first step
PARENTS
gstate
CHILDREN
astruct_set_symmetries,symmetry_set_n_sym,wrtout
SOURCE
519 subroutine wvl_descr_atoms_set_sym(wvl, efield, irrzon, nsppol, nsym, phnons, & 520 & symafm, symrel, tnons, tolsym) 521 522 use defs_wvltypes 523 use m_ab7_symmetry 524 #if defined HAVE_BIGDFT 525 use BigDFT_API, only: astruct_set_symmetries 526 #endif 527 528 !This section has been created automatically by the script Abilint (TD). 529 !Do not modify the following lines by hand. 530 #undef ABI_FUNC 531 #define ABI_FUNC 'wvl_descr_atoms_set_sym' 532 !End of the abilint section 533 534 implicit none 535 536 !Arguments ------------------------------------ 537 !scalars 538 integer, intent(in) :: nsppol,nsym 539 real(dp), intent(in) :: tolsym 540 type(wvl_internal_type), intent(inout) :: wvl 541 !arrays 542 integer, intent(in) :: symafm(3,3,nsym), symrel(3,3,nsym) 543 integer, target, intent(in) :: irrzon(:,:,:) 544 real(dp), target, intent(in) :: phnons(:,:,:) 545 real(dp), intent(in) :: efield(3), tnons(3,nsym) 546 547 !Local variables------------------------------- 548 !scalars 549 #if defined HAVE_BIGDFT 550 integer :: errno 551 character(len=500) :: message 552 #endif 553 554 ! ********************************************************************* 555 556 #if defined HAVE_BIGDFT 557 558 write(message, '(a,a)' ) ch10,& 559 & ' wvl_descr_atoms_set_sym: Create symmetries for the wvl%atoms object.' 560 call wrtout(std_out,message,'COLL') 561 562 wvl%atoms%astruct%sym%symObj = -1 563 nullify(wvl%atoms%astruct%sym%irrzon) 564 nullify(wvl%atoms%astruct%sym%phnons) 565 call astruct_set_symmetries(wvl%atoms%astruct, (nsym <= 1), tolsym, efield, nsppol) 566 if (nsym > 1) then 567 call symmetry_set_n_sym(wvl%atoms%astruct%sym%symObj, nsym, symrel, tnons, symafm, errno) 568 end if 569 wvl%atoms%astruct%sym%irrzon => irrzon 570 wvl%atoms%astruct%sym%phnons => phnons 571 572 #else 573 BIGDFT_NOTENABLED_ERROR() 574 if (.false.) write(std_out,*) nsppol,nsym,tolsym,wvl%h(1),symafm(1,1,1),symrel(1,1,1),& 575 & irrzon(1,1,1),phnons(1,1,1),efield(1),tnons(1,1) 576 #endif 577 578 end subroutine wvl_descr_atoms_set_sym
defs_wvltypes/wvl_descr_free [ Functions ]
[ Top ] [ defs_wvltypes ] [ Functions ]
NAME
wvl_descr_free
FUNCTION
Free the wvl%atoms% datastructure (deallocate or nullify)
INPUTS
wvl <type(wvl_internal_type)>=internal variables for wavelets
OUTPUT
wvl <type(wvl_internal_type)>=internal variables for wavelets
PARENTS
gstate,wvl_memory
CHILDREN
deallocate_atoms_data,f_free_ptr
SOURCE
339 subroutine wvl_descr_free(wvl) 340 341 use defs_wvltypes 342 #if defined HAVE_BIGDFT 343 use BigDFT_API, only : deallocate_atoms_data 344 use dynamic_memory 345 #endif 346 347 !This section has been created automatically by the script Abilint (TD). 348 !Do not modify the following lines by hand. 349 #undef ABI_FUNC 350 #define ABI_FUNC 'wvl_descr_free' 351 !End of the abilint section 352 353 implicit none 354 355 !Arguments ------------------------------------ 356 !scalars 357 type(wvl_internal_type), intent(inout) :: wvl 358 !arrays 359 360 !Local variables------------------------------- 361 !scalars 362 363 ! ********************************************************************* 364 365 #if defined HAVE_BIGDFT 366 !These arrays are pointers on memory handled by ABINIT. 367 nullify(wvl%atoms%astruct%sym%irrzon) 368 nullify(wvl%atoms%astruct%sym%phnons) 369 if (associated(wvl%atoms%nlccpar)) then 370 call f_free_ptr(wvl%atoms%nlccpar) 371 end if 372 call deallocate_atoms_data(wvl%atoms) 373 #endif 374 if(allocated(wvl%npspcode_paw_init_guess)) then 375 ABI_DEALLOCATE(wvl%npspcode_paw_init_guess) 376 end if 377 end subroutine wvl_descr_free
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
184 subroutine wvl_descr_psp_fill(gth_params, ipsp, ixc, nelpsp, nzatom, pspunit) 185 186 use defs_datatypes 187 #if defined HAVE_BIGDFT 188 use BigDFT_API, only: atomic_info, UNINITIALIZED, psp_from_data 189 #endif 190 191 !This section has been created automatically by the script Abilint (TD). 192 !Do not modify the following lines by hand. 193 #undef ABI_FUNC 194 #define ABI_FUNC 'wvl_descr_psp_fill' 195 !End of the abilint section 196 197 implicit none 198 199 !Arguments ------------------------------------ 200 integer, intent(in) :: ipsp, pspunit, nzatom, nelpsp, ixc 201 type(pseudopotential_gth_type), intent(inout) :: gth_params 202 !Local variables------------------------------- 203 #if defined HAVE_BIGDFT 204 integer :: ios, ii, nzatom_, nelpsp_, npspcode_, ixc_ 205 real(dp) :: ehomo, radfine 206 logical :: exists 207 character(len = 2) :: symbol 208 character(len=100) :: line 209 character(len=500) :: message 210 #endif 211 212 ! *************************************************************************** 213 214 #if defined HAVE_BIGDFT 215 216 ! check if gth_params%psppar have been set 217 if (any(gth_params%psppar == UNINITIALIZED(1._dp))) then 218 call atomic_info(nzatom, nelpsp, symbol = symbol) 219 ixc_ = ixc 220 if (ixc_>0.and.ixc_< 10) ixc_=1 221 if (ixc_>0.and.ixc_>=10) ixc_=11 222 call psp_from_data(symbol, nzatom_, nelpsp_, npspcode_, ixc_, & 223 & gth_params%psppar(:,:,ipsp), exists) 224 if(.not. exists) then 225 write(message,'(a,a,a,a)')ch10,& 226 & "wvl_descr_psp_fill : bug, chemical element not found in BigDFT table",ch10,& 227 & "Action: upgrade BigDFT table" 228 call wrtout(ab_out,message,'COLL') 229 MSG_BUG(message) 230 end if 231 gth_params%set(ipsp) = .true. 232 end if 233 234 ! Try to read radii from pspunit 235 if (pspunit /= 0) then 236 read (pspunit, '(a100)', iostat = ios) line 237 if (ios /= 0) then 238 line='' 239 end if 240 ! We try to read the values from the pseudo. 241 read (line, *, iostat = ios) gth_params%radii_cf(ipsp, 1), gth_params%radii_cf(ipsp, 2), & 242 & gth_params%radii_cf(ipsp, 3) 243 if (ios /= 0 .or. gth_params%radii_cf(ipsp, 3) < zero) then 244 read (line, *, iostat = ios) gth_params%radii_cf(ipsp, 1), gth_params%radii_cf(ipsp, 2) 245 gth_params%radii_cf(ipsp, 3) = UNINITIALIZED(gth_params%radii_cf(ipsp, 3)) 246 end if 247 if (ios /= 0 .or. gth_params%radii_cf(ipsp, 1) < zero .or. gth_params%radii_cf(ipsp, 2) < zero) then 248 gth_params%radii_cf(ipsp, 1) = UNINITIALIZED(gth_params%radii_cf(ipsp, 1)) 249 gth_params%radii_cf(ipsp, 2) = UNINITIALIZED(gth_params%radii_cf(ipsp, 2)) 250 write(message, '(a,a,a,a,a,a,a)' ) '-', ch10,& 251 & '- wvl_descr_psp_fill : COMMENT -',ch10,& 252 & "- the pseudo-potential does not include geometric informations,",ch10,& 253 & '- values have been computed.' 254 call wrtout(ab_out,message,'COLL') 255 call wrtout(std_out, message,'COLL') 256 end if 257 end if 258 259 ! Update radii. 260 if (gth_params%radii_cf(ipsp, 1) == UNINITIALIZED(gth_params%radii_cf(ipsp, 1))) then 261 call atomic_info(nzatom, nelpsp, ehomo = ehomo) 262 !assigning the radii by calculating physical parameters 263 gth_params%radii_cf(ipsp, 1)=1._dp/sqrt(abs(2._dp*ehomo)) 264 end if 265 if (gth_params%radii_cf(ipsp, 2) == UNINITIALIZED(gth_params%radii_cf(ipsp, 2))) then 266 radfine = gth_params%psppar(0, 0, ipsp) 267 do ii = 1, 4, 1 268 if (gth_params%psppar(ii, 0, ipsp) /= zero) then 269 radfine = min(radfine, gth_params%psppar(ii, 0, ipsp)) 270 end if 271 end do 272 gth_params%radii_cf(ipsp, 2)=radfine 273 end if 274 if (gth_params%radii_cf(ipsp, 3) == UNINITIALIZED(gth_params%radii_cf(ipsp, 3))) then 275 gth_params%radii_cf(ipsp, 3)=gth_params%radii_cf(ipsp, 2) 276 end if 277 278 ! Set flag. 279 if (gth_params%radii_cf(ipsp, 1) >= 0.d0 .and. gth_params%radii_cf(ipsp, 2) >= 0.d0) then 280 gth_params%hasGeometry(ipsp) = .true. 281 else 282 gth_params%hasGeometry(ipsp) = .false. 283 end if 284 285 write(message, '(a,f12.7,a,f12.7,a,f12.7)' )& 286 & ' radii_cf(1)=', gth_params%radii_cf(ipsp, 1),& 287 & '; radii_cf(2)=', gth_params%radii_cf(ipsp, 2),& 288 & '; rad_cov=', gth_params%radii_cf(ipsp, 3) 289 call wrtout(ab_out,message,'COLL') 290 call wrtout(std_out, message,'COLL') 291 292 ! Some consistency checks on radii, to be moved earlier, but need to update test refs. 293 ! gth_params%radii_cf(ipsp, 3) = min(gth_params%radii_cf(ipsp, 3), gth_params%radii_cf(ipsp, 2)) 294 295 ! This was one before 296 ! maxrad=zero 297 ! do ii=0,2,1 298 ! if (ii==1) maxrad=zero 299 ! if (gth_params%psppar(ii,0,ipsp)/=zero) maxrad=max(maxrad,gth_params%psppar(ii,0,ipsp)) 300 ! end do 301 ! if (maxrad== zero) then 302 ! gth_params%radii_cf(ipsp,3)=zero 303 ! else 304 ! gth_params%radii_cf(ipsp,3)=max( & 305 !& min(dtset%wvl_crmult*psps%gth_params%radii_cf(ipsp,1), & 306 !& 15._dp*maxrad)/dtset%wvl_frmult,gth_params%radii_cf(ipsp,2)) 307 ! end if 308 309 #else 310 BIGDFT_NOTENABLED_ERROR() 311 if (.false.) write(std_out,*) ipsp,pspunit,nzatom,nelpsp,ixc,gth_params%psppar 312 #endif 313 314 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
78 subroutine wvl_descr_psp_set(filoccup, nsppol, psps, spinat, wvl) 79 80 use defs_datatypes 81 use defs_wvltypes 82 #if defined HAVE_BIGDFT 83 use BigDFT_API, only: aoig_set,UNINITIALIZED,dict_init,dict_free,dictionary, & 84 & operator(//), bigdft_mpi, dict_set 85 use BigDFT_API, only: psp_data_merge_to_dict, psp_dict_fill_all, atomic_info, & 86 & psp_dict_analyse, atomic_data_set_from_dict, merge_input_file_to_dict 87 #endif 88 89 !This section has been created automatically by the script Abilint (TD). 90 !Do not modify the following lines by hand. 91 #undef ABI_FUNC 92 #define ABI_FUNC 'wvl_descr_psp_set' 93 !End of the abilint section 94 95 implicit none 96 97 !Arguments ------------------------------------ 98 !scalars 99 integer, intent(in) :: nsppol 100 type(wvl_internal_type), intent(inout) :: wvl 101 type(pseudopotential_type), intent(in) :: psps 102 character(len = *), intent(in) :: filoccup 103 !arrays 104 real(dp),intent(in) :: spinat(:,:) 105 106 !Local variables------------------------------- 107 #if defined HAVE_BIGDFT 108 integer :: ityp,pspcod 109 logical :: exists 110 real(dp) :: radii_cf(3) 111 character(len=2) :: symbol 112 character(len=27) :: filename 113 type(dictionary), pointer :: dict 114 #endif 115 116 ! ********************************************************************* 117 118 #if defined HAVE_BIGDFT 119 120 !We create the atoms_data structure, the part that is dependent from psp. 121 do ityp=1,size(psps%pspcod) 122 wvl%atoms%psppar(:,:,ityp)= psps%gth_params%psppar(:,:,ityp) 123 wvl%atoms%npspcode(ityp) = merge(psps%pspcod(ityp),7,psps%pspcod(ityp)/=17) 124 wvl%atoms%ixcpsp(ityp) = psps%pspxc(ityp) 125 wvl%atoms%nzatom(ityp) = int(psps%znucltypat(ityp)) 126 wvl%atoms%nelpsp(ityp) = int(psps%ziontypat(ityp)) 127 end do 128 wvl%atoms%donlcc = .false. 129 130 call dict_init(dict) 131 radii_cf = UNINITIALIZED(1._dp) 132 do ityp = 1, wvl%atoms%astruct%ntypes, 1 133 call atomic_info(wvl%atoms%nzatom(ityp), wvl%atoms%nelpsp(ityp), & 134 & symbol = symbol) 135 write(wvl%atoms%astruct%atomnames(ityp), "(A)") symbol 136 filename = 'psppar.' // trim(symbol) 137 pspcod=psps%pspcod(ityp);if (pspcod==17) pspcod=7 ! PAW specificity 138 call psp_data_merge_to_dict(dict // filename, int(psps%znucltypat(ityp)), & 139 & int(psps%ziontypat(ityp)), pspcod, psps%pspxc(ityp), & 140 & psps%gth_params%psppar(0:4,0:6,ityp), radii_cf, & 141 & UNINITIALIZED(1._dp), UNINITIALIZED(1._dp)) 142 call psp_dict_fill_all(dict, trim(symbol), psps%pspxc(ityp)) 143 end do 144 145 call psp_dict_analyse(dict, wvl%atoms) 146 inquire(file = filoccup, exist = exists) 147 if (exists) then 148 call merge_input_file_to_dict(dict//"Atomic occupation",filoccup,bigdft_mpi) 149 end if 150 ! Need to add spinat in the dictionary, later 151 if (spinat(1,1) == zero .and. .false.) then 152 call dict_set(dict // "spinat", spinat(1,1)) 153 end if 154 call atomic_data_set_from_dict(dict, "Atomic occupation", wvl%atoms, nsppol) 155 call dict_free(dict) 156 157 #else 158 BIGDFT_NOTENABLED_ERROR() 159 if (.false.) write(std_out,*) nsppol,wvl%h(1),psps%npsp,filoccup,spinat(1,1) 160 #endif 161 162 end subroutine wvl_descr_psp_set