TABLE OF CONTENTS
ABINIT/m_paw2wvl [ Modules ]
NAME
paw2wvl
FUNCTION
COPYRIGHT
Copyright (C) 2011-2024 ABINIT group (T. Rangel, MT) 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
16 #if defined HAVE_CONFIG_H 17 #include "config.h" 18 #endif 19 20 #include "abi_common.h" 21 22 module m_paw2wvl 23 24 use defs_basis 25 use defs_wvltypes 26 use m_abicore 27 use m_errors 28 29 use m_pawtab, only : pawtab_type 30 use m_paw_ij, only : paw_ij_type 31 32 implicit none 33 34 private
ABINIT/paw2wvl [ Functions ]
NAME
paw2wvl
FUNCTION
Points WVL objects to PAW objects
INPUTS
argin(sizein)=description
OUTPUT
argout(sizeout)=description
SIDE EFFECTS
NOTES
SOURCE
65 subroutine paw2wvl(pawtab,proj,wvl) 66 67 implicit none 68 69 !Arguments ------------------------------------ 70 type(pawtab_type),intent(in)::pawtab(:) 71 type(wvl_internal_type), intent(inout) :: wvl 72 type(wvl_projectors_type),intent(inout)::proj 73 74 !Local variables------------------------------- 75 #if defined HAVE_BIGDFT 76 integer :: ib,ig 77 integer :: itypat,jb,ll,lmnmax,lmnsz,nn,ntypat,ll_,nn_ 78 integer ::maxmsz,msz1,ptotgau,max_lmn2_size 79 logical :: test_wvl 80 real(dp) :: a1 81 integer,allocatable :: msz(:) 82 #endif 83 84 !extra variables, use to debug 85 ! integer::nr,unitp,i_shell,ii,ir,ng 86 ! real(dp)::step,rmax 87 ! real(dp),allocatable::r(:), y(:) 88 ! complex::fac,arg 89 ! complex(dp),allocatable::f(:),g(:,:) 90 91 ! ********************************************************************* 92 93 !DEBUG 94 !write (std_out,*) ' paw2wvl : enter' 95 !ENDDEBUG 96 97 #if defined HAVE_BIGDFT 98 99 ntypat=size(pawtab) 100 101 !wvl object must be allocated in pawtab 102 if (ntypat>0) then 103 test_wvl=.true. 104 do itypat=1,ntypat 105 if (pawtab(itypat)%has_wvl==0.or.(.not.associated(pawtab(itypat)%wvl))) test_wvl=.false. 106 end do 107 if (.not.test_wvl) then 108 ABI_BUG('pawtab%wvl must be allocated!') 109 end if 110 end if 111 112 !Find max mesh size 113 ABI_MALLOC(msz,(ntypat)) 114 do itypat=1,ntypat 115 msz(itypat)=pawtab(itypat)%wvl%rholoc%msz 116 end do 117 maxmsz=maxval(msz(1:ntypat)) 118 119 ABI_MALLOC(wvl%rholoc%msz,(ntypat)) 120 ABI_MALLOC(wvl%rholoc%d,(maxmsz,4,ntypat)) 121 ABI_MALLOC(wvl%rholoc%rad,(maxmsz,ntypat)) 122 ABI_MALLOC(wvl%rholoc%radius,(ntypat)) 123 124 do itypat=1,ntypat 125 msz1=pawtab(itypat)%wvl%rholoc%msz 126 msz(itypat)=msz1 127 wvl%rholoc%msz(itypat)=msz1 128 wvl%rholoc%d(1:msz1,:,itypat)=pawtab(itypat)%wvl%rholoc%d(1:msz1,:) 129 wvl%rholoc%rad(1:msz1,itypat)=pawtab(itypat)%wvl%rholoc%rad(1:msz1) 130 wvl%rholoc%radius(itypat)=pawtab(itypat)%rpaw 131 end do 132 ABI_FREE(msz) 133 ! 134 !Now fill projectors type: 135 ! 136 ABI_MALLOC(proj%G,(ntypat)) 137 ! 138 !nullify all for security: 139 !do itypat=1,ntypat 140 !call nullify_gaussian_basis(proj%G(itypat)) 141 !end do 142 proj%G(:)%nat=1 !not used 143 proj%G(:)%ncplx=2 !Complex gaussians 144 ! 145 !Obtain dimensions: 146 !jb=0 147 !ptotgau=0 148 do itypat=1,ntypat 149 ! do ib=1,pawtab(itypat)%basis_size 150 ! jb=jb+1 151 ! end do 152 ! ptotgau=ptotgau+pawtab(itypat)%wvl%ptotgau 153 ptotgau=pawtab(itypat)%wvl%ptotgau 154 proj%G(itypat)%nexpo=ptotgau 155 proj%G(itypat)%nshltot=pawtab(itypat)%basis_size 156 end do 157 !proj%G%nexpo=ptotgau 158 !proj%G%nshltot=jb 159 ! 160 !Allocations 161 do itypat=1,ntypat 162 ABI_MALLOC(proj%G(itypat)%ndoc ,(proj%G(itypat)%nshltot)) 163 ABI_MALLOC(proj%G(itypat)%nam ,(proj%G(itypat)%nshltot)) 164 ABI_MALLOC(proj%G(itypat)%xp ,(proj%G(itypat)%ncplx,proj%G(itypat)%nexpo)) 165 ABI_MALLOC(proj%G(itypat)%psiat,(proj%G(itypat)%ncplx,proj%G(itypat)%nexpo)) 166 end do 167 168 !jb=0 169 do itypat=1,ntypat 170 proj%G(itypat)%ndoc(:)=pawtab(itypat)%wvl%pngau(:) 171 end do 172 ! 173 do itypat=1,ntypat 174 proj%G(itypat)%xp(:,:)=pawtab(itypat)%wvl%parg(:,:) 175 proj%G(itypat)%psiat(:,:)=pawtab(itypat)%wvl%pfac(:,:) 176 end do 177 178 !Change the real part of psiat to the form adopted in BigDFT: 179 !Here we use exp{(a+ib)x^2}, where a is a negative number. 180 !In BigDFT: exp{-0.5(x/c)^2}exp{i(bx)^2}, and c is positive 181 !Hence c=sqrt(0.5/abs(a)) 182 do itypat=1,ntypat 183 do ig=1,proj%G(itypat)%nexpo 184 a1=proj%G(itypat)%xp(1,ig) 185 a1=(sqrt(0.5/abs(a1))) 186 proj%G(itypat)%xp(1,ig)=a1 187 end do 188 end do 189 190 !debug 191 !write(*,*)'paw2wvl 178: erase me set gaussians real equal to hgh (for Li)' 192 !ABI_FREE(proj%G(1)%xp) 193 !ABI_FREE(proj%G(1)%psiat) 194 !ABI_MALLOC(proj%G(1)%psiat,(2,2)) 195 !ABI_MALLOC(proj%G(1)%xp,(2,2)) 196 !proj%G(1)%ndoc(:)=1 !1 gaussian per shell 197 !proj%G(1)%nexpo=2 ! two gaussians in total 198 !proj%G(1)%xp=zero 199 !proj%G(1)%xp(1,1)=0.666375d0 200 !proj%G(1)%xp(1,2)=1.079306d0 201 !proj%G(1)%psiat(1,:)=1.d0 202 !proj%G(1)%psiat(2,:)=zero 203 204 ! 205 !begin debug 206 ! 207 !write(*,*)'paw2wvl, comment me' 208 !and comment out variables 209 !rmax=2.d0 210 !nr=rmax/(0.0001d0) 211 !ABI_MALLOC(r,(nr)) 212 !ABI_MALLOC(f,(nr)) 213 !ABI_MALLOC(y,(nr)) 214 !step=rmax/real(nr-1,dp) 215 !do ir=1,nr 216 !r(ir)=real(ir-1,dp)*step 217 !end do 218 !! 219 !unitp=400 220 !do itypat=1,ntypat 221 !ig=0 222 !do i_shell=1,proj%G(itypat)%nshltot 223 !unitp=unitp+1 224 !f(:)=czero 225 !ng=proj%G(itypat)%ndoc(i_shell) 226 !ABI_MALLOC(g,(nr,ng)) 227 !do ii=1,ng 228 !ig=ig+1 229 !fac=cmplx(proj%G(itypat)%psiat(1,ig),proj%G(itypat)%psiat(2,ig)) 230 !a1=-0.5d0/(proj%G(itypat)%xp(1,ig)**2) 231 !arg=cmplx(a1,proj%G(itypat)%xp(2,ig)) 232 !g(:,ii)=fac*exp(arg*r(:)**2) 233 !f(:)=f(:)+g(:,ii) 234 !end do 235 !do ir=1,nr 236 !write(unitp,'(9999f16.7)')r(ir),real(f(ir)),(real(g(ir,ii)),& 237 !& ii=1,ng) 238 !end do 239 !ABI_FREE(g) 240 !end do 241 !end do 242 !ABI_FREE(r) 243 !ABI_FREE(f) 244 !ABI_FREE(y) 245 !stop 246 !end debug 247 248 249 !jb=0 250 do itypat=1,ntypat 251 jb=0 252 ll_ = 0 253 nn_ = 0 254 do ib=1,pawtab(itypat)%lmn_size 255 ll=pawtab(itypat)%indlmn(1,ib) 256 nn=pawtab(itypat)%indlmn(3,ib) 257 ! write(*,*)ll,pawtab(itypat)%indlmn(2,ib),nn 258 if(ib>1 .and. ll == ll_ .and. nn==nn_) cycle 259 jb=jb+1 260 ! proj%G%nam(jb)= pawtab(itypat)%indlmn(1,ib)+1!l quantum number 261 ! write(*,*)jb,proj%G%nshltot,proj%G%nam(jb) 262 proj%G(itypat)%nam(jb)= pawtab(itypat)%indlmn(1,ib)+1 !l quantum number 263 ! 1 is added due to BigDFT convention 264 ! write(*,*)jb,pawtab(itypat)%indlmn(1:3,ib) 265 ll_=ll 266 nn_=nn 267 end do 268 end do 269 !Nullify remaining objects 270 do itypat=1,ntypat 271 nullify(proj%G(itypat)%rxyz) 272 nullify(proj%G(itypat)%nshell) 273 end do 274 !nullify(proj%G%rxyz) 275 276 !now the index l,m,n objects: 277 lmnmax=maxval(pawtab(:)%lmn_size) 278 ABI_MALLOC(wvl%paw%indlmn,(6,lmnmax,ntypat)) 279 wvl%paw%lmnmax=lmnmax 280 wvl%paw%ntypes=ntypat 281 wvl%paw%indlmn=0 282 do itypat=1,ntypat 283 lmnsz=pawtab(itypat)%lmn_size 284 wvl%paw%indlmn(1:6,1:lmnsz,itypat)=pawtab(itypat)%indlmn(1:6,1:lmnsz) 285 end do 286 287 !allocate and copy sij 288 !max_lmn2_size=max(pawtab(:)%lmn2_size,1) 289 max_lmn2_size=lmnmax*(lmnmax+1)/2 290 ABI_MALLOC(wvl%paw%sij,(max_lmn2_size,ntypat)) 291 !sij is not yet calculated here. 292 !We copy this in gstate after call to pawinit 293 !do itypat=1,ntypat 294 !wvl%paw%sij(1:pawtab(itypat)%lmn2_size,itypat)=pawtab(itypat)%sij(:) 295 !end do 296 297 ABI_MALLOC(wvl%paw%rpaw,(ntypat)) 298 do itypat=1,ntypat 299 wvl%paw%rpaw(itypat)= pawtab(itypat)%rpaw 300 end do 301 302 if(allocated(wvl%npspcode_paw_init_guess)) then 303 ABI_FREE(wvl%npspcode_paw_init_guess) 304 end if 305 ABI_MALLOC(wvl%npspcode_paw_init_guess,(ntypat)) 306 do itypat=1,ntypat 307 wvl%npspcode_paw_init_guess(itypat)=pawtab(itypat)%wvl%npspcode_init_guess 308 end do 309 310 #else 311 if (.false.) write(std_out) pawtab(1)%mesh_size,wvl%h(1),proj%nlpsp 312 #endif 313 314 !DEBUG 315 !write (std_out,*) ' paw2wvl : exit' 316 !stop 317 !ENDDEBUG 318 319 end subroutine paw2wvl
ABINIT/paw2wvl_ij [ Functions ]
NAME
paw2wvl_ij
FUNCTION
FIXME: add description.
INPUTS
argin(sizein)=description
OUTPUT
argout(sizeout)=description
SIDE EFFECTS
NOTES
SOURCE
341 subroutine paw2wvl_ij(option,paw_ij,wvl) 342 343 #if defined HAVE_BIGDFT 344 use BigDFT_API, only : nullify_paw_ij_objects 345 #endif 346 implicit none 347 348 !Arguments ------------------------------------ 349 integer,intent(in)::option 350 type(wvl_internal_type), intent(inout)::wvl 351 type(paw_ij_type),intent(in) :: paw_ij(:) 352 !Local variables------------------------------- 353 #if defined HAVE_BIGDFT 354 integer :: iatom,iaux,my_natom 355 character(len=500) :: message 356 #endif 357 358 ! ************************************************************************* 359 360 DBG_ENTER("COLL") 361 362 #if defined HAVE_BIGDFT 363 my_natom=size(paw_ij) 364 365 !Option==1: allocate and copy 366 if(option==1) then 367 ABI_MALLOC(wvl%paw%paw_ij,(my_natom)) 368 do iatom=1,my_natom 369 call nullify_paw_ij_objects(wvl%paw%paw_ij(iatom)) 370 wvl%paw%paw_ij(iatom)%cplex =paw_ij(iatom)%qphase 371 wvl%paw%paw_ij(iatom)%cplex_dij =paw_ij(iatom)%cplex_dij 372 wvl%paw%paw_ij(iatom)%has_dij =paw_ij(iatom)%has_dij 373 wvl%paw%paw_ij(iatom)%has_dijfr =0 374 wvl%paw%paw_ij(iatom)%has_dijhartree =0 375 wvl%paw%paw_ij(iatom)%has_dijhat =0 376 wvl%paw%paw_ij(iatom)%has_dijso =0 377 wvl%paw%paw_ij(iatom)%has_dijU =0 378 wvl%paw%paw_ij(iatom)%has_dijxc =0 379 wvl%paw%paw_ij(iatom)%has_dijxc_val =0 380 wvl%paw%paw_ij(iatom)%has_exexch_pot =0 381 wvl%paw%paw_ij(iatom)%has_pawu_occ =0 382 wvl%paw%paw_ij(iatom)%lmn_size =paw_ij(iatom)%lmn_size 383 wvl%paw%paw_ij(iatom)%lmn2_size =paw_ij(iatom)%lmn2_size 384 wvl%paw%paw_ij(iatom)%ndij =paw_ij(iatom)%ndij 385 wvl%paw%paw_ij(iatom)%nspden =paw_ij(iatom)%nspden 386 wvl%paw%paw_ij(iatom)%nsppol =paw_ij(iatom)%nsppol 387 if (paw_ij(iatom)%has_dij/=0) then 388 iaux=paw_ij(iatom)%cplex_dij*paw_ij(iatom)%lmn2_size 389 ABI_MALLOC(wvl%paw%paw_ij(iatom)%dij,(iaux,paw_ij(iatom)%ndij)) 390 wvl%paw%paw_ij(iatom)%dij(:,:)=paw_ij(iatom)%dij(:,:) 391 end if 392 end do 393 394 ! Option==2: deallocate 395 elseif(option==2) then 396 do iatom=1,my_natom 397 wvl%paw%paw_ij(iatom)%has_dij=0 398 if (associated(wvl%paw%paw_ij(iatom)%dij)) then 399 ABI_FREE(wvl%paw%paw_ij(iatom)%dij) 400 end if 401 end do 402 ABI_FREE(wvl%paw%paw_ij) 403 404 ! Option==3: only copy 405 elseif(option==3) then 406 do iatom=1,my_natom 407 wvl%paw%paw_ij(iatom)%cplex =paw_ij(iatom)%qphase 408 wvl%paw%paw_ij(iatom)%cplex_dij =paw_ij(iatom)%cplex_dij 409 wvl%paw%paw_ij(iatom)%lmn_size =paw_ij(iatom)%lmn_size 410 wvl%paw%paw_ij(iatom)%lmn2_size =paw_ij(iatom)%lmn2_size 411 wvl%paw%paw_ij(iatom)%ndij =paw_ij(iatom)%ndij 412 wvl%paw%paw_ij(iatom)%nspden =paw_ij(iatom)%nspden 413 wvl%paw%paw_ij(iatom)%nsppol =paw_ij(iatom)%nsppol 414 wvl%paw%paw_ij(iatom)%dij(:,:) =paw_ij(iatom)%dij(:,:) 415 end do 416 417 else 418 message = 'paw2wvl_ij: option should be equal to 1, 2 or 3' 419 ABI_ERROR(message) 420 end if 421 422 #else 423 if (.false.) write(std_out,*) option,wvl%h(1),paw_ij(1)%ndij 424 #endif 425 426 DBG_EXIT("COLL") 427 428 end subroutine paw2wvl_ij
ABINIT/wvl_cprjreorder [ Functions ]
NAME
wvl_cprjreorder
FUNCTION
Change the order of a wvl-cprj datastructure From unsorted cprj to atom-sorted cprj (atm_indx=atindx) From atom-sorted cprj to unsorted cprj (atm_indx=atindx1)
INPUTS
atm_indx(natom)=index table for atoms From unsorted wvl%paw%cprj to atom-sorted wvl%paw%cprj (atm_indx=atindx) From atom-sorted wvl%paw%cprj to unsorted wvl%paw%cprj (atm_indx=atindx1)
OUTPUT
SOURCE
521 subroutine wvl_cprjreorder(wvl,atm_indx) 522 523 #if defined HAVE_BIGDFT 524 use BigDFT_API,only : cprj_objects,cprj_paw_alloc,cprj_clean 525 use dynamic_memory 526 #endif 527 implicit none 528 529 !Arguments ------------------------------------ 530 !scalars 531 !arrays 532 integer,intent(in) :: atm_indx(:) 533 type(wvl_internal_type),intent(inout),target :: wvl 534 535 !Local variables------------------------------- 536 #if defined HAVE_BIGDFT 537 !scalars 538 integer :: iexit,ii,jj,kk,n1atindx,n1cprj,n2cprj,ncpgr 539 character(len=100) :: msg 540 !arrays 541 integer,allocatable :: nlmn(:) 542 type(cprj_objects),pointer :: cprj(:,:) 543 type(cprj_objects),allocatable :: cprj_tmp(:,:) 544 #endif 545 546 ! ************************************************************************* 547 548 DBG_ENTER("COLL") 549 550 #if defined HAVE_BIGDFT 551 cprj => wvl%paw%cprj 552 553 n1cprj=size(cprj,dim=1);n2cprj=size(cprj,dim=2) 554 n1atindx=size(atm_indx,dim=1) 555 if (n1cprj==0.or.n2cprj==0.or.n1atindx<=1) return 556 if (n1cprj/=n1atindx) then 557 msg='wrong sizes!' 558 ABI_BUG(msg) 559 end if 560 561 !Nothing to do when the atoms are already sorted 562 iexit=1;ii=0 563 do while (iexit==1.and.ii<n1atindx) 564 ii=ii+1 565 if (atm_indx(ii)/=ii) iexit=0 566 end do 567 if (iexit==1) return 568 569 ABI_MALLOC(nlmn,(n1cprj)) 570 do ii=1,n1cprj 571 nlmn(ii)=cprj(ii,1)%nlmn 572 end do 573 ncpgr=cprj(1,1)%ncpgr 574 575 ABI_MALLOC(cprj_tmp,(n1cprj,n2cprj)) 576 call cprj_paw_alloc(cprj_tmp,ncpgr,nlmn) 577 do jj=1,n2cprj 578 do ii=1,n1cprj 579 cprj_tmp(ii,jj)%nlmn=nlmn(ii) 580 cprj_tmp(ii,jj)%ncpgr=ncpgr 581 cprj_tmp(ii,jj)%cp(:,:)=cprj(ii,jj)%cp(:,:) 582 if (ncpgr>0) cprj_tmp(ii,jj)%dcp(:,:,:)=cprj(ii,jj)%dcp(:,:,:) 583 end do 584 end do 585 586 call cprj_clean(cprj) 587 588 do jj=1,n2cprj 589 do ii=1,n1cprj 590 kk=atm_indx(ii) 591 cprj(kk,jj)%nlmn=nlmn(ii) 592 cprj(kk,jj)%ncpgr=ncpgr 593 cprj(kk,jj)%cp=f_malloc_ptr((/2,nlmn(ii)/),id='cprj%cp') 594 cprj(kk,jj)%cp(:,:)=cprj_tmp(ii,jj)%cp(:,:) 595 if (ncpgr>0) then 596 cprj(kk,jj)%dcp=f_malloc_ptr((/2,ncpgr,nlmn(ii)/),id='cprj%dcp') 597 cprj(kk,jj)%dcp(:,:,:)=cprj_tmp(kk,jj)%dcp(:,:,:) 598 end if 599 end do 600 end do 601 602 call cprj_clean(cprj_tmp) 603 ABI_FREE(cprj_tmp) 604 ABI_FREE(nlmn) 605 606 #else 607 if (.false.) write(std_out,*) atm_indx(1),wvl%h(1) 608 #endif 609 610 DBG_EXIT("COLL") 611 612 end subroutine wvl_cprjreorder
ABINIT/wvl_paw_free [ Functions ]
NAME
wvl_paw_free
FUNCTION
Frees memory for WVL+PAW implementation
INPUTS
ntypat = number of atom types wvl= wvl type wvl_proj= wvl projector type
OUTPUT
SIDE EFFECTS
NOTES
SOURCE
452 subroutine wvl_paw_free(wvl) 453 454 implicit none 455 456 !Arguments ------------------------------------ 457 type(wvl_internal_type),intent(inout) :: wvl 458 459 !Local variables------------------------------- 460 461 ! ************************************************************************* 462 463 #if defined HAVE_BIGDFT 464 465 !PAW objects 466 if( associated(wvl%paw%spsi)) then 467 ABI_FREE(wvl%paw%spsi) 468 end if 469 if( associated(wvl%paw%indlmn)) then 470 ABI_FREE(wvl%paw%indlmn) 471 end if 472 if( associated(wvl%paw%sij)) then 473 ABI_FREE(wvl%paw%sij) 474 end if 475 if( associated(wvl%paw%rpaw)) then 476 ABI_FREE(wvl%paw%rpaw) 477 end if 478 479 !rholoc 480 if( associated(wvl%rholoc%msz )) then 481 ABI_FREE(wvl%rholoc%msz) 482 end if 483 if( associated(wvl%rholoc%d )) then 484 ABI_FREE(wvl%rholoc%d) 485 end if 486 if( associated(wvl%rholoc%rad)) then 487 ABI_FREE(wvl%rholoc%rad) 488 end if 489 if( associated(wvl%rholoc%radius)) then 490 ABI_FREE(wvl%rholoc%radius) 491 end if 492 493 #else 494 if (.false.) write(std_out,*) wvl%h(1) 495 #endif 496 497 !paw%paw_ij and paw%cprj are allocated and deallocated inside vtorho 498 499 end subroutine wvl_paw_free