TABLE OF CONTENTS


ABINIT/fillcell [ Functions ]

[ Top ] [ Functions ]

NAME

 fillcell

FUNCTION

 Computes the atomic position of all the atoms in the unit cell starting
 with the symmetry operations and the atoms from the asymetric unit cell.

COPYRIGHT

 Copyright (C) 1999-2018 ABINIT group (RC,XG)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .

INPUTS

  natrd = number of atoms in the assymetric unit cell
  natom = total number of atoms (to be checked)
  nsym = number of symmetry operations
  symafm(nsym)=(anti)ferromagnetic part of symmetry operations
  symrel(3,3,nsym)=symmetry operations in real space in terms
   of primitive translations
  tnons(3,nsym)=nonsymmorphic translations for symmetry operations
  tolsym=tolerance on symmetries
  typat(1:natrd)=type integer for each atom in cell
  xred(3,1:natrd)=reduced dimensionless atomic coordinates

OUTPUT

SIDE EFFECTS

  At input, for the assymetric unit cell
  nucdipmom(3,1:natrd)=nuclear magnetic dipole moments of the atoms
  spinat(3,1:natrd)=spin-magnetization of the atoms
  typat(1:natrd)=type integer for each atom in cell
  xred(3,1:natrd)=reduced dimensionless atomic coordinates

  At output, for the complete unit cell
  nucdipmom(3,1:natom)=nuclear magnetic dipole moments of the atoms
  spinat(3,1:natom)=spin-magnetization of the atoms
  typat(1:natom)=type integer for each atom in cell
  xred(3,1:natom)=reduced dimensionless atomic coordinates

PARENTS

      ingeo

CHILDREN

SOURCE

 50 #if defined HAVE_CONFIG_H
 51 #include "config.h"
 52 #endif
 53 
 54 #include "abi_common.h"
 55 
 56 
 57 subroutine fillcell(natom,natrd,nsym,nucdipmom,spinat,symafm,symrel,tnons,tolsym,typat,xred)
 58 
 59  use defs_basis
 60  use m_errors
 61  use m_profiling_abi
 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 'fillcell'
 67 !End of the abilint section
 68 
 69  implicit none
 70 
 71 !Arguments ------------------------------------
 72 !scalars
 73  integer,intent(in) :: natom,natrd,nsym
 74 !arrays
 75  integer,intent(in) :: symafm(nsym),symrel(3,3,nsym)
 76  integer,intent(inout) :: typat(natom)
 77  real(dp),intent(in) :: tolsym
 78  real(dp),intent(in) :: tnons(3,nsym)
 79  real(dp),intent(inout) :: nucdipmom(3,natom),spinat(3,natom),xred(3,natom)
 80 
 81 !Local variables ------------------------------
 82 !scalars
 83  integer :: curat,flagch,flageq,ii,iij,jj,kk
 84  character(len=500) :: message
 85 !arrays
 86  integer :: bcktypat(nsym*natrd)
 87  real(dp) :: bckat(3),bcknucdipmom(3,nsym*natrd)
 88  real(dp) :: bckspinat(3,nsym*natrd),bckxred(3,nsym*natrd)
 89 
 90 ! *************************************************************************
 91 
 92 !DEBUG
 93 !write(std_out,*)' fillcell : enter with nsym, natrd= ',nsym,natrd
 94 !write(std_out,*)' Describe the different symmetry operations (index,symrel,tnons,symafm)'
 95 !do ii=1,nsym
 96 !write(std_out,'(i3,2x,9i3,3es12.2,i3)')ii,symrel(:,:,ii),tnons(:,ii),symafm(ii)
 97 !end do
 98 !write(std_out,*)' Describe the input atoms (index,typat,xred,spinat)'
 99 !do jj=1,natrd
100 !write(std_out,'(i3,2x,i3,6es12.2)')jj,typat(jj),xred(:,jj),spinat(:,jj)
101 !end do
102 !ENDDEBUG
103 
104  curat=0
105 
106 !Cycle over all the symmetry operations
107  do ii=1,nsym
108 
109 !  Cycle over all the atoms in the assymetric unit cell
110    do jj=1,natrd
111 
112 !    Symmetry operation application
113      bckat(:)=matmul(symrel(:,:,ii),xred(:,jj))+tnons(:,ii)
114 
115 !    Normalization of the coordinates in [0,1)
116      do iij=1,3
117        do while (bckat(iij)<-tolsym)
118          bckat(iij)=bckat(iij)+1.0d0
119        end do
120        do while (bckat(iij)>=1.0d0-tolsym)
121          bckat(iij)=bckat(iij)-1.0d0
122        end do
123      end do
124 
125 !    Check for duplicate atoms
126      flagch=0
127      do kk=1,curat
128        flageq=0
129        if ( abs(bckxred(1,kk)-bckat(1))<tolsym  .and. &
130 &       abs(bckxred(2,kk)-bckat(2))<tolsym  .and. &
131 &       abs(bckxred(3,kk)-bckat(3))<tolsym       ) exit
132        flagch=flagch+1
133      end do
134 
135      if (flagch==curat) then
136 !      Add the obtained atom to the bckxred list
137        curat=curat+1
138        bckxred(:,curat)=bckat
139        bcktypat(curat)=typat(jj)
140        bcknucdipmom(:,curat)=nucdipmom(:,jj)
141        bckspinat(:,curat)=spinat(:,jj)*symafm(ii)
142      end if
143 
144    end do
145 
146  end do
147 
148 !DEBUG
149 !write(std_out,*)' fillcell : Proposed coordinates ='
150 !do ii=1,curat
151 !write(std_out,'(i4,3es16.6)' )ii,bckxred(:,ii)
152 !end do
153 !ENDDEBUG
154 
155  if (curat>natom) then
156    write(message, '(a,i3,a,a,i7,a,a,a,a)' )&
157 &   'The number of atoms obtained from symmetries, ',curat,ch10,&
158 &   'is greater than the input number of atoms, natom=',natom,ch10,&
159 &   'This is not allowed.  ',ch10,&
160 &   'Action: modify natom or the symmetry data in the input file.'
161    MSG_ERROR(message)
162  end if
163 
164  if (curat<natom) then
165    write(message, '(a,i3,a,a,i7,a,a,a,a)' )&
166 &   'fillcell : The number of atoms obtained from symmetries, ',curat,ch10,&
167 &   'is lower than the input number of atoms, natom=',natom,ch10,&
168 &   'This is not allowed.  ',ch10,&
169 &   'Action: modify natom or the symmetry data in the input file.'
170    MSG_ERROR(message)
171  end if
172 
173 !Assignment of symmetry to xred
174  xred(:,1:natom)=bckxred(:,1:natom)
175  typat(1:natom)=bcktypat(1:natom)
176  nucdipmom(1:3,1:natom)=bcknucdipmom(1:3,1:natom)
177  spinat(1:3,1:natom)=bckspinat(1:3,1:natom)
178 
179 !DEBUG
180 !write(std_out,*)' fillcell : exit with natom=',natom
181 !write(std_out,*)' Describe the output atoms (index,typat,xred,spinat)'
182 !do jj=1,natom
183 !write(std_out,'(i3,2x,i3,6es12.2)')jj,typat(jj),xred(:,jj),spinat(:,jj)
184 !end do
185 !ENDDEBUG
186 
187 end subroutine fillcell