TABLE OF CONTENTS


ABINIT/xcart2xred [ Functions ]

[ Top ] [ Functions ]

NAME

 xcart2xred

FUNCTION

 Convert from cartesian coordinates xcart(3,natom) in bohr to
 dimensionless reduced coordinates xred(3,natom) by using
 xred(mu,ia)=gprimd(1,mu)*xcart(1,ia)
            +gprimd(2,mu)*xcart(2,ia)
            +gprimd(3,mu)*xcart(3,ia)
 where gprimd is the inverse of rprimd
 Note that the reverse operation is deon by xred2xcart

COPYRIGHT

 Copyright (C) 1998-2017 ABINIT group (DCA, XG, GMR)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

  natom=number of atoms in unit cell
  rprimd(3,3)=dimensional real space primitive translations (bohr)
  xcart(3,natom)=cartesian coordinates of atoms (bohr)

OUTPUT

  xred(3,natom)=dimensionless reduced coordinates of atoms

SIDE EFFECTS

PARENTS

      driver,evdw_wannier,ingeo,m_cut3d,m_dens,m_effective_potential
      m_effective_potential_file,m_fit_polynomial_coeff,m_mep
      m_paw_pwaves_lmn,m_pred_lotf,mkcore_paw,mkcore_wvl,mover_effpot
      pawmkaewf,pimd_langevin_npt,pimd_langevin_nvt,pimd_nosehoover_npt
      pimd_nosehoover_nvt,prcref,prcref_PMA,pred_delocint,pred_diisrelax
      pred_isokinetic,pred_isothermal,pred_langevin,pred_moldyn,pred_nose
      pred_srkna14,pred_steepdesc,pred_velverlet,pred_verlet,relaxpol
      wrt_moldyn_netcdf,wvl_setboxgeometry

CHILDREN

      matr3inv

SOURCE

47 #if defined HAVE_CONFIG_H
48 #include "config.h"
49 #endif
50 
51 #include "abi_common.h"
52 
53 
54 subroutine xcart2xred(natom,rprimd,xcart,xred)
55 
56  use defs_basis
57  use m_errors
58  use m_profiling_abi
59 
60 !This section has been created automatically by the script Abilint (TD).
61 !Do not modify the following lines by hand.
62 #undef ABI_FUNC
63 #define ABI_FUNC 'xcart2xred'
64  use interfaces_32_util
65 !End of the abilint section
66 
67  implicit none
68 
69 !Arguments ------------------------------------
70 !scalars
71  integer,intent(in) :: natom
72 !arrays
73  real(dp),intent(in) :: rprimd(3,3),xcart(3,natom)
74  real(dp),intent(out) :: xred(3,natom)
75 
76 !Local variables-------------------------------
77 !scalars
78  integer :: iatom,mu
79 !arrays
80  real(dp) :: gprimd(3,3)
81 
82 ! *************************************************************************
83 
84  call matr3inv(rprimd,gprimd)
85  do iatom=1,natom
86    do mu=1,3
87      xred(mu,iatom)= gprimd(1,mu)*xcart(1,iatom)+gprimd(2,mu)*xcart(2,iatom)+&
88 &     gprimd(3,mu)*xcart(3,iatom)
89    end do
90  end do
91 
92 end subroutine xcart2xred