Cleanup. Repeated software moved to separate routines.

git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@9817 c92efa57-630b-4861-b058-cf58834340f0
This commit is contained in:
dalcorso 2013-01-19 17:20:23 +00:00
parent 545f2fa036
commit 9513510809
9 changed files with 153 additions and 106 deletions

View File

@ -129,6 +129,7 @@ random_matrix.o \
read_wfc_rspace_and_fwfft.o \
rotate_dvscf_star.o \
rotate_and_add_dyn.o \
rotate_pattern_add.o \
run_pwscf.o \
save_ph_input.o \
set_asr_c.o \

View File

@ -232,18 +232,9 @@ subroutine d2ionq (nat, ntyp, ityp, zv, tau, alat, omega, q, at, &
! The dynamical matrix was computed in cartesian axis and now we put
! it on the basis of the modes
!
do nu_i = 1, nmodes
do nu_j = 1, nmodes
work = (0.d0, 0.d0)
do nb_jcart = 1, 3 * nat
do na_icart = 1, 3 * nat
work = work + CONJG(u (na_icart, nu_i) ) * &
dy3 (na_icart, nb_jcart) * u (nb_jcart, nu_j)
enddo
enddo
dyn (nu_i, nu_j) = dyn (nu_i, nu_j) - work
enddo
enddo
dy3 = -dy3
!
CALL rotate_pattern_add(nat, u, dyn, dy3)
!
call stop_clock ('d2ionq')
return

View File

@ -69,20 +69,13 @@ subroutine dynmat0_new
! rotate again in the pattern basis
!
call zcopy (9 * nat * nat, dyn, 1, dynwrk, 1)
do nu_i = 1, 3 * nat
do nu_j = 1, 3 * nat
wrk = (0.d0, 0.d0)
do nb_jcart = 1, 3 * nat
do na_icart = 1, 3 * nat
wrk = wrk + CONJG(u (na_icart, nu_i) ) * &
dynwrk (na_icart, nb_jcart) * &
u (nb_jcart, nu_j)
enddo
enddo
dyn (nu_i, nu_j) = wrk
enddo
enddo
dyn=(0.d0, 0.d0)
CALL rotate_pattern_add(nat, u, dyn, dynwrk)
endif
! call tra_write_matrix('dynmat0 dyn',dyn,u,nat)
dyn_rec(:,:)=dyn(:,:)
done_irr(0) = 1

View File

@ -253,20 +253,9 @@ SUBROUTINE dynmat_us()
!
! We rotate the dynamical matrix on the basis of patterns
!
DO nu_i = 1, 3 * nat
DO nu_j = 1, 3 * nat
work = (0.0d0, 0.0d0)
DO na_jcart = 1, 3 * nat
DO na_icart = 1, 3 * nat
work = work + CONJG (u (na_icart, nu_i) ) * &
dynwrk (na_icart, na_jcart) * &
u (na_jcart, nu_j)
ENDDO
ENDDO
dyn (nu_i, nu_j) = dyn (nu_i, nu_j) + work
ENDDO
ENDDO
CALL rotate_pattern_add(nat, u, dyn, dynwrk)
IF (noncolin) THEN
DEALLOCATE (deff_nc)
ELSE

View File

@ -108,18 +108,7 @@ subroutine dynmatcc
!
! rotate in the pattern basis and add to dynmat
!
do nu_i = 1, 3 * nat
do nu_j = 1, 3 * nat
wrk = (0.d0, 0.d0)
do nb_jcart = 1, 3 * nat
do na_icart = 1, 3 * nat
wrk = wrk + CONJG(u (na_icart, nu_i) ) * dynwrk (na_icart, &
nb_jcart) * u (nb_jcart, nu_j)
enddo
enddo
dyn (nu_i, nu_j) = dyn (nu_i, nu_j) + wrk
enddo
enddo
CALL rotate_pattern_add(nat, u, dyn, dynwrk)
!
deallocate (work)
deallocate (vxc)

View File

@ -54,17 +54,9 @@ subroutine q2qstar_ph (dyn, at, bg, nat, nsym, s, invs, irt, rtau, &
nsq = nsym / nq
if (nsq * nq /= nsym) call errore ('q2star_ph', 'wrong degeneracy', 1)
!
! Writes dyn.mat. dyn(3*nat,3*nat) on the 4-index array phi(3,3,nat,nta)
! Writes dyn.mat. dyn(3*nat,3*nat) on the 4-index array phi(3,3,nat,nat)
!
do i = 1, 3 * nat
na = (i - 1) / 3 + 1
icar = i - 3 * (na - 1)
do j = 1, 3 * nat
nb = (j - 1) / 3 + 1
jcar = j - 3 * (nb - 1)
phi (icar, jcar, na, nb) = dyn (i, j)
enddo
enddo
CALL scompact_dyn(nat, dyn, phi)
!
! Go to crystal coordinates
!
@ -105,17 +97,9 @@ subroutine q2qstar_ph (dyn, at, bg, nat, nsym, s, invs, irt, rtau, &
enddo
enddo
!
! Saves 4-index array phi(3,3,nat,nta) on the dyn.mat. dyn(3*nat,3*nat)
! Saves 4-index array phi2(3,3,nat,nat) on the dyn.mat. dyn(3*nat,3*nat)
!
do i = 1, 3 * nat
na = (i - 1) / 3 + 1
icar = i - 3 * (na - 1)
do j = 1, 3 * nat
nb = (j - 1) / 3 + 1
jcar = j - 3 * (nb - 1)
dyn (i, j) = phi2 (icar, jcar, na, nb)
enddo
enddo
CALL compact_dyn(nat, dyn, phi2)
endif
!
! For each q of the star rotates phi with the appropriate sym.op. -> phi

View File

@ -0,0 +1,132 @@
!
! Copyright (C) 2012 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 rotate_pattern_add(nat, u, dyn, dynwrk)
! This routine rotates the dynamical matrix dynwork written
! in cartesian basis to the basis of the patterns u and adds it to
! the dynamical matrix dyn that is supposed to be in the basis of the
! patterns.
USE kinds, ONLY : DP
IMPLICIT NONE
INTEGER, INTENT(IN) :: nat
COMPLEX(DP), INTENT(IN) :: u(3*nat, 3*nat)
COMPLEX(DP), INTENT(INOUT) :: dyn(3*nat, 3*nat)
COMPLEX(DP), INTENT(IN) :: dynwrk(3*nat, 3*nat)
COMPLEX(DP) :: work
INTEGER :: nu_i, nu_j, na_icart, na_jcart
DO nu_i = 1, 3 * nat
DO nu_j = 1, 3 * nat
work = (0.0d0, 0.0d0)
DO na_jcart = 1, 3 * nat
DO na_icart = 1, 3 * nat
work = work + CONJG (u (na_icart, nu_i) ) * &
dynwrk (na_icart, na_jcart) * &
u (na_jcart, nu_j)
ENDDO
ENDDO
dyn (nu_i, nu_j) = dyn (nu_i, nu_j) + work
ENDDO
ENDDO
RETURN
END SUBROUTINE rotate_pattern_add
!
!----------------------------------------------------------------------
SUBROUTINE dyn_pattern_to_cart(nat, u, dyn, phi)
! This routine transforms the dynamical matrix dyn, written in the basis
! of the pattern, in the dynamical matrix phi, in the cartesian basis.
USE kinds, ONLY : DP
IMPLICIT NONE
INTEGER, INTENT(IN) :: nat
COMPLEX(DP), INTENT(IN) :: u(3*nat, 3*nat)
COMPLEX(DP), INTENT(IN) :: dyn(3*nat, 3*nat)
COMPLEX(DP), INTENT(OUT) :: phi(3, 3, nat, nat)
COMPLEX(DP) :: work
INTEGER :: i, j, icart, jcart, na, nb, mu, nu
DO i = 1, 3 * nat
na = (i - 1) / 3 + 1
icart = i - 3 * (na - 1)
DO j = 1, 3 * nat
nb = (j - 1) / 3 + 1
jcart = j - 3 * (nb - 1)
work = (0.d0, 0.d0)
DO mu = 1, 3 * nat
DO nu = 1, 3 * nat
work = work + u (i, mu) * dyn (mu, nu) * CONJG(u (j, nu) )
ENDDO
ENDDO
phi (icart, jcart, na, nb) = work
ENDDO
ENDDO
RETURN
END SUBROUTINE dyn_pattern_to_cart
!
!-----------------------------------------------------------------------
SUBROUTINE compact_dyn(nat, dyn, phi)
!-----------------------------------------------------------------------
!
! This routine writes the dynamical matrix from a 3,3,nat,nat array
! to a 3*nat,3*nat array.
!
USE kinds, ONLY : DP
IMPLICIT NONE
INTEGER, INTENT(IN) :: nat
COMPLEX(DP), INTENT(IN) :: phi(3,3,nat,nat)
COMPLEX(DP), INTENT(OUT) :: dyn(3*nat, 3*nat)
INTEGER :: na, nb, icart, jcart, imode, jmode
DO na = 1, nat
DO icart = 1, 3
imode = 3 * ( na - 1 ) + icart
DO nb = 1, nat
DO jcart = 1, 3
jmode = 3 * ( nb - 1 ) + jcart
dyn (imode, jmode) = phi (icart, jcart, na, nb)
END DO
END DO
END DO
END DO
RETURN
END SUBROUTINE compact_dyn
!-----------------------------------------------------------------------
SUBROUTINE scompact_dyn(nat, dyn, phi)
!-----------------------------------------------------------------------
!
! This routine writes the dynamical matrix from a 3*nat,3*nat array
! to a 3,3,nat,nat array.
!
USE kinds, ONLY : DP
IMPLICIT NONE
INTEGER, INTENT(IN) :: nat
COMPLEX(DP), INTENT(OUT) :: phi(3,3,nat,nat)
COMPLEX(DP), INTENT(IN) :: dyn(3*nat, 3*nat)
INTEGER :: na, nb, icart, jcart, imode, jmode
DO na = 1, nat
DO icart = 1, 3
imode = 3 * ( na - 1 ) + icart
DO nb = 1, nat
DO jcart = 1, 3
jmode = 3 * ( nb - 1 ) + jcart
phi (icart, jcart, na, nb) = dyn (imode, jmode)
END DO
END DO
END DO
END DO
RETURN
END SUBROUTINE scompact_dyn

View File

@ -110,17 +110,7 @@ subroutine set_irr_new (xq, u, npert, nirr, eigen)
!
! We copy the dynamical matrix in a bidimensional array
!
do na = 1, nat
do nb = 1, nat
do ipol = 1, 3
imode = ipol + 3 * (na - 1)
do jpol = 1, 3
jmode = jpol + 3 * (nb - 1)
phi (imode, jmode) = wdyn (ipol, jpol, na, nb)
enddo
enddo
enddo
enddo
CALL compact_dyn(nat, phi, wdyn)
!
! Diagonalize
!

View File

@ -55,21 +55,7 @@ subroutine symdyn_munu_new (dyn, u, xq, s, invs, rtau, irt, at, &
!
! First we transform in the cartesian coordinates
!
do i = 1, 3 * nat
na = (i - 1) / 3 + 1
icart = i - 3 * (na - 1)
do j = 1, 3 * nat
nb = (j - 1) / 3 + 1
jcart = j - 3 * (nb - 1)
work = (0.d0, 0.d0)
do mu = 1, 3 * nat
do nu = 1, 3 * nat
work = work + u (i, mu) * dyn (mu, nu) * CONJG(u (j, nu) )
enddo
enddo
phi (icart, jcart, na, nb) = work
enddo
enddo
CALL dyn_pattern_to_cart(nat, u, dyn, phi)
!
! Then we transform to the crystal axis
!
@ -92,17 +78,9 @@ subroutine symdyn_munu_new (dyn, u, xq, s, invs, rtau, irt, at, &
enddo
enddo
!
! rewrite the dynamical matrix on the array dyn with dimension 3nat x 3
! rewrite the dynamical matrix on the array dyn with dimension 3nat x 3nat
!
do i = 1, 3 * nat
na = (i - 1) / 3 + 1
icart = i - 3 * (na - 1)
do j = 1, 3 * nat
nb = (j - 1) / 3 + 1
jcart = j - 3 * (nb - 1)
dyn (i, j) = phi (icart, jcart, na, nb)
enddo
CALL compact_dyn(nat, dyn, phi)
enddo
return
end subroutine symdyn_munu_new