TABLE OF CONTENTS
ABINIT/m_rwwf [ Modules ]
NAME
m_rwwf
FUNCTION
Read/Write wavefunctions.
COPYRIGHT
Copyright (C) 1998-2024 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 .
SOURCE
16 #if defined HAVE_CONFIG_H 17 #include "config.h" 18 #endif 19 20 #include "abi_common.h" 21 22 module m_rwwf 23 24 use defs_basis 25 use m_errors 26 use m_wffile 27 use m_abicore 28 use m_xmpi 29 #if defined HAVE_MPI2 30 use mpi 31 #endif 32 use m_nctk 33 #ifdef HAVE_NETCDF 34 use netcdf 35 #endif 36 37 use defs_abitypes, only : mpi_type 38 use m_time, only : timab 39 40 implicit none 41 42 private 43 44 #if defined HAVE_MPI1 45 include 'mpif.h' 46 #endif
m_rwwf/readwf [ Functions ]
[ Top ] [ m_rwwf ] [ 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.
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
SOURCE
292 subroutine readwf(cg,eigen,formeig,headform,icg,ikpt,isppol,kg_k,mband,mcg,mpi_enreg,& 293 & nband,nband_disk,npw,nspinor,occ,option,optkg,wff) 294 295 !Arguments ------------------------------------ 296 integer,intent(in) :: formeig,headform,icg,ikpt,isppol,mband,mcg,nband,npw 297 integer,intent(in) :: nspinor,option,optkg 298 integer,intent(inout) :: nband_disk 299 integer,intent(inout),target :: kg_k(3,optkg*npw) 300 real(dp),intent(inout),target :: cg(2,mcg),eigen((2*mband)**formeig*mband),occ(mband) 301 type(MPI_type),intent(in) :: mpi_enreg 302 type(wffile_type),intent(inout) :: wff 303 304 !Local variables------------------------------- 305 integer :: iband,indxx,ios,ipw,ispinor_index,nband1,ncid_hdr 306 integer :: npw1,npwso,npwso1,npwsotot,npwtot,nrec,nspinor1,nspinortot,unitwf,use_f90 307 integer :: band1,band2,ierr 308 character(len=500) :: msg 309 character(len=fnlen) :: fname 310 integer :: ikpt_this_proc,ispinor 311 integer,allocatable :: ind_cg_mpi_to_seq(:) 312 #ifdef HAVE_NETCDF 313 integer :: kg_varid,eig_varid,occ_varid,cg_varid,h1_varid,mband_varid,ncerr,ii,mband_file 314 integer,allocatable :: gall(:,:) 315 real(dp),allocatable :: cg_all(:,:,:),h1mat(:,:,:) 316 #endif 317 318 ! ********************************************************************* 319 320 !write(std_out,*)"rwwf with ikpt, isppol, option, etsf",ikpt, isppol, option, wff%iomode == IO_MODE_ETSF 321 322 !Check the options 323 if ( ALL(option /= [1,3,-1,-2,-4])) then 324 write(msg,'(a,i0,a)')'The argument option should be -4, -2, -1, 1 or 3. However, option=',option,'.' 325 ABI_BUG(msg) 326 end if 327 328 npwtot=npw; npwso=npw*nspinor 329 npwsotot=npwso 330 nspinortot=min(2,(1+mpi_enreg%paral_spinor)*nspinor) 331 332 unitwf=wff%unwff 333 ncid_hdr=unitwf 334 335 use_f90=0 336 if (wff%iomode==IO_MODE_FORTRAN.or.(wff%iomode==IO_MODE_FORTRAN_MASTER.and.wff%master==wff%me)) use_f90=1 337 338 if (option==1.or.option==-2) then 339 340 ! Compute mapping my_gtable --> sequential gtable. 341 if (any(wff%iomode==[IO_MODE_MPI, IO_MODE_ETSF]).and.nband>0) then 342 call xmpi_sum(npwsotot,wff%spaceComm_mpiio,ios) 343 npwtot=npwsotot/nspinortot 344 ABI_MALLOC(ind_cg_mpi_to_seq,(npwso)) 345 if (allocated(mpi_enreg%my_kgtab)) then 346 ikpt_this_proc=mpi_enreg%my_kpttab(ikpt) 347 if ( ikpt_this_proc <= 0 ) then 348 ABI_BUG("rwwf: ikpt_this_proc <= 0") 349 end if 350 do ispinor=1,nspinor 351 ispinor_index=ispinor 352 if (mpi_enreg%nproc_spinor>1) ispinor_index=mpi_enreg%me_spinor + 1 353 ind_cg_mpi_to_seq(1+npw*(ispinor-1):npw*ispinor)=npwtot*(ispinor_index-1) & 354 & + mpi_enreg%my_kgtab(1:npw,ikpt_this_proc) 355 end do 356 else 357 ind_cg_mpi_to_seq(1:npwso) = (/(ipw,ipw=1,npwso)/) 358 end if 359 end if 360 end if 361 362 !--------------------------------------------------------------------------- 363 !Read the first record: npw, nspinor, nband_disk 364 !--------------------------------------------------------------------------- 365 366 if (headform>=40.or.headform==0) then ! headform==0 refers to the current headform 367 call WffReadNpwRec(ios,ikpt,isppol,nband_disk,npw1,nspinor1,wff) 368 npwso1=npw1*nspinor1 369 if(ios/=0)then 370 inquire (unit=unitwf, NAME=fname) 371 write(msg,'(3a,i4,2a,i4,4a)') & 372 & 'Reading option of rwwf. Trying to read',ch10,& 373 & 'the (npw,nspinor,nband) record of a wf file, unit=',unitwf,ch10,& 374 & 'gave iostat=',ios,'. Your file is likely not correct.',ch10,& 375 & 'Action: check your input wf file:',trim(fname) 376 ABI_ERROR(msg) 377 end if 378 379 else 380 ! Old format 381 if(use_f90==1)then 382 read (unitwf,iostat=ios) npwso1,nband_disk 383 else if(wff%iomode==IO_MODE_MPI)then 384 call xderiveRRecInit(wff,ios) 385 call xderiveRead(wff,npwso1,ios) 386 call xderiveRead(wff,nband_disk,ios) 387 call xderiveRRecEnd(wff,ios) 388 end if 389 if(ios/=0)then 390 inquire (unit=unitwf, NAME=fname) 391 write(msg,'(3a,i4,2a,i4,4a)') & 392 & 'Reading option of rwwf. Trying to read',ch10,& 393 & 'the (npw,nband) record of a wf file with headform <40 , unit=',unitwf,ch10,& 394 & 'gave iostat=',ios,'. Your file is likely not correct.',ch10,& 395 & 'Action: check your input wf file:',trim(fname) 396 ABI_ERROR(msg) 397 end if 398 end if ! headform 399 400 if (option==1.or.option==-2) then ! Will read the wavefunction and/or kg data, so check npw and nspinor 401 402 if (headform>=40.or.headform==0) then ! New format. headform==0 refers to the current headform 403 if (npwtot/=npw1) then 404 write(msg,'(3a,i0,a,i0,a)') & 405 & 'Reading option of rwwf. One should have npwtot=npw1',ch10,& 406 & 'However, npwtot= ',npwtot,', and npw1= ',npw1,'.' 407 ABI_BUG(msg) 408 end if 409 if(nspinortot/=nspinor1)then 410 write(msg,'(3a,i0,a,i0,a)') & 411 & 'Reading option of rwwf. One should have nspinor=nspinor1',ch10,& 412 & 'However, nspinortot= ',nspinortot,', and nspinor1= ',nspinor1,'.' 413 ABI_BUG(msg) 414 end if 415 else ! Treat the Old format. 416 if(npwsotot/=npwso1)then 417 write(msg,'(3a,i0,a,i0,a)') & 418 & 'Reading option of rwwf. One should have npwso=npwso1',ch10,& 419 & 'However, npwsotot= ',npwsotot,', and npwso1= ',npwso1,'.' 420 ABI_BUG(msg) 421 end if 422 end if ! headform 423 ! 424 end if ! option==1.or.option==2 425 426 !--------------------------------------------------------------------------- 427 !Read the second record: (k+G) vectors 428 !--------------------------------------------------------------------------- 429 430 if (headform>=40.or.headform==0) then ! headform==0 refers to the current headform 431 432 if ((option==1.or.option==-2.or.option==3).and.optkg/=0 )then 433 434 if(use_f90==1)then 435 read(unitwf,iostat=ios) kg_k(1:3,1:npw) 436 437 else if(wff%iomode==IO_MODE_MPI)then 438 call xderiveRRecInit(wff,ios) 439 if (allocated(mpi_enreg%my_kgtab)) then 440 ikpt_this_proc=mpi_enreg%my_kpttab(ikpt) 441 call xderiveRead(wff,kg_k(1:3,1:npw),3,npw,wff%spaceComm_mpiio,& 442 & mpi_enreg%my_kgtab(1:npw,ikpt_this_proc),ios) 443 else 444 ! MG The call below uses MPI_SCAN but here we want to read the full set of G. 445 call xderiveRead(wff,kg_k(1:3,1:npw),3,npw,wff%spaceComm_mpiio,ios) 446 call xmpi_read_int2d(wff,kg_k(1:3,1:npw),wff%spaceComm_mpiio,xmpio_collective,ios) 447 end if 448 call xderiveRRecEnd(wff,ios) 449 450 #ifdef HAVE_NETCDF 451 else if (wff%iomode == IO_MODE_ETSF) then 452 ! Read reduced_coordinates_of_plane_waves for this k point (npw1 is npw_disk). 453 ! TODO: spinor parallelism 454 NCF_CHECK(nf90_inq_varid(wff%unwff, "reduced_coordinates_of_plane_waves", kg_varid)) 455 if (npw == npw1) then 456 ncerr = nf90_get_var(wff%unwff, kg_varid, kg_k, start=[1,1,ikpt], count=[3,npw1,1]) 457 NCF_CHECK_MSG(ncerr, "getting kg_k") 458 else 459 write(std_out,*)"ETSF Reading distributed kg_k" 460 ABI_MALLOC(gall, (3, npw1)) 461 ncerr = nf90_get_var(wff%unwff, kg_varid, gall, start=[1,1,ikpt], count=[3,npw1,1]) 462 NCF_CHECK_MSG(ncerr, "getting kg_k") 463 do ipw=1,npw 464 kg_k(:,ipw) = gall(:,ind_cg_mpi_to_seq(ipw)) 465 end do 466 ABI_FREE(gall) 467 end if 468 #endif 469 end if 470 471 else ! option 472 call WffReadSkipRec(ios,1,wff) ! Skip the record 473 end if 474 475 if(ios/=0)then 476 inquire (unit=unitwf, NAME=fname) 477 write(msg,'(3a,i4,2a,i4,4a)') & 478 & 'Reading option of rwwf. Trying to read',ch10,& 479 & 'the k+g record of a wf file, unit=',unitwf,ch10,& 480 & 'gave iostat=',ios,'. Your file is likely not correct.',ch10,& 481 & 'Action: check your input wf file:',trim(fname) 482 ABI_ERROR(msg) 483 end if 484 end if ! headform 485 486 !--------------------------------------------------------------------------- 487 !Read the third record: eigenvalues 488 !--------------------------------------------------------------------------- 489 !The reading of occ should be enabled, BUT taking into account 490 !of headform of the disk file : occ was NOT present in the disk files with headform=22 491 492 nband1 = min(nband,nband_disk) 493 494 !===== Case formeig=0: read eigenvalues ===== 495 if (formeig==0) then 496 497 if (option==1.or.option==3.or.option==-4) then 498 if(use_f90==1)then 499 read (unitwf,iostat=ios) eigen(1:nband1) 500 else if(wff%iomode==IO_MODE_MPI)then 501 call xderiveRRecInit(wff,ios) 502 call xderiveRead(wff,eigen,nband1,xmpi_comm_self,ios) 503 call xderiveRRecEnd(wff,ios) 504 end if 505 else 506 call WffReadSkipRec(ios,1,wff) 507 end if ! option 508 509 if(ios/=0)then 510 inquire (unit=unitwf, NAME=fname) 511 write(msg,'(3a,i4,2a,i4,4a)') & 512 & 'Reading option of rwwf. Trying to read',ch10,& 513 & 'an eigenvalue record of a wf file, unit=',unitwf,ch10,& 514 & 'gave iostat=',ios,'. Your file is likely not correct.',ch10,& 515 & 'Action: check your input wf file.',trim(fname) 516 ABI_ERROR(msg) 517 end if 518 519 #ifdef HAVE_NETCDF 520 if (wff%iomode == IO_MODE_ETSF) then 521 ! get eigenvalues and occupations 522 NCF_CHECK(nf90_inq_varid(wff%unwff, "eigenvalues", eig_varid)) 523 ncerr = nf90_get_var(wff%unwff, eig_varid, eigen, start=[1,ikpt,isppol], count=[nband1,1,1]) 524 NCF_CHECK_MSG(ncerr, "getting eig_k") 525 526 NCF_CHECK(nf90_inq_varid(wff%unwff, "occupations", occ_varid)) 527 ncerr = nf90_get_var(wff%unwff, occ_varid, occ, start=[1,ikpt,isppol], count=[nband1,1,1]) 528 NCF_CHECK_MSG(ncerr, "getting occ_k") 529 end if 530 #endif 531 532 ! ===== Case formeig=1: read matrix of eigenvalues ===== 533 ! Will be written later (together with wave-functions) 534 else if(formeig==1)then 535 end if ! formeig 536 537 !--------------------------------------------------------------------------- 538 !Read the wave-function coefficients 539 !--------------------------------------------------------------------------- 540 541 !Select bands 542 nband1=min(nband,nband_disk) 543 if(nband1>0.and.option/=-1)then 544 545 ! ===== Case formeig=0: read only wave-functions ===== 546 if (formeig==0) then 547 548 if (option==1.or.option==-2) then 549 550 if (use_f90==1) then 551 do iband=1,nband1 552 ipw=(iband-1)*npwso+icg 553 read(unitwf,iostat=ios) cg(1:2,ipw+1:ipw+npwso) 554 if (ios/=0) exit 555 end do 556 else if (wff%iomode==IO_MODE_MPI) then 557 call WffReadWrite_mpio(wff,1,cg,mcg,icg,nband1,npwso,npwsotot,ind_cg_mpi_to_seq,ios) 558 end if 559 560 else if (option/=-4) then 561 do iband=1,nband1 562 call WffReadSkipRec(ios,1,wff) ! Skip the record 563 end do 564 end if ! option 565 566 if(ios/=0)then 567 inquire (unit=unitwf, NAME=fname) 568 write(msg,'(3a,i4,2a,i4,4a)') & 569 & 'Reading option of rwwf. Trying to read',ch10,& 570 & 'a RF wf record of a wf file, unit=',unitwf,ch10,& 571 & 'gave iostat=',ios,'. Your file is likely not correct.',ch10,& 572 & 'Action: check your input wf file.',trim(fname) 573 ABI_ERROR(msg) 574 end if 575 576 #ifdef HAVE_NETCDF 577 if (wff%iomode == IO_MODE_ETSF) then 578 ! The coefficients_of_wavefunctions on file have shape [cplex, mpw, nspinor, mband, nkpt, nsppol] 579 NCF_CHECK(nf90_inq_varid(wff%unwff, "coefficients_of_wavefunctions", cg_varid)) 580 if (npw == npw1) then 581 ncerr = nf90_get_var(wff%unwff, cg_varid, cg(:, icg+1:), start=[1,1,1,1,ikpt,isppol], & 582 count=[2,npw,nspinor,nband1,1,1]) 583 NCF_CHECK_MSG(ncerr, "getting cg_k") 584 else 585 write(std_out,*)"ETSF Reading distributed cg" 586 ABI_MALLOC_OR_DIE(cg_all, (2, npw1*nspinor, nband1), ierr) 587 ncerr = nf90_get_var(wff%unwff, cg_varid, cg_all, start=[1,1,1,1,ikpt,isppol], & 588 count=[2,npw1,nspinor,nband1,1,1]) 589 NCF_CHECK_MSG(ncerr, "getting cg_k") 590 ii = icg 591 do iband=1,nband1 592 do ipw=1,npw 593 ii = ii + 1 594 cg(:,ii) = cg_all(:,ind_cg_mpi_to_seq(ipw),iband) 595 end do 596 end do 597 ABI_FREE(cg_all) 598 end if 599 end if 600 #endif 601 602 ! ===== Case formeig=1: read eigenvalues, occupations and wave-functions ===== 603 else if(formeig==1)then 604 ! write(std_out,*)"nband1",nband1 605 ! ABI_CHECK(nband1==nband_disk,"nband != nband_disk") 606 607 if (wff%iomode /= IO_MODE_ETSF) then 608 indxx=0 609 do iband=1,nband1 610 611 if (option==1.or.option==3.or.option==-4) then 612 if(use_f90==1)then 613 read (unitwf,iostat=ios) eigen(1+indxx:2*nband1+indxx) 614 else if(wff%iomode==IO_MODE_MPI)then 615 ! Should use an access with a "view" 616 call xderiveRRecInit(wff,ios) 617 call xderiveRead(wff,eigen(1+indxx:2*nband1+indxx),2*nband1,xmpi_comm_self,ios) 618 call xderiveRRecEnd(wff,ios) 619 end if 620 indxx=indxx+2*nband1 621 else 622 call WffReadSkipRec(ios,1,wff) ! Skip the record 623 end if 624 625 if(ios/=0)then 626 inquire (unit=unitwf, NAME=fname) 627 write(msg,'(3a,i4,2a,i4,4a)') & 628 & 'Reading option of rwwf. Trying to read',ch10,& 629 & 'a RF eigenvalue record of a wf file, unit=',unitwf,ch10,& 630 & 'gave iostat=',ios,'. Your file is likely not correct.',ch10,& 631 & 'Action: check your input wf file.',trim(fname) 632 ABI_ERROR(msg) 633 end if 634 635 if(option==1.or.option==-2)then 636 ipw=(iband-1)*npwso+icg 637 if(use_f90==1)then 638 ipw=(iband-1)*npwso+icg 639 read(unitwf,iostat=ios) cg(1:2,ipw+1:ipw+npwso) 640 else if(wff%iomode==IO_MODE_MPI)then 641 ! Should use an access with a "view" 642 call xderiveRRecInit(wff,ios) 643 call xderiveRead(wff,cg(1:2,ipw+1:ipw+npwso),2,npwso,wff%spaceComm_mpiio,ind_cg_mpi_to_seq,ios) 644 call xderiveRRecEnd(wff,ios) 645 end if 646 else if (option/=-4) then 647 call WffReadSkipRec(ios,1,wff) ! Skip the record 648 end if ! option 649 650 if(ios/=0)then 651 inquire (unit=unitwf, NAME=fname) 652 write(msg,'(3a,i4,2a,i4,4a)') & 653 & 'Reading option of rwwf. Trying to read',ch10,& 654 & 'a RF wf record of a wf file, unit=',unitwf,ch10,& 655 & 'gave iostat=',ios,'. Your file is likely not correct.',ch10,& 656 & 'Action: check your input wf file.',trim(fname) 657 ABI_ERROR(msg) 658 end if 659 end do ! iband 660 661 else 662 #ifdef HAVE_NETCDF 663 ! ETSF-IO 664 if (any(option == [1, 3, -4])) then 665 ! Read eigen. Remember that the matrix on file has shape [2, mband, mband, nkpt, nspin] 666 ! whereas the eigen array used in Abinit is packed. Read the full matrix first, then pack data. 667 NCF_CHECK(nf90_inq_dimid(wff%unwff, "max_number_of_states", mband_varid)) 668 NCF_CHECK(nf90_inquire_dimension(wff%unwff, mband_varid, len=mband_file)) 669 h1_varid = nctk_idname(wff%unwff, "h1_matrix_elements") 670 671 ABI_MALLOC(h1mat, (2, mband_file, mband_file)) 672 ncerr = nf90_get_var(wff%unwff, h1_varid, h1mat, start=[1,1,1,ikpt,isppol]) 673 NCF_CHECK_MSG(ncerr, "getting h1_matrix_elements") 674 675 indxx=1 676 do band2=1,nband1 677 do band1=1,nband1 678 eigen(indxx:indxx+1) = h1mat(:,band1,band2) 679 indxx = indxx + 2 680 end do 681 end do 682 ABI_FREE(h1mat) 683 end if 684 685 if (any(option == [1, -2])) then 686 ! Read wavefunctions. 687 ! The coefficients_of_wavefunctions on file have shape [cplex, mpw, nspinor, mband, nkpt, nsppol] 688 ABI_CHECK(npw == npw1, "npw != npw1 not coded") 689 690 cg_varid = nctk_idname(wff%unwff, "coefficients_of_wavefunctions") 691 ncerr = nf90_get_var(wff%unwff, cg_varid, cg(:, icg+1:), start=[1,1,1,1,ikpt,isppol], & 692 count=[2,npw,nspinor,nband1,1,1]) 693 NCF_CHECK_MSG(ncerr, "getting cg_k") 694 end if 695 #endif 696 end if 697 698 end if ! formeig == 1 699 700 end if ! nband >0 701 702 !If fewer than all bands were read wind disk file forward to end of bands for this k point. 703 !Will have to fill the non-filled bands outside of this routine ... 704 if (nband<nband_disk .or. option==-1) then 705 nrec=(formeig+1)*(nband_disk-nband) 706 if(option==-1)nrec=(formeig+1)*nband_disk 707 call WffReadSkipRec(ios,nrec,wff) 708 end if 709 710 !--------------------------------------------------------------------------- 711 ! Free memory 712 !--------------------------------------------------------------------------- 713 if (allocated(ind_cg_mpi_to_seq)) then 714 ABI_FREE(ind_cg_mpi_to_seq) 715 end if 716 717 end subroutine readwf
m_rwwf/rwwf [ Functions ]
[ Top ] [ m_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.
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.
SOURCE
127 subroutine rwwf(cg,eigen,formeig,headform,icg,ikpt,isppol,kg_k,mband,mcg,mpi_enreg,& 128 & nband,nband_disk,npw,nspinor,occ,option,optkg,tim_rwwf,wff) 129 130 !Arguments ------------------------------------ 131 integer,intent(in) :: formeig,headform,icg,ikpt,isppol,mband,mcg,nband,npw 132 integer,intent(inout) :: nband_disk 133 integer,intent(in) :: nspinor,option,optkg,tim_rwwf 134 integer,intent(inout),target :: kg_k(3,optkg*npw) 135 real(dp),intent(inout),target :: cg(2,mcg),eigen((2*mband)**formeig*mband),occ(mband) 136 type(wffile_type),intent(inout) :: wff 137 type(MPI_type), intent(in) :: mpi_enreg 138 139 !Local variables------------------------------- 140 !scalars 141 character(len=500) :: msg 142 !arrays 143 real(dp) :: tsec(2) 144 145 ! ************************************************************************* 146 147 call timab(270+tim_rwwf,1,tsec) 148 149 !Might check that icg+npw*nband*nspinor is smaller than mcg 150 151 !Check that nband is smaller than mband, if one will not skip the records. 152 if (nband>mband .and. option/=-1)then 153 write(msg,'(a,i0,a,i0,a)')' One should have nband<=mband. However, nband=',nband,', and mband=',mband,'.' 154 ABI_BUG(msg) 155 end if 156 157 !Check that formeig is 0 or 1. 158 if ( ALL(formeig/=(/0,1/)) ) then 159 write(msg,'(a,i0,a)')' The argument formeig should be 0 or 1. However, formeig=',formeig,'.' 160 ABI_BUG(msg) 161 end if 162 163 !Check the value of option 164 if ( ALL( option /= (/1,2,3,4,5,-1,-2,-4/) )) then 165 write(msg,'(a,i0,a)')' The argument option should be between 1 and 5, or -1, -2, -4. However, option=',option,'.' 166 ABI_BUG(msg) 167 end if 168 169 if (option/=2.and.option/=4 .and. option/=5 ) then ! read 170 call readwf(cg,eigen,formeig,headform,icg,ikpt,isppol,kg_k,mband,mcg,mpi_enreg,& 171 & nband,nband_disk,npw,nspinor,occ,option,optkg,wff) 172 else ! write 173 call writewf(cg,eigen,formeig,icg,ikpt,isppol,kg_k,mband,mcg,mpi_enreg,& 174 & nband,nband_disk,npw,nspinor,occ,option,optkg,wff) 175 end if 176 177 call timab(270+tim_rwwf,2,tsec) 178 179 end subroutine rwwf
m_rwwf/WffReadSkipK [ Functions ]
[ Top ] [ m_rwwf ] [ Functions ]
NAME
WffReadSkipK
FUNCTION
(Wavefunction file, read action : skip one k-point blok) This subroutine skips the block of records related to one k point, and one spin-polarization, that contains the wavefunctions as well as the eigenvalues and occupations, in a wavefunction file that has been already initialized.
INPUTS
formeig=format of the eigenvalues 0 => vector of eigenvalues (for Ground-State files) 1 => hermitian matrix of eigenvalues (for Response-Function files) 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 ikpt=index of current k point (only needed for error message) isppol=spin polarization currently treated (only needed for error message) mpi_enreg=information about MPI parallelization wff=structured info for wavefunction file
OUTPUT
NOTES
SOURCE
210 subroutine WffReadSkipK(formeig,headform,ikpt,isppol,mpi_enreg,wff) 211 212 !Arguments ------------------------------------ 213 !scalars 214 integer,intent(in) :: formeig,headform,ikpt,isppol 215 type(MPI_type),intent(in) :: mpi_enreg 216 type(wffile_type),intent(inout) :: wff 217 218 !Local variables------------------------------- 219 !scalars 220 integer :: icg,mband,mcg,nband,nband_disk,option,optkg,tim_rwwf 221 integer,parameter :: nspinor1=1,npw1=1 222 !arrays 223 integer,allocatable :: kg_dum(:,:) 224 real(dp) :: cg_dum(2,1),occ_dum(1) 225 real(dp),allocatable :: eig_dum(:) 226 227 ! ************************************************************************* 228 229 ! No need to skip if netcdf 230 if (wff%iomode == IO_MODE_ETSF) return 231 232 option=-1 233 tim_rwwf=0 ; mcg=1 ; mband=1 ; icg=0 ; optkg=0 ; nband=0 234 ABI_MALLOC(eig_dum,(2**formeig)) 235 ABI_MALLOC(kg_dum,(3,optkg*npw1)) 236 237 call rwwf(cg_dum,eig_dum,formeig,headform,icg,ikpt,isppol,kg_dum,mband,mcg,mpi_enreg,nband,& 238 & nband_disk,npw1,nspinor1,occ_dum,option,optkg,tim_rwwf,wff) 239 240 ABI_FREE(eig_dum) 241 ABI_FREE(kg_dum) 242 243 end subroutine WffReadSkipK
m_rwwf/writewf [ Functions ]
[ Top ] [ m_rwwf ] [ 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.
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
SOURCE
767 subroutine writewf(cg,eigen,formeig,icg,ikpt,isppol,kg_k,mband,mcg,mpi_enreg,& 768 & nband,nband_disk,npw,nspinor,occ,option,optkg,wff) 769 770 !Arguments ------------------------------------ 771 integer,intent(in) :: formeig,icg,ikpt,isppol,mband,mcg,nband,nband_disk,npw,nspinor,option,optkg 772 integer,intent(in),target :: kg_k(3,optkg*npw) 773 real(dp),intent(in),target :: cg(2,mcg),eigen((2*mband)**formeig*mband),occ(mband) 774 type(MPI_type),intent(in) :: mpi_enreg 775 type(wffile_type),intent(inout) :: wff 776 777 !Local variables------------------------------- 778 integer :: iband,ii,ios,ipw,ispinor_index,nband2,ncid_hdr,npwso,npwsotot,npwtot,nspinortot 779 integer :: unitwf,use_f90 780 character(len=500) :: msg 781 integer :: ikpt_this_proc,ispinor,me_cart_3d 782 integer,allocatable :: ind_cg_mpi_to_seq(:) 783 real(dp),ABI_CONTIGUOUS pointer :: cg_ptr(:,:) 784 #ifdef HAVE_NETCDF 785 integer :: kg_varid,eig_varid,occ_varid,cg_varid,ncerr 786 character(len=nctk_slen) :: kdep 787 #endif 788 #ifdef HAVE_MPI 789 integer(kind=MPI_OFFSET_KIND) :: off(1) 790 #endif 791 792 ! ********************************************************************* 793 794 !Check the options 795 if ( ALL(option /= (/2,4,5/)) ) then 796 write(msg,'(a,i0)')' The argument option should be 2, 4 or 5. However, option=',option 797 ABI_BUG(msg) 798 end if 799 800 if(wff%iomode==IO_MODE_MPI)then 801 if(option==5)then 802 write(msg,'(a,i0)')' With MPI-IO activated, the argument option should be 2 or 4. However, option=',option 803 ABI_BUG(msg) 804 end if 805 end if 806 807 if (wff%iomode==IO_MODE_ETSF) then 808 ! outkss still calls this routine! 809 ABI_WARNING("You should use outwff or ncwrite_cg to write data in netcdf format!") 810 end if 811 812 !Check that nband_disk is not larger than nband (only for writing) 813 if (nband<nband_disk) then 814 write(msg,'(3a,i5,a,i5,a)') & 815 & 'Writing option of rwwf. One should have nband<=nband_disk',ch10,& 816 & 'However, nband= ',nband,', and nband_disk= ',nband_disk,'.' 817 ABI_BUG(msg) 818 end if 819 820 npwtot=npw; npwso=npw*nspinor 821 unitwf=wff%unwff; ncid_hdr=unitwf 822 npwsotot=npwso 823 nspinortot=min(2,(1+mpi_enreg%paral_spinor)*nspinor) 824 825 use_f90=0 826 if (wff%iomode==IO_MODE_FORTRAN.or.(wff%iomode ==IO_MODE_FORTRAN_MASTER.and.wff%master==wff%me)) use_f90=1 827 828 if (wff%iomode==IO_MODE_MPI) then 829 830 call xmpi_sum(npwsotot,wff%spaceComm_mpiio,ios) 831 npwtot=npwsotot/nspinortot 832 833 if (option/=4) then 834 ABI_MALLOC(ind_cg_mpi_to_seq,(npwso)) 835 if (allocated(mpi_enreg%my_kgtab)) then 836 ikpt_this_proc=mpi_enreg%my_kpttab(ikpt) 837 do ispinor=1,nspinor 838 ispinor_index=ispinor 839 if (mpi_enreg%nproc_spinor > 1) ispinor_index = mpi_enreg%me_spinor + 1 840 ind_cg_mpi_to_seq(1+npw*(ispinor-1):npw*ispinor)=npwtot*(ispinor_index-1) & 841 & + mpi_enreg%my_kgtab(1:npw,ikpt_this_proc) 842 end do 843 else 844 ind_cg_mpi_to_seq(1:npwso) = (/(ipw,ipw=1,npwso)/) 845 end if 846 !write(std_out,*)"MPI-IO ind_cg_mpi_to_seq", ind_cg_mpi_to_seq(1:5) 847 end if 848 end if 849 850 !--------------------------------------------------------------------------- 851 !Write the first record: npw, nspinor, nband_disk 852 !--------------------------------------------------------------------------- 853 !Not modified for netCDF: no need to add writing of nband_disk,npw,nspinor 854 855 call WffWriteNpwRec(ios,nband_disk,npwtot,nspinortot,wff,opt_paral=2) 856 857 !--------------------------------------------------------------------------- 858 !Write the second record: (k+G) vectors 859 !--------------------------------------------------------------------------- 860 861 if (optkg/=0.and.option/=4) then 862 if(use_f90==1)then 863 write(unitwf) kg_k(1:3,1:optkg*npw) 864 else if (wff%iomode==IO_MODE_MPI) then 865 866 if (allocated(mpi_enreg%my_kgtab)) then 867 me_cart_3d=xmpi_comm_rank(mpi_enreg%comm_bandspinorfft) 868 ikpt_this_proc=mpi_enreg%my_kpttab(ikpt) 869 call xderiveWRecInit(wff,ios,me_cart_3d) 870 if (mpi_enreg%me_spinor==0) then 871 call xderiveWrite(wff,kg_k,3,npw,mpi_enreg%comm_bandfft, & 872 & mpi_enreg%my_kgtab(1:npw,ikpt_this_proc),ios) 873 end if 874 call xderiveWRecEnd(wff,ios,me_cart_3d) 875 else 876 ! MG does it work if we are not using FFT distribution ? 877 call xderiveWRecInit(wff,ios ) 878 if (mpi_enreg%me_spinor==0) then 879 call xderiveWrite(wff,kg_k,3,optkg*npw,Wff%spaceComm_mpiio,ios) 880 end if 881 call xderiveWRecEnd(wff,ios) 882 end if 883 884 #ifdef HAVE_NETCDF 885 else if (wff%iomode == IO_MODE_ETSF) then 886 ! Write the reduced_coordinates_of_plane_waves for this k point. 887 NCF_CHECK(nf90_inq_varid(wff%unwff, "reduced_coordinates_of_plane_waves", kg_varid)) 888 NCF_CHECK(nf90_get_att(wff%unwff, kg_varid, "k_dependent", kdep)) 889 if (kdep == "no") then 890 ncerr = nf90_put_var(wff%unwff, kg_varid, kg_k, start=[1,1], count=[3,npw]) 891 else 892 ncerr = nf90_put_var(wff%unwff, kg_varid, kg_k, start=[1,1,ikpt], count=[3,npw,1]) 893 end if 894 NCF_CHECK_MSG(ncerr, "putting kg_k") 895 #endif 896 end if ! end if wff%iomode 897 else ! Still skip the record 898 if (use_f90==1) then 899 write(unitwf) 900 else if (wff%iomode==IO_MODE_MPI) then 901 call xderiveWRecInit(wff,wff%spaceComm_mpiio,ios) 902 call xderiveWRecEnd(wff,wff%spaceComm_mpiio,ios) 903 end if 904 end if 905 906 !--------------------------------------------------------------------------- 907 !Write the third record: eigenvalues and occupations 908 !--------------------------------------------------------------------------- 909 910 !===== Case formeig=0: write eigenvalues and occupations ===== 911 if (formeig==0) then 912 if (use_f90==1) then 913 write(unitwf) (eigen(iband),iband=1,nband_disk),(occ(iband),iband=1,nband_disk) 914 else if(wff%iomode==IO_MODE_MPI) then 915 if (wff%me_mpiio==0) then 916 call xderiveWRecInit(wff,ios) 917 call xderiveWrite(wff,eigen,nband_disk,xmpi_comm_self,ios) 918 call xderiveWrite(wff,occ,nband_disk,xmpi_comm_self,ios) 919 call xderiveWRecEnd(wff,ios) 920 end if 921 #ifdef HAVE_MPI 922 off(1)=wff%offwff 923 call MPI_BCAST(off,1,wff%offset_mpi_type,0,wff%spaceComm_mpiio,ios) 924 wff%offwff=off(1) 925 #endif 926 #ifdef HAVE_NETCDF 927 else if (wff%iomode == IO_MODE_ETSF) then 928 ! Write eigenvalues and occupation factors. 929 NCF_CHECK(nf90_inq_varid(wff%unwff, "eigenvalues", eig_varid)) 930 ncerr = nf90_put_var(wff%unwff, eig_varid, eigen, start=[1,ikpt,isppol], count=[mband,1,1]) 931 NCF_CHECK_MSG(ncerr, "putting eig_k") 932 933 NCF_CHECK(nf90_inq_varid(wff%unwff, "occupations", occ_varid)) 934 ncerr = nf90_put_var(wff%unwff, occ_varid, occ, start=[1,ikpt,isppol], count=[mband,1,1]) 935 NCF_CHECK_MSG(ncerr, "putting occ_k") 936 #endif 937 end if 938 939 ! ===== Case formeig=1: write matrix of eigenvalues ===== 940 ! Will be written later (together with wave-functions) 941 else if(formeig==1)then 942 end if ! formeig 943 944 !--------------------------------------------------------------------------- 945 !Write the wave-function coefficients 946 !--------------------------------------------------------------------------- 947 948 !===== Case formeig=0: write only wave-functions ===== 949 if (formeig==0) then 950 ! If option=4, do not write wave functions 951 if (option/=4) then 952 if (use_f90==1) then 953 do iband=1,nband_disk 954 ipw=(iband-1)*npwso+icg 955 if(option/=5)then 956 write(unitwf) cg(1:2,ipw+1:ipw+npwso) ! VALGRIND complains some elements of cg are not initialized, but written 957 else 958 write(unitwf) 959 end if 960 end do 961 else if(wff%iomode==IO_MODE_MPI)then 962 cg_ptr => cg ! Need pointer to bypass "inout" intent attribute 963 call WffReadWrite_mpio(wff,2,cg_ptr,mcg,icg,nband_disk,npwso,npwsotot,ind_cg_mpi_to_seq,ios) 964 nullify(cg_ptr) 965 #ifdef HAVE_NETCDF 966 else if (wff%iomode == IO_MODE_ETSF .and. option/=5) then 967 968 NCF_CHECK(nf90_inq_varid(wff%unwff, "reduced_coordinates_of_plane_waves", kg_varid)) 969 NCF_CHECK(nf90_get_att(wff%unwff, kg_varid, "k_dependent", kdep)) 970 971 ! The coefficients_of_wavefunctions on file have shape [cplex, mpw, nspinor, mband, nkpt, nsppol] 972 NCF_CHECK(nf90_inq_varid(wff%unwff, "coefficients_of_wavefunctions", cg_varid)) 973 !write(std_out,*)"put cg, count: ",[2,npw,nspinor,nband,1,1] 974 ncerr = nf90_put_var(wff%unwff, cg_varid, cg(:, icg+1:), start=[1,1,1,1,ikpt,isppol], & 975 count=[2,npw,nspinor,nband,1,1]) 976 NCF_CHECK_MSG(ncerr, "putting cg_k") 977 !write(std_out,*)"after cg" 978 #endif 979 end if 980 end if ! option/=4 981 982 ! ===== Case formeig=1: write eigenvalues and wave-functions ===== 983 else if(formeig==1)then 984 985 ! Not available for NETCDF and ETSF_IO 986 ABI_CHECK(wff%iomode /= IO_MODE_ETSF, "ETSF-write-eigen1 not coded!") 987 988 ! ABI_CHECK(nband_disk==nband,"nband_disk!=nband") 989 990 nband2=2*nband_disk 991 do iband=1,nband_disk 992 ipw=(iband-1)*npwso+icg 993 ii=(iband-1)*nband2 994 if(use_f90==1)then 995 write(unitwf) eigen(1+ii:nband2+ii) 996 if (option/=5) then 997 write(unitwf) cg(1:2,1+ipw:npwso+ipw) 998 else 999 write(unitwf) 1000 end if 1001 else if(wff%iomode==IO_MODE_MPI)then 1002 ! Should use an access with a "view" 1003 call xderiveWRecInit(wff,ios) 1004 call xderiveWrite(wff,eigen(1+ii:ii+nband2),nband2,wff%spaceComm_mpiio,ios) 1005 call xderiveWRecEnd(wff,ios) 1006 if (option/=4) then 1007 call xderiveWRecInit(wff,ios) 1008 call xderiveWrite(wff,cg(1:2,ipw+1:ipw+npwso),2,npwso,wff%spaceComm_mpiio,ios) 1009 call xderiveWRecEnd(wff,ios) 1010 end if 1011 end if 1012 end do 1013 1014 end if ! formeig 1015 1016 !--------------------------------------------------------------------------- 1017 !Final statements 1018 !--------------------------------------------------------------------------- 1019 1020 if (allocated(ind_cg_mpi_to_seq)) then 1021 ABI_FREE(ind_cg_mpi_to_seq) 1022 end if 1023 1024 RETURN 1025 1026 !Silence compiler warning. 1027 ABI_UNUSED((/ii,mpi_enreg%me/)) 1028 1029 end subroutine writewf