! ! Copyright (C) 2001-2005 Quantum ESPRESSO group ! This file is distributed under the terms of the ! GNU General Public License. See the file `License' ! in the root directory of the present distribution, ! or http://www.gnu.org/copyleft/gpl.txt . ! ! !---------------------------------------------------------------------------- SUBROUTINE read_file_xspectra(xread_wf) !---------------------------------------------------------------------------- ! ! ... This routine allocates space for all quantities already computed ! ... in the pwscf program and reads them from the data file. ! USE kinds, ONLY : DP USE ions_base, ONLY : nat, nsp, ityp, tau, if_pos USE basis, ONLY : natomwfc USE cell_base, ONLY : tpiba2, alat,omega, at, bg, ibrav USE force_mod, ONLY : force USE klist, ONLY : nkstot, nks, xk, wk, npk USE lsda_mod, ONLY : lsda, nspin, current_spin, isk USE wvfct, ONLY : nbnd, nbndx, et, wg, npwx USE symm_base, ONLY : irt, nsym, ftau, s, d1, d2, d3 USE ktetra, ONLY : tetra, ntetra USE extfield, ONLY : forcefield, tefield USE cellmd, ONLY : cell_factor, lmovecell USE gvect, ONLY : gg, ecutwfc, ngm, g, nr1, nr2, nr3, nrxx,& nrx1, nrx2, nrx3, eigts1, eigts2, eigts3, & nl, gstart USE gsmooth, ONLY : ngms, nls, nrx1s, nr1s, nr2s, nr3s ! USE spin_orb, ONLY : so, lspinorb ! CG USE spin_orb, ONLY : lspinorb USE scf, ONLY : rho, rho_core, rhog_core, v USE wavefunctions_module, ONLY : psic USE vlocal, ONLY : strf USE io_files, ONLY : tmp_dir, prefix, iunpun, nwordwfc, iunwfc USE buffers, ONLY : open_buffer, close_buffer USE uspp_param, ONLY : upf USE noncollin_module, ONLY : noncolin, npol USE mp_global, ONLY : kunit USE pw_restart, ONLY : pw_readfile USE xml_io_base, ONLY : pp_check_file USE uspp, ONLY : okvan, becsum USE paw_variables, ONLY : okpaw, ddd_PAW USE paw_onecenter, ONLY : paw_potential USE paw_init, ONLY : paw_init_onecenter, allocate_paw_internals USE ldaU, ONLY : eth USE dfunct, ONLY : newd ! USE paw_gipaw, ONLY : set_paw_upf ! ! IMPLICIT NONE ! INTEGER :: i, is, ik, ibnd, nb, nt, ios, isym, ierr REAL(DP) :: rdum(1,1), ehart, etxc, vtxc, etotefield, epaw, charge REAL(DP) :: sr(3,3,48) LOGICAL :: exst ! LOGICAL :: xread_wf ! ! ! ! ... first we check if the file can be used for post-processing ! IF ( .NOT. pp_check_file() ) & CALL infomsg( 'read_file', 'file ' // TRIM( tmp_dir ) // TRIM( prefix ) & & // '.save not guaranteed to be safe for post-processing' ) ! ! ... here we read the variables that dimension the system ! ... in parallel execution, only root proc reads the file ! ... and then broadcasts the values to all other procs ! ! ... a reset of the internal flags is necessary because some codes call ! ... read_file() more than once ! CALL pw_readfile( 'reset', ierr ) CALL pw_readfile( 'dim', ierr ) ! CALL errore( 'read_file ', 'problem reading file ' // & & TRIM( tmp_dir ) // TRIM( prefix ) // '.save', ierr ) ! ! ... allocate space for atomic positions, symmetries, forces, tetrahedra ! IF ( nat < 0 ) & CALL errore( 'read_file', 'wrong number of atoms', 1 ) ! ! ... allocation ! ALLOCATE( ityp( nat ) ) ! ALLOCATE( tau( 3, nat ) ) ALLOCATE( if_pos( 3, nat ) ) ALLOCATE( force( 3, nat ) ) ! IF ( tefield ) ALLOCATE( forcefield( 3, nat ) ) ! ALLOCATE( irt( 48, nat ) ) ALLOCATE( tetra( 4, MAX( ntetra, 1 ) ) ) ! ! ... here we read all the variables defining the system ! ... in parallel execution, only root proc read the file ! ... and then broadcast the values to all other procs ! !------------------------------------------------------------------------------- ! ... XML punch-file !------------------------------------------------------------------------------- ! CALL set_dimensions() ! ! ... check whether LSDA ! IF ( lsda ) THEN ! nspin = 2 npol = 1 ! ELSE IF ( noncolin ) THEN ! nspin = 4 npol = 2 current_spin = 1 ! ELSE ! nspin = 1 npol = 1 current_spin = 1 ! END IF ! cell_factor = 1.D0 lmovecell = .FALSE. ! ! ... allocate memory for eigenvalues and weights (read from file) ! nbndx = nbnd ALLOCATE( et( nbnd, nkstot ) , wg( nbnd, nkstot ) ) ! CALL pw_readfile( 'nowave', ierr ) ! ! ... distribute across pools k-points and related variables. ! ... nks is defined by the following routine as the number ! ... of k-points in the current pool ! CALL divide_et_impera( xk, wk, isk, lsda, nkstot, nks ) ! CALL poolscatter( nbnd, nkstot, et, nks, et ) CALL poolscatter( nbnd, nkstot, wg, nks, wg ) ! ! ... read pseudopotentials ! CALL pw_readfile( 'pseudo', ierr ) ! CALL readpp() ! ! DO nt = 1, nsp call set_paw_upf(nt, upf(nt)) ENDDO ! ! okvan = ANY ( upf(:)%tvanp ) okpaw = ANY ( upf(1:nsp)%tpawp ) ! ! ... check for spin-orbit pseudopotentials ! ! this seems to be obsolete ! DO nt = 1, nsp ! ! ! so(nt) = upf(nt)%has_so ! ! ! END DO ! ! ! ! IF ( .NOT. xread_wf ) THEN CALL read_k_points() ! nkstot=nks ! CALL divide_et_impera( xk, wk, isk, lsda, nkstot, nks ) DEALLOCATE ( et, wg ) ALLOCATE( et( nbnd, nkstot ) , wg( nbnd, nkstot ) ) END IF ! IF ( .NOT. lspinorb ) CALL average_pp ( nsp ) ! ! ... allocate memory for G- and R-space fft arrays ! CALL pre_init() CALL allocate_fft() CALL ggen() ! ! ... allocate the potential and wavefunctions ! CALL allocate_locpot() CALL allocate_nlpot() IF (okpaw) THEN CALL allocate_paw_internals() CALL paw_init_onecenter() CALL d_matrix(d1,d2,d3) ENDIF CALL allocate_wfc() ! ! ... read the charge density ! CALL pw_readfile( 'rho', ierr ) ! ! ... re-calculate the local part of the pseudopotential vltot ! ... and the core correction charge (if any) - This is done here ! ... for compatibility with the previous version of read_file ! CALL init_vloc() ! CALL struc_fact( nat, tau, nsp, ityp, ngm, g, bg, & nr1, nr2, nr3, strf, eigts1, eigts2, eigts3 ) ! CALL setlocal() ! CALL set_rhoc() ! ! ... bring rho to G-space ! DO is = 1, nspin ! psic(:) = rho%of_r(:,is) ! CALL cft3( psic, nr1, nr2, nr3, nrx1, nrx2, nrx3, -1 ) ! rho%of_g(:,is) = psic(nl(:)) ! END DO ! ! ... recalculate the potential ! CALL v_of_rho( rho, rho_core, rhog_core, & ehart, etxc, vtxc, eth, etotefield, charge, v ) ! ! ... reads the wavefunctions and writes them in 'distributed' form ! ... to unit iunwfc (for compatibility) ! ! if(xread_wf) then nwordwfc = nbnd*npwx*npol ! CALL open_buffer ( iunwfc, 'wfc', nwordwfc, nks, exst ) ! CALL pw_readfile( 'wave', ierr ) ! nks=nkstot if(lsda) then nks=nkstot/2 CALL set_kup_and_kdw( xk, wk, isk, nks, npk ) endif CALL divide_et_impera( xk, wk, isk, lsda, nkstot, nks ) endif ! CALL init_us_1() IF (okpaw) then CALL compute_becsum(1) CALL PAW_potential(rho%bec, ddd_PAW, epaw) ENDIF CALL newd() ! if(xread_wf) CALL close_buffer ( iunwfc, 'KEEP' ) ! ! RETURN ! CONTAINS ! !------------------------------------------------------------------------ SUBROUTINE set_dimensions() !------------------------------------------------------------------------ ! USE constants, ONLY : pi USE cell_base, ONLY : alat, tpiba, tpiba2 USE gvect, ONLY : ecutwfc, dual, gcutm USE gsmooth, ONLY : gcutms, doublegrid ! ! ! ... Set the units in real and reciprocal space ! tpiba = 2.D0 * pi / alat tpiba2 = tpiba**2 ! ! ... Compute the cut-off of the G vectors ! gcutm = dual * ecutwfc / tpiba2 ! doublegrid = ( dual > 4.D0 ) ! IF ( doublegrid ) THEN ! gcutms = 4.D0 * ecutwfc / tpiba2 ! ELSE ! gcutms = gcutm ! END IF ! END SUBROUTINE set_dimensions ! ! subroutine read_k_points USE ktetra, ONLY : nk1, nk2, nk3, k1, k2, k3 USE io_global, ONLY : ionode_id USE klist, ONLY : npk USE constants, ONLY : degspin USE parser, ONLY : read_line use mp, ONLY : mp_bcast implicit none INTEGER :: npool, nkl, nkr, nkbl, iks, ike CHARACTER(LEN=256) :: input_line INTEGER i,j,k,n ! Define new k-point mesh CALL read_line( input_line ) READ(input_line, *) nk1, nk2, nk3, k1, k2 ,k3 IF ( k1 < 0 .OR. k1 > 1 .OR. & k2 < 0 .OR. k2 > 1 .OR. & k3 < 0 .OR. k3 > 1 ) CALL errore & ('card_kpoints', 'invalid offsets: must be 0 or 1', 1) IF ( nk1 <= 0 .OR. nk2 <= 0 .OR. nk3 <= 0 ) CALL errore & ('card_kpoints', 'invalid values for nk1, nk2, nk3', 1) CALL mp_bcast( k1, ionode_id ) CALL mp_bcast( k2, ionode_id ) CALL mp_bcast( k3, ionode_id ) CALL mp_bcast( nk1, ionode_id ) CALL mp_bcast( nk2, ionode_id ) CALL mp_bcast( nk3, ionode_id ) ! if(nsym.eq.1) then nks=nk1*nk2*nk3 do i=1,nk1 do j=1,nk2 do k=1,nk3 ! this is nothing but consecutive ordering n = (k-1) + (j-1)*nk3 + (i-1)*nk2*nk3 + 1 ! xkg are the components of the complete grid in crystal axis xk(1,n) = DBLE(i-1)/nk1 + DBLE(k1)/2/nk1 xk(2,n) = DBLE(j-1)/nk2 + DBLE(k2)/2/nk2 xk(3,n) = DBLE(k-1)/nk3 + DBLE(k3)/2/nk3 end do end do end do wk(1:nks)=1.0/DBLE(nks) call cryst_to_cart(nks,xk,bg,1) ! else ! CALL kpoint_grid( nsym, s, bg, npk, k1, k2, k3, & ! nk1, nk2, nk3, nks, xk, wk ) ! endif if(lsda) then CALL set_kup_and_kdw( xk, wk, isk, nks, npk ) ELSE IF ( noncolin ) THEN call errore('define_and_distribute_k_points', & 'noncolinear not implemented', 1 ) else isk(1:nks)=1 wk(1:nks)=2.d0/dfloat(nks) endif end subroutine read_k_points ! END SUBROUTINE read_file_xspectra