mirror of https://gitlab.com/QEF/q-e.git
Missing variables initialization.
In some cases this can cause issues if the compiler does not set to 0 but to random values.
This commit is contained in:
parent
3af2945d8e
commit
2ad2c31a36
|
@ -50,7 +50,7 @@
|
|||
!! Set of symmetry operations
|
||||
INTEGER, INTENT(in) :: irt(48, nat)
|
||||
!! For each atoms give the rotated atoms
|
||||
REAL(KIND = DP), INTENT(inout) :: sxq(3, 48)
|
||||
REAL(KIND = DP), INTENT(in) :: sxq(3, 48)
|
||||
!! Symmetry matrix
|
||||
REAL(KIND = DP), INTENT(inout) :: rtau(3, 48, nat)
|
||||
!! the relative position of the rotated atom to the original one
|
||||
|
@ -132,8 +132,6 @@
|
|||
!!
|
||||
INTEGER :: nqs
|
||||
!!
|
||||
INTEGER :: axis
|
||||
!!
|
||||
INTEGER :: nrws
|
||||
!!
|
||||
INTEGER :: ierr
|
||||
|
@ -205,8 +203,8 @@
|
|||
COMPLEX(KIND = DP), ALLOCATABLE :: dyn(:, :, :, :) ! 3,3,nat,nat
|
||||
!! Dynamical matrix
|
||||
!
|
||||
axis = 3
|
||||
!
|
||||
q(:, :) = zero
|
||||
!
|
||||
! the call to set_ndnmbr is just a trick to get quickly
|
||||
! a file label by exploiting an existing subroutine
|
||||
! (if you look at the sub you will find that the original
|
||||
|
|
|
@ -415,6 +415,7 @@
|
|||
bmat(:, :, :, :) = czero
|
||||
cu(:, :, :) = czero
|
||||
cuq(:, :, :) = czero
|
||||
sxq(:, :) = zero
|
||||
!
|
||||
! read interatomic force constat matrix from q2r
|
||||
IF (lifc) THEN
|
||||
|
@ -496,7 +497,7 @@
|
|||
INQUIRE(FILE = TRIM(filename), EXIST = exst)
|
||||
IF (.NOT. exst) CALL errore('elphon_shuffle_wrap', &
|
||||
'cannot open file for reading or writing', ierr)
|
||||
CALL read_disp_pattern_only (iunpattern, filename, iq_irr, ierr)
|
||||
CALL read_disp_pattern_only(iunpattern, filename, iq_irr, ierr)
|
||||
IF (ierr /= 0) CALL errore('elphon_shuffle_wrap', ' Problem with modes file', 1)
|
||||
ENDIF
|
||||
!
|
||||
|
@ -533,7 +534,7 @@
|
|||
! ######################### star of q #########################
|
||||
!
|
||||
sym_smallq(:) = 0
|
||||
CALL star_q2(xq, at, bg, nsym, s, invs, nq, sxq, isq, imq, .TRUE., sym_smallq)
|
||||
CALL star_q2(xq, at, bg, nsym, s, invs, t_rev, nq, sxq, isq, imq, .TRUE., sym_smallq)
|
||||
IF (fixsym) THEN
|
||||
IF (epw_noinv) imq = 1 ! Any non-zero integer is ok.
|
||||
ENDIF
|
||||
|
@ -719,7 +720,7 @@
|
|||
! are equal to 5+ digits).
|
||||
! For any volunteers, please write to giustino@civet.berkeley.edu
|
||||
!
|
||||
CALL elphon_shuffle(iq_irr, nqc_irr, nqc, gmapsym(:,isym1), eigv(:,isym1), isym, xq0, .FALSE.)
|
||||
CALL elphon_shuffle(iq_irr, nqc_irr, nqc, gmapsym(:, isym1), eigv(:, isym1), isym, xq0, .FALSE.)
|
||||
!
|
||||
! bring epmatq in the mode representation of iq_first,
|
||||
! and then in the cartesian representation of iq
|
||||
|
@ -763,7 +764,7 @@
|
|||
!
|
||||
xq0 = -xq0
|
||||
!
|
||||
CALL elphon_shuffle(iq_irr, nqc_irr, nqc, gmapsym(:,isym1), eigv(:,isym1), isym, xq0, .TRUE.)
|
||||
CALL elphon_shuffle(iq_irr, nqc_irr, nqc, gmapsym(:, isym1), eigv(:, isym1), isym, xq0, .TRUE.)
|
||||
! bring epmatq in the mode representation of iq_first,
|
||||
! and then in the cartesian representation of iq
|
||||
!
|
||||
|
|
|
@ -436,7 +436,7 @@
|
|||
!---------------------------------------------------------------------------
|
||||
!
|
||||
!-----------------------------------------------------------------------
|
||||
SUBROUTINE star_q2(xq, at, bg, nsym, s, invs, nq, sxq, isq, imq, verbosity, sym_smallq)
|
||||
SUBROUTINE star_q2(xq, at, bg, nsym, s, invs, t_rev, nq, sxq, isq, imq, verbosity, sym_smallq)
|
||||
!-----------------------------------------------------------------------
|
||||
!!
|
||||
!! Generate the star of q vectors that are equivalent to the input one
|
||||
|
@ -457,6 +457,8 @@
|
|||
!! Symmetry operations
|
||||
INTEGER, INTENT(in) :: invs(48)
|
||||
!! list of inverse operation indices
|
||||
INTEGER, INTENT(in) :: t_rev(48)
|
||||
!! Time-reveral sym
|
||||
INTEGER, INTENT(out) :: nq
|
||||
!! degeneracy of the star of q
|
||||
INTEGER, INTENT(out) :: isq(48)
|
||||
|
@ -469,7 +471,7 @@
|
|||
!! direct lattice vectors
|
||||
REAL(KIND = DP), INTENT(in) :: bg(3, 3)
|
||||
!! reciprocal lattice vectors
|
||||
REAL(KIND = DP), INTENT(out) :: sxq(3, 48)
|
||||
REAL(KIND = DP), INTENT(inout) :: sxq(3, 48)
|
||||
!! list of vectors in the star of q
|
||||
!
|
||||
! Local variables
|
||||
|
@ -499,6 +501,7 @@
|
|||
!! Tolerence
|
||||
!
|
||||
zero(:) = 0.d0
|
||||
saq(:, :) = 0.0d0
|
||||
!
|
||||
! go to crystal coordinates
|
||||
!
|
||||
|
@ -525,6 +528,8 @@
|
|||
+ bg(i, 2) * raq(2) &
|
||||
+ bg(i, 3) * raq(3)
|
||||
ENDDO
|
||||
IF (t_rev(isym) == 1) raq = -raq
|
||||
!
|
||||
DO iq = 1, nq
|
||||
IF (eqvect(raq, saq(1, iq), zero, accep)) THEN
|
||||
isq(isym) = iq
|
||||
|
|
|
@ -57,6 +57,7 @@ subroutine star_q1(xq, at, bg, nsym, s, invs, nq, sxq, isq, imq, verbosity,&
|
|||
!
|
||||
!
|
||||
zero(:) = 0.d0
|
||||
saq(:, :) = 0.0d0
|
||||
!
|
||||
! go to crystal coordinates
|
||||
!
|
||||
|
@ -176,6 +177,7 @@ subroutine star_q (xq, at, bg, nsym, s, invs, nq, sxq, isq, imq, verbosity)
|
|||
!
|
||||
!
|
||||
zero(:) = 0.d0
|
||||
saq(:, :) = 0.0d0
|
||||
!
|
||||
! go to crystal coordinates
|
||||
!
|
||||
|
|
Loading…
Reference in New Issue