TABLE OF CONTENTS
ABINIT/prep_getghc [ Functions ]
NAME
prep_getghc
FUNCTION
this routine prepares the data to the call of getghc.
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (FBottin,MT,GZ,MD,FDahm) 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
blocksize= size of block for FFT cpopt=flag defining the status of cprjin%cp(:)=<Proj_i|Cnk> scalars (see below, side effects) cwavef(2,npw*my_nspinor*blocksize)=planewave coefficients of wavefunction. gs_hamk <type(gs_hamiltonian_type)>=all data for the hamiltonian at k gvnlc=matrix elements <G|Vnonlocal|C> lambda=factor to be used when computing <G|H-lambda.S|C> - only for sij_opt=-1 Typically lambda is the eigenvalue (or its guess) mpi_enreg=informations about mpi parallelization prtvol=control print volume and debugging output sij_opt= -PAW ONLY- if 0, only matrix elements <G|H|C> have to be computed (S=overlap) if 1, matrix elements <G|S|C> have to be computed in gsc in addition to ghc if -1, matrix elements <G|H-lambda.S|C> have to be computed in ghc (gsc not used)
OUTPUT
gwavef=(2,npw*my_nspinor*blocksize)=matrix elements <G|H|C> (if sij_opt>=0) or <G|H-lambda.S|C> (if sij_opt=-1). swavef=(2,npw*my_nspinor*blocksize)=matrix elements <G|S|C>.
SIDE EFFECTS
====== if gs_hamk%usepaw==1 cwaveprj(natom,my_nspinor*bandpp)= wave functions at k projected with nl projectors
PARENTS
chebfi,lobpcgwf,m_lobpcgwf,mkresi
CHILDREN
dcopy,multithreaded_getghc,prep_index_wavef_bandpp,prep_sort_wavef_spin prep_wavef_sym_do,prep_wavef_sym_undo,timab,xmpi_alltoallv
SOURCE
48 #if defined HAVE_CONFIG_H 49 #include "config.h" 50 #endif 51 52 #include "abi_common.h" 53 54 subroutine prep_getghc(cwavef,gs_hamk,gvnlc,gwavef,swavef,lambda,blocksize,& 55 & mpi_enreg,prtvol,sij_opt,cpopt,cwaveprj,& 56 & already_transposed) ! optional argument 57 58 use defs_basis 59 use defs_abitypes 60 use m_profiling_abi 61 use m_errors 62 use m_xmpi 63 use m_bandfft_kpt, only : bandfft_kpt,bandfft_kpt_get_ikpt 64 65 use m_pawcprj, only : pawcprj_type 66 use m_hamiltonian, only : gs_hamiltonian_type 67 68 !This section has been created automatically by the script Abilint (TD). 69 !Do not modify the following lines by hand. 70 #undef ABI_FUNC 71 #define ABI_FUNC 'prep_getghc' 72 use interfaces_18_timing 73 use interfaces_66_wfs, except_this_one => prep_getghc 74 !End of the abilint section 75 76 implicit none 77 78 !Arguments ------------------------------------ 79 !scalars 80 integer,intent(in) :: blocksize,cpopt,prtvol,sij_opt 81 logical, intent(in),optional :: already_transposed 82 real(dp),intent(in) :: lambda 83 type(gs_hamiltonian_type),intent(inout) :: gs_hamk 84 type(mpi_type),intent(inout) :: mpi_enreg 85 !arrays 86 real(dp),intent(inout) :: cwavef(:,:),gvnlc (:,:),gwavef(:,:),swavef(:,:) 87 type(pawcprj_type), intent(inout) :: cwaveprj(:,:) 88 89 !Local variables------------------------------- 90 !scalars 91 integer,parameter :: tim_getghc=6 92 integer :: bandpp,bandpp_sym,idatarecv0,ier,ikpt_this_proc,iscalc,mcg,my_nspinor 93 integer :: nbval,ndatarecv,ndatarecv_tot,ndatasend_sym,nproc_band,nproc_fft 94 integer :: old_me_g0,spaceComm=0 95 logical :: flag_inv_sym, do_transpose 96 character(len=100) :: msg 97 !arrays 98 integer,allocatable :: index_wavef_band(:),index_wavef_send(:),index_wavef_spband(:) 99 integer,allocatable :: rdisplsloc(:),recvcountsloc(:),sdisplsloc(:),sendcountsloc(:) 100 integer,ABI_CONTIGUOUS pointer :: kg_k_gather_sym(:,:) 101 integer,ABI_CONTIGUOUS pointer :: rdispls(:),rdispls_sym(:) 102 integer,ABI_CONTIGUOUS pointer :: recvcounts(:),recvcounts_sym(:),recvcounts_sym_tot(:) 103 integer,ABI_CONTIGUOUS pointer :: sdispls(:),sdispls_sym(:) 104 integer,ABI_CONTIGUOUS pointer :: sendcounts(:),sendcounts_sym(:),sendcounts_sym_all(:) 105 integer,ABI_CONTIGUOUS pointer :: tab_proc(:) 106 real(dp) :: tsec(2) 107 real(dp),allocatable,target :: cwavef_alltoall1(:,:),gvnlc_alltoall1(:,:) 108 real(dp),allocatable,target :: gwavef_alltoall1(:,:),swavef_alltoall1(:,:) 109 real(dp),allocatable,target :: cwavef_alltoall2(:,:),gvnlc_alltoall2(:,:) 110 real(dp),allocatable,target :: gwavef_alltoall2(:,:),swavef_alltoall2(:,:) 111 real(dp),pointer :: ewavef_alltoall_sym(:,:),gvnlc_alltoall_sym(:,:) 112 real(dp),pointer :: gwavef_alltoall_sym(:,:) 113 real(dp),pointer :: swavef_alltoall_sym(:,:) 114 115 ! ************************************************************************* 116 117 call timab(630,1,tsec) 118 call timab(631,3,tsec) 119 120 !Some inits 121 nproc_band = mpi_enreg%nproc_band 122 nproc_fft = mpi_enreg%nproc_fft 123 bandpp = mpi_enreg%bandpp 124 my_nspinor = max(1,gs_hamk%nspinor/mpi_enreg%nproc_spinor) 125 126 do_transpose = .true. 127 if(present(already_transposed)) then 128 if(already_transposed) do_transpose = .false. 129 end if 130 131 flag_inv_sym = (gs_hamk%istwf_k==2 .and. any(gs_hamk%ngfft(7) == [401,402,312])) 132 if (flag_inv_sym) then 133 gs_hamk%istwf_k = 1 134 if (modulo(bandpp,2)==0) bandpp_sym = bandpp/2 135 if (modulo(bandpp,2)/=0) bandpp_sym = bandpp 136 end if 137 138 !Check sizes 139 mcg=2*gs_hamk%npw_fft_k*my_nspinor*bandpp 140 if (do_transpose) mcg=2*gs_hamk%npw_k*my_nspinor*blocksize 141 if (size(cwavef)<mcg) then 142 msg='wrong size for cwavef!' 143 MSG_BUG(msg) 144 end if 145 if (size(gwavef)<mcg) then 146 msg='wrong size for gwavef!' 147 MSG_BUG(msg) 148 end if 149 if (size(gvnlc)<mcg) then 150 msg='wrong size for gvnlc!' 151 MSG_BUG(msg) 152 end if 153 if (sij_opt==1) then 154 if (size(swavef)<mcg) then 155 msg='wrong size for swavef!' 156 MSG_BUG(msg) 157 end if 158 end if 159 if (gs_hamk%usepaw==1.and.cpopt>=0) then 160 if (size(cwaveprj)<gs_hamk%natom*my_nspinor*bandpp) then 161 msg='wrong size for cwaveprj!' 162 MSG_BUG(msg) 163 end if 164 end if 165 166 !==================================================================================== 167 168 spaceComm=mpi_enreg%comm_fft 169 if(mpi_enreg%paral_kgb==1) spaceComm=mpi_enreg%comm_band 170 171 ikpt_this_proc=bandfft_kpt_get_ikpt() 172 173 ABI_ALLOCATE(sendcountsloc,(nproc_band)) 174 ABI_ALLOCATE(sdisplsloc ,(nproc_band)) 175 ABI_ALLOCATE(recvcountsloc,(nproc_band)) 176 ABI_ALLOCATE(rdisplsloc ,(nproc_band)) 177 178 recvcounts =>bandfft_kpt(ikpt_this_proc)%recvcounts(:) 179 sendcounts =>bandfft_kpt(ikpt_this_proc)%sendcounts(:) 180 rdispls =>bandfft_kpt(ikpt_this_proc)%rdispls (:) 181 sdispls =>bandfft_kpt(ikpt_this_proc)%sdispls (:) 182 ndatarecv = bandfft_kpt(ikpt_this_proc)%ndatarecv 183 184 if (flag_inv_sym ) then 185 idatarecv0 = bandfft_kpt(ikpt_this_proc)%idatarecv0 186 ndatarecv_tot = bandfft_kpt(ikpt_this_proc)%ndatarecv_tot 187 ndatasend_sym = bandfft_kpt(ikpt_this_proc)%ndatasend_sym 188 kg_k_gather_sym =>bandfft_kpt(ikpt_this_proc)%kg_k_gather_sym(:,:) 189 rdispls_sym =>bandfft_kpt(ikpt_this_proc)%rdispls_sym(:) 190 recvcounts_sym =>bandfft_kpt(ikpt_this_proc)%recvcounts_sym(:) 191 recvcounts_sym_tot =>bandfft_kpt(ikpt_this_proc)%recvcounts_sym_tot(:) 192 sdispls_sym =>bandfft_kpt(ikpt_this_proc)%sdispls_sym(:) 193 sendcounts_sym =>bandfft_kpt(ikpt_this_proc)%sendcounts_sym(:) 194 sendcounts_sym_all =>bandfft_kpt(ikpt_this_proc)%sendcounts_sym_all(:) 195 tab_proc =>bandfft_kpt(ikpt_this_proc)%tab_proc(:) 196 end if 197 iscalc=(sij_opt+1)/2 ! 0 if S not calculated, 1 otherwise 198 nbval=(ndatarecv*my_nspinor*bandpp)*iscalc 199 200 if ( ((.not.flag_inv_sym) .and. bandpp==1 .and. mpi_enreg%paral_spinor==0 .and. my_nspinor==2 ).or. & 201 & ((.not.flag_inv_sym) .and. bandpp>1) .or. flag_inv_sym ) then 202 ABI_ALLOCATE(cwavef_alltoall1,(2,ndatarecv*my_nspinor*bandpp)) 203 ABI_ALLOCATE(gwavef_alltoall1,(2,ndatarecv*my_nspinor*bandpp)) 204 ABI_ALLOCATE(swavef_alltoall1,(2,ndatarecv*my_nspinor*bandpp)) 205 ABI_ALLOCATE(gvnlc_alltoall1,(2,ndatarecv*my_nspinor*bandpp)) 206 swavef_alltoall1(:,:)=zero 207 gvnlc_alltoall1(:,:)=zero 208 cwavef_alltoall1(:,:)=zero 209 gwavef_alltoall1(:,:)=zero 210 end if 211 ABI_ALLOCATE(cwavef_alltoall2,(2,ndatarecv*my_nspinor*bandpp)) 212 ABI_ALLOCATE(gwavef_alltoall2,(2,ndatarecv*my_nspinor*bandpp)) 213 ABI_ALLOCATE(swavef_alltoall2,(2,ndatarecv*my_nspinor*bandpp)) 214 ABI_ALLOCATE(gvnlc_alltoall2,(2,ndatarecv*my_nspinor*bandpp)) 215 swavef_alltoall2(:,:)=zero 216 gvnlc_alltoall2(:,:)=zero 217 cwavef_alltoall2(:,:)=zero 218 gwavef_alltoall2(:,:)=zero 219 220 recvcountsloc(:)=recvcounts(:)*2*my_nspinor*bandpp 221 rdisplsloc(:)=rdispls(:)*2*my_nspinor*bandpp 222 sendcountsloc(:)=sendcounts(:)*2*my_nspinor 223 sdisplsloc(:)=sdispls(:)*2*my_nspinor 224 call timab(631,2,tsec) 225 226 if(do_transpose) then 227 call timab(545,3,tsec) 228 if ( ((.not.flag_inv_sym) .and. bandpp==1 .and. mpi_enreg%paral_spinor==0 .and. my_nspinor==2 ).or. & 229 & ((.not.flag_inv_sym) .and. bandpp>1) .or. flag_inv_sym ) then 230 call xmpi_alltoallv(cwavef,sendcountsloc,sdisplsloc,cwavef_alltoall1,& 231 & recvcountsloc,rdisplsloc,spaceComm,ier) 232 else 233 call xmpi_alltoallv(cwavef,sendcountsloc,sdisplsloc,cwavef_alltoall2,& 234 & recvcountsloc,rdisplsloc,spaceComm,ier) 235 end if 236 call timab(545,2,tsec) 237 else 238 ! Here, we cheat, and use DCOPY to bypass some compiler's overzealous bound-checking 239 ! (ndatarecv*my_nspinor*bandpp might be greater than the declared size of cwavef) 240 call DCOPY(2*ndatarecv*my_nspinor*bandpp, cwavef, 1, cwavef_alltoall2, 1) 241 end if 242 243 if(gs_hamk%istwf_k==2) then 244 old_me_g0=mpi_enreg%me_g0 245 if (mpi_enreg%me_fft==0) then 246 mpi_enreg%me_g0=1 247 else 248 mpi_enreg%me_g0=0 249 end if 250 end if 251 252 !==================================================================== 253 if ((.not.(flag_inv_sym)) .and. (bandpp==1)) then 254 if (do_transpose .and. mpi_enreg%paral_spinor==0.and.my_nspinor==2)then 255 call timab(632,3,tsec) 256 ! Sort to have all ispinor=1 first, then all ispinor=2 257 call prep_sort_wavef_spin(nproc_band,my_nspinor,ndatarecv,recvcounts,rdispls,index_wavef_spband) 258 cwavef_alltoall2(:,:)=cwavef_alltoall1(:,index_wavef_spband) 259 call timab(632,2,tsec) 260 end if 261 262 call timab(635,3,tsec) 263 call multithreaded_getghc(cpopt,cwavef_alltoall2,cwaveprj,gwavef_alltoall2,swavef_alltoall2(:,1:nbval),& 264 & gs_hamk,gvnlc_alltoall2,lambda,mpi_enreg,1,prtvol,sij_opt,tim_getghc,0) 265 call timab(635,2,tsec) 266 267 if (do_transpose .and. mpi_enreg%paral_spinor==0.and.my_nspinor==2)then 268 call timab(634,3,tsec) 269 gwavef_alltoall1(:,index_wavef_spband)=gwavef_alltoall2(:,:) 270 if (sij_opt==1) swavef_alltoall1(:,index_wavef_spband)=swavef_alltoall2(:,:) 271 gvnlc_alltoall1(:,index_wavef_spband)=gvnlc_alltoall2(:,:) 272 ABI_DEALLOCATE(index_wavef_spband) 273 call timab(634,2,tsec) 274 end if 275 276 else if ((.not.(flag_inv_sym)) .and. (bandpp>1)) then 277 ! ------------------------------------------------------------- 278 ! Computation of the index to class the waves functions below bandpp 279 ! ------------------------------------------------------------- 280 281 if(do_transpose) then 282 call timab(632,3,tsec) 283 call prep_index_wavef_bandpp(nproc_band,bandpp,& 284 & my_nspinor,ndatarecv, recvcounts,rdispls, index_wavef_band) 285 ! ------------------------------------------------------- 286 ! Sorting of the waves functions below bandpp 287 ! ------------------------------------------------------- 288 cwavef_alltoall2(:,:) = cwavef_alltoall1(:,index_wavef_band) 289 call timab(632,2,tsec) 290 end if 291 292 ! ---------------------- 293 ! Fourier transformation 294 ! ---------------------- 295 call timab(636,3,tsec) 296 call multithreaded_getghc(cpopt,cwavef_alltoall2,cwaveprj,gwavef_alltoall2,swavef_alltoall2,gs_hamk,& 297 & gvnlc_alltoall2,lambda,mpi_enreg,bandpp,prtvol,sij_opt,tim_getghc,0) 298 call timab(636,2,tsec) 299 300 ! ----------------------------------------------------- 301 ! Sorting of waves functions below the processors 302 ! ----------------------------------------------------- 303 if(do_transpose) then 304 call timab(634,3,tsec) 305 gwavef_alltoall1(:,index_wavef_band) = gwavef_alltoall2(:,:) 306 if (sij_opt==1) swavef_alltoall1(:,index_wavef_band) = swavef_alltoall2(:,:) 307 gvnlc_alltoall1(:,index_wavef_band) = gvnlc_alltoall2(:,:) 308 ABI_DEALLOCATE(index_wavef_band) 309 call timab(634,2,tsec) 310 end if 311 312 313 else if (flag_inv_sym) then 314 315 ! ------------------------------------------------------------- 316 ! Computation of the index to class the waves functions below bandpp 317 ! ------------------------------------------------------------- 318 if(do_transpose) then 319 call timab(632,3,tsec) 320 call prep_index_wavef_bandpp(nproc_band,bandpp,& 321 & my_nspinor,ndatarecv,& 322 & recvcounts,rdispls,& 323 & index_wavef_band) 324 325 ! ------------------------------------------------------- 326 ! Sorting the wave functions below bandpp 327 ! ------------------------------------------------------- 328 cwavef_alltoall2(:,:) = cwavef_alltoall1(:,index_wavef_band) 329 end if 330 331 ! ------------------------------------------------------------ 332 ! We associate the waves functions by two 333 ! ------------------------------------------------------------ 334 call prep_wavef_sym_do(mpi_enreg,bandpp,my_nspinor,& 335 & ndatarecv,& 336 & ndatarecv_tot,ndatasend_sym,tab_proc,& 337 & cwavef_alltoall2,& 338 & sendcounts_sym,sdispls_sym,& 339 & recvcounts_sym,rdispls_sym,& 340 & ewavef_alltoall_sym,& 341 & index_wavef_send) 342 343 ! ------------------------------------------------------------ 344 ! Allocation 345 ! ------------------------------------------------------------ 346 ABI_ALLOCATE(gwavef_alltoall_sym,(2,ndatarecv_tot*bandpp_sym)) 347 ABI_ALLOCATE(swavef_alltoall_sym,(2,(ndatarecv_tot*bandpp_sym)*iscalc)) 348 ABI_ALLOCATE(gvnlc_alltoall_sym ,(2,ndatarecv_tot*bandpp_sym)) 349 350 gwavef_alltoall_sym(:,:)=zero 351 swavef_alltoall_sym(:,:)=zero 352 gvnlc_alltoall_sym(:,:)=zero 353 354 call timab(632,2,tsec) 355 356 ! ------------------------------------------------------------ 357 ! Fourier calculation 358 ! ------------------------------------------------------------ 359 call timab(637,3,tsec) 360 call multithreaded_getghc(cpopt,ewavef_alltoall_sym,cwaveprj,gwavef_alltoall_sym,swavef_alltoall_sym,gs_hamk,& 361 & gvnlc_alltoall_sym,lambda,mpi_enreg,bandpp_sym,prtvol,sij_opt,tim_getghc,1,& 362 & kg_fft_k=kg_k_gather_sym) 363 call timab(637,2,tsec) 364 365 call timab(633,3,tsec) 366 367 ! ------------------------------------------------------------ 368 ! We dissociate each wave function in two waves functions 369 ! gwavef is classed below of bandpp 370 ! ------------------------------------------------------------ 371 call prep_wavef_sym_undo(mpi_enreg,bandpp,my_nspinor,& 372 & ndatarecv,& 373 & ndatarecv_tot,ndatasend_sym,idatarecv0,& 374 & gwavef_alltoall2,& 375 & sendcounts_sym,sdispls_sym,& 376 & recvcounts_sym,rdispls_sym,& 377 & gwavef_alltoall_sym,& 378 & index_wavef_send) 379 if (sij_opt==1)then 380 call prep_wavef_sym_undo(mpi_enreg,bandpp,my_nspinor,& 381 & ndatarecv,& 382 & ndatarecv_tot,ndatasend_sym,idatarecv0,& 383 & swavef_alltoall2,& 384 & sendcounts_sym,sdispls_sym,& 385 & recvcounts_sym,rdispls_sym,& 386 & swavef_alltoall_sym,& 387 & index_wavef_send) 388 end if 389 call prep_wavef_sym_undo(mpi_enreg,bandpp,my_nspinor,& 390 & ndatarecv,& 391 & ndatarecv_tot,ndatasend_sym,idatarecv0,& 392 & gvnlc_alltoall2,& 393 & sendcounts_sym,sdispls_sym,& 394 & recvcounts_sym,rdispls_sym,& 395 & gvnlc_alltoall_sym,& 396 & index_wavef_send) 397 398 ABI_DEALLOCATE(ewavef_alltoall_sym) 399 ABI_DEALLOCATE(index_wavef_send) 400 ABI_DEALLOCATE(gwavef_alltoall_sym) 401 ABI_DEALLOCATE(swavef_alltoall_sym) 402 ABI_DEALLOCATE(gvnlc_alltoall_sym) 403 404 ! ------------------------------------------- 405 ! We call getghc to calculate the nl matrix elements. 406 ! -------------------------------------------- 407 gs_hamk%istwf_k=2 408 !!write(std_out,*)"Setting iswfk_k to 2" 409 410 old_me_g0=mpi_enreg%me_g0 411 if (mpi_enreg%me_fft==0) then 412 mpi_enreg%me_g0=1 413 else 414 mpi_enreg%me_g0=0 415 end if 416 call timab(633,2,tsec) 417 418 call timab(638,3,tsec) 419 call multithreaded_getghc(cpopt,cwavef_alltoall2,cwaveprj,gwavef_alltoall2,swavef_alltoall2,gs_hamk,& 420 & gvnlc_alltoall2,lambda,mpi_enreg,bandpp,prtvol,sij_opt,tim_getghc,2) 421 call timab(638,2,tsec) 422 423 call timab(634,3,tsec) 424 mpi_enreg%me_g0=old_me_g0 425 426 gs_hamk%istwf_k=1 427 428 ! ------------------------------------------------------- 429 ! Sorting the wave functions below the processors 430 ! ------------------------------------------------------- 431 if(do_transpose) then 432 ! cwavef_alltoall(:,index_wavef_band) = cwavef_alltoall(:,:) ! NOT NEEDED 433 gwavef_alltoall1(:,index_wavef_band) = gwavef_alltoall2(:,:) 434 if (sij_opt==1) swavef_alltoall1(:,index_wavef_band) = swavef_alltoall2(:,:) 435 gvnlc_alltoall1(:,index_wavef_band) = gvnlc_alltoall2(:,:) 436 ABI_DEALLOCATE(index_wavef_band) 437 call timab(634,2,tsec) 438 end if 439 440 end if 441 !==================================================================== 442 443 if (gs_hamk%istwf_k==2) mpi_enreg%me_g0=old_me_g0 444 445 if(do_transpose) then 446 447 call timab(545,3,tsec) 448 if ( ((.not.flag_inv_sym) .and. bandpp==1 .and. mpi_enreg%paral_spinor==0 .and. my_nspinor==2 ).or. & 449 & ((.not.flag_inv_sym) .and. bandpp>1) .or. flag_inv_sym ) then 450 if (sij_opt==1) then 451 call xmpi_alltoallv(swavef_alltoall1,recvcountsloc,rdisplsloc,swavef,& 452 & sendcountsloc,sdisplsloc,spaceComm,ier) 453 end if 454 call xmpi_alltoallv(gvnlc_alltoall1,recvcountsloc,rdisplsloc,gvnlc,& 455 & sendcountsloc,sdisplsloc,spaceComm,ier) 456 call xmpi_alltoallv(gwavef_alltoall1,recvcountsloc,rdisplsloc,gwavef,& 457 & sendcountsloc,sdisplsloc,spaceComm,ier) 458 else 459 if (sij_opt==1) then 460 call xmpi_alltoallv(swavef_alltoall2,recvcountsloc,rdisplsloc,swavef,& 461 & sendcountsloc,sdisplsloc,spaceComm,ier) 462 end if 463 call xmpi_alltoallv(gvnlc_alltoall2,recvcountsloc,rdisplsloc,gvnlc,& 464 & sendcountsloc,sdisplsloc,spaceComm,ier) 465 call xmpi_alltoallv(gwavef_alltoall2,recvcountsloc,rdisplsloc,gwavef,& 466 & sendcountsloc,sdisplsloc,spaceComm,ier) 467 end if 468 469 call timab(545,2,tsec) 470 else 471 if(sij_opt == 1) then 472 call DCOPY(2*ndatarecv*my_nspinor*bandpp, swavef_alltoall2, 1, swavef, 1) 473 end if 474 call DCOPY(2*ndatarecv*my_nspinor*bandpp, gvnlc_alltoall2, 1, gvnlc, 1) 475 call DCOPY(2*ndatarecv*my_nspinor*bandpp, gwavef_alltoall2, 1, gwavef, 1) 476 end if 477 478 !==================================================================== 479 if (flag_inv_sym) then 480 gs_hamk%istwf_k = 2 481 end if 482 !==================================================================== 483 ABI_DEALLOCATE(sendcountsloc) 484 ABI_DEALLOCATE(sdisplsloc) 485 ABI_DEALLOCATE(recvcountsloc) 486 ABI_DEALLOCATE(rdisplsloc) 487 ABI_DEALLOCATE(cwavef_alltoall2) 488 ABI_DEALLOCATE(gwavef_alltoall2) 489 ABI_DEALLOCATE(gvnlc_alltoall2) 490 ABI_DEALLOCATE(swavef_alltoall2) 491 492 if ( ((.not.flag_inv_sym) .and. bandpp==1 .and. mpi_enreg%paral_spinor==0 .and. my_nspinor==2 ).or. & 493 & ((.not.flag_inv_sym) .and. bandpp>1) .or. flag_inv_sym ) then 494 ABI_DEALLOCATE(cwavef_alltoall1) 495 ABI_DEALLOCATE(gwavef_alltoall1) 496 ABI_DEALLOCATE(gvnlc_alltoall1) 497 ABI_DEALLOCATE(swavef_alltoall1) 498 end if 499 500 call timab(630,2,tsec) 501 502 end subroutine prep_getghc