TABLE OF CONTENTS


ABINIT/gensymshub4 [ Functions ]

[ Top ] [ Functions ]

NAME

 gensymshub4

FUNCTION

 Assigns the Bravais magnetic translations for shubnikov type IV
 symmetry groups starting from the Fedorov space group
 and the translation, generator of the anti-ferromagnetic
 operations (as input). It will double nsym and tnons. In the end both
 symrel and tnons are completely determined.

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 .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

 genafm(3) = translation, generator of the anti-ferromagnetic symmatry operations
 msym = maximum number of symmetry operations

OUTPUT

 symafm(msym)= (anti)ferromagnetic part of symmetry operations

SIDE EFFECTS

 nsym=number of symmetry operations, without magnetic operations at input,
  and with magnetic operations at output
 symrel(3,3,msym)=symmetry operations in real space in terms
  of primitive translations, without magnetic operations at input,
  and with magnetic operations at output
 tnons(3,msym)=nonsymmorphic translations for symmetry operations
  without magnetic operations at input,
  and with magnetic operations at output

PARENTS

      ingeo

CHILDREN

SOURCE

44 #if defined HAVE_CONFIG_H
45 #include "config.h"
46 #endif
47 
48 #include "abi_common.h"
49 
50 
51 subroutine gensymshub4(genafm,msym,nsym,symafm,symrel,tnons)
52 
53  use defs_basis
54  use m_profiling_abi
55  use m_errors
56 
57 !This section has been created automatically by the script Abilint (TD).
58 !Do not modify the following lines by hand.
59 #undef ABI_FUNC
60 #define ABI_FUNC 'gensymshub4'
61 !End of the abilint section
62 
63  implicit none
64 
65 !Arguments ------------------------------------
66 !scalars
67  integer,intent(in) :: msym
68  integer,intent(inout) :: nsym
69 !arrays
70  integer,intent(inout) :: symafm(msym),symrel(3,3,msym)
71  real(dp),intent(in) :: genafm(3)
72  real(dp),intent(inout) :: tnons(3,msym)
73 
74 !Local variables ------------------------------
75 !scalars
76  integer :: ii
77  character(len=500) :: message
78 
79 ! *************************************************************************
80 
81  if(msym<2*nsym)then
82    write(message, '(3a)' )&
83 &   'The number of symmetries in the Shubnikov type IV space group',ch10,&
84 &   'is larger than the maximal allowed number of symmetries.'
85    MSG_ERROR(message)
86  end if
87 
88  do ii=1,nsym
89    tnons(:,nsym+ii)=tnons(:,ii)+genafm(:)
90    symrel(:,:,nsym+ii)=symrel(:,:,ii)
91    symafm(ii)=1
92    symafm(nsym+ii)=-1
93  end do
94  nsym=nsym*2
95 
96 end subroutine gensymshub4