TABLE OF CONTENTS


ABINIT/setsym [ Functions ]

[ Top ] [ Functions ]

NAME

 setsym

FUNCTION

 Set up irreducible zone in  G space by direct calculation.
 Do not call this routine if nsym=1 (only identity symmetry).
 Only indsym and symrec get returned if iscf=0.
 symrec needed to symmetrize coordinate gradients in sygrad.
 (symrec is redundant and could be removed later in favor of symrel)

COPYRIGHT

 Copyright (C) 1998-2018 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

 iscf=(<= 0  =>non-SCF), >0 => SCF
 natom=number of atoms in unit cell
 nfft=(effective) number of FFT grid points (for this processor)
 ngfft(18)=contain all needed information about 3D FFT, see ~abinit/doc/variables/vargs.htm#ngfft
 nspden=number of spin-density components
 nsppol=1 for unpolarized, 2 for spin-polarized
 nsym=number of symmetries in space group (at least 1)
 symafm(nsym)=(anti)ferromagnetic part of symmetry operations
 symrel(3,3,nsym)=symmetry operations in terms of real space
 primitive translations
 tnons(3,nsym)=nonsymmorphic translations of space group in terms
 of real space primitive translations (may be 0)
 typat(natom)=atom type (integer) for each atom
 xred(3,natom)=atomic coordinates in terms of real space primitive translations

OUTPUT

 indsym(4,nsym,natom)=indirect indexing of atom labels--see subroutine
   symatm for definition (if nsym>1)
 irrzon(nfft,2,(nspden/nsppol)-3*(nspden/4))=irreducible zone data
 phnons(2,nfft,(nspden/nsppol)-3*(nspden/4))=nonsymmorphic translation phases
 symrec(3,3,nsym)=symmetry operations in terms of reciprocal
   space primitive translations (if nsym>1)

NOTES

 nsppol and nspden are needed in case of (anti)ferromagnetic symmetry operations

PARENTS

      dfpt_looppert,gstate,m_dvdb,nonlinear,respfn,scfcv

CHILDREN

      chkgrp,irrzg,mati3inv,symatm,symdet,timab

SOURCE

 55 #if defined HAVE_CONFIG_H
 56 #include "config.h"
 57 #endif
 58 
 59 #include "abi_common.h"
 60 
 61 subroutine setsym(indsym,irrzon,iscf,natom,nfft,ngfft,nspden,nsppol,nsym,phnons,&
 62 & symafm,symrec,symrel,tnons,typat,xred)
 63 
 64  use defs_basis
 65  use m_profiling_abi
 66  use m_errors
 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 'setsym'
 72  use interfaces_18_timing
 73  use interfaces_32_util
 74  use interfaces_41_geometry
 75  use interfaces_56_recipspace, except_this_one => setsym
 76 !End of the abilint section
 77 
 78  implicit none
 79 
 80 !Arguments ------------------------------------
 81 !scalars
 82  integer,intent(in) :: iscf,natom,nfft,nspden,nsppol,nsym
 83 !arrays
 84  integer,intent(in) :: ngfft(18),symafm(nsym),symrel(3,3,nsym),typat(natom)
 85  integer,intent(out) :: indsym(4,nsym,natom)
 86  integer,intent(inout) :: irrzon(nfft,2,(nspden/nsppol)-3*(nspden/4)) !vz_i
 87  integer,intent(out) :: symrec(3,3,nsym)
 88  real(dp),intent(in) :: tnons(3,nsym),xred(3,natom)
 89  real(dp),intent(out) :: phnons(2,nfft,(nspden/nsppol)-3*(nspden/4))
 90 
 91 !Local variables-------------------------------
 92 !scalars
 93  integer :: isym,ierr
 94  real(dp) :: tolsym8
 95 !arrays
 96  integer,allocatable :: determinant(:)
 97  real(dp) :: tsec(2)
 98 
 99 ! *************************************************************************
100 
101  call timab(6,1,tsec)
102 
103 !Check that symmetries have unity determinant
104  ABI_ALLOCATE(determinant,(nsym))
105  call symdet(determinant,nsym,symrel)
106  ABI_DEALLOCATE(determinant)
107 
108 
109 !Get the symmetry matrices in terms of reciprocal basis
110  do isym=1,nsym
111    call mati3inv(symrel(:,:,isym),symrec(:,:,isym))
112  end do
113 
114 !Check for group closure
115  call chkgrp(nsym,symafm,symrel,ierr)
116  ABI_CHECK(ierr==0,"Error in group closure")
117 
118  call chkgrp(nsym,symafm,symrec,ierr)
119  ABI_CHECK(ierr==0,"Error in group closure")
120 
121 !Obtain a list of rotated atom labels:
122  tolsym8=tol8
123  call symatm(indsym,natom,nsym,symrec,tnons,tolsym8,typat,xred)
124 
125 !If non-SCF calculation, or nsym==1, do not need IBZ data
126  if ( (iscf>0 .or. iscf==-3) .and. nsym>1 ) then
127 !  Locate irreducible zone in reciprocal space for symmetrization:
128    call irrzg(irrzon,nspden,nsppol,nsym,ngfft(1),ngfft(2),ngfft(3),phnons,symafm,symrel,tnons)
129  end if
130 
131  call timab(6,2,tsec)
132 
133 end subroutine setsym