TABLE OF CONTENTS
- ABINIT/m_abi2big
- ABINIT/wvl_setBoxGeometry
- ABINIT/wvl_setngfft
- m_abi2big/wvl_eigen_abi2big
- m_abi2big/wvl_occ_abi2big
- m_abi2big/wvl_occopt_abi2big
- m_abi2big/wvl_rho_abi2big
- m_abi2big/wvl_rhov_abi2big_1D_1D
- m_abi2big/wvl_rhov_abi2big_1D_2D
- m_abi2big/wvl_rhov_abi2big_1D_4D
- m_abi2big/wvl_rhov_abi2big_2D_1D
- m_abi2big/wvl_rhov_abi2big_2D_2D
- m_abi2big/wvl_rhov_abi2big_2D_4D
- m_abi2big/wvl_rhov_abi2big_gen
- m_abi2big/wvl_vhartr_abi2big
- m_abi2big/wvl_vtrial_abi2big
- m_abi2big/wvl_vxc_abi2big
ABINIT/m_abi2big [ Modules ]
NAME
m_abi2big
FUNCTION
Module to copy objects from ABINIT to BigDFT and viceversa.
COPYRIGHT
Copyright (C) 2012-2024 ABINIT group (TR,DC,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_abi2big 23 24 use defs_basis 25 use defs_wvltypes 26 use m_errors 27 use m_xmpi 28 29 use m_geometry, only : mkrdim, xcart2xred, xred2xcart 30 31 implicit none 32 33 private 34 35 public :: wvl_vtrial_abi2big 36 !to copy vtrial to wvl_den%rhov and viceversa. 37 38 public :: wvl_rho_abi2big 39 !to copy a density from ABINIT to BigDFT or viceversa 40 41 public :: wvl_vhartr_abi2big 42 ! to copy V_hartree from ABINIT to BigDFT and viceversa 43 44 public :: wvl_vxc_abi2big 45 ! to copy Vxc from ABINIT to BigDFT and viceversa 46 47 public :: wvl_occ_abi2big 48 ! to copy occupations from/to ABINIT to/from BigDFT 49 50 public :: wvl_eigen_abi2big 51 ! to copy eigenvalues from/to ABINIT to/from BigDFT 52 53 public :: wvl_occopt_abi2big 54 ! maps occupation method in ABINIT and BigDFT 55 56 public :: wvl_rhov_abi2big 57 !generic routine to copy a density or potential from/to 58 !ABINIT to/from BigDFT. 59 interface wvl_rhov_abi2big 60 module procedure wvl_rhov_abi2big_2D_4D 61 module procedure wvl_rhov_abi2big_1D_4D 62 module procedure wvl_rhov_abi2big_2D_2D 63 module procedure wvl_rhov_abi2big_1D_2D 64 module procedure wvl_rhov_abi2big_2D_1D 65 module procedure wvl_rhov_abi2big_1D_1D 66 end interface wvl_rhov_abi2big 67 68 logical,parameter :: hmem=.false. !high memory 69 !! Set hmem=.false. if memory is limited. It will copy element by element. 70 !! If hmem=.true. all elements are copied at once. 71 72 public :: wvl_setngfft 73 74 public :: wvl_setBoxGeometry 75 76 contains
ABINIT/wvl_setBoxGeometry [ Functions ]
NAME
wvl_setBoxGeometry
FUNCTION
When wavelets are used, the box definition needs to be changed. The box size is recomputed knowing some psp information 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.
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
SOURCE
1218 subroutine wvl_setBoxGeometry(prtvol, radii, rprimd, xred, wvl, wvl_crmult, wvl_frmult) 1219 1220 #if defined HAVE_BIGDFT 1221 use BigDFT_API, only: system_size,nullify_locreg_descriptors 1222 #endif 1223 implicit none 1224 1225 !Arguments ------------------------------------ 1226 !scalars 1227 integer,intent(in) :: prtvol 1228 real(dp), intent(in) :: wvl_crmult, wvl_frmult 1229 type(wvl_internal_type), intent(inout) :: wvl 1230 !arrays 1231 real(dp),intent(in) :: radii(:,:) 1232 real(dp),intent(inout) :: rprimd(3,3),xred(:,:) 1233 1234 !Local variables------------------------------- 1235 #if defined HAVE_BIGDFT 1236 !scalars 1237 integer :: ii 1238 logical,parameter :: OCLconv=.false. 1239 character(len=500) :: message 1240 !arrays 1241 real(dp) :: rprim(3,3),acell(3) 1242 real(dp),allocatable :: xcart(:,:) 1243 #endif 1244 1245 ! ********************************************************************* 1246 1247 #if defined HAVE_BIGDFT 1248 if (prtvol == 0) then 1249 write(message, '(a,a,a,a)' ) ch10,& 1250 & ' wvl_setBoxGeometry : Changing the box for wavelets computation.' 1251 call wrtout(std_out,message,'COLL') 1252 end if 1253 1254 !Store xcart for each atom 1255 ABI_MALLOC(xcart,(3, wvl%atoms%astruct%nat)) 1256 call xred2xcart(wvl%atoms%astruct%nat, rprimd, xcart, xred) 1257 1258 call nullify_locreg_descriptors(wvl%Glr) 1259 call system_size(wvl%atoms, xcart, radii, wvl_crmult, & 1260 & wvl_frmult, wvl%h(1), wvl%h(2), wvl%h(3), OCLconv, wvl%Glr, wvl%shift) 1261 1262 acell(:) = wvl%atoms%astruct%cell_dim(:) 1263 1264 if (prtvol == 0) then 1265 write(message, '(a,3F12.6)' ) & 1266 & ' | acell is now: ', acell 1267 call wrtout(std_out,message,'COLL') 1268 write(message, '(a,2I5,a,a,2I5,a,a,2I5)' ) & 1269 & ' | nfl1, nfu1: ', wvl%Glr%d%nfl1, wvl%Glr%d%nfu1, ch10, & 1270 & ' | nfl2, nfu2: ', wvl%Glr%d%nfl2, wvl%Glr%d%nfu2, ch10, & 1271 & ' | nfl3, nfu3: ', wvl%Glr%d%nfl3, wvl%Glr%d%nfu3 1272 call wrtout(std_out,message,'COLL') 1273 end if 1274 1275 !Change the metric to orthogonal one 1276 rprim(:, :) = real(0., dp) 1277 do ii = 1, 3, 1 1278 rprim(ii,ii) = real(1., dp) 1279 end do 1280 call mkrdim(acell, rprim, rprimd) 1281 1282 !Save shifted atom positions into xred 1283 call xcart2xred(wvl%atoms%astruct%nat, rprimd, xcart, xred) 1284 ABI_FREE(xcart) 1285 1286 if (prtvol == 0) then 1287 write(message, '(a,3I12)' ) & 1288 & ' | box size for datas: ', wvl%Glr%d%n1i, wvl%Glr%d%n2i, wvl%Glr%d%n3i 1289 call wrtout(std_out,message,'COLL') 1290 write(message, '(a,3I12)' ) & 1291 & ' | box size for wavelets:', wvl%Glr%d%n1, wvl%Glr%d%n2, wvl%Glr%d%n3 1292 call wrtout(std_out,message,'COLL') 1293 end if 1294 1295 #else 1296 BIGDFT_NOTENABLED_ERROR() 1297 if (.false.) write(std_out,*) prtvol,wvl_crmult,wvl_frmult,wvl%h(1),& 1298 & radii(1,1),rprimd(1,1),xred(1,1) 1299 #endif 1300 1301 end subroutine wvl_setBoxGeometry
ABINIT/wvl_setngfft [ Functions ]
NAME
wvl_setngfft
FUNCTION
When wavelets are used, the FFT grid is used to store potentials and density. The size of the grid takes into account the two resolution in wavelet description and also the distribution over processor in the parallel case. The FFT grid is not in strict terms an FFT grid but rather a real space grid. Its dimensions are not directly compatible with FFTs. This is not relevant when using the wavelet part of the code and in the Poisson solver the arrays are extended to match FFT dimensions internally. But for other parts of the code, this must be taken into account. see doc/variables/vargs.html#ngfft for details about ngfft
SIDE EFFECTS
mpi_enreg=information about MPI parallelization (description of the density and potentials scatterring is allocated and updated). dtset <type(dataset_type)>=the FFT grid is changed.
SOURCE
1122 subroutine wvl_setngfft(me_wvl, mgfft, nfft, ngfft, nproc_wvl, n1i, n2i, n3i,n3d) 1123 1124 implicit none 1125 1126 !Arguments ------------------------------------ 1127 !scalars 1128 integer, intent(inout) :: mgfft, nfft 1129 integer, intent(in) :: n1i, n2i, n3i,n3d, nproc_wvl, me_wvl 1130 !arrays 1131 integer, intent(inout) :: ngfft(18) 1132 1133 !Local variables------------------------------- 1134 !scalars 1135 #if defined HAVE_BIGDFT 1136 character(len=500) :: message 1137 #endif 1138 1139 ! ************************************************************************* 1140 1141 #if defined HAVE_BIGDFT 1142 write(message, '(a,a,a,a)' ) ch10,& 1143 & ' wvl_setngfft : Changing the FFT grid definition.' 1144 call wrtout(std_out,message,'COLL') 1145 1146 !Change nfft and ngfft 1147 !Now ngfft will use the density definition (since the potential size 1148 !is always smaller than the density one). ???? 1149 ngfft(1) = n1i 1150 ngfft(2) = n2i 1151 ngfft(3) = n3i 1152 1153 nfft = n1i*n2i*n3d 1154 !Set up fft array dimensions ngfft(4,5,6) to avoid cache conflicts 1155 !Code paste from getng() 1156 ngfft(4) = 2 * (ngfft(1) / 2) + 1 1157 ngfft(5) = 2 * (ngfft(2) / 2) + 1 1158 ngfft(6) = ngfft(3) 1159 if (nproc_wvl == 0) then 1160 ngfft(9) = 0 ! paral_fft 1161 ngfft(10) = 1 ! nproc_fft 1162 ngfft(11) = 0 ! me_fft 1163 ngfft(12) = 0 ! n2proc 1164 ngfft(13) = 0 ! n3proc 1165 else 1166 ngfft(9) = 1 ! paral_fft 1167 ngfft(10) = nproc_wvl 1168 ngfft(11) = me_wvl 1169 ngfft(12) = ngfft(2) 1170 ngfft(13) = n3d 1171 end if 1172 1173 write(message, '(a,3I12)' ) & 1174 & ' | ngfft(1:3) is now: ', ngfft(1:3) 1175 call wrtout(std_out,message,'COLL') 1176 write(message, '(a,3I12)' ) & 1177 & ' | ngfft(4:6) is now: ', ngfft(4:6) 1178 call wrtout(std_out,message,'COLL') 1179 1180 !Set mgfft 1181 mgfft= max(ngfft(1), ngfft(2), ngfft(3)) 1182 1183 #else 1184 BIGDFT_NOTENABLED_ERROR() 1185 if (.false.) write(std_out,*) mgfft,nfft,n1i,n2i,n3i,n3d,nproc_wvl,me_wvl,ngfft(1) 1186 #endif 1187 1188 end subroutine wvl_setngfft
m_abi2big/wvl_eigen_abi2big [ Functions ]
[ Top ] [ m_abi2big ] [ Functions ]
NAME
wvl_eigen_abi2big
FUNCTION
Copies eigenvalues in ABINIT to BigDFT objects and viceversa.
INPUTS
opt= 1) copy from ABINIT to BigDFT 2) copy from BigDFT to ABINIT nsppol= number of spin polarization
OUTPUT
SIDE EFFECTS
occ is copied to wfs%ks%orbs%occup, or viceversa, depending on "opt" (see above).
SOURCE
478 subroutine wvl_eigen_abi2big(mband,nkpt,nsppol,eigen,opt,wvl_wfs) 479 480 implicit none 481 482 !Arguments ------------------------------------ 483 integer , intent(in) :: mband,nkpt,nsppol,opt 484 real(dp) , intent(inout) :: eigen(mband*nkpt*nsppol) 485 type(wvl_wf_type), intent(inout) :: wvl_wfs 486 487 !Local variables------------------------------- 488 #if defined HAVE_BIGDFT 489 integer :: ii,norb,norbd,norbu 490 character(len=100) :: message 491 #endif 492 493 ! ************************************************************************* 494 495 DBG_ENTER("COLL") 496 497 #if defined HAVE_BIGDFT 498 !PENDING: I am not sure this will work for nsppol==2 499 !check also the parallel case. 500 501 norbu=wvl_wfs%ks%orbs%norbu 502 norbd=wvl_wfs%ks%orbs%norbd 503 norb =wvl_wfs%ks%orbs%norb 504 if(opt==1) then !ABINIT -> BIGDFT 505 if (nsppol == 1) then 506 wvl_wfs%ks%orbs%eval=eigen 507 else 508 wvl_wfs%ks%orbs%eval(1:norbu)=eigen(1:norbu) 509 wvl_wfs%ks%orbs%eval(norbu + 1:norb)= & 510 & eigen(mband + 1:mband + norbd) 511 end if 512 elseif(opt==2) then !BigDFT -> ABINIT 513 if (nsppol == 1) then 514 do ii=1,norb 515 eigen(ii)=wvl_wfs%ks%orbs%eval(ii) 516 end do 517 else 518 eigen(1:norbu) = wvl_wfs%ks%orbs%eval(1:norbu) 519 eigen(mband + 1:mband + norbd) = & 520 & wvl_wfs%ks%orbs%eval(norbu + 1:norb) 521 end if 522 else 523 message='wvl_eigen_abi2big: wrong option' 524 ABI_BUG(message) 525 end if 526 527 #else 528 BIGDFT_NOTENABLED_ERROR() 529 if (.false.) write(std_out,*) mband,nkpt,nsppol,opt,eigen(1),wvl_wfs%ks 530 #endif 531 532 DBG_EXIT("COLL") 533 534 end subroutine wvl_eigen_abi2big
m_abi2big/wvl_occ_abi2big [ Functions ]
[ Top ] [ m_abi2big ] [ Functions ]
NAME
wvl_occ_abi2big
FUNCTION
Copies occupations in ABINIT to BigDFT objects and viceversa.
INPUTS
opt= 1) copy from ABINIT to BigDFT 2) copy from BigDFT to ABINIT nsppol= number of spin polarization
OUTPUT
SIDE EFFECTS
occ is copied to wfs%ks%orbs%occup, or viceversa, depending on "opt" (see above).
SOURCE
398 subroutine wvl_occ_abi2big(mband,nkpt,nsppol,occ,opt,wvl_wfs) 399 400 implicit none 401 402 !Arguments ------------------------------------ 403 integer , intent(in) :: mband,nkpt,nsppol,opt 404 real(dp) , intent(inout) :: occ(mband*nkpt*nsppol) 405 type(wvl_wf_type), intent(inout) :: wvl_wfs 406 407 !Local variables------------------------------- 408 #if defined HAVE_BIGDFT 409 integer :: norb,norbd,norbu,ii 410 character(len=100) :: message 411 #endif 412 413 ! ************************************************************************* 414 415 DBG_ENTER("COLL") 416 417 #if defined HAVE_BIGDFT 418 !PENDING: I am not sure this will work for nsppol==2 419 !check also the parallel case. 420 421 norbu=wvl_wfs%ks%orbs%norbu 422 norbd=wvl_wfs%ks%orbs%norbd 423 norb =wvl_wfs%ks%orbs%norb 424 if(opt==1) then !ABINIT -> BIGDFT 425 if (nsppol == 1) then 426 do ii=1,norb 427 wvl_wfs%ks%orbs%occup(ii)=occ(ii) 428 end do 429 else 430 wvl_wfs%ks%orbs%occup(1:norbu)=occ(1:norbu) 431 wvl_wfs%ks%orbs%occup(norbu + 1:norb)= & 432 & occ(mband + 1:mband + norbd) 433 end if 434 elseif(opt==2) then !BigDFT -> ABINIT 435 if (nsppol == 1) then 436 do ii=1,norb 437 occ=wvl_wfs%ks%orbs%occup 438 end do 439 else 440 occ(1:norbu) = wvl_wfs%ks%orbs%occup(1:norbu) 441 occ(mband + 1:mband + norbd) = & 442 & wvl_wfs%ks%orbs%occup(norbu + 1:norb) 443 end if 444 else 445 message='wvl_occ_abi2big: wrong option' 446 ABI_BUG(message) 447 end if 448 449 #else 450 BIGDFT_NOTENABLED_ERROR() 451 if (.false.) write(std_out,*) mband,nkpt,nsppol,opt,occ(1),wvl_wfs%ks 452 #endif 453 454 DBG_EXIT("COLL") 455 456 end subroutine wvl_occ_abi2big
m_abi2big/wvl_occopt_abi2big [ Functions ]
[ Top ] [ m_abi2big ] [ Functions ]
NAME
wvl_occopt_abi2big
FUNCTION
Copies occopt in ABINIT to BigDFT objects and viceversa.
INPUTS
opt= 1) copy from ABINIT to BigDFT 2) copy from BigDFT to ABINIT
OUTPUT
SIDE EFFECTS
NOTES
Several smearing schemes do not exists in both codes such as the SMEARING_DIST_ERF in BigDFT.
SOURCE
558 subroutine wvl_occopt_abi2big(occopt_abi,occopt_big,opt) 559 560 #if defined HAVE_BIGDFT 561 use BigDFT_API, only : & 562 & SMEARING_DIST_FERMI, SMEARING_DIST_COLD1, SMEARING_DIST_COLD2,& 563 & SMEARING_DIST_METPX 564 #endif 565 implicit none 566 567 !Arguments ------------------------------------ 568 integer , intent(inout) :: occopt_abi,occopt_big 569 integer , intent(in) :: opt 570 571 !Local variables------------------------------- 572 #if defined HAVE_BIGDFT 573 character(len=500) :: message 574 #endif 575 576 ! ************************************************************************* 577 578 DBG_ENTER("COLL") 579 580 #if defined HAVE_BIGDFT 581 582 if(opt==1) then !ABINIT -> BIGDFT 583 if(occopt_abi==3) then 584 occopt_big=SMEARING_DIST_FERMI 585 elseif(occopt_abi==4) then 586 occopt_big=SMEARING_DIST_COLD1 587 elseif(occopt_abi==5) then 588 occopt_big=SMEARING_DIST_COLD2 589 elseif(occopt_abi==6) then 590 occopt_big=SMEARING_DIST_METPX 591 else 592 write(message,'(4a)') ch10,& 593 & ' wvl_occopt_abi2big: occopt does not have a corresponding option in BigDFT.',ch10,& 594 & ' Action: change the value of occopt to a number between 3 and 6' 595 ABI_ERROR(message) 596 end if 597 elseif(opt==2) then !BigDFT -> ABINIT 598 if(occopt_big==SMEARING_DIST_FERMI) then 599 occopt_abi=3 600 elseif(occopt_big==SMEARING_DIST_COLD1) then 601 occopt_abi=4 602 elseif(occopt_big==SMEARING_DIST_COLD2) then 603 occopt_abi=5 604 elseif(occopt_big==SMEARING_DIST_METPX) then 605 occopt_abi=6 606 else 607 ! One should never get here. 608 write(message,'(4a)') ch10,& 609 & ' wvl_occopt_abi2big: occopt in BigDFT does not have a corresponding option in ABINIT.',ch10,& 610 & ' Action: contact the ABINIT group' 611 ABI_ERROR(message) 612 end if 613 else 614 message='wvl_occopt_abi2big: wrong option' 615 ABI_BUG(message) 616 end if 617 618 #else 619 BIGDFT_NOTENABLED_ERROR() 620 if (.false.) write(std_out,*) occopt_abi,occopt_big,opt 621 #endif 622 623 DBG_EXIT("COLL") 624 625 end subroutine wvl_occopt_abi2big
m_abi2big/wvl_rho_abi2big [ Functions ]
[ Top ] [ m_abi2big ] [ Functions ]
NAME
wvl_rho_abi2big
FUNCTION
Copies the density from ABINIT to BigDFT, or viceversa.
INPUTS
opt= 1) copy from ABINIT to BigDFT 2) copy from BigDFT to ABINIT rhor(nfft,nspden)= trial potential (ABINIT array) wvl_den= density-potential BigDFT object
OUTPUT
SIDE EFFECTS
Density copied from ABINIT to BigDFT or viceversa. It verifies that (or sets) wvl_den%rhov_is= ELECTRONIC_DENSITY.
NOTES
It uses the generic routine wvl_rhov_abi2big.
SOURCE
185 subroutine wvl_rho_abi2big(opt,rhor,wvl_den) 186 187 #if defined HAVE_BIGDFT 188 use BigDFT_API, only : ELECTRONIC_DENSITY 189 #endif 190 implicit none 191 192 !Arguments ------------------------------------ 193 integer , intent(in) :: opt 194 real(dp) , intent(inout) :: rhor(:,:) 195 type(wvl_denspot_type), intent(inout) :: wvl_den 196 197 !Local variables------------------------------- 198 #if defined HAVE_BIGDFT 199 character(len=100) :: message 200 #endif 201 202 ! ************************************************************************* 203 204 DBG_ENTER("COLL") 205 206 #if defined HAVE_BIGDFT 207 ! write(message,'(2a)') ch10,'wvl_rho_abi2big: but why are you copying me :..o(' 208 ! call wrtout(std_out,message,'COLL') 209 210 if(opt==1) then !ABINIT -> BIGDFT 211 212 call wvl_rhov_abi2big(opt,rhor,wvl_den%denspot%rhov) 213 wvl_den%denspot%rhov_is = ELECTRONIC_DENSITY 214 215 elseif(opt==2) then !BigDFT -> ABINIT 216 217 if(wvl_den%denspot%rhov_is .ne. ELECTRONIC_DENSITY) then 218 message='wvl_rho_abi2big: rhov should contain the ELECTRONIC_DENSITY' 219 ABI_BUG(message) 220 end if 221 call wvl_rhov_abi2big(opt,rhor,wvl_den%denspot%rhov) 222 223 else 224 message='wvl_rho_abi2big: wrong option' 225 ABI_BUG(message) 226 end if 227 228 #else 229 BIGDFT_NOTENABLED_ERROR() 230 if (.false.) write(std_out,*) opt,rhor(1,1),wvl_den%symObj 231 #endif 232 233 DBG_EXIT("COLL") 234 235 end subroutine wvl_rho_abi2big
m_abi2big/wvl_rhov_abi2big_1D_1D [ Functions ]
[ Top ] [ m_abi2big ] [ Functions ]
NAME
wvl_rhov_abi2big_1D_1D
FUNCTION
Copies a density/potential from ABINIT to BigDFT or viceversa. Target : ABINIT 2D arrays (without spin), BigDFT 1D arrays (with or without spin)
INPUTS
opt= 1: copy from ABINIT to BigDFT, 2: copy from BigDFT to ABINIT rhov_abi(:) = density/potential array in ABINIT rhov_big(:) = density/potential array in BigDFT [shift] = shift to be applied in rhov_abi array (parallelism)
OUTPUT
SIDE EFFECTS
At output rhov_abi is copied to rhov_abinit, or viceversa
SOURCE
917 subroutine wvl_rhov_abi2big_1D_1D(opt,rhov_abi,rhov_big,shift) 918 919 implicit none 920 921 !Arguments ------------------------------------ 922 integer,intent(in) :: opt 923 integer,intent(in),optional :: shift 924 real(dp) :: rhov_abi(:),rhov_big(:) 925 926 !Local variables------------------------------- 927 #if defined HAVE_BIGDFT 928 integer :: nfft_abi,nfft_big,shift_ 929 #endif 930 931 ! ************************************************************************* 932 933 #if defined HAVE_BIGDFT 934 nfft_abi=size(rhov_abi) ; nfft_big=size(rhov_big) 935 shift_=0;if (present(shift)) shift_=shift 936 call wvl_rhov_abi2big_gen(nfft_abi,nfft_big,1,opt,rhov_abi,rhov_big,shift_) 937 #else 938 BIGDFT_NOTENABLED_ERROR() 939 if (.false. .and. present(shift)) write(std_out,*) opt,rhov_abi(1),rhov_big(1) 940 #endif 941 942 end subroutine wvl_rhov_abi2big_1D_1D
m_abi2big/wvl_rhov_abi2big_1D_2D [ Functions ]
[ Top ] [ m_abi2big ] [ Functions ]
NAME
wvl_rhov_abi2big_1D_2D
FUNCTION
Copies a density/potential from ABINIT to BigDFT or viceversa. Target : ABINIT 1D arrays (without spin), BigDFT 2D arrays (with spin)
INPUTS
opt= 1: copy from ABINIT to BigDFT, 2: copy from BigDFT to ABINIT rhov_abi(:) = density/potential array in ABINIT rhov_big(:,:) = density/potential array in BigDFT [shift] = shift to be applied in rhov_abi array (parallelism)
OUTPUT
SIDE EFFECTS
At output rhov_abi is copied to rhov_abinit, or viceversa
SOURCE
813 subroutine wvl_rhov_abi2big_1D_2D(opt,rhov_abi,rhov_big,shift) 814 815 implicit none 816 817 !Arguments ------------------------------------ 818 integer,intent(in) :: opt 819 integer,intent(in),optional :: shift 820 real(dp) :: rhov_abi(:),rhov_big(:,:) 821 822 !Local variables------------------------------- 823 #if defined HAVE_BIGDFT 824 integer :: nfft_abi,nfft_big,shift_ 825 character(len=100) :: message 826 #endif 827 828 ! ************************************************************************* 829 830 #if defined HAVE_BIGDFT 831 if (size(rhov_big,2)/=1) then 832 message='wvl_rhov_abi2big: ABINIT and BigDFT objects do not have the same nspden' 833 ABI_BUG(message) 834 end if 835 nfft_abi=size(rhov_abi) ; nfft_big=size(rhov_big) 836 shift_=0;if (present(shift)) shift_=shift 837 call wvl_rhov_abi2big_gen(nfft_abi,nfft_big,1,opt,rhov_abi,rhov_big,shift_) 838 #else 839 BIGDFT_NOTENABLED_ERROR() 840 if (.false. .and. present(shift)) write(std_out,*) opt,rhov_abi(1),rhov_big(1,1) 841 #endif 842 843 end subroutine wvl_rhov_abi2big_1D_2D
m_abi2big/wvl_rhov_abi2big_1D_4D [ Functions ]
[ Top ] [ m_abi2big ] [ Functions ]
NAME
wvl_rhov_abi2big_1D_4D
FUNCTION
Copies a density/potential from ABINIT to BigDFT or viceversa. Target : ABINIT 1D arrays (without spin), BigDFT 4D arrays (with spin)
INPUTS
opt= 1: copy from ABINIT to BigDFT, 2: copy from BigDFT to ABINIT rhov_abi(:) = density/potential array in ABINIT rhov_big(:,:,:,:) = density/potential array in BigDFT [shift] = shift to be applied in rhov_abi array (parallelism)
OUTPUT
SIDE EFFECTS
At output rhov_abi is copied to rhov_abinit, or viceversa
SOURCE
704 subroutine wvl_rhov_abi2big_1D_4D(opt,rhov_abi,rhov_big,shift) 705 706 implicit none 707 708 !Arguments ------------------------------------ 709 integer,intent(in) :: opt 710 integer,intent(in),optional :: shift 711 real(dp) :: rhov_abi(:),rhov_big(:,:,:,:) 712 713 !Local variables------------------------------- 714 #if defined HAVE_BIGDFT 715 integer :: nfft_abi,nfft_big,shift_ 716 character(len=100) :: message 717 #endif 718 719 ! ************************************************************************* 720 721 #if defined HAVE_BIGDFT 722 if (size(rhov_big,4)/=1) then 723 message='wvl_rhov_abi2big: ABINIT and BigDFT objects do not have the same nspden' 724 ABI_BUG(message) 725 end if 726 nfft_abi=size(rhov_abi) ; nfft_big=size(rhov_big) 727 shift_=0;if (present(shift)) shift_=shift 728 call wvl_rhov_abi2big_gen(nfft_abi,nfft_big,1,opt,rhov_abi,rhov_big,shift_) 729 #else 730 BIGDFT_NOTENABLED_ERROR() 731 if (.false. .and. present(shift)) write(std_out,*) opt,rhov_abi(1),rhov_big(1,1,1,1) 732 #endif 733 734 end subroutine wvl_rhov_abi2big_1D_4D
m_abi2big/wvl_rhov_abi2big_2D_1D [ Functions ]
[ Top ] [ m_abi2big ] [ Functions ]
NAME
wvl_rhov_abi2big_2D_1D
FUNCTION
Copies a density/potential from ABINIT to BigDFT or viceversa. Target : ABINIT 2D arrays (with spin), BigDFT 1D arrays (with or without spin)
INPUTS
opt= 1: copy from ABINIT to BigDFT, 2: copy from BigDFT to ABINIT rhov_abi(:,:) = density/potential array in ABINIT rhov_big(:) = density/potential array in BigDFT [shift] = shift to be applied in rhov_abi array (parallelism)
OUTPUT
SIDE EFFECTS
At output rhov_abi is copied to rhov_abinit, or viceversa
SOURCE
867 subroutine wvl_rhov_abi2big_2D_1D(opt,rhov_abi,rhov_big,shift) 868 869 implicit none 870 871 !Arguments ------------------------------------ 872 integer,intent(in) :: opt 873 integer,intent(in),optional :: shift 874 real(dp) :: rhov_abi(:,:),rhov_big(:) 875 876 !Local variables------------------------------- 877 #if defined HAVE_BIGDFT 878 integer :: nfft_abi,nfft_big,nspden,shift_ 879 #endif 880 881 ! ************************************************************************* 882 883 #if defined HAVE_BIGDFT 884 nspden=size(rhov_abi,2) 885 nfft_abi=size(rhov_abi)/nspden ; nfft_big=size(rhov_big)/nspden 886 shift_=0;if (present(shift)) shift_=shift 887 call wvl_rhov_abi2big_gen(nfft_abi,nfft_big,nspden,opt,rhov_abi,rhov_big,shift_) 888 #else 889 BIGDFT_NOTENABLED_ERROR() 890 if (.false. .and. present(shift)) write(std_out,*) opt,rhov_abi(1,1),rhov_big(1) 891 #endif 892 893 end subroutine wvl_rhov_abi2big_2D_1D
m_abi2big/wvl_rhov_abi2big_2D_2D [ Functions ]
[ Top ] [ m_abi2big ] [ Functions ]
NAME
wvl_rhov_abi2big_2D_2D
FUNCTION
Copies a density/potential from ABINIT to BigDFT or viceversa. Target : ABINIT 2D arrays (with spin), BigDFT 2D arrays (with spin)
INPUTS
opt= 1: copy from ABINIT to BigDFT, 2: copy from BigDFT to ABINIT rhov_abi(:,:) = density/potential array in ABINIT rhov_big(:,:) = density/potential array in BigDFT [shift] = shift to be applied in rhov_abi array (parallelism)
OUTPUT
SIDE EFFECTS
At output rhov_abi is copied to rhov_abinit, or viceversa
SOURCE
758 subroutine wvl_rhov_abi2big_2D_2D(opt,rhov_abi,rhov_big,shift) 759 760 implicit none 761 762 !Arguments ------------------------------------ 763 integer,intent(in) :: opt 764 integer,intent(in),optional :: shift 765 real(dp) :: rhov_abi(:,:),rhov_big(:,:) 766 767 !Local variables------------------------------- 768 #if defined HAVE_BIGDFT 769 integer :: nfft_abi,nfft_big,nspden,shift_ 770 character(len=100) :: message 771 #endif 772 773 ! ************************************************************************* 774 775 #if defined HAVE_BIGDFT 776 nspden=size(rhov_abi,2) 777 if (size(rhov_big,2)/=nspden) then 778 message='wvl_rhov_abi2big: ABINIT and BigDFT objects do not have the same nspden' 779 ABI_BUG(message) 780 end if 781 nfft_abi=size(rhov_abi)/nspden ; nfft_big=size(rhov_big)/nspden 782 shift_=0;if (present(shift)) shift_=shift 783 call wvl_rhov_abi2big_gen(nfft_abi,nfft_big,nspden,opt,rhov_abi,rhov_big,shift_) 784 #else 785 BIGDFT_NOTENABLED_ERROR() 786 if (.false. .and. present(shift)) write(std_out,*) opt,rhov_abi(1,1),rhov_big(1,1) 787 #endif 788 789 end subroutine wvl_rhov_abi2big_2D_2D
m_abi2big/wvl_rhov_abi2big_2D_4D [ Functions ]
[ Top ] [ m_abi2big ] [ Functions ]
NAME
wvl_rhov_abi2big_2D_4D
FUNCTION
Copies a density/potential from ABINIT to BigDFT or viceversa. Target : ABINIT 2D arrays (with spin), BigDFT 4D arrays (with spin)
INPUTS
opt= 1: copy from ABINIT to BigDFT, 2: copy from BigDFT to ABINIT rhov_abi(:,:) = density/potential array in ABINIT rhov_big(:,:,:,:) = density/potential array in BigDFT [shift] = shift to be applied in rhov_abi array (parallelism)
OUTPUT
SIDE EFFECTS
At output rhov_abi is copied to rhov_abinit, or viceversa
SOURCE
649 subroutine wvl_rhov_abi2big_2D_4D(opt,rhov_abi,rhov_big,shift) 650 651 implicit none 652 653 !Arguments ------------------------------------ 654 integer,intent(in) :: opt 655 integer,intent(in),optional :: shift 656 real(dp) :: rhov_abi(:,:),rhov_big(:,:,:,:) 657 658 !Local variables------------------------------- 659 #if defined HAVE_BIGDFT 660 integer :: nfft_abi,nfft_big,nspden,shift_ 661 character(len=100) :: message 662 #endif 663 664 ! ************************************************************************* 665 666 #if defined HAVE_BIGDFT 667 nspden=size(rhov_abi,2) 668 if (size(rhov_big,4)/=nspden) then 669 message='wvl_rhov_abi2big: ABINIT and BigDFT objects do not have the same nspden' 670 ABI_BUG(message) 671 end if 672 nfft_abi=size(rhov_abi)/nspden ; nfft_big=size(rhov_big)/nspden 673 shift_=0;if (present(shift)) shift_=shift 674 call wvl_rhov_abi2big_gen(nfft_abi,nfft_big,nspden,opt,rhov_abi,rhov_big,shift_) 675 #else 676 BIGDFT_NOTENABLED_ERROR() 677 if (.false. .and. present(shift)) write(std_out,*) opt,rhov_abi(1,1),rhov_big(1,1,1,1) 678 #endif 679 680 end subroutine wvl_rhov_abi2big_2D_4D
m_abi2big/wvl_rhov_abi2big_gen [ Functions ]
[ Top ] [ m_abi2big ] [ Functions ]
NAME
wvl_rhov_abi2big_gen
FUNCTION
Copies a density/potential from ABINIT to BigDFT or viceversa. This is a generic routine to copy objects.
INPUTS
nfft_abi = size of rhov_abi nfft_big = size of rhov_big nspden = number of spin components opt= 1) copy from ABINIT to BigDFT 2) copy from BigDFT to ABINIT rhov_abi(nfft_abi,nspden) = density/potential array in ABINIT rhov_big(nfft_big,nspden) = density/potential array in BigDFT shift = shift to be applied in rhov_abi array (parallelism)
OUTPUT
SIDE EFFECTS
At output rhov_abi is copied to rhov_abinit, or viceversa
NOTES
This routine is duplicated: This option is faster but it requires more memory. Notice that we cannot point the variables since the spin convention is not the same in BigDFT and ABINIT. In ABINIT: index 1 is for the total spin (spin up + spin down) and index 2 is for spin up. In BigDFT: indices 1 and 2 are for spin up and down, respectively.
SOURCE
978 subroutine wvl_rhov_abi2big_gen(nfft_abi,nfft_big,nspden,opt,rhov_abi,rhov_big,shift) 979 980 implicit none 981 982 !Arguments ------------------------------------ 983 integer,intent(in) :: nfft_abi,nfft_big,nspden,opt,shift 984 real(dp) :: rhov_abi(nfft_abi,nspden),rhov_big(nfft_big,nspden) 985 986 !Local variables------------------------------- 987 #if defined HAVE_BIGDFT 988 integer :: ifft,jfft,nfft 989 real(dp) :: tmpUp,tmpDown,tmpTot 990 character(len=100) :: message 991 real(dp),allocatable :: rhoup(:),rhodn(:),rhotot(:) 992 #endif 993 994 ! ************************************************************************* 995 996 DBG_ENTER("COLL") 997 998 #if defined HAVE_BIGDFT 999 !No objects to copy; in BigDFT by default they have size of 1! 1000 if(size(rhov_big)==1.and.size(rhov_abi)==0) return 1001 1002 nfft=nfft_big;if (nfft_big+shift>nfft_abi) nfft=nfft-shift 1003 1004 if (nfft_abi<nfft+shift) then 1005 message='wvl_rhov_abi2big: cannot handle nfft(abi)<nfft(big)+shift case' 1006 ABI_BUG(message) 1007 end if 1008 if(nspden==4) then 1009 message='wvl_rhov_abi2big: nspden=4 not yet supported' 1010 ABI_ERROR(message) 1011 end if 1012 1013 if (hmem.and.nspden==2) then 1014 if (opt==1) then 1015 ABI_MALLOC(rhoup,(nfft)) 1016 ABI_MALLOC(rhodn,(nfft)) 1017 else if (opt==2) then 1018 ABI_MALLOC(rhotot,(nfft)) 1019 end if 1020 end if 1021 1022 if (opt==1) then !ABINIT -> BIGDFT 1023 if (nspden==2) then 1024 if (hmem) then 1025 rhoup(1:nfft)=rhov_abi(shift+1:shift+nfft,2) 1026 rhodn(1:nfft)=rhov_abi(shift+1:shift+nfft,1)-rhoup(1:nfft) 1027 rhov_big(:,1)=rhoup(:) 1028 rhov_big(:,2)=rhodn(:) 1029 else 1030 do ifft=1,nfft 1031 jfft=shift+ifft 1032 ! We change convention for BigDFT 1033 tmpDown=rhov_abi(jfft,1)-rhov_abi(jfft,2) 1034 tmpUp =rhov_abi(jfft,2) 1035 rhov_big(ifft,1)=tmpUp 1036 rhov_big(ifft,2)=tmpDown 1037 end do 1038 end if !hmem 1039 else !nspden==1 1040 if (hmem) then 1041 rhov_big(1:nfft,1)=rhov_abi(shift+1:shift+nfft,1) 1042 else 1043 do ifft=1,nfft 1044 rhov_big(ifft,1)=rhov_abi(shift+ifft,1) 1045 end do 1046 end if!hmem 1047 end if !nspden 1048 1049 else if (opt==2) then !BigDFT -> ABINIT 1050 if (nspden==2) then 1051 if (hmem) then 1052 rhotot(:)=rhov_big(:,1)+rhov_big(:,2) 1053 rhov_abi(shift+1:shift+nfft,1)=rhotot(1:nfft) 1054 rhov_abi(shift+1:shift+nfft,2)=rhov_big(1:nfft,1) 1055 else 1056 do ifft=1,nfft 1057 jfft=shift+ifft 1058 ! We change convention for BigDFT 1059 tmpTot=rhov_big(ifft,1)+rhov_big(ifft,2) 1060 rhov_abi(jfft,1)=tmpTot 1061 rhov_abi(jfft,2)=rhov_big(ifft,1) !Spin Up 1062 end do 1063 end if !hmem 1064 else if (nspden==1) then 1065 if (hmem) then 1066 rhov_abi(shift+1:shift+nfft,1)=rhov_big(1:nfft,1) 1067 else 1068 do ifft=1,nfft 1069 rhov_abi(shift+ifft,1)=rhov_big(ifft,1) 1070 end do 1071 end if !hmem 1072 end if !nspden 1073 1074 else 1075 message='wvl_rhov_abi2big_gen: wrong option' 1076 ABI_BUG(message) 1077 end if 1078 1079 if (hmem.and.nspden==2) then 1080 if (opt==1) then 1081 ABI_FREE(rhoup) 1082 ABI_FREE(rhodn) 1083 else if (opt==2) then 1084 ABI_FREE(rhotot) 1085 end if !opt 1086 end if 1087 1088 #else 1089 BIGDFT_NOTENABLED_ERROR() 1090 if (.false.) write(std_out,*) nfft_abi,nfft_big,nspden,opt,shift,rhov_big(1,1),rhov_abi(1,1) 1091 #endif 1092 1093 DBG_EXIT("COLL") 1094 1095 end subroutine wvl_rhov_abi2big_gen
m_abi2big/wvl_vhartr_abi2big [ Functions ]
[ Top ] [ m_abi2big ] [ Functions ]
NAME
wvl_vhartr_abi2big
FUNCTION
Copies vhartree in ABINIT to BigDFT objects and viceversa.
INPUTS
opt= 1) copy from ABINIT to BigDFT 2) copy from BigDFT to ABINIT vhartr(nfft)= Hartree potential (ABINIT array) wvl_den= density-potential BigDFT object
OUTPUT
SIDE EFFECTS
vhartr is copied to wvl_den, or viceversa, depending on "opt" (see above). It verifies that (or sets) wvl_den%rhov_is = HARTREE_POTENTIAL
NOTES
It uses the generic routine wvl_rhov_abi2big.
SOURCE
262 subroutine wvl_vhartr_abi2big(opt,vhartr,wvl_den) 263 264 #if defined HAVE_BIGDFT 265 use BigDFT_API, only : HARTREE_POTENTIAL 266 #endif 267 implicit none 268 269 !Arguments ------------------------------------ 270 integer , intent(in) :: opt 271 real(dp) , intent(inout) :: vhartr(:) 272 type(wvl_denspot_type), intent(inout) :: wvl_den 273 274 !Local variables------------------------------- 275 #if defined HAVE_BIGDFT 276 integer :: shiftV 277 character(len=100) :: message 278 #endif 279 280 ! ************************************************************************* 281 282 DBG_ENTER("COLL") 283 284 #if defined HAVE_BIGDFT 285 ! write(message,'(2a)') ch10, 'wvl_vhartr_abi2big: but why are you copying me :..o(' 286 ! call wrtout(std_out,message,'COLL') 287 288 shiftV=wvl_den%denspot%dpbox%ndims(1)*wvl_den%denspot%dpbox%ndims(2) & 289 & *wvl_den%denspot%dpbox%i3xcsh 290 291 if(opt==1) then !ABINIT -> BIGDFT 292 293 call wvl_rhov_abi2big(opt,vhartr,wvl_den%denspot%rhov,shift=shiftV) 294 wvl_den%denspot%rhov_is = HARTREE_POTENTIAL 295 296 elseif(opt==2) then !BigDFT -> ABINIT 297 298 if(wvl_den%denspot%rhov_is .ne. HARTREE_POTENTIAL) then 299 message='wvl_vhartr_abi2big: rhov should contain the HARTREE_POTENTIAL' 300 ABI_BUG(message) 301 end if 302 call wvl_rhov_abi2big(opt,vhartr,wvl_den%denspot%rhov,shift=shiftV) 303 304 else 305 message='wvl_vhartr_abi2big: wrong option' 306 ABI_BUG(message) 307 end if 308 309 #else 310 BIGDFT_NOTENABLED_ERROR() 311 if (.false.) write(std_out,*) opt,vhartr(1),wvl_den%symObj 312 #endif 313 314 DBG_EXIT("COLL") 315 316 end subroutine wvl_vhartr_abi2big
m_abi2big/wvl_vtrial_abi2big [ Functions ]
[ Top ] [ m_abi2big ] [ Functions ]
NAME
wvl_vtrial_abi2big
FUNCTION
Copies vtrial in ABINIT to BigDFT objects and viceversa.
INPUTS
opt= 1) copy from ABINIT to BigDFT 2) copy from BigDFT to ABINIT vtrial(nfft,nspden)= trial potential (ABINIT array) wvl_den= density-potential BigDFT object
OUTPUT
SIDE EFFECTS
vtrial is copied to wvl_den, or viceversa, depending on "opt" (see above). It verifies that (or sets) wvl_den%rhov_is = KS_POTENTIAL.
NOTES
It uses the generic routine wvl_rhov_abi2big.
SOURCE
103 subroutine wvl_vtrial_abi2big(opt,vtrial,wvl_den) 104 105 #if defined HAVE_BIGDFT 106 use BigDFT_API, only : KS_POTENTIAL 107 #endif 108 implicit none 109 110 !Arguments ------------------------------------ 111 integer,intent(in) :: opt 112 real(dp),intent(inout) :: vtrial(:,:) 113 type(wvl_denspot_type),intent(inout) :: wvl_den 114 115 !Local variables------------------------------- 116 #if defined HAVE_BIGDFT 117 integer :: shiftV 118 character(len=100) :: message 119 #endif 120 121 ! ************************************************************************* 122 123 DBG_ENTER("COLL") 124 125 #if defined HAVE_BIGDFT 126 ! write(message,'(2a)') ch10,' wvl_vtrial_abi2big: but why are you copying me :..o(' 127 ! call wrtout(std_out,message,'COLL') 128 129 shiftV=wvl_den%denspot%dpbox%ndims(1)*wvl_den%denspot%dpbox%ndims(2) & 130 & *wvl_den%denspot%dpbox%i3xcsh 131 132 if(opt==1) then !ABINIT -> BIGDFT 133 134 call wvl_rhov_abi2big(opt,vtrial,wvl_den%denspot%rhov,shift=shiftV) 135 wvl_den%denspot%rhov_is = KS_POTENTIAL 136 137 elseif(opt==2) then !BigDFT -> ABINIT 138 139 if(wvl_den%denspot%rhov_is .ne. KS_POTENTIAL) then 140 message='wvl_vtrial_abi2big: rhov should contain the KS_POTENTIAL' 141 ABI_BUG(message) 142 end if 143 144 call wvl_rhov_abi2big(opt,vtrial,wvl_den%denspot%rhov,shift=shiftV) 145 146 else 147 message='wvl_vtrial_abi2big: wrong option' 148 ABI_BUG(message) 149 end if 150 151 #else 152 BIGDFT_NOTENABLED_ERROR() 153 if (.false.) write(std_out,*) opt,vtrial(1,1),wvl_den%symObj 154 #endif 155 156 DBG_EXIT("COLL") 157 158 end subroutine wvl_vtrial_abi2big
m_abi2big/wvl_vxc_abi2big [ Functions ]
[ Top ] [ m_abi2big ] [ Functions ]
NAME
wvl_vxc_abi2big
FUNCTION
It copies the Vxc potential from ABINIT to BigDFT or viceversa.
INPUTS
opt= 1) copy from ABINIT to BigDFT 2) copy from BigDFT to ABINIT vxc(nfft,nspden)= trial potential (ABINIT array) wvl_den= density-potential BigDFT object
OUTPUT
SIDE EFFECTS
vxc is copied to wvl_den, or viceversa, depending on "opt" (see above).
NOTES
It uses the generic routine wvl_rhov_abi2big.
SOURCE
342 subroutine wvl_vxc_abi2big(opt,vxc,wvl_den) 343 344 implicit none 345 346 !Arguments ------------------------------------ 347 integer,intent(in) :: opt 348 real(dp),intent(inout) :: vxc(:,:) 349 type(wvl_denspot_type), intent(inout) :: wvl_den 350 351 !Local variables------------------------------- 352 #if defined HAVE_BIGDFT 353 integer :: shiftV 354 #endif 355 356 ! ************************************************************************* 357 358 DBG_ENTER("COLL") 359 360 #if defined HAVE_BIGDFT 361 ! write(message,'(2a)') ch10, 'wvl_vxc_abi2big: but why are you copying me :..o(' 362 ! call wrtout(std_out,message,'COLL') 363 364 shiftV=wvl_den%denspot%dpbox%ndims(1)*wvl_den%denspot%dpbox%ndims(2) & 365 & *wvl_den%denspot%dpbox%i3xcsh 366 367 call wvl_rhov_abi2big(opt,vxc,wvl_den%denspot%v_xc,shift=shiftV) 368 369 #else 370 BIGDFT_NOTENABLED_ERROR() 371 if (.false.) write(std_out,*) opt,vxc(1,1),wvl_den%symObj 372 #endif 373 374 DBG_EXIT("COLL") 375 376 end subroutine wvl_vxc_abi2big