TABLE OF CONTENTS


ABINIT/freeze_displ_allmodes [ Functions ]

[ Top ] [ Functions ]

NAME

 freeze_displ_allmodes

FUNCTION

  From a given set of phonon modes, generate and output supercells and
  displaced configurations of atoms.
  Typically useful to follow soft modes and see distorsions of crystal structures

COPYRIGHT

 Copyright (C) 2010-2018 ABINIT group (MJV)
 This file is distributed under the terms of the
 GNU General Public Licence, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

 amu(ntypat) = mass of the atoms (atomic mass unit)
 displ(2,3*natom,3*natom) = phonon mode displacements (complex)
 freeze_displ = amplitude of the displacement to freeze into the supercell
 natom = number of atoms in the unit cell
 ntypat = number of atom types
 phfrq(3*natom) = phonon frequencies
 qphnrm = norm of phonon q vector (should be 1 or 0)
 qphon = phonon wavevector
 rprimd(3,3) = dimensionfull primitive translations in real space
 typat(natom) = integer label of each type of atom (1,2,...)
 xcart(3,natom) = cartesian coords of atoms in unit cell (bohr)

OUTPUT

 for the moment only prints to file, but could also return pointer to supercell object, with
 rprimd and atomic positions, for further use

NOTES

 freeze_displ could be determined automatically from a temperature and the phonon frequency,
 as the average displacement of the mode with a Bose distribution.

PARENTS

      m_phonons

CHILDREN

      destroy_supercell,freeze_displ_supercell,init_supercell_for_qpt
      prt_supercell_for_qpt

SOURCE

 49 #if defined HAVE_CONFIG_H
 50 #include "config.h"
 51 #endif
 52 
 53 #include "abi_common.h"
 54 
 55 subroutine freeze_displ_allmodes(displ, freeze_displ, natom, outfile_radix, phfreq,  &
 56 &         qphon, rprimd, typat, xcart, znucl)
 57 
 58 
 59  use defs_basis
 60  use m_profiling_abi
 61  use m_supercell
 62 
 63 !This section has been created automatically by the script Abilint (TD).
 64 !Do not modify the following lines by hand.
 65 #undef ABI_FUNC
 66 #define ABI_FUNC 'freeze_displ_allmodes'
 67 !End of the abilint section
 68 
 69  implicit none
 70 
 71 ! arguments
 72 ! scalar
 73  integer,intent(in) :: natom
 74  character(len=*),intent(in) :: outfile_radix
 75  real(dp), intent(in) :: freeze_displ
 76 
 77 !arrays
 78  integer,intent(in) :: typat(natom)
 79 
 80  real(dp),intent(in) :: displ(2,3*natom,3*natom)
 81  real(dp),intent(in) :: rprimd(3,3)
 82  real(dp),intent(in) :: phfreq(3*natom)
 83  real(dp),intent(in) :: qphon(3)
 84  real(dp),intent(in) :: xcart(3,natom)
 85  real(dp),intent(in) :: znucl(:) 
 86 
 87 ! local vars
 88  integer :: jmode
 89  type(supercell_type) :: scell
 90 
 91 ! *************************************************************************
 92 
 93 !determine supercell needed to freeze phonon
 94  call init_supercell_for_qpt(natom, qphon, rprimd, typat, xcart, znucl, scell)
 95 
 96  do jmode = 1, 3*natom
 97 ! reset positions
 98    scell%xcart = scell%xcart_ref
 99 
100 !  displace atoms according to phonon jmode
101    call freeze_displ_supercell(displ(:,:,jmode), freeze_displ, scell)
102 
103 !  print out everything for this wavevector and mode
104    call prt_supercell_for_qpt (phfreq(jmode), jmode, outfile_radix, scell)
105  end do
106 
107  call destroy_supercell (scell)
108 
109 end subroutine freeze_displ_allmodes