2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! Copyright (C) 2001 PWSCF 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 .
|
|
|
|
!
|
2004-10-26 17:32:48 +08:00
|
|
|
#include "f_defs.h"
|
2005-05-25 10:58:35 +08:00
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2003-02-08 00:04:36 +08:00
|
|
|
subroutine d3_summary
|
2003-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! This routine writes on output the quantities which have been read
|
|
|
|
! from the punch file, and the quantities computed in the d3_setup
|
|
|
|
! file.
|
|
|
|
!
|
|
|
|
! if iverbosity = 0 only a partial summary is done.
|
|
|
|
!
|
2006-01-16 04:18:53 +08:00
|
|
|
USE ions_base, ONLY : nat, ityp, ntyp => nsp, atm, tau, amass
|
2003-11-06 21:06:49 +08:00
|
|
|
USE io_global, ONLY : stdout
|
2004-02-04 19:06:45 +08:00
|
|
|
USE kinds, only : DP
|
2003-01-20 05:58:50 +08:00
|
|
|
use pwcom
|
2004-06-01 01:55:33 +08:00
|
|
|
USE uspp_param, ONLY : rinner, nqlc, nqf, lll, nbeta, psd, iver, tvanp
|
2004-04-28 18:25:36 +08:00
|
|
|
USE atom, ONLY: mesh, numeric, xmin, dx, nlcc
|
2004-04-02 21:01:15 +08:00
|
|
|
USE control_flags, ONLY : iverbosity
|
2003-01-20 05:58:50 +08:00
|
|
|
use phcom
|
|
|
|
use d3com
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
integer :: i, l, nt, mu, nu, ipol, apol, na, isymq, isym, nsymtot, &
|
|
|
|
ik, ib, irr, imode0, iaux
|
|
|
|
! generic counter
|
|
|
|
! counter on angular momenta
|
|
|
|
! counter on atomic types
|
|
|
|
! counter on modes
|
|
|
|
! counter on modes
|
|
|
|
! counter on polarizations
|
|
|
|
! counter on polarizations
|
|
|
|
! counter on atoms
|
|
|
|
! counter on symmetries
|
|
|
|
! counter on symmetries
|
|
|
|
! counter on symmetries
|
|
|
|
! counter on k points
|
|
|
|
! counter on beta functions
|
|
|
|
! counter on irreducible representation
|
|
|
|
! the first mode
|
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
real (DP) :: ft1, ft2, ft3, sr (3, 3), xkg (3)
|
2003-01-20 05:58:50 +08:00
|
|
|
! fractionary translation
|
|
|
|
! fractionary translation
|
|
|
|
! fractionary translation
|
|
|
|
! the symmetry matrix in cartesian coord
|
|
|
|
! k point in crystal coordinates
|
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
character :: ps * 5
|
2003-01-20 05:58:50 +08:00
|
|
|
! the name of the pseudo
|
2006-04-06 20:13:54 +08:00
|
|
|
WRITE( stdout, 100) title, crystal, ibrav, alat, omega, nat, ntyp, &
|
2003-01-20 05:58:50 +08:00
|
|
|
ecutwfc, ecutwfc * dual
|
|
|
|
100 format (/,5x,a75,/,/,5x, 'crystal is ',a20,/,/,5x, &
|
|
|
|
& 'bravais-lattice index = ',i12,/,5x, &
|
|
|
|
& 'lattice parameter (a_0) = ',f12.4,' a.u.',/,5x, &
|
|
|
|
& 'unit-cell volume = ',f12.4,' (a.u.)^3',/,5x, &
|
|
|
|
& 'number of atoms/cell = ',i12,/,5x, &
|
|
|
|
& 'number of atomic types = ',i12,/,5x, &
|
|
|
|
& 'kinetic-energy cut-off = ',f12.4,' Ry',/,5x, &
|
|
|
|
& 'charge densisty cut-off = ',f12.4,' Ry',/,5x,/)
|
|
|
|
!
|
2003-06-14 00:55:38 +08:00
|
|
|
! and here more detailed information. Description of the unit cell
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(2(3x,3(2x,"celldm(",i1,")=",f11.5),/))') (i, &
|
2003-01-20 05:58:50 +08:00
|
|
|
celldm (i) , i = 1, 6)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x, &
|
2003-01-20 05:58:50 +08:00
|
|
|
& "crystal axes: (cart. coord. in units of a_0)",/, &
|
|
|
|
& 3(15x,"a(",i1,") = (",3f8.4," ) ",/ ) )') (apol, &
|
|
|
|
& (at (ipol, apol) , ipol = 1, 3) , apol = 1, 3)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x, &
|
2003-01-20 05:58:50 +08:00
|
|
|
&"reciprocal axes: (cart. coord. in units 2 pi/a_0)",/, &
|
|
|
|
& 3(15x,"b(",i1,") = (",3f8.4," ) ",/ ) )') (apol, &
|
|
|
|
& (bg (ipol, apol) , ipol = 1, 3) , apol = 1, 3)
|
|
|
|
!
|
|
|
|
! description of the atoms inside the unit cell
|
|
|
|
!
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/, 5x,"Atoms inside the unit cell: ")')
|
|
|
|
WRITE( stdout, '(/,3x,"Cartesian axes")')
|
|
|
|
WRITE( stdout, '(/,5x,"site n. atom mass ", &
|
2003-01-20 05:58:50 +08:00
|
|
|
& " positions (a_0 units)")')
|
|
|
|
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(7x,i2,5x,a6,f8.4," tau(",i2, ") = (",3f11.5," )")') &
|
2003-01-20 05:58:50 +08:00
|
|
|
(na, atm (ityp (na) ) , amass (ityp (na) ) / amconv, na, &
|
|
|
|
(tau (ipol, na ) , ipol = 1, 3) , na = 1, nat)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/,5x,"Computing dynamical matrix for ")')
|
|
|
|
WRITE( stdout, '(20x,"q = (",3f11.5," )")') (xq (ipol) , ipol = 1, 3)
|
2003-02-08 00:04:36 +08:00
|
|
|
if (q0mode_todo (1) .le.0) then
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/,5x,"Computing all the modes ")')
|
2003-02-08 00:04:36 +08:00
|
|
|
else
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/,5x,"Computing only selected modes: ")')
|
2003-02-08 00:04:36 +08:00
|
|
|
do i = 1, 3 * nat
|
2003-11-06 21:06:49 +08:00
|
|
|
if (q0mode (i) ) WRITE( stdout, '(5x,"Mode to be computed: ",i5)') i
|
2003-01-20 05:58:50 +08:00
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
!
|
|
|
|
! description of symmetries
|
|
|
|
!
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, * )
|
2003-02-08 00:04:36 +08:00
|
|
|
if (nsymg0.le.1) then
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,"No symmetry! for q=0 ")')
|
2003-02-08 00:04:36 +08:00
|
|
|
else
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,i2," + 1 = ",i2," q=0 Sym.Ops. ",/)') &
|
2003-01-20 05:58:50 +08:00
|
|
|
nsymg0, nsymg0 + 1
|
|
|
|
endif
|
2003-02-08 00:04:36 +08:00
|
|
|
if (.not.lgamma) then
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, * )
|
2003-02-08 00:04:36 +08:00
|
|
|
if (nsymq.le.1.and..not.minus_q) then
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,"No symmetry!")')
|
2003-02-08 00:04:36 +08:00
|
|
|
else
|
|
|
|
if (minus_q) then
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,i2," Sym.Ops. (with q -> -q+G )",/)') &
|
2003-01-20 05:58:50 +08:00
|
|
|
nsymq + 1
|
2003-02-08 00:04:36 +08:00
|
|
|
else
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,i2," Sym.Ops. (no q -> -q+G )",/)') &
|
2003-01-20 05:58:50 +08:00
|
|
|
nsymq
|
|
|
|
endif
|
|
|
|
endif
|
|
|
|
endif
|
2003-02-08 00:04:36 +08:00
|
|
|
if (iverbosity.eq.1) then
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(36x,"s",24x,"frac. trans.")')
|
2003-02-08 00:04:36 +08:00
|
|
|
if (minus_q) then
|
|
|
|
iaux = 0
|
|
|
|
else
|
|
|
|
iaux = 1
|
2003-01-20 05:58:50 +08:00
|
|
|
endif
|
2003-02-08 00:04:36 +08:00
|
|
|
do isymq = iaux, nsymg0
|
|
|
|
if (isymq.eq.0) then
|
|
|
|
isym = irotmq
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/,5x,"This transformation sends q -> -q+G")')
|
2003-02-08 00:04:36 +08:00
|
|
|
else
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
! It seems to me variable irgq is useless !
|
|
|
|
! isym = irgq(isymq)
|
2003-02-08 00:04:36 +08:00
|
|
|
isym = isymq
|
2003-01-20 05:58:50 +08:00
|
|
|
endif
|
2003-02-08 00:04:36 +08:00
|
|
|
if (isymq.eq.nsymq + 1) then
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(//,5x,&
|
2003-01-20 05:58:50 +08:00
|
|
|
&"In the following are listed symmetries of the crystal")')
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,"not belonging to the small group of q")')
|
2003-01-20 05:58:50 +08:00
|
|
|
endif
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/6x,"isym = ",i2,5x,a45/)') isymq, sname (isym)
|
2003-02-08 00:04:36 +08:00
|
|
|
call s_axis_to_cart (s (1, 1, isym), sr, at, bg)
|
2003-01-20 05:58:50 +08:00
|
|
|
if (ftau (1, isym) .ne.0.or.ftau (2, isym) .ne.0.or.ftau (3, &
|
|
|
|
isym) .ne.0) then
|
|
|
|
ft1 = at (1, 1) * ftau (1, isym) / nr1 + at (1, 2) * ftau ( &
|
|
|
|
2, isym) / nr2 + at (1, 3) * ftau (3, isym) / nr3
|
|
|
|
ft2 = at (2, 1) * ftau (1, isym) / nr1 + at (2, 2) * ftau ( &
|
|
|
|
2, isym) / nr2 + at (2, 3) * ftau (3, isym) / nr3
|
|
|
|
ft3 = at (3, 1) * ftau (1, isym) / nr1 + at (3, 2) * ftau ( &
|
|
|
|
2, isym) / nr2 + at (3, 3) * ftau (3, isym) / nr3
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(1x,"cryst.",3x,"s(",i2,") = (",3(i6,5x), &
|
2003-01-20 05:58:50 +08:00
|
|
|
&" ) f =( ",f10.7," )")') isymq, (s (1, ipol, isym),&
|
General cleanup of intrinsic functions:
conversion to real => DBLE
(including real part of a complex number)
conversion to complex => CMPLX
complex conjugate => CONJG
imaginary part => AIMAG
All functions are uppercase.
CMPLX is preprocessed by f_defs.h and performs an explicit cast:
#define CMPLX(a,b) cmplx(a,b,kind=DP)
This implies that 1) f_defs.h must be included whenever a CMPLX is present,
2) CMPLX should stay in a single line, 3) DP must be defined.
All occurrences of real, float, dreal, dfloat, dconjg, dimag, dcmplx
removed - please do not reintroduce any of them.
Tested only with ifc7 and g95 - beware unintended side effects
Maybe not the best solution (explicit casts everywhere would be better)
but it can be easily changed with a script if the need arises.
The following code might be used to test for possible trouble:
program test_intrinsic
implicit none
integer, parameter :: dp = selected_real_kind(14,200)
real (kind=dp) :: a = 0.123456789012345_dp
real (kind=dp) :: b = 0.987654321098765_dp
complex (kind=dp) :: c = ( 0.123456789012345_dp, 0.987654321098765_dp)
print *, ' A = ', a
print *, ' DBLE(A)= ', DBLE(a)
print *, ' C = ', c
print *, 'CONJG(C)= ', CONJG(c)
print *, 'DBLE(c),AIMAG(C) = ', DBLE(c), AIMAG(c)
print *, 'CMPLX(A,B,kind=dp)= ', CMPLX( a, b, kind=dp)
end program test_intrinsic
Note that CMPLX and REAL without a cast yield single precision numbers on
ifc7 and g95 !!!
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@2133 c92efa57-630b-4861-b058-cf58834340f0
2005-08-27 01:44:42 +08:00
|
|
|
ipol = 1, 3) , DBLE (ftau (1, isym) ) / DBLE (nr1)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(17x," (",3(i6,5x), &
|
2003-01-20 05:58:50 +08:00
|
|
|
& " ) ( ",f10.7," )")') &
|
|
|
|
(s (2, ipol, &
|
General cleanup of intrinsic functions:
conversion to real => DBLE
(including real part of a complex number)
conversion to complex => CMPLX
complex conjugate => CONJG
imaginary part => AIMAG
All functions are uppercase.
CMPLX is preprocessed by f_defs.h and performs an explicit cast:
#define CMPLX(a,b) cmplx(a,b,kind=DP)
This implies that 1) f_defs.h must be included whenever a CMPLX is present,
2) CMPLX should stay in a single line, 3) DP must be defined.
All occurrences of real, float, dreal, dfloat, dconjg, dimag, dcmplx
removed - please do not reintroduce any of them.
Tested only with ifc7 and g95 - beware unintended side effects
Maybe not the best solution (explicit casts everywhere would be better)
but it can be easily changed with a script if the need arises.
The following code might be used to test for possible trouble:
program test_intrinsic
implicit none
integer, parameter :: dp = selected_real_kind(14,200)
real (kind=dp) :: a = 0.123456789012345_dp
real (kind=dp) :: b = 0.987654321098765_dp
complex (kind=dp) :: c = ( 0.123456789012345_dp, 0.987654321098765_dp)
print *, ' A = ', a
print *, ' DBLE(A)= ', DBLE(a)
print *, ' C = ', c
print *, 'CONJG(C)= ', CONJG(c)
print *, 'DBLE(c),AIMAG(C) = ', DBLE(c), AIMAG(c)
print *, 'CMPLX(A,B,kind=dp)= ', CMPLX( a, b, kind=dp)
end program test_intrinsic
Note that CMPLX and REAL without a cast yield single precision numbers on
ifc7 and g95 !!!
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@2133 c92efa57-630b-4861-b058-cf58834340f0
2005-08-27 01:44:42 +08:00
|
|
|
&isym) , ipol = 1, 3) , DBLE (ftau (2, isym) ) / DBLE (nr2)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(17x," (",3(i6,5x), &
|
2003-01-20 05:58:50 +08:00
|
|
|
& " ) ( ",f10.7," )"/)') (s (3, ipol, &
|
General cleanup of intrinsic functions:
conversion to real => DBLE
(including real part of a complex number)
conversion to complex => CMPLX
complex conjugate => CONJG
imaginary part => AIMAG
All functions are uppercase.
CMPLX is preprocessed by f_defs.h and performs an explicit cast:
#define CMPLX(a,b) cmplx(a,b,kind=DP)
This implies that 1) f_defs.h must be included whenever a CMPLX is present,
2) CMPLX should stay in a single line, 3) DP must be defined.
All occurrences of real, float, dreal, dfloat, dconjg, dimag, dcmplx
removed - please do not reintroduce any of them.
Tested only with ifc7 and g95 - beware unintended side effects
Maybe not the best solution (explicit casts everywhere would be better)
but it can be easily changed with a script if the need arises.
The following code might be used to test for possible trouble:
program test_intrinsic
implicit none
integer, parameter :: dp = selected_real_kind(14,200)
real (kind=dp) :: a = 0.123456789012345_dp
real (kind=dp) :: b = 0.987654321098765_dp
complex (kind=dp) :: c = ( 0.123456789012345_dp, 0.987654321098765_dp)
print *, ' A = ', a
print *, ' DBLE(A)= ', DBLE(a)
print *, ' C = ', c
print *, 'CONJG(C)= ', CONJG(c)
print *, 'DBLE(c),AIMAG(C) = ', DBLE(c), AIMAG(c)
print *, 'CMPLX(A,B,kind=dp)= ', CMPLX( a, b, kind=dp)
end program test_intrinsic
Note that CMPLX and REAL without a cast yield single precision numbers on
ifc7 and g95 !!!
git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@2133 c92efa57-630b-4861-b058-cf58834340f0
2005-08-27 01:44:42 +08:00
|
|
|
& isym) , ipol = 1, 3) , DBLE (ftau (3, isym) ) / DBLE (nr3)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(1x,"cart.",3x,"s(",i2,") = (",3f11.7, &
|
2003-01-20 05:58:50 +08:00
|
|
|
& " ) f =( ",f10.7," )")') isymq, (sr (1 &
|
|
|
|
&, ipol) , ipol = 1, 3) , ft1
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(17x," (",3f11.7, &
|
2003-01-20 05:58:50 +08:00
|
|
|
& " ) ( ",f10.7," )")') (sr (2, ipol) &
|
|
|
|
& , ipol = 1, 3) , ft2
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(17x," (",3f11.7, &
|
2003-01-20 05:58:50 +08:00
|
|
|
& " ) ( ",f10.7," )"/)') (sr (3, ipol &
|
|
|
|
&) , ipol = 1, 3) , ft3
|
2003-02-08 00:04:36 +08:00
|
|
|
else
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(1x,"cryst.",3x,"s(",i2,") = (",3(i6,5x), &
|
2003-01-20 05:58:50 +08:00
|
|
|
& " )")') isymq, (s (1, ipol, isym) , ipol = &
|
|
|
|
&1, 3)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(17x," (",3(i6,5x)," )")') (s (2, ipol, isym) &
|
2003-01-20 05:58:50 +08:00
|
|
|
, ipol = 1, 3)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(17x," (",3(i6,5x)," )"/)') (s (3, ipol, &
|
2003-01-20 05:58:50 +08:00
|
|
|
isym) , ipol = 1, 3)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(1x,"cart.",3x,"s(",i2,") = (",3f11.7, &
|
2003-01-20 05:58:50 +08:00
|
|
|
& " )")') isymq, (sr (1, ipol) , ipol = 1, 3)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(17x," (",3f11.7," )")') (sr (2, ipol) , &
|
2003-01-20 05:58:50 +08:00
|
|
|
ipol = 1, 3)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(17x," (",3f11.7," )"/)') (sr (3, ipol) , &
|
2003-01-20 05:58:50 +08:00
|
|
|
ipol = 1, 3)
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
!
|
|
|
|
! Description of the reciprocal lattice vectors
|
|
|
|
!
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/5x,"G cutoff =",f10.4," (", &
|
2003-01-20 05:58:50 +08:00
|
|
|
& i7," G-vectors)"," FFT grid: (",i3, &
|
|
|
|
& ",",i3,",",i3,")")') gcutm, ngm, nr1, nr2, nr3
|
|
|
|
|
2003-11-06 21:06:49 +08:00
|
|
|
if (doublegrid) WRITE( stdout, '(5x,"G cutoff =",f10.4," (", &
|
2003-01-20 05:58:50 +08:00
|
|
|
& i7," G-vectors)"," smooth grid: (",i3, &
|
|
|
|
& ",",i3,",",i3,")")') gcutms, ngms, &
|
|
|
|
&nr1s, nr2s, nr3s
|
2003-02-08 00:04:36 +08:00
|
|
|
if (degauss.eq.0.d0) then
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,"number of k points=",i5)') nkstot
|
2003-02-08 00:04:36 +08:00
|
|
|
else
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,"number of k points=",i5, &
|
2006-10-23 20:32:54 +08:00
|
|
|
& " gaussian broad. (Ry)=",f8.4,5x, &
|
2003-01-20 05:58:50 +08:00
|
|
|
& "ngauss = ",i3)') nkstot, degauss, ngauss
|
|
|
|
endif
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(23x,"cart. coord. in units 2pi/a_0")')
|
2003-02-08 00:04:36 +08:00
|
|
|
do ik = 1, nkstot
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(8x,"k(",i5,") = (",3f12.7,"), wk =",f12.7)') ik, &
|
2003-01-20 05:58:50 +08:00
|
|
|
(xk (ipol, ik) , ipol = 1, 3) , wk (ik)
|
|
|
|
enddo
|
2003-02-08 00:04:36 +08:00
|
|
|
if (iverbosity.eq.1) then
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/23x,"cryst. coord.")')
|
2003-02-08 00:04:36 +08:00
|
|
|
do ik = 1, nkstot
|
|
|
|
do ipol = 1, 3
|
2003-01-20 05:58:50 +08:00
|
|
|
! xkg are the compone
|
|
|
|
xkg (ipol) = at (1, ipol) * xk (1, ik) + at (2, ipol) * xk (2, &
|
|
|
|
ik) + at (3, ipol) * xk (3, ik)
|
|
|
|
! of xk in the crysta
|
|
|
|
! rec. lattice basis
|
|
|
|
enddo
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(8x,"k(",i5,") = (",3f12.7,"), wk =",f12.7)') &
|
2003-01-20 05:58:50 +08:00
|
|
|
ik, (xkg (ipol) , ipol = 1, 3) , wk (ik)
|
|
|
|
enddo
|
|
|
|
|
|
|
|
endif
|
2003-02-08 00:04:36 +08:00
|
|
|
do nt = 1, ntyp
|
|
|
|
if (tvanp (nt) ) then
|
|
|
|
ps = '(US)'
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/5x,"pseudo",i2," is ",a2, &
|
2003-01-20 05:58:50 +08:00
|
|
|
& 1x,a5," zval =",f5.1," lmax=",i2, &
|
|
|
|
& " lloc=",i2)') nt, psd (nt) , ps, zp (nt) , lmax (nt) &
|
|
|
|
&, lloc (nt)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,"Version ", 3i3, " of US pseudo code")') &
|
2003-01-20 05:58:50 +08:00
|
|
|
(iver (i, nt) , i = 1, 3)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/,5x,"Using log mesh of ", i3, " points",/)') &
|
2003-01-20 05:58:50 +08:00
|
|
|
mesh (nt)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,"The pseudopotential has ",i2, &
|
2003-01-20 05:58:50 +08:00
|
|
|
& " beta functions with: ",/)') nbeta (nt)
|
2003-02-08 00:04:36 +08:00
|
|
|
do ib = 1, nbeta (nt)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(15x," l(",i1,") = ",i3)') ib, lll (ib, nt)
|
2003-01-20 05:58:50 +08:00
|
|
|
|
|
|
|
enddo
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/,5x,"Q(r) pseudized with ", &
|
2003-01-20 05:58:50 +08:00
|
|
|
& i2," coefficients, rinner = ",3f8.3, /, &
|
|
|
|
& 58x,2f8.3)') nqf (nt) , (rinner (i, nt) , i = 1, nqlc ( &
|
|
|
|
&nt) )
|
2003-02-08 00:04:36 +08:00
|
|
|
else
|
|
|
|
if (nlc (nt) .eq.1.and.nnl (nt) .eq.1) then
|
|
|
|
ps = '(vbc)'
|
|
|
|
elseif (nlc (nt) .eq.2.and.nnl (nt) .eq.3) then
|
|
|
|
ps = '(bhs)'
|
|
|
|
elseif (nlc (nt) .eq.1.and.nnl (nt) .eq.3) then
|
|
|
|
ps = '(our)'
|
|
|
|
else
|
|
|
|
ps = ' '
|
2003-01-20 05:58:50 +08:00
|
|
|
|
|
|
|
endif
|
|
|
|
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/5x,"pseudo",i2," is ",a2, &
|
2003-01-20 05:58:50 +08:00
|
|
|
& 1x,a5," zval =",f5.1," lmax=",i2, &
|
|
|
|
& " lloc=",i2)') nt, psd (nt) , ps, zp (nt) , lmax (nt) &
|
|
|
|
&, lloc (nt)
|
2003-02-08 00:04:36 +08:00
|
|
|
if (numeric (nt) ) then
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,"(in numerical form: ",i3, &
|
2003-01-20 05:58:50 +08:00
|
|
|
&" grid points",", xmin = ",f5.2,", dx = ", &
|
|
|
|
&f6.4,")")') mesh (nt) , xmin (nt) , dx (nt)
|
2003-02-08 00:04:36 +08:00
|
|
|
else
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/14x,"i=",7x,"1",13x,"2",10x,"3")')
|
|
|
|
WRITE( stdout, '(/5x,"core")')
|
|
|
|
WRITE( stdout, '(5x,"alpha =",4x,3g13.5)') (alpc (i, nt) , i = &
|
2003-01-20 05:58:50 +08:00
|
|
|
1, 2)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,"a(i) =",4x,3g13.5)') (cc (i, nt) , i = 1, 2)
|
2003-02-08 00:04:36 +08:00
|
|
|
do l = 0, lmax (nt)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/5x,"l = ",i2)') l
|
|
|
|
WRITE( stdout, '(5x,"alpha =",4x,3g13.5)') (alps (i, l, nt) , &
|
2003-01-20 05:58:50 +08:00
|
|
|
i = 1, 3)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,"a(i) =",4x,3g13.5)') (aps (i, l, nt) , i = 1, &
|
2003-01-20 05:58:50 +08:00
|
|
|
&3)
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,"a(i+3)=",4x,3g13.5)') (aps (i, l, nt) , i &
|
2003-01-20 05:58:50 +08:00
|
|
|
= 4, 6)
|
|
|
|
enddo
|
2003-11-06 21:06:49 +08:00
|
|
|
if (nlcc (nt) ) WRITE( stdout, 200) a_nlcc (nt), b_nlcc (nt), &
|
2003-01-20 05:58:50 +08:00
|
|
|
alpha_nlcc (nt)
|
|
|
|
200 format(/5x,'nonlinear core correction: ', &
|
|
|
|
& 'rho(r) = ( a + b r^2) exp(-alpha r^2)', &
|
|
|
|
& /,5x,'a =',4x,g11.5, &
|
|
|
|
& /,5x,'b =',4x,g11.5, &
|
|
|
|
& /,5x,'alpha=',4x,g11.5)
|
|
|
|
endif
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
!
|
|
|
|
! Representation for q=0
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
if (.not.lgamma) then
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(//5x,"Atomic displacements (q=0 Repr):")')
|
|
|
|
WRITE( stdout, '(5x,"There are ",i5, &
|
2003-01-20 05:58:50 +08:00
|
|
|
& " irreducible representations")') nirrg0
|
2003-02-08 00:04:36 +08:00
|
|
|
imode0 = 0
|
|
|
|
do irr = 1, nirrg0
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/, 5x,"Representation ",i5,i7, &
|
2003-01-20 05:58:50 +08:00
|
|
|
& " modes - To be done")') irr, npertg0 (irr)
|
2003-02-08 00:04:36 +08:00
|
|
|
if (iverbosity.eq.1) then
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,"Phonon polarizations are as follows:",/)')
|
2003-02-08 00:04:36 +08:00
|
|
|
if (npertg0 (irr) .eq.1) then
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(20x," mode # ",i3)') imode0 + 1
|
|
|
|
WRITE( stdout, '(20x," (",2f10.5," ) ")') ( (ug0 (mu, nu) , nu = &
|
2003-01-20 05:58:50 +08:00
|
|
|
& imode0 + 1, imode0 + npertg0 (irr) ) , mu = 1, 3 * nat)
|
2003-02-08 00:04:36 +08:00
|
|
|
elseif (npertg0 (irr) .eq.2) then
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(2(10x," mode # ",i3,16x))') imode0 + 1, &
|
2003-01-20 05:58:50 +08:00
|
|
|
imode0 + 2
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(2(10x," (",2f10.5," ) "))') ( (ug0 (mu, nu),&
|
2003-01-20 05:58:50 +08:00
|
|
|
nu = imode0 + 1, imode0 + npertg0 (irr) ) , mu = 1, 3 * nat)
|
2003-02-08 00:04:36 +08:00
|
|
|
else
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(4x,3(" mode # ",i3,13x))') imode0 + 1, &
|
2003-01-20 05:58:50 +08:00
|
|
|
imode0 + 2, imode0 + 3
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '((5x,3("(",2f10.5," ) ")))') ( (ug0 (mu, &
|
2003-01-20 05:58:50 +08:00
|
|
|
nu) , nu = imode0 + 1, imode0 + npertg0 (irr) ) , mu = 1, &
|
|
|
|
3 * nat)
|
|
|
|
endif
|
2003-02-08 00:04:36 +08:00
|
|
|
imode0 = imode0 + npertg0 (irr)
|
2003-01-20 05:58:50 +08:00
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
!
|
|
|
|
! Representation for a generic q
|
|
|
|
!
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(//5x,"Atomic displacements:")')
|
|
|
|
WRITE( stdout, '(5x,"There are ",i5," irreducible representations") &
|
2003-01-20 05:58:50 +08:00
|
|
|
&') nirr
|
2003-02-08 00:04:36 +08:00
|
|
|
imode0 = 0
|
|
|
|
do irr = 1, nirr
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/, 5x,"Representation ",i5,i7, &
|
2003-01-20 05:58:50 +08:00
|
|
|
& " modes - To be done")') irr, npert (irr)
|
2003-02-08 00:04:36 +08:00
|
|
|
if (iverbosity.eq.1) then
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(5x,"Phonon polarizations are as follows:",/)')
|
2003-02-08 00:04:36 +08:00
|
|
|
if (npert (irr) .eq.1) then
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(20x," mode # ",i3)') imode0 + 1
|
|
|
|
WRITE( stdout, '(20x," (",2f10.5," ) ")') ( (u (mu, nu) , nu = &
|
2003-01-20 05:58:50 +08:00
|
|
|
imode0 + 1, imode0 + npert (irr) ) , mu = 1, 3 * nat)
|
2003-02-08 00:04:36 +08:00
|
|
|
elseif (npert (irr) .eq.2) then
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(2(10x," mode # ",i3,16x))') imode0 + 1, &
|
2003-01-20 05:58:50 +08:00
|
|
|
imode0 + 2
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(2(10x," (",2f10.5," ) "))') ( (u (mu, nu) , &
|
2003-01-20 05:58:50 +08:00
|
|
|
nu = imode0 + 1, imode0 + npert (irr) ) , mu = 1, 3 * nat)
|
2003-02-08 00:04:36 +08:00
|
|
|
else
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(4x,3(" mode # ",i3,13x))') imode0 + 1, imode0 &
|
2003-01-20 05:58:50 +08:00
|
|
|
+ 2, imode0 + 3
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '((5x,3("(",2f10.5," ) ")))') ( (u (mu, nu) , &
|
2003-01-20 05:58:50 +08:00
|
|
|
nu = imode0 + 1, imode0 + npert (irr) ) , mu = 1, 3 * nat)
|
|
|
|
endif
|
2003-02-08 00:04:36 +08:00
|
|
|
imode0 = imode0 + npert (irr)
|
2003-01-20 05:58:50 +08:00
|
|
|
endif
|
|
|
|
|
|
|
|
enddo
|
2003-11-06 21:06:49 +08:00
|
|
|
WRITE( stdout, '(/20x,"** Complex Version **")')
|
2005-05-25 10:58:35 +08:00
|
|
|
!
|
|
|
|
CALL flush_unit( stdout )
|
|
|
|
!
|
2003-02-08 00:04:36 +08:00
|
|
|
return
|
2003-01-20 05:58:50 +08:00
|
|
|
end subroutine d3_summary
|