TABLE OF CONTENTS
ABINIT/readwf [ Functions ]
NAME
readwf
FUNCTION
This subroutine reads the block of records related to one k point, and one spin-polarization, that contains the wavefunctions (as well as the eigenvalues). The disk file unitwf should have been prepared outside of this routine, in order to read the correct records.
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (DCA,XG,GMR,MVer,MB,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 .
INPUTS
formeig=format of the eigenvalues 0 => vector of eigenvalues 1 => hermitian matrix of eigenvalues icg=shift to be given to the location of the cg array ikpt=index of current k point (only needed for error message) isppol=spin polarization currently treated (only needed for error message) kg_k(3,optkg*npw)=k+g data (only if optkg==1) mband=maximum number of bands (dimension of cg, eigen and occ) mcg=dimension of cg nband=number of bands actually in cg, eigen and occ nband_disk=number of bands on the disk file npw=number of plane waves (on current proc) nspinor=number of spinorial components of the wavefunctions (on current proc) option= 1 for reading cg, eigen and occ, -1 for reading/skipping, -2 for reading cg only 3 for reading the eigenvalues only -4 for reading a file written with 4 optkg= if 1 , read or write kg_k ; if 0, do not care about kg_k wff=struct info for wavefunction
SIDE EFFECTS
Current kpt and spin updated cg(2,npw*nspinor*mband)=planewave coefficients of wavefunctions, eigen((2*mband)**formeig *mband)=array for holding eigenvalues (hartree) occ(mband)=array for holding electronic occupations
NOTES
WARNING : occ is not read in the present status of this routine WARNING : skipping k-blocks is also done in the randac subroutine WARNING : reading the two first records is also done in the rdnpw routine
PARENTS
rwwf
CHILDREN
mpi_bcast,wffreadwrite_mpio,wffwritenpwrec,xderivewrecend xderivewrecinit,xderivewrite,xmpi_sum
SOURCE
231 #if defined HAVE_CONFIG_H 232 #include "config.h" 233 #endif 234 235 #include "abi_common.h" 236 237 subroutine readwf(cg,eigen,formeig,headform,icg,ikpt,isppol,kg_k,mband,mcg,mpi_enreg,& 238 & nband,nband_disk,npw,nspinor,occ,option,optkg,wff) 239 240 use defs_basis 241 use defs_abitypes 242 use m_profiling_abi 243 use m_xmpi 244 use m_errors 245 use m_wffile 246 use m_nctk 247 #ifdef HAVE_MPI2 248 use mpi 249 #endif 250 #ifdef HAVE_NETCDF 251 use netcdf 252 #endif 253 254 !This section has been created automatically by the script Abilint (TD). 255 !Do not modify the following lines by hand. 256 #undef ABI_FUNC 257 #define ABI_FUNC 'readwf' 258 !End of the abilint section 259 260 implicit none 261 262 #if defined HAVE_MPI1 263 include 'mpif.h' 264 #endif 265 266 !Arguments ------------------------------------ 267 integer,intent(in) :: formeig,headform,icg,ikpt,isppol,mband,mcg,nband,npw 268 integer,intent(in) :: nspinor,option,optkg 269 integer,intent(inout) :: nband_disk 270 integer,intent(inout),target :: kg_k(3,optkg*npw) 271 real(dp),intent(inout),target :: cg(2,mcg),eigen((2*mband)**formeig*mband),occ(mband) 272 type(MPI_type),intent(in) :: mpi_enreg 273 type(wffile_type),intent(inout) :: wff 274 275 !Local variables------------------------------- 276 integer :: iband,indxx,ios,ipw,ispinor_index,nband1,ncid_hdr 277 integer :: npw1,npwso,npwso1,npwsotot,npwtot,nrec,nspinor1,nspinortot,unitwf,use_f90 278 integer :: band1,band2,ierr 279 character(len=500) :: msg 280 character(len=fnlen) :: fname 281 integer :: ikpt_this_proc,ispinor 282 integer,allocatable :: ind_cg_mpi_to_seq(:) 283 #ifdef HAVE_NETCDF 284 integer :: kg_varid,eig_varid,occ_varid,cg_varid,h1_varid,mband_varid,ncerr,ii,mband_file 285 integer,allocatable :: gall(:,:) 286 real(dp),allocatable :: cg_all(:,:,:),h1mat(:,:,:) 287 #endif 288 289 ! ********************************************************************* 290 291 !write(std_out,*)"rwwf with ikpt, isppol, option, etsf",ikpt, isppol, option, wff%iomode == IO_MODE_ETSF 292 293 !Check the options 294 if ( ALL(option /= [1,3,-1,-2,-4])) then 295 write(msg,'(a,i0,a)')'The argument option should be -4, -2, -1, 1 or 3. However, option=',option,'.' 296 MSG_BUG(msg) 297 end if 298 299 npwtot=npw; npwso=npw*nspinor 300 npwsotot=npwso 301 nspinortot=min(2,(1+mpi_enreg%paral_spinor)*nspinor) 302 303 unitwf=wff%unwff 304 ncid_hdr=unitwf 305 306 use_f90=0 307 if (wff%iomode==IO_MODE_FORTRAN.or.(wff%iomode==IO_MODE_FORTRAN_MASTER.and.wff%master==wff%me)) use_f90=1 308 309 if (option==1.or.option==-2) then 310 311 ! Compute mapping my_gtable --> sequential gtable. 312 if (any(wff%iomode==[IO_MODE_MPI, IO_MODE_ETSF]).and.nband>0) then 313 call xmpi_sum(npwsotot,wff%spaceComm_mpiio,ios) 314 npwtot=npwsotot/nspinortot 315 ABI_ALLOCATE(ind_cg_mpi_to_seq,(npwso)) 316 if (allocated(mpi_enreg%my_kgtab)) then 317 ikpt_this_proc=mpi_enreg%my_kpttab(ikpt) 318 if ( ikpt_this_proc <= 0 ) then 319 MSG_BUG("rwwf: ikpt_this_proc <= 0") 320 end if 321 do ispinor=1,nspinor 322 ispinor_index=ispinor 323 if (mpi_enreg%nproc_spinor>1) ispinor_index=mpi_enreg%me_spinor + 1 324 ind_cg_mpi_to_seq(1+npw*(ispinor-1):npw*ispinor)=npwtot*(ispinor_index-1) & 325 & + mpi_enreg%my_kgtab(1:npw,ikpt_this_proc) 326 end do 327 else 328 ind_cg_mpi_to_seq(1:npwso) = (/(ipw,ipw=1,npwso)/) 329 end if 330 end if 331 end if 332 333 !--------------------------------------------------------------------------- 334 !Read the first record: npw, nspinor, nband_disk 335 !--------------------------------------------------------------------------- 336 337 if (headform>=40.or.headform==0) then ! headform==0 refers to the current headform 338 call WffReadNpwRec(ios,ikpt,isppol,nband_disk,npw1,nspinor1,wff) 339 npwso1=npw1*nspinor1 340 if(ios/=0)then 341 inquire (unit=unitwf, NAME=fname) 342 write(msg,'(3a,i4,2a,i4,4a)') & 343 & 'Reading option of rwwf. Trying to read',ch10,& 344 & 'the (npw,nspinor,nband) record of a wf file, unit=',unitwf,ch10,& 345 & 'gave iostat=',ios,'. Your file is likely not correct.',ch10,& 346 & 'Action: check your input wf file:',trim(fname) 347 MSG_ERROR(msg) 348 end if 349 350 else 351 ! Old format 352 if(use_f90==1)then 353 read (unitwf,iostat=ios) npwso1,nband_disk 354 else if(wff%iomode==IO_MODE_MPI)then 355 call xderiveRRecInit(wff,ios) 356 call xderiveRead(wff,npwso1,ios) 357 call xderiveRead(wff,nband_disk,ios) 358 call xderiveRRecEnd(wff,ios) 359 end if 360 if(ios/=0)then 361 inquire (unit=unitwf, NAME=fname) 362 write(msg,'(3a,i4,2a,i4,4a)') & 363 & 'Reading option of rwwf. Trying to read',ch10,& 364 & 'the (npw,nband) record of a wf file with headform <40 , unit=',unitwf,ch10,& 365 & 'gave iostat=',ios,'. Your file is likely not correct.',ch10,& 366 & 'Action: check your input wf file:',trim(fname) 367 MSG_ERROR(msg) 368 end if 369 end if ! headform 370 371 if (option==1.or.option==-2) then ! Will read the wavefunction and/or kg data, so check npw and nspinor 372 373 if (headform>=40.or.headform==0) then ! New format. headform==0 refers to the current headform 374 if (npwtot/=npw1) then 375 write(msg,'(3a,i0,a,i0,a)') & 376 & 'Reading option of rwwf. One should have npwtot=npw1',ch10,& 377 & 'However, npwtot= ',npwtot,', and npw1= ',npw1,'.' 378 MSG_BUG(msg) 379 end if 380 if(nspinortot/=nspinor1)then 381 write(msg,'(3a,i0,a,i0,a)') & 382 & 'Reading option of rwwf. One should have nspinor=nspinor1',ch10,& 383 & 'However, nspinortot= ',nspinortot,', and nspinor1= ',nspinor1,'.' 384 MSG_BUG(msg) 385 end if 386 else ! Treat the Old format. 387 if(npwsotot/=npwso1)then 388 write(msg,'(3a,i0,a,i0,a)') & 389 & 'Reading option of rwwf. One should have npwso=npwso1',ch10,& 390 & 'However, npwsotot= ',npwsotot,', and npwso1= ',npwso1,'.' 391 MSG_BUG(msg) 392 end if 393 end if ! headform 394 ! 395 end if ! option==1.or.option==2 396 397 !--------------------------------------------------------------------------- 398 !Read the second record: (k+G) vectors 399 !--------------------------------------------------------------------------- 400 401 if (headform>=40.or.headform==0) then ! headform==0 refers to the current headform 402 403 if ((option==1.or.option==-2.or.option==3).and.optkg/=0 )then 404 405 if(use_f90==1)then 406 read(unitwf,iostat=ios) kg_k(1:3,1:npw) 407 408 else if(wff%iomode==IO_MODE_MPI)then 409 call xderiveRRecInit(wff,ios) 410 if (allocated(mpi_enreg%my_kgtab)) then 411 ikpt_this_proc=mpi_enreg%my_kpttab(ikpt) 412 call xderiveRead(wff,kg_k(1:3,1:npw),3,npw,wff%spaceComm_mpiio,& 413 & mpi_enreg%my_kgtab(1:npw,ikpt_this_proc),ios) 414 else 415 ! MG The call below uses MPI_SCAN but here we want to read the full set of G. 416 call xderiveRead(wff,kg_k(1:3,1:npw),3,npw,wff%spaceComm_mpiio,ios) 417 call xmpi_read_int2d(wff,kg_k(1:3,1:npw),wff%spaceComm_mpiio,xmpio_collective,ios) 418 end if 419 call xderiveRRecEnd(wff,ios) 420 421 #ifdef HAVE_NETCDF 422 else if (wff%iomode == IO_MODE_ETSF) then 423 ! Read reduced_coordinates_of_plane_waves for this k point (npw1 is npw_disk). 424 ! TODO: spinor parallelism 425 NCF_CHECK(nf90_inq_varid(wff%unwff, "reduced_coordinates_of_plane_waves", kg_varid)) 426 if (npw == npw1) then 427 ncerr = nf90_get_var(wff%unwff, kg_varid, kg_k, start=[1,1,ikpt], count=[3,npw1,1]) 428 NCF_CHECK_MSG(ncerr, "getting kg_k") 429 else 430 write(std_out,*)"ETSF Reading distributed kg_k" 431 ABI_MALLOC(gall, (3, npw1)) 432 ncerr = nf90_get_var(wff%unwff, kg_varid, gall, start=[1,1,ikpt], count=[3,npw1,1]) 433 NCF_CHECK_MSG(ncerr, "getting kg_k") 434 do ipw=1,npw 435 kg_k(:,ipw) = gall(:,ind_cg_mpi_to_seq(ipw)) 436 end do 437 ABI_FREE(gall) 438 end if 439 #endif 440 end if 441 442 else ! option 443 call WffReadSkipRec(ios,1,wff) ! Skip the record 444 end if 445 446 if(ios/=0)then 447 inquire (unit=unitwf, NAME=fname) 448 write(msg,'(3a,i4,2a,i4,4a)') & 449 & 'Reading option of rwwf. Trying to read',ch10,& 450 & 'the k+g record of a wf file, unit=',unitwf,ch10,& 451 & 'gave iostat=',ios,'. Your file is likely not correct.',ch10,& 452 & 'Action: check your input wf file:',trim(fname) 453 MSG_ERROR(msg) 454 end if 455 end if ! headform 456 457 !--------------------------------------------------------------------------- 458 !Read the third record: eigenvalues 459 !--------------------------------------------------------------------------- 460 !The reading of occ should be enabled, BUT taking into account 461 !of headform of the disk file : occ was NOT present in the disk files with headform=22 462 463 nband1 = min(nband,nband_disk) 464 465 !===== Case formeig=0: read eigenvalues ===== 466 if (formeig==0) then 467 468 if (option==1.or.option==3.or.option==-4) then 469 if(use_f90==1)then 470 read (unitwf,iostat=ios) eigen(1:nband1) 471 else if(wff%iomode==IO_MODE_MPI)then 472 call xderiveRRecInit(wff,ios) 473 call xderiveRead(wff,eigen,nband1,xmpi_comm_self,ios) 474 call xderiveRRecEnd(wff,ios) 475 end if 476 else 477 call WffReadSkipRec(ios,1,wff) 478 end if ! option 479 480 if(ios/=0)then 481 inquire (unit=unitwf, NAME=fname) 482 write(msg,'(3a,i4,2a,i4,4a)') & 483 & 'Reading option of rwwf. Trying to read',ch10,& 484 & 'an eigenvalue record of a wf file, unit=',unitwf,ch10,& 485 & 'gave iostat=',ios,'. Your file is likely not correct.',ch10,& 486 & 'Action: check your input wf file.',trim(fname) 487 MSG_ERROR(msg) 488 end if 489 490 #ifdef HAVE_NETCDF 491 if (wff%iomode == IO_MODE_ETSF) then 492 ! get eigenvalues and occupations 493 NCF_CHECK(nf90_inq_varid(wff%unwff, "eigenvalues", eig_varid)) 494 ncerr = nf90_get_var(wff%unwff, eig_varid, eigen, start=[1,ikpt,isppol], count=[nband1,1,1]) 495 NCF_CHECK_MSG(ncerr, "getting eig_k") 496 497 NCF_CHECK(nf90_inq_varid(wff%unwff, "occupations", occ_varid)) 498 ncerr = nf90_get_var(wff%unwff, occ_varid, occ, start=[1,ikpt,isppol], count=[nband1,1,1]) 499 NCF_CHECK_MSG(ncerr, "getting occ_k") 500 end if 501 #endif 502 503 ! ===== Case formeig=1: read matrix of eigenvalues ===== 504 ! Will be written later (together with wave-functions) 505 else if(formeig==1)then 506 end if ! formeig 507 508 !--------------------------------------------------------------------------- 509 !Read the wave-function coefficients 510 !--------------------------------------------------------------------------- 511 512 !Select bands 513 nband1=min(nband,nband_disk) 514 if(nband1>0.and.option/=-1)then 515 516 ! ===== Case formeig=0: read only wave-functions ===== 517 if (formeig==0) then 518 519 if (option==1.or.option==-2) then 520 521 if (use_f90==1) then 522 do iband=1,nband1 523 ipw=(iband-1)*npwso+icg 524 read(unitwf,iostat=ios) cg(1:2,ipw+1:ipw+npwso) 525 if (ios/=0) exit 526 end do 527 else if (wff%iomode==IO_MODE_MPI) then 528 call WffReadWrite_mpio(wff,1,cg,mcg,icg,nband1,npwso,npwsotot,ind_cg_mpi_to_seq,ios) 529 end if 530 531 else if (option/=-4) then 532 do iband=1,nband1 533 call WffReadSkipRec(ios,1,wff) ! Skip the record 534 end do 535 end if ! option 536 537 if(ios/=0)then 538 inquire (unit=unitwf, NAME=fname) 539 write(msg,'(3a,i4,2a,i4,4a)') & 540 & 'Reading option of rwwf. Trying to read',ch10,& 541 & 'a RF wf record of a wf file, unit=',unitwf,ch10,& 542 & 'gave iostat=',ios,'. Your file is likely not correct.',ch10,& 543 & 'Action: check your input wf file.',trim(fname) 544 MSG_ERROR(msg) 545 end if 546 547 #ifdef HAVE_NETCDF 548 if (wff%iomode == IO_MODE_ETSF) then 549 ! The coefficients_of_wavefunctions on file have shape [cplex, mpw, nspinor, mband, nkpt, nsppol] 550 NCF_CHECK(nf90_inq_varid(wff%unwff, "coefficients_of_wavefunctions", cg_varid)) 551 if (npw == npw1) then 552 ncerr = nf90_get_var(wff%unwff, cg_varid, cg(:, icg+1:), start=[1,1,1,1,ikpt,isppol], & 553 count=[2,npw,nspinor,nband1,1,1]) 554 NCF_CHECK_MSG(ncerr, "getting cg_k") 555 else 556 write(std_out,*)"ETSF Reading distributed cg" 557 ABI_STAT_MALLOC(cg_all, (2, npw1*nspinor, nband1), ierr) 558 ABI_CHECK(ierr==0, "oom cg_all") 559 ncerr = nf90_get_var(wff%unwff, cg_varid, cg_all, start=[1,1,1,1,ikpt,isppol], & 560 count=[2,npw1,nspinor,nband1,1,1]) 561 NCF_CHECK_MSG(ncerr, "getting cg_k") 562 ii = icg 563 do iband=1,nband1 564 do ipw=1,npw 565 ii = ii + 1 566 cg(:,ii) = cg_all(:,ind_cg_mpi_to_seq(ipw),iband) 567 end do 568 end do 569 ABI_FREE(cg_all) 570 end if 571 end if 572 #endif 573 574 ! ===== Case formeig=1: read eigenvalues, occupations and wave-functions ===== 575 else if(formeig==1)then 576 ! write(std_out,*)"nband1",nband1 577 ! ABI_CHECK(nband1==nband_disk,"nband != nband_disk") 578 579 if (wff%iomode /= IO_MODE_ETSF) then 580 indxx=0 581 do iband=1,nband1 582 583 if (option==1.or.option==3.or.option==-4) then 584 if(use_f90==1)then 585 read (unitwf,iostat=ios) eigen(1+indxx:2*nband1+indxx) 586 else if(wff%iomode==IO_MODE_MPI)then 587 ! Should use an access with a "view" 588 call xderiveRRecInit(wff,ios) 589 call xderiveRead(wff,eigen(1+indxx:2*nband1+indxx),2*nband1,xmpi_comm_self,ios) 590 call xderiveRRecEnd(wff,ios) 591 end if 592 indxx=indxx+2*nband1 593 else 594 call WffReadSkipRec(ios,1,wff) ! Skip the record 595 end if 596 597 if(ios/=0)then 598 inquire (unit=unitwf, NAME=fname) 599 write(msg,'(3a,i4,2a,i4,4a)') & 600 & 'Reading option of rwwf. Trying to read',ch10,& 601 & 'a RF eigenvalue record of a wf file, unit=',unitwf,ch10,& 602 & 'gave iostat=',ios,'. Your file is likely not correct.',ch10,& 603 & 'Action: check your input wf file.',trim(fname) 604 MSG_ERROR(msg) 605 end if 606 607 if(option==1.or.option==-2)then 608 ipw=(iband-1)*npwso+icg 609 if(use_f90==1)then 610 ipw=(iband-1)*npwso+icg 611 read(unitwf,iostat=ios) cg(1:2,ipw+1:ipw+npwso) 612 else if(wff%iomode==IO_MODE_MPI)then 613 ! Should use an access with a "view" 614 call xderiveRRecInit(wff,ios) 615 call xderiveRead(wff,cg(1:2,ipw+1:ipw+npwso),2,npwso,wff%spaceComm_mpiio,ind_cg_mpi_to_seq,ios) 616 call xderiveRRecEnd(wff,ios) 617 end if 618 else if (option/=-4) then 619 call WffReadSkipRec(ios,1,wff) ! Skip the record 620 end if ! option 621 622 if(ios/=0)then 623 inquire (unit=unitwf, NAME=fname) 624 write(msg,'(3a,i4,2a,i4,4a)') & 625 & 'Reading option of rwwf. Trying to read',ch10,& 626 & 'a RF wf record of a wf file, unit=',unitwf,ch10,& 627 & 'gave iostat=',ios,'. Your file is likely not correct.',ch10,& 628 & 'Action: check your input wf file.',trim(fname) 629 MSG_ERROR(msg) 630 end if 631 end do ! iband 632 633 else 634 #ifdef HAVE_NETCDF 635 ! ETSF-IO 636 if (any(option == [1, 3, -4])) then 637 ! Read eigen. Remember that the matrix on file has shape [2, mband, mband, nkpt, nspin] 638 ! whereas the eigen array used in Abinit is packed. Read the full matrix first, then pack data. 639 NCF_CHECK(nf90_inq_dimid(wff%unwff, "max_number_of_states", mband_varid)) 640 NCF_CHECK(nf90_inquire_dimension(wff%unwff, mband_varid, len=mband_file)) 641 h1_varid = nctk_idname(wff%unwff, "h1_matrix_elements") 642 643 ABI_MALLOC(h1mat, (2, mband_file, mband_file)) 644 ncerr = nf90_get_var(wff%unwff, h1_varid, h1mat, start=[1,1,1,ikpt,isppol]) 645 NCF_CHECK_MSG(ncerr, "getting h1_matrix_elements") 646 647 indxx=1 648 do band2=1,nband1 649 do band1=1,nband1 650 eigen(indxx:indxx+1) = h1mat(:,band1,band2) 651 indxx = indxx + 2 652 end do 653 end do 654 ABI_FREE(h1mat) 655 end if 656 657 if (any(option == [1, -2])) then 658 ! Read wavefunctions. 659 ! The coefficients_of_wavefunctions on file have shape [cplex, mpw, nspinor, mband, nkpt, nsppol] 660 ABI_CHECK(npw == npw1, "npw != npw1 not coded") 661 662 cg_varid = nctk_idname(wff%unwff, "coefficients_of_wavefunctions") 663 ncerr = nf90_get_var(wff%unwff, cg_varid, cg(:, icg+1:), start=[1,1,1,1,ikpt,isppol], & 664 count=[2,npw,nspinor,nband1,1,1]) 665 NCF_CHECK_MSG(ncerr, "getting cg_k") 666 end if 667 #endif 668 end if 669 670 end if ! formeig == 1 671 672 end if ! nband >0 673 674 !If fewer than all bands were read wind disk file forward to end of bands for this k point. 675 !Will have to fill the non-filled bands outside of this routine ... 676 if (nband<nband_disk .or. option==-1) then 677 nrec=(formeig+1)*(nband_disk-nband) 678 if(option==-1)nrec=(formeig+1)*nband_disk 679 call WffReadSkipRec(ios,nrec,wff) 680 end if 681 682 !--------------------------------------------------------------------------- 683 ! Free memory 684 !--------------------------------------------------------------------------- 685 if (allocated(ind_cg_mpi_to_seq)) then 686 ABI_DEALLOCATE(ind_cg_mpi_to_seq) 687 end if 688 689 end subroutine readwf
ABINIT/rwwf [ Functions ]
NAME
rwwf
FUNCTION
This subroutine reads (different options) or write (option=2) the block of records related to one k point, and one spin-polarization, that contains the wavefunctions (as well as the eigenvalues and occupations). If called with option -1, the records will be skipped. If called with option -2, only the wavefunctions are read. The disk file unitwf should have been prepared outside of this routine, in order to read or write the correct records.
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (DCA,XG,GMR,MVer,MB,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 .
INPUTS
formeig=format of the eigenvalues 0 => vector of eigenvalues 1 => hermitian matrix of eigenvalues headform=format of the header of the wf file, also governing the k block format in case headform=0, use the default (current) format and headform icg=shift to be given to the location of the cg array ikpt=index of current k point (only needed for error message) isppol=spin polarization currently treated (only needed for error message) mband=maximum number of bands (dimension of cg, eigen and occ) mcg=dimention of cg nband=number of bands actually in cg, eigen and occ (if writing mode : must be larger or equal to nband_disk, only nband_disk bands are written ; if reading mode : can be equal, larger or smaller than nband_disk, but cg, eigen and occ will not be completely filled if nband>nband_disk) nband_disk=number of bands on the disk file npw=number of plane waves nspinor=number of spinorial components of the wavefunctions (on current proc) option= 2 for writing cg, eigen and occ, 1 for reading cg and eigen, -1 for reading/skipping, -2 for reading cg only 3 for reading the eigenvalues only 4 for writing a file containing only eigenvalues and occupations (need to be read with option 4) 5 for writing a file containing only eigenvalues and occupations, that can however be read with option 3) -4 for reading a file written with 4 (different from 3 which reads a normal option 2 file) optkg= if 1 , read or write kg_k ; if 0, do not care about kg_k tim_rwwf=timing code of the calling routine (set to 0 if not attributed) wff=struct info for wavefunction | unitwf=unit number for wavefunction
SIDE EFFECTS
cg(2,npw*nspinor*mband)=planewave coefficients of wavefunctions, input if option=2; output if option=1 or -2 eigen((2*mband)**formeig *mband)=array for holding eigenvalues (hartree) input if option=2 or 4 or 5; output if option=1 kg_k(3,optkg*npw)=k+g data (only if optkg==1) input if option=2; output if option=1 or -2 nband_disk=number of bands on disk input if option=2 or 4 or 5; output in the other cases occ(mband)=array for holding eigenvalues (hartree) input if option=2 or 4 or 5; output if option=1 no meaning if frmeig/=0
NOTES
WARNING : occ is not read in the present status of this routine WARNING : skipping k-blocks is also done in the randac subroutine WARNING : reading the two first records is also done in the rdnpw routine WARNING : writing the two first records is also done in the dfpt_vtowfk routine TODO! if (mpi_enreg%flag_ind_kg_mpi_to_seq==1 ) then Some arguments are contained in the wff datastructure, and should be eliminated. option 3 should be called -3 (reading -> negative option) and others (-1,1) re-shuffled.
PARENTS
WffReadEigK,WffReadSkipK,initwf,m_iowf,m_wfk,newkpt,uderiv
CHILDREN
mpi_bcast,wffreadwrite_mpio,wffwritenpwrec,xderivewrecend xderivewrecinit,xderivewrite,xmpi_sum
SOURCE
85 #if defined HAVE_CONFIG_H 86 #include "config.h" 87 #endif 88 89 #include "abi_common.h" 90 91 subroutine rwwf(cg,eigen,formeig,headform,icg,ikpt,isppol,kg_k,mband,mcg,mpi_enreg,& 92 & nband,nband_disk,npw,nspinor,occ,option,optkg,tim_rwwf,wff) 93 94 use defs_basis 95 use defs_abitypes 96 use m_errors 97 use m_wffile 98 use m_profiling_abi 99 #if defined HAVE_MPI2 && !defined FC_G95 100 use mpi 101 #endif 102 #ifdef HAVE_NETCDF 103 use netcdf 104 #endif 105 106 !This section has been created automatically by the script Abilint (TD). 107 !Do not modify the following lines by hand. 108 #undef ABI_FUNC 109 #define ABI_FUNC 'rwwf' 110 use interfaces_18_timing 111 use interfaces_56_io_mpi, except_this_one => rwwf 112 !End of the abilint section 113 114 implicit none 115 116 #if defined HAVE_MPI1 || (defined HAVE_MPI && defined FC_G95) 117 include 'mpif.h' 118 #endif 119 120 !Arguments ------------------------------------ 121 integer,intent(in) :: formeig,headform,icg,ikpt,isppol,mband,mcg,nband,npw 122 integer,intent(inout) :: nband_disk 123 integer,intent(in) :: nspinor,option,optkg,tim_rwwf 124 integer,intent(inout),target :: kg_k(3,optkg*npw) 125 real(dp),intent(inout),target :: cg(2,mcg),eigen((2*mband)**formeig*mband),occ(mband) 126 type(wffile_type),intent(inout) :: wff 127 type(MPI_type), intent(in) :: mpi_enreg 128 129 !Local variables------------------------------- 130 !scalars 131 character(len=500) :: msg 132 !arrays 133 real(dp) :: tsec(2) 134 135 ! ************************************************************************* 136 137 call timab(270+tim_rwwf,1,tsec) 138 139 !Might check that icg+npw*nband*nspinor is smaller than mcg 140 141 !Check that nband is smaller than mband, if one will not skip the records. 142 if (nband>mband .and. option/=-1)then 143 write(msg,'(a,i0,a,i0,a)')' One should have nband<=mband. However, nband=',nband,', and mband=',mband,'.' 144 MSG_BUG(msg) 145 end if 146 147 !Check that formeig is 0 or 1. 148 if ( ALL(formeig/=(/0,1/)) ) then 149 write(msg,'(a,i0,a)')' The argument formeig should be 0 or 1. However, formeig=',formeig,'.' 150 MSG_BUG(msg) 151 end if 152 153 !Check the value of option 154 if ( ALL( option /= (/1,2,3,4,5,-1,-2,-4/) )) then 155 write(msg,'(a,i0,a)')' The argument option should be between 1 and 5, or -1, -2, -4. However, option=',option,'.' 156 MSG_BUG(msg) 157 end if 158 159 if (option/=2.and.option/=4 .and. option/=5 ) then ! read 160 call readwf(cg,eigen,formeig,headform,icg,ikpt,isppol,kg_k,mband,mcg,mpi_enreg,& 161 & nband,nband_disk,npw,nspinor,occ,option,optkg,wff) 162 else ! write 163 call writewf(cg,eigen,formeig,icg,ikpt,isppol,kg_k,mband,mcg,mpi_enreg,& 164 & nband,nband_disk,npw,nspinor,occ,option,optkg,wff) 165 end if 166 167 call timab(270+tim_rwwf,2,tsec) 168 169 end subroutine rwwf
ABINIT/writewf [ Functions ]
NAME
writewf
FUNCTION
This subroutine writes the block of records related to one k point, and one spin-polarization, that contains the wavefunctions (as well as the eigenvalues and occupations). The disk file unitwf should have been prepared outside of this routine, in order to write the correct records.
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (DCA,XG,GMR,MVer,MB,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 .
INPUTS
cg(2,npw*nspinor*mband)=planewave coefficients of wavefunctions, eigen((2*mband)**formeig *mband)=array for holding eigenvalues (hartree) formeig=format of the eigenvalues 0 => vector of eigenvalues 1 => hermitian matrix of eigenvalues icg=shift to be given to the location of the cg array ikpt=index of current k point (only needed for error message) isppol=spin polarization currently treated (only needed for error message) kg_k(3,optkg*npw)=k+g data (only if optkg==1) mband=maximum number of bands (dimension of cg, eigen and occ) mcg=dimension of cg nband=number of bands actually in cg, eigen and occ (must be larger or equal to nband_disk, only nband_disk bands are written) nband_disk=number of bands on the disk file npw=number of plane waves nspinor=number of spinorial components of the wavefunctions (on current proc) occ(mband)=array for holding electronic occupations option= 2 for writing cg, eigen and occ, 4 for writing a file containing only eigenvalues and occupations optkg= if 1 , read or write kg_k ; if 0, do not care about kg_k wff=struct info for wavefunction
OUTPUT
(none, only writing)
SIDE EFFECTS
NOTES
WARNING : skipping k-blocks is also done in the randac subroutine WARNING : writing the two first records is also done in the dfpt_vtowfk routine
PARENTS
rwwf
CHILDREN
mpi_bcast,wffreadwrite_mpio,wffwritenpwrec,xderivewrecend xderivewrecinit,xderivewrite,xmpi_sum
SOURCE
752 #if defined HAVE_CONFIG_H 753 #include "config.h" 754 #endif 755 756 #include "abi_common.h" 757 758 subroutine writewf(cg,eigen,formeig,icg,ikpt,isppol,kg_k,mband,mcg,mpi_enreg,& 759 & nband,nband_disk,npw,nspinor,occ,option,optkg,wff) 760 761 use defs_basis 762 use defs_abitypes 763 use m_xmpi 764 use m_errors 765 use m_wffile 766 use m_profiling_abi 767 use m_nctk 768 #ifdef HAVE_MPI2 769 use mpi 770 #endif 771 #ifdef HAVE_NETCDF 772 use netcdf 773 #endif 774 775 !This section has been created automatically by the script Abilint (TD). 776 !Do not modify the following lines by hand. 777 #undef ABI_FUNC 778 #define ABI_FUNC 'writewf' 779 !End of the abilint section 780 781 implicit none 782 783 #if defined HAVE_MPI1 784 include 'mpif.h' 785 #endif 786 787 !Arguments ------------------------------------ 788 integer,intent(in) :: formeig,icg,ikpt,isppol,mband,mcg,nband,nband_disk,npw,nspinor,option,optkg 789 integer,intent(in),target :: kg_k(3,optkg*npw) 790 real(dp),intent(in),target :: cg(2,mcg),eigen((2*mband)**formeig*mband),occ(mband) 791 type(MPI_type),intent(in) :: mpi_enreg 792 type(wffile_type),intent(inout) :: wff 793 794 !Local variables------------------------------- 795 integer :: iband,ii,ios,ipw,ispinor_index,nband2,ncid_hdr,npwso,npwsotot,npwtot,nspinortot 796 integer :: unitwf,use_f90 797 character(len=500) :: msg 798 integer :: ikpt_this_proc,ispinor,me_cart_3d 799 integer,allocatable :: ind_cg_mpi_to_seq(:) 800 real(dp),ABI_CONTIGUOUS pointer :: cg_ptr(:,:) 801 #ifdef HAVE_NETCDF 802 integer :: kg_varid,eig_varid,occ_varid,cg_varid,ncerr 803 character(len=nctk_slen) :: kdep 804 #endif 805 806 ! ********************************************************************* 807 808 !Check the options 809 if ( ALL(option /= (/2,4,5/)) ) then 810 write(msg,'(a,i0)')' The argument option should be 2, 4 or 5. However, option=',option 811 MSG_BUG(msg) 812 end if 813 814 if(wff%iomode==IO_MODE_MPI)then 815 if(option==5)then 816 write(msg,'(a,i0)')' With MPI-IO activated, the argument option should be 2 or 4. However, option=',option 817 MSG_BUG(msg) 818 end if 819 end if 820 821 if (wff%iomode==IO_MODE_ETSF) then 822 ! outkss still calls this routine! 823 MSG_WARNING("You should use outwff or ncwrite_cg to write data in netcdf format!") 824 end if 825 826 !Check that nband_disk is not larger than nband (only for writing) 827 if (nband<nband_disk) then 828 write(msg,'(3a,i5,a,i5,a)') & 829 & 'Writing option of rwwf. One should have nband<=nband_disk',ch10,& 830 & 'However, nband= ',nband,', and nband_disk= ',nband_disk,'.' 831 MSG_BUG(msg) 832 end if 833 834 npwtot=npw; npwso=npw*nspinor 835 unitwf=wff%unwff; ncid_hdr=unitwf 836 npwsotot=npwso 837 nspinortot=min(2,(1+mpi_enreg%paral_spinor)*nspinor) 838 839 use_f90=0 840 if (wff%iomode==IO_MODE_FORTRAN.or.(wff%iomode ==IO_MODE_FORTRAN_MASTER.and.wff%master==wff%me)) use_f90=1 841 842 if (wff%iomode==IO_MODE_MPI) then 843 844 call xmpi_sum(npwsotot,wff%spaceComm_mpiio,ios) 845 npwtot=npwsotot/nspinortot 846 847 if (option/=4) then 848 ABI_ALLOCATE(ind_cg_mpi_to_seq,(npwso)) 849 if (allocated(mpi_enreg%my_kgtab)) then 850 ikpt_this_proc=mpi_enreg%my_kpttab(ikpt) 851 do ispinor=1,nspinor 852 ispinor_index=ispinor 853 if (mpi_enreg%nproc_spinor > 1) ispinor_index = mpi_enreg%me_spinor + 1 854 ind_cg_mpi_to_seq(1+npw*(ispinor-1):npw*ispinor)=npwtot*(ispinor_index-1) & 855 & + mpi_enreg%my_kgtab(1:npw,ikpt_this_proc) 856 end do 857 else 858 ind_cg_mpi_to_seq(1:npwso) = (/(ipw,ipw=1,npwso)/) 859 end if 860 !write(std_out,*)"MPI-IO ind_cg_mpi_to_seq", ind_cg_mpi_to_seq(1:5) 861 end if 862 end if 863 864 !--------------------------------------------------------------------------- 865 !Write the first record: npw, nspinor, nband_disk 866 !--------------------------------------------------------------------------- 867 !Not modified for netCDF: no need to add writing of nband_disk,npw,nspinor 868 869 call WffWriteNpwRec(ios,nband_disk,npwtot,nspinortot,wff,opt_paral=2) 870 871 !--------------------------------------------------------------------------- 872 !Write the second record: (k+G) vectors 873 !--------------------------------------------------------------------------- 874 875 if (optkg/=0.and.option/=4) then 876 if(use_f90==1)then 877 write(unitwf) kg_k(1:3,1:optkg*npw) 878 else if (wff%iomode==IO_MODE_MPI) then 879 880 if (allocated(mpi_enreg%my_kgtab)) then 881 me_cart_3d=xmpi_comm_rank(mpi_enreg%comm_bandspinorfft) 882 ikpt_this_proc=mpi_enreg%my_kpttab(ikpt) 883 call xderiveWRecInit(wff,ios,me_cart_3d) 884 if (mpi_enreg%me_spinor==0) then 885 call xderiveWrite(wff,kg_k,3,npw,mpi_enreg%comm_bandfft, & 886 & mpi_enreg%my_kgtab(1:npw,ikpt_this_proc),ios) 887 end if 888 call xderiveWRecEnd(wff,ios,me_cart_3d) 889 else 890 ! MG does it work if we are not using FFT distribution ? 891 call xderiveWRecInit(wff,ios ) 892 if (mpi_enreg%me_spinor==0) then 893 call xderiveWrite(wff,kg_k,3,optkg*npw,Wff%spaceComm_mpiio,ios) 894 end if 895 call xderiveWRecEnd(wff,ios) 896 end if 897 898 #ifdef HAVE_NETCDF 899 else if (wff%iomode == IO_MODE_ETSF) then 900 ! Write the reduced_coordinates_of_plane_waves for this k point. 901 NCF_CHECK(nf90_inq_varid(wff%unwff, "reduced_coordinates_of_plane_waves", kg_varid)) 902 NCF_CHECK(nf90_get_att(wff%unwff, kg_varid, "k_dependent", kdep)) 903 if (kdep == "no") then 904 ncerr = nf90_put_var(wff%unwff, kg_varid, kg_k, start=[1,1], count=[3,npw]) 905 else 906 ncerr = nf90_put_var(wff%unwff, kg_varid, kg_k, start=[1,1,ikpt], count=[3,npw,1]) 907 end if 908 NCF_CHECK_MSG(ncerr, "putting kg_k") 909 #endif 910 end if ! end if wff%iomode 911 else ! Still skip the record 912 if (use_f90==1) then 913 write(unitwf) 914 else if (wff%iomode==IO_MODE_MPI) then 915 call xderiveWRecInit(wff,wff%spaceComm_mpiio,ios) 916 call xderiveWRecEnd(wff,wff%spaceComm_mpiio,ios) 917 end if 918 end if 919 920 !--------------------------------------------------------------------------- 921 !Write the third record: eigenvalues and occupations 922 !--------------------------------------------------------------------------- 923 924 !===== Case formeig=0: write eigenvalues and occupations ===== 925 if (formeig==0) then 926 if (use_f90==1) then 927 write(unitwf) (eigen(iband),iband=1,nband_disk),(occ(iband),iband=1,nband_disk) 928 else if(wff%iomode==IO_MODE_MPI) then 929 if (wff%me_mpiio==0) then 930 call xderiveWRecInit(wff,ios) 931 call xderiveWrite(wff,eigen,nband_disk,xmpi_comm_self,ios) 932 call xderiveWrite(wff,occ,nband_disk,xmpi_comm_self,ios) 933 call xderiveWRecEnd(wff,ios) 934 end if 935 #ifdef HAVE_MPI 936 call MPI_BCAST(wff%offwff,1,wff%offset_mpi_type,0,wff%spaceComm_mpiio,ios) 937 #endif 938 #ifdef HAVE_NETCDF 939 else if (wff%iomode == IO_MODE_ETSF) then 940 ! Write eigenvalues and occupation factors. 941 NCF_CHECK(nf90_inq_varid(wff%unwff, "eigenvalues", eig_varid)) 942 ncerr = nf90_put_var(wff%unwff, eig_varid, eigen, start=[1,ikpt,isppol], count=[mband,1,1]) 943 NCF_CHECK_MSG(ncerr, "putting eig_k") 944 945 NCF_CHECK(nf90_inq_varid(wff%unwff, "occupations", occ_varid)) 946 ncerr = nf90_put_var(wff%unwff, occ_varid, occ, start=[1,ikpt,isppol], count=[mband,1,1]) 947 NCF_CHECK_MSG(ncerr, "putting occ_k") 948 #endif 949 end if 950 951 ! ===== Case formeig=1: write matrix of eigenvalues ===== 952 ! Will be written later (together with wave-functions) 953 else if(formeig==1)then 954 end if ! formeig 955 956 !--------------------------------------------------------------------------- 957 !Write the wave-function coefficients 958 !--------------------------------------------------------------------------- 959 960 !===== Case formeig=0: write only wave-functions ===== 961 if (formeig==0) then 962 ! If option=4, do not write wave functions 963 if (option/=4) then 964 if (use_f90==1) then 965 do iband=1,nband_disk 966 ipw=(iband-1)*npwso+icg 967 if(option/=5)then 968 write(unitwf) cg(1:2,ipw+1:ipw+npwso) ! VALGRIND complains some elements of cg are not initialized, but written 969 else 970 write(unitwf) 971 end if 972 end do 973 else if(wff%iomode==IO_MODE_MPI)then 974 cg_ptr => cg ! Need pointer to bypass "inout" intent attribute 975 call WffReadWrite_mpio(wff,2,cg_ptr,mcg,icg,nband_disk,npwso,npwsotot,ind_cg_mpi_to_seq,ios) 976 nullify(cg_ptr) 977 #ifdef HAVE_NETCDF 978 else if (wff%iomode == IO_MODE_ETSF .and. option/=5) then 979 980 NCF_CHECK(nf90_inq_varid(wff%unwff, "reduced_coordinates_of_plane_waves", kg_varid)) 981 NCF_CHECK(nf90_get_att(wff%unwff, kg_varid, "k_dependent", kdep)) 982 983 ! The coefficients_of_wavefunctions on file have shape [cplex, mpw, nspinor, mband, nkpt, nsppol] 984 NCF_CHECK(nf90_inq_varid(wff%unwff, "coefficients_of_wavefunctions", cg_varid)) 985 !write(std_out,*)"put cg, count: ",[2,npw,nspinor,nband,1,1] 986 ncerr = nf90_put_var(wff%unwff, cg_varid, cg(:, icg+1:), start=[1,1,1,1,ikpt,isppol], & 987 count=[2,npw,nspinor,nband,1,1]) 988 NCF_CHECK_MSG(ncerr, "putting cg_k") 989 !write(std_out,*)"after cg" 990 #endif 991 end if 992 end if ! option/=4 993 994 ! ===== Case formeig=1: write eigenvalues and wave-functions ===== 995 else if(formeig==1)then 996 997 ! Not available for NETCDF and ETSF_IO 998 ABI_CHECK(wff%iomode /= IO_MODE_ETSF, "ETSF-write-eigen1 not coded!") 999 1000 ! ABI_CHECK(nband_disk==nband,"nband_disk!=nband") 1001 1002 nband2=2*nband_disk 1003 do iband=1,nband_disk 1004 ipw=(iband-1)*npwso+icg 1005 ii=(iband-1)*nband2 1006 if(use_f90==1)then 1007 write(unitwf) eigen(1+ii:nband2+ii) 1008 if (option/=5) then 1009 write(unitwf) cg(1:2,1+ipw:npwso+ipw) 1010 else 1011 write(unitwf) 1012 end if 1013 else if(wff%iomode==IO_MODE_MPI)then 1014 ! Should use an access with a "view" 1015 call xderiveWRecInit(wff,ios) 1016 call xderiveWrite(wff,eigen(1+ii:ii+nband2),nband2,wff%spaceComm_mpiio,ios) 1017 call xderiveWRecEnd(wff,ios) 1018 if (option/=4) then 1019 call xderiveWRecInit(wff,ios) 1020 call xderiveWrite(wff,cg(1:2,ipw+1:ipw+npwso),2,npwso,wff%spaceComm_mpiio,ios) 1021 call xderiveWRecEnd(wff,ios) 1022 end if 1023 end if 1024 end do 1025 1026 end if ! formeig 1027 1028 !--------------------------------------------------------------------------- 1029 !Final statements 1030 !--------------------------------------------------------------------------- 1031 1032 if (allocated(ind_cg_mpi_to_seq)) then 1033 ABI_DEALLOCATE(ind_cg_mpi_to_seq) 1034 end if 1035 1036 RETURN 1037 1038 !Silence compiler warning. 1039 ABI_UNUSED((/ii,mpi_enreg%me/)) 1040 1041 end subroutine writewf