TABLE OF CONTENTS


ABINIT/xfpack_x2vin [ Functions ]

[ Top ] [ Functions ]

NAME

 xfpack_x2vin

FUNCTION

 Old option=1, transfer xred, acell, and rprim to vin

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (XG,MJV)
 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

 acell0(3)=reference length scales of primitive translations (bohr), needed
   for some values of optcell.
 natom=number of atoms in cell
 ndim=dimension of vin arrays
 nsym=order of group.
 rprimd0(3,3)=reference real space primitive translations,
   needed for some values of optcell.
 optcell=option for the optimisation of the unit cell. Described in abinit_help.
  Depending on its value, different part of acell and rprim
  are contained in vin.
 symrel(3,3,nsym)=symmetry operators in terms of action on primitive translations
 ucvol=unit cell volume (bohr^3), needed for some values of optcell.
 ucvol0=reference unit cell volume (bohr^3), needed for some values of optcell.

OUTPUT

  (see side effects)

SIDE EFFECTS

 Input/Output variables
 acell(3)=length scales of primitive translations (bohr)
 rprim(3,3)=dimensionless real space primitive translations
 vin(ndim)=vector that contains xred and some quantity derived
   from acell and rprim, depending on the value of optcell.
 xred(3,natom)=reduced dimensionless atomic coordinates

PARENTS

      pred_bfgs,pred_delocint,pred_lbfgs,pred_verlet,xfh_recover_deloc
      xfh_recover_new

CHILDREN

      matr3inv,metric,mkrdim,strainsym

SOURCE

 52 #if defined HAVE_CONFIG_H
 53 #include "config.h"
 54 #endif
 55 
 56 #include "abi_common.h"
 57 
 58 
 59 subroutine xfpack_x2vin(acell,acell0,natom,ndim,nsym,optcell,&
 60   & rprim,rprimd0,symrel,ucvol,ucvol0,vin,xred)
 61 
 62  use defs_basis
 63  use m_errors
 64  use m_profiling_abi
 65 
 66  use m_geometry,   only : mkrdim, metric
 67 
 68 !This section has been created automatically by the script Abilint (TD).
 69 !Do not modify the following lines by hand.
 70 #undef ABI_FUNC
 71 #define ABI_FUNC 'xfpack_x2vin'
 72  use interfaces_32_util
 73  use interfaces_45_geomoptim, except_this_one => xfpack_x2vin
 74 !End of the abilint section
 75 
 76  implicit none
 77 
 78 !Arguments ------------------------------------
 79 !scalars
 80  integer,intent(in) :: natom,ndim,nsym,optcell
 81  real(dp),intent(in) :: ucvol0
 82  real(dp),intent(inout) :: ucvol !vz_i
 83 !arrays
 84  integer,intent(in) :: symrel(3,3,nsym)
 85  real(dp),intent(in) :: acell0(3),rprimd0(3,3)
 86  real(dp),intent(in) :: acell(3),rprim(3,3)
 87  real(dp),intent(in) :: xred(3,natom)
 88  real(dp),intent(out) :: vin(ndim)
 89 
 90 !Local variables-------------------------------
 91 !scalars
 92  integer :: ii,jj,kk
 93  real(dp) :: scale
 94  character(len=500) :: message
 95 !arrays
 96  real(dp) :: gmet(3,3),gprimd(3,3),gprimd0(3,3),rmet(3,3),rprimd(3,3)
 97  real(dp) :: rprimd_symm(3,3),scaling(3,3)
 98 
 99 ! *************************************************************************
100 
101 !!DEBUG
102 !write(ab_out,*) ''
103 !write(ab_out,*) 'xfpack_x2vin'
104 !write(ab_out,*) 'natom=',natom
105 !write(ab_out,*) 'ndim=',ndim
106 !write(ab_out,*) 'nsym=',nsym
107 !write(ab_out,*) 'optcell=',optcell
108 !write(ab_out,*) 'ucvol=',ucvol
109 !write(ab_out,*) 'xred='
110 !do kk=1,natom
111 !write(ab_out,*) xred(:,kk)
112 !end do
113 !write(ab_out,*) 'VECTOR INPUT (vin) xfpack_x2vin INPUT'
114 !do ii=1,ndim,3
115 !if (ii+2<=ndim)then
116 !write(ab_out,*) ii,vin(ii:ii+2)
117 !else
118 !write(ab_out,*) ii,vin(ii:ndim)
119 !end if
120 !end do
121 !!DEBUG
122 
123 
124 !##########################################################
125 !### 1. Test for compatible ndim
126 
127  if(optcell==0 .and. ndim/=3*natom)then
128    write(message,'(a,a,a,i4,a,i4,a)' )&
129 &   '  When optcell=0, ndim MUST be equal to 3*natom,',ch10,&
130 &   '  while ndim=',ndim,' and 3*natom=',3*natom,'.'
131    MSG_BUG(message)
132  end if
133 
134  if( (optcell==1 .or. optcell==4 .or. optcell==5 .or. optcell==6) &
135 & .and. ndim/=3*natom+1)then
136    write(message,'(a,a,a,i4,a,i4,a)' )&
137 &   '  When optcell=1,4,5 or 6, ndim MUST be equal to 3*natom+1,',ch10,&
138 &   '  while ndim=',ndim,' and 3*natom+1=',3*natom+1,'.'
139    MSG_BUG(message)
140  end if
141 
142  if( (optcell==2 .or. optcell==3) &
143 & .and. ndim/=3*natom+6)then
144    write(message,'(a,a,a,i4,a,i4,a)' )&
145 &   '  When optcell=2 or 3, ndim MUST be equal to 3*natom+6,',ch10,&
146 &   '  while ndim=',ndim,' and 3*natom+6=',3*natom+6,'.'
147    MSG_BUG(message)
148  end if
149 
150  if( optcell>=7 .and. ndim/=3*natom+3)then
151    write(message,'(a,a,a,i4,a,i4,a)' )&
152 &   '  When optcell=7,8 or 9, ndim MUST be equal to 3*natom+3,',ch10,&
153 &   '  while ndim=',ndim,' and 3*natom+3=',3*natom+3,'.'
154    MSG_BUG(message)
155  end if
156 
157 !##########################################################
158 !### 2. option=1, transfer xred, acell, and rprim to vin
159 
160 !Get vin from xred, acell, and rprim
161  vin(1:3*natom)= reshape(xred(:,:), (/3*natom/) )
162 
163  if(optcell/=0)then
164    call mkrdim(acell,rprim,rprimd)
165    call strainsym(nsym,rprimd0,rprimd,rprimd_symm,symrel)
166    call metric(gmet,gprimd,-1,rmet,rprimd_symm,ucvol)
167 
168    if(optcell==1)then
169 
170 !    vin(3*natom+1)=ucvol**third
171      vin(3*natom+1)=(ucvol/ucvol0)**third
172 
173    else if(optcell==2 .or. optcell==3 .or. optcell>=7)then
174 
175 !    Generates gprimd0
176      call matr3inv(rprimd0,gprimd0)
177      do ii=1,3
178        do jj=1,3
179          scaling(ii,jj)=0.0_dp
180          do kk=1,3
181            scaling(ii,jj)=scaling(ii,jj)+rprimd_symm(ii,kk)*gprimd0(jj,kk)
182          end do
183        end do
184      end do
185 !    Rescale if the volume must be preserved
186      if(optcell==3)then
187        scale=(ucvol0/ucvol)**third
188        scaling(:,:)=scale*scaling(:,:)
189      end if
190      if(optcell==2 .or. optcell==3)then
191        vin(3*natom+1)=scaling(1,1) ; vin(3*natom+4)=(scaling(2,3)+scaling(3,2))*0.5_dp
192        vin(3*natom+2)=scaling(2,2) ; vin(3*natom+5)=(scaling(1,3)+scaling(3,1))*0.5_dp
193        vin(3*natom+3)=scaling(3,3) ; vin(3*natom+6)=(scaling(1,2)+scaling(2,1))*0.5_dp
194      else if(optcell>=7)then
195        vin(3*natom+1)=scaling(1,1)
196        vin(3*natom+2)=scaling(2,2)
197        vin(3*natom+3)=scaling(3,3)
198        if(optcell==7)vin(3*natom+1)=(scaling(2,3)+scaling(3,2))*0.5_dp
199        if(optcell==8)vin(3*natom+2)=(scaling(1,3)+scaling(3,1))*0.5_dp
200        if(optcell==9)vin(3*natom+3)=(scaling(1,2)+scaling(2,1))*0.5_dp
201      end if
202 
203    else if(optcell==4 .or. optcell==5 .or. optcell==6)then
204 
205      vin(3*natom+1)=acell(optcell-3)/acell0(optcell-3)
206 
207    end if
208 
209  end if
210 
211 
212 end subroutine xfpack_x2vin