TABLE OF CONTENTS


ABINIT/symsgtetra [ Functions ]

[ Top ] [ Functions ]

NAME

 symsgtetra

FUNCTION

 Yields all the TETRAGONAL symmetry operations starting from the space group symbol.
 according to the International Tables of Crystallography, 1983.

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

 msym = default number of symmetries
 nsym = the number of symmetry operations
 shubnikov= magnetic type of the space group to be generated
 spgorig = the origin choice (1 or 2) for the axes system
 spgroup = the numeric symbol of the space groups
 spgroupma= number of the magnetic space group

OUTPUT

 symafm(nsym)=(anti)ferromagnetic part of symmetry operations
 symrel(3,3,nsym) = 3D matrix containg symmetry operations
 tnons(3,nsym) = 2D matrix containing translations associated

PARENTS

      gensymspgr

CHILDREN

      bldgrp,spgdata

SOURCE

 38 #if defined HAVE_CONFIG_H
 39 #include "config.h"
 40 #endif
 41 
 42 #include "abi_common.h"
 43 
 44 
 45 subroutine symsgtetra(msym,nsym,shubnikov,spgaxor,spgorig,spgroup,&
 46 &   spgroupma,symafm,symrel,tnons)
 47 
 48  use defs_basis
 49  use m_profiling_abi
 50 
 51 !This section has been created automatically by the script Abilint (TD).
 52 !Do not modify the following lines by hand.
 53 #undef ABI_FUNC
 54 #define ABI_FUNC 'symsgtetra'
 55  use interfaces_41_geometry, except_this_one => symsgtetra
 56 !End of the abilint section
 57 
 58  implicit none
 59 
 60 !Arguments ------------------------------------
 61 !scalars
 62  integer,intent(in) :: msym,nsym,shubnikov,spgaxor,spgorig,spgroup
 63  integer,intent(in) :: spgroupma
 64 !arrays
 65  integer,intent(inout) :: symafm(msym),symrel(3,3,msym) !vz_i
 66  real(dp),intent(inout) :: tnons(3,msym) !vz_i
 67 
 68 !Local variables ------------------------------
 69 !integer :: isym
 70 !scalars
 71  integer :: nogen,sporder
 72  character(len=1) :: brvsb
 73  character(len=15) :: intsb,ptintsb,ptschsb,schsb
 74  character(len=35) :: intsbl
 75 !arrays
 76  integer :: gen4m(3,3),genmmm(3,3),genmmp(3,3),genmpm(3,3),genmpp(3,3)
 77  integer :: genpmm(3,3),genpmp(3,3),genppm(3,3)
 78 
 79 ! *************************************************************************
 80 !DEBUG
 81 !write(std_out,*) ' symsgtetra: ',spgroup,shubnikov,spgroupma
 82 !ENDDEBUG
 83  nogen=0
 84 
 85  tnons(:,:)=zero
 86 !The identity operation belongs to all space groups
 87  symrel(:,:,1)=0 ; symrel(1,1,1)=1 ; symrel(2,2,1)=1 ; symrel(3,3,1)=1
 88 
 89 !The next operation belongs to most TETRAGONAL space groups, except 81,82,111,121.
 90  symrel(:,:,2)=0 ; symrel(1,2,2)=-1 ; symrel(2,1,2)=1 ; symrel(3,3,2)=1
 91 
 92 !Predefine some generators
 93  genmpp(:,:)=0 ; genmpp(1,1)=-1 ; genmpp(2,2)= 1 ; genmpp(3,3)= 1
 94  genpmp(:,:)=0 ; genpmp(1,1)= 1 ; genpmp(2,2)=-1 ; genpmp(3,3)= 1
 95  genppm(:,:)=0 ; genppm(1,1)= 1 ; genppm(2,2)= 1 ; genppm(3,3)=-1
 96  genpmm(:,:)=0 ; genpmm(1,1)= 1 ; genpmm(2,2)=-1 ; genpmm(3,3)=-1
 97  genmpm(:,:)=0 ; genmpm(1,1)=-1 ; genmpm(2,2)= 1 ; genmpm(3,3)=-1
 98  genmmp(:,:)=0 ; genmmp(1,1)=-1 ; genmmp(2,2)=-1 ; genmmp(3,3)= 1
 99  genmmm(:,:)=0 ; genmmm(1,1)=-1 ; genmmm(2,2)=-1 ; genmmm(3,3)=-1
100  gen4m(:,:)=0  ; gen4m(1,2)=-1 ; gen4m(2,1)=1 ; gen4m(3,3)=-1
101 
102 !Default non-magnetic behaviour
103  symafm(1:nsym)=1
104 
105 
106 !assigns the generators to each space group
107  select case (spgroup)
108 !  TETRAGONAL space groups
109  case (75,79,83,87)        !P4, I4, P4/m, I4/m
110    nogen=2
111  case (76,80)                !P41, I41
112    tnons(:,2)=(/0.d0,0.d0,0.25d0/)
113    nogen=2
114  case (77,84)                !P42, P42/m
115    tnons(:,2)=(/0.d0,0.d0,0.5d0/)
116    nogen=2
117  case (78)                !P43
118    tnons(:,2)=(/0.d0,0.d0,0.75d0/)
119    nogen=2
120  case (81,82)                !PB4, IB4
121    symrel(:,:,1)=0 ; symrel(1,1,1)=1 ; symrel(2,2,1)=1 ; symrel(3,3,1)=1
122    symrel(:,:,2)=0 ; symrel(1,1,2)=-1; symrel(2,2,2)=-1; symrel(3,3,2)=1
123    symrel(:,:,3)=0 ; symrel(1,2,3)=1 ; symrel(2,1,3)=-1 ; symrel(3,3,3)=-1
124    symrel(:,:,4)=0 ; symrel(1,2,4)=-1 ; symrel(2,1,4)=1 ; symrel(3,3,4)=-1
125    nogen=0
126    if (shubnikov==3) then
127      symafm(3:4)=-1
128    end if
129  case (85,86,88)                !P4/n
130    symrel(:,:,1)=0 ; symrel(1,1,1)=1 ; symrel(2,2,1)=1 ; symrel(3,3,1)=1 ; symrel(:,:,5) = -1*symrel(:,:,1)
131    symrel(:,:,2)=0 ; symrel(1,1,2)=-1; symrel(2,2,2)=-1; symrel(3,3,2)=1 ; symrel(:,:,6) = -1*symrel(:,:,2)
132    symrel(:,:,3)=0 ; symrel(1,2,3)=-1 ; symrel(2,1,3)=1 ; symrel(3,3,3)=1 ; symrel(:,:,7) = -1*symrel(:,:,3)
133    symrel(:,:,4)=0 ; symrel(1,2,4)=1 ; symrel(2,1,4)=-1 ; symrel(3,3,4)=1 ; symrel(:,:,8) = -1*symrel(:,:,4)
134    nogen=0
135    if(spgorig==1) then
136      select case (spgroup)
137      case (85)                !P4/n
138        tnons(:,3)=(/0.5d0,0.5d0,0.d0/) ; tnons(:,4)=(/0.5d0,0.5d0,0.d0/)
139        tnons(:,5)=(/0.5d0,0.5d0,0.d0/) ; tnons(:,6)=(/0.5d0,0.5d0,0.d0/)
140      case (86)                !P42/n
141        tnons(:,3)=(/0.5d0,0.5d0,0.5d0/) ; tnons(:,4)=(/0.5d0,0.5d0,0.5d0/)
142        tnons(:,5)=(/0.5d0,0.5d0,0.5d0/) ; tnons(:,6)=(/0.5d0,0.5d0,0.5d0/)
143      case (88)                 !I41/a
144        tnons(:,2)=(/0.5d0,0.5d0,0.5d0/)
145        tnons(:,3)=(/0.0d0,0.5d0,0.25d0/)
146        tnons(:,4)=(/0.5d0,0.0d0,0.75d0/)
147        tnons(:,5)=(/0.0d0,0.5d0,0.25d0/)
148        tnons(:,6)=(/0.5d0,0.0d0,0.75d0/)
149        tnons(:,7)=(/0.5d0,0.5d0,0.5d0/)
150      end select
151    else
152      select case (spgroup)
153      case (85)          !P4/n
154        tnons(:,2)=(/0.5d0,0.5d0,0.d0/) ;  tnons(:,6)=(/0.5d0,0.5d0,0.d0/)
155        tnons(:,3)=(/0.5d0,0.0d0,0.d0/) ;  tnons(:,7)=(/0.5d0,0.0d0,0.d0/)
156        tnons(:,4)=(/0.0d0,0.5d0,0.d0/) ;  tnons(:,8)=(/0.0d0,0.5d0,0.d0/)
157      case (86)          !P42/n
158        tnons(:,2)=(/0.5d0,0.5d0,0.0d0/) ; tnons(:,6)=(/0.5d0,0.5d0,0.0d0/)
159        tnons(:,3)=(/0.0d0,0.5d0,0.5d0/) ; tnons(:,7)=(/0.0d0,0.5d0,0.5d0/)
160        tnons(:,4)=(/0.5d0,0.0d0,0.5d0/) ; tnons(:,8)=(/0.5d0,0.0d0,0.5d0/)
161      case (88)          !I41/a
162        tnons(:,2)=(/0.5d0,0.0d0,0.5d0/) ;  tnons(:,6)=(/0.5d0,0.0d0,0.5d0/)
163        tnons(:,3)=(/0.75d0,0.25d0,0.25d0/)
164        tnons(:,7)=(/0.25d0,0.75d0,0.75d0/)
165        tnons(:,4)=(/0.75d0,0.75d0,0.75d0/)
166        tnons(:,8)=(/0.25d0,0.25d0,0.25d0/)
167      end select
168    end if
169    if (shubnikov==3) then
170      select case (spgroupma)
171      case(61,69,83)
172        symafm(3)=-1;symafm(4)=-1;symafm(7)=-1; symafm(8)=-1
173      case(62,70,84)
174        symafm(5)=-1;symafm(6)=-1;symafm(7)=-1; symafm(8)=-1
175      case(63,71,85)
176        symafm(3)=-1;symafm(5)=-1;symafm(4)=-1; symafm(6)=-1
177      end select
178    end if
179  case (89,97)                !P422, I422
180    symrel(:,:,3) = genmpm(:,:)
181    nogen=3
182  case (90)                !P4212
183    tnons(:,2)=(/0.5d0,0.5d0,0.0d0/)
184    symrel(:,:,3) = genmpm(:,:)
185    tnons(:,3)=(/0.5d0,0.5d0,0.d0/)
186    nogen=3
187  case (91)                !P4122
188    tnons(:,2)=(/0.d0,0.d0,0.25d0/)
189    symrel(:,:,3) = genmpm(:,:)
190    nogen=3
191  case (92)                !P41212
192    tnons(:,2)=(/0.5d0,0.5d0,0.25d0/)
193    symrel(:,:,3) = genmpm(:,:)
194    tnons(:,3)=(/0.5d0,0.5d0,0.25d0/)
195    nogen=3
196  case (93)                !P4222
197    tnons(:,2)=(/0.d0,0.d0,0.5d0/)
198    symrel(:,:,3) = genmpm(:,:)
199    nogen=3
200  case (94)                !P42212
201    tnons(:,2)=(/0.5d0,0.5d0,0.5d0/)
202    symrel(:,:,3) = genmpm(:,:)
203    tnons(:,3)=(/0.5d0,0.5d0,0.5d0/)
204    nogen=3
205  case (95)                !P4322
206    tnons(:,2)=(/0.d0,0.d0,0.75d0/)
207    symrel(:,:,3) = genmpm(:,:)
208    nogen=3
209  case (96)                !P43212
210    tnons(:,2)=(/0.5d0,0.5d0,0.75d0/)
211    symrel(:,:,3) = genmpm(:,:)
212    tnons(:,3)=(/0.5d0,0.5d0,0.75d0/)
213    nogen=3
214  case (98)                !I4122
215    tnons(:,2)=(/0.d0,0.5d0,0.25d0/)
216    symrel(:,:,3) = genmpm(:,:)
217    tnons(:,3)=(/0.5d0,0.0d0,0.75d0/)
218    nogen=3
219  case (99,107,123,139)        !P4mm, I4mm, P4/mmm, I4/mmm
220    symrel(:,:,3) = genpmp(:,:)
221    nogen=3
222  case (100)                !P4bm
223    symrel(:,:,3) = genpmp(:,:)
224    tnons(:,3)=(/0.5d0,0.5d0,0.d0/)
225    nogen=3
226  case (101,132)                !P42cm, P42/mcm
227    tnons(:,2)=(/0.d0,0.d0,0.5d0/)
228    symrel(:,:,3) = genpmp(:,:)
229    tnons(:,3)=(/0.d0,0.d0,0.5d0/)
230    nogen=3
231  case (102)                !P42nm
232    tnons(:,2)=(/0.5d0,0.5d0,0.5d0/)
233    symrel(:,:,3) = genpmp(:,:)
234    tnons(:,3)=(/0.5d0,0.5d0,0.5d0/)
235    nogen=3
236  case (103,124)                !P4cc, P4/mcc
237    symrel(:,:,3) = genpmp(:,:)
238    tnons(:,3)=(/0.d0,0.d0,0.5d0/)
239    nogen=3
240  case (104)                !P4nc
241    symrel(:,:,3) = genpmp(:,:)
242    tnons(:,3)=(/0.5d0,0.5d0,0.5d0/)
243    nogen=3
244  case (105,131)                !P42mc, P42/mmc
245    tnons(:,2)=(/0.d0,0.d0,0.5d0/)
246    symrel(:,:,3) = genpmp(:,:)
247    nogen=3
248  case (106,135)                !P42bc, P42/mbc
249    tnons(:,2)=(/0.d0,0.d0,0.5d0/)
250    symrel(:,:,3) = genpmp(:,:)
251    tnons(:,3)=(/0.5d0,0.5d0,0.d0/)
252    nogen=3
253  case (108)                !I4cm
254    symrel(:,:,3) = genpmp(:,:)
255    tnons(:,3)=(/0.d0,0.d0,0.5d0/)
256    nogen=3
257  case (109)                !I41md
258    tnons(:,2)=(/0.d0,0.5d0,0.25d0/)
259    symrel(:,:,3) = genpmp(:,:)
260    symrel(:,:,4) = genmmp(:,:)
261    tnons(:,4)=(/0.5d0,0.5d0,0.5d0/)
262    nogen=4
263  case (110)                !I41cd
264    tnons(:,2)=(/0.d0,0.5d0,0.25d0/)
265    symrel(:,:,3) = genpmp(:,:)
266    tnons(:,3)=(/0.d0,0.d0,0.5d0/)
267    symrel(:,:,4) = genmmp(:,:)
268    tnons(:,4)=(/0.5d0,0.5d0,0.5d0/)
269    nogen=4
270  case (111,121)                !PB42m, IB42m
271    symrel(:,:,2) = gen4m(:,:)
272    symrel(:,:,3) = genmpm(:,:)
273    nogen=3
274  case (112)                !PB42c
275    symrel(:,:,2) = gen4m(:,:)
276    symrel(:,:,3) = genmpm(:,:)
277    tnons(:,3)=(/0.d0,0.d0,0.5d0/)
278    nogen=3
279  case (113)                !PB421m
280    symrel(:,:,2) = gen4m(:,:)
281    symrel(:,:,3) = genmpm(:,:)
282    tnons(:,3)=(/0.5d0,0.5d0,0.d0/)
283    nogen=3
284  case (114)                !PB421c
285    symrel(:,:,2) = gen4m(:,:)
286    symrel(:,:,3) = genmpm(:,:)
287    tnons(:,3)=(/0.5d0,0.5d0,0.5d0/)
288    nogen=3
289  case (115,119)                !PB4m2,IB4m2
290    symrel(:,:,2) = gen4m(:,:)
291    symrel(:,:,3) = genpmp(:,:)
292    nogen=3
293  case (116,120)                !PB4c2, IB4c2
294    symrel(:,:,2) = gen4m(:,:)
295    symrel(:,:,3) = genpmp(:,:)
296    tnons(:,3)=(/0.d0,0.d0,0.5d0/)
297    nogen=3
298  case (117)                !PB4b2
299    symrel(:,:,2) = gen4m(:,:)
300    symrel(:,:,3) = genpmp(:,:)
301    tnons(:,3)=(/0.5d0,0.5d0,0.d0/)
302    nogen=3
303  case (118)                !PB4n2
304    symrel(:,:,2) = gen4m(:,:)
305    symrel(:,:,3) = genpmp(:,:)
306    tnons(:,3)=(/0.5d0,0.5d0,0.5d0/)
307    nogen=3
308  case (122)                !IB42d
309    symrel(:,:,2)=0 ; symrel(1,2,2)=-1 ; symrel(2,1,2)=1 ; symrel(3,3,2)=-1
310    symrel(:,:,3) = genmpm(:,:)
311    tnons(:,3)=(/0.5d0,0.d0,0.75d0/)
312    nogen=3
313  case (125,126,129,130,133,134,137,138,141,142)
314    symrel(:,:,1)=0 ; symrel(1,1,1)=1 ; symrel(2,2,1)=1 ; symrel(3,3,1)=1
315    symrel(:,:,2)=0 ; symrel(1,1,2)=-1; symrel(2,2,2)=-1; symrel(3,3,2)=1
316    symrel(:,:,3)=0 ; symrel(1,2,3)=-1 ; symrel(2,1,3)=1 ; symrel(3,3,3)=1
317    symrel(:,:,4)=0 ; symrel(1,2,4)=1 ; symrel(2,1,4)=-1 ; symrel(3,3,4)=1
318    symrel(:,:,5) = genmpm(:,:)
319    symrel(:,:,6) = genmmm(:,:)
320    if (spgorig==1) then
321      select case (spgroup)
322      case (125)                !P4/nbm
323        tnons(:,6)=(/0.5d0,0.5d0,0.d0/)
324      case (126)         !P4/nnc
325        tnons(:,6)=(/0.5d0,0.5d0,0.5d0/)
326      case (129)                !P4/nmm
327        tnons(:,3)=(/0.5d0,0.5d0,0.d0/) ; tnons(:,4)=(/0.5d0,0.5d0,0.d0/)
328        tnons(:,5)=(/0.5d0,0.5d0,0.d0/) ; tnons(:,6)=(/0.5d0,0.5d0,0.d0/)
329      case (130)                !P4/ncc
330        tnons(:,3)=(/0.5d0,0.5d0,0.d0/) ; tnons(:,4)=(/0.5d0,0.5d0,0.d0/)
331        tnons(:,5)=(/0.5d0,0.5d0,0.5d0/)
332        tnons(:,6)=(/0.5d0,0.5d0,0.d0/)
333      case (133)         !P42/nbc
334        tnons(:,3)=(/0.5d0,0.5d0,0.5d0/) ; tnons(:,4)=(/0.5d0,0.5d0,0.5d0/)
335        tnons(:,5)=(/0.d0,0.d0,0.5d0/)
336        tnons(:,6)=(/0.5d0,0.5d0,0.5d0/)
337      case (134)         !P42/nnm
338        tnons(:,3)=(/0.5d0,0.5d0,0.5d0/) ; tnons(:,4)=(/0.5d0,0.5d0,0.5d0/)
339        tnons(:,6)=(/0.5d0,0.5d0,0.5d0/)
340      case (137)         !P42/nmc
341        tnons(:,3)=(/0.5d0,0.5d0,0.5d0/) ; tnons(:,4)=(/0.5d0,0.5d0,0.5d0/)
342        tnons(:,5)=(/0.5d0,0.5d0,0.5d0/) ; tnons(:,6)=(/0.5d0,0.5d0,0.5d0/)
343      case (138)         !P42/ncm
344        tnons(:,3)=(/0.5d0,0.5d0,0.5d0/) ; tnons(:,4)=(/0.5d0,0.5d0,0.5d0/)
345        tnons(:,5)=(/0.5d0,0.5d0,0.d0/) ; tnons(:,6)=(/0.5d0,0.5d0,0.5d0/)
346      case (141)         !I41/amd
347        tnons(:,2)=(/0.5d0,0.5d0,0.5d0/)
348        tnons(:,3)=(/0.d0,0.5d0,0.25d0/)
349        tnons(:,4)=(/0.5d0,0.d0,0.75d0/) ; tnons(:,5)=(/0.5d0,0.d0,0.75d0/)
350        tnons(:,6)=(/0.d0,0.5d0,0.25d0/)
351      case (142)         !I41/acd
352        tnons(:,2)=(/0.5d0,0.5d0,0.5d0/)
353        tnons(:,3)=(/0.d0,0.5d0,0.25d0/)
354        tnons(:,4)=(/0.5d0,0.d0,0.75d0/)
355        tnons(:,5)=(/0.5d0,0.d0,0.25d0/)
356        tnons(:,6)=(/0.d0,0.5d0,0.25d0/)
357      end select
358    else
359      select case (spgroup)
360      case (125)         !P4/nbm
361        tnons(:,2)=(/0.5d0,0.5d0,0.d0/) ; tnons(:,3)=(/0.5d0,0.d0,0.d0/)
362        tnons(:,4)=(/0.0d0,0.5d0,0.d0/) ; tnons(:,5)=(/0.5d0,0.d0,0.d0/)
363      case (126)                !P4/nnc
364        tnons(:,2)=(/0.5d0,0.5d0,0.d0/) ; tnons(:,3)=(/0.5d0,0.d0,0.d0/)
365        tnons(:,4)=(/0.0d0,0.5d0,0.d0/) ; tnons(:,5)=(/0.5d0,0.d0,0.5d0/)
366      case (129)                !P4/nmm
367        tnons(:,2)=(/0.5d0,0.5d0,0.d0/) ; tnons(:,3)=(/0.5d0,0.d0,0.d0/)
368        tnons(:,4)=(/0.0d0,0.5d0,0.d0/) ;  tnons(:,5)=(/0.0d0,0.5d0,0.d0/)
369      case (130)         !P4/ncc
370        tnons(:,2)=(/0.5d0,0.5d0,0.d0/) ;  tnons(:,3)=(/0.5d0,0.d0,0.d0/)
371        tnons(:,4)=(/0.0d0,0.5d0,0.d0/) ; tnons(:,5)=(/0.0d0,0.5d0,0.5d0/)
372      case (133)        !P42/nbc
373        tnons(:,2)=(/0.5d0,0.5d0,0.d0/) ; tnons(:,3)=(/0.5d0,0.d0,0.5d0/)
374        tnons(:,4)=(/0.0d0,0.5d0,0.5d0/) ; tnons(:,5)=(/0.5d0,0.d0,0.d0/)
375      case (134)                !P42/nnm
376        tnons(:,2)=(/0.5d0,0.5d0,0.d0/) ; tnons(:,3)=(/0.5d0,0.d0,0.5d0/)
377        tnons(:,4)=(/0.0d0,0.5d0,0.5d0/) ; tnons(:,5)=(/0.5d0,0.d0,0.5d0/)
378      case (137)         !P42/nmc
379        tnons(:,2)=(/0.5d0,0.5d0,0.d0/) ; tnons(:,3)=(/0.5d0,0.d0,0.5d0/)
380        tnons(:,4)=(/0.0d0,0.5d0,0.5d0/) ; tnons(:,5)=(/0.d0,0.5d0,0.d0/)
381      case (138)         !P42/ncm
382        tnons(:,2)=(/0.5d0,0.5d0,0.d0/) ; tnons(:,3)=(/0.5d0,0.d0,0.5d0/)
383        tnons(:,4)=(/0.0d0,0.5d0,0.5d0/) ; tnons(:,5)=(/0.d0,0.5d0,0.5d0/)
384      case (141)         !I41/amd
385        tnons(:,2)=(/0.5d0,0.d0,0.5d0/) ; tnons(:,3)=(/0.25d0,0.75d0,0.25d0/)
386        tnons(:,4)=(/0.25d0,0.25d0,0.75d0/) ;  tnons(:,5)=(/0.5d0,0.d0,0.5d0/)
387      case (142)         !I41/acd
388        tnons(:,2)=(/0.5d0,0.d0,0.5d0/) ; tnons(:,3)=(/0.25d0,0.75d0,0.25d0/)
389        tnons(:,4)=(/0.25d0,0.25d0,0.75d0/) ; tnons(:,5)=(/0.5d0,0.d0,0.d0/)
390      end select
391    end if
392    nogen=6
393  case (127)                !P4/mbm
394    symrel(:,:,3) = genmpm(:,:)
395    tnons(:,3)=(/0.5d0,0.5d0,0.d0/)
396    nogen=3
397  case (128)                !P4/mnc
398    symrel(:,:,3) = genmpm(:,:)
399    tnons(:,3)=(/0.5d0,0.5d0,0.5d0/)
400    nogen=3
401  case (136)                !P42/mnm
402    tnons(:,2)=(/0.5d0,0.5d0,0.5d0/)
403    symrel(:,:,3) = genmpm(:,:)
404    tnons(:,3)=(/0.5d0,0.5d0,0.5d0/)
405    nogen=3
406  case (140)                !I4/mcm
407    symrel(:,:,3) = genmpm(:,:)
408    tnons(:,3)=(/0.d0,0.d0,0.5d0/)
409    nogen=3
410  end select
411 
412  if(shubnikov==3)then
413    select case(spgroupma)
414    case(3,9,15,21,27,31,45,47,53,55,77,79,&
415 &     89,97,105,113,121,129,137,145,153,159,166,174,182,190,198,206,&
416 &     214,222,230,236,242,248,254,262,270,278,285,293,301,309,317,&
417 &     323,330,336,343,346,355,358,391,392,403,404,&
418 &     439,442,451,454,487,490,499,500,535,&
419 &     538,545,546)
420      symafm(2)=-1
421    case(90,98,106,114,122,130,138,146,154,160,167,175,183,191,199,&
422 &     207,215,223,231,237,243,249,255,263,271,279,287,295,303,311,319,&
423 &     325,331,337,345,347,357,359,389,393,401,405,441,443,453,455,&
424 &     489,491,497,501,537,539,543,547)
425      symafm(3)=-1
426    case(91,99,107,115,123,131,139,147,155,161,165,173,181,189,&
427 &     197,205,213,221,229,235,241,247,253,261,269,277,286,294,302,310,&
428 &     318,324,329,335,342,344,354,356,390,394,402,406,&
429 &     438,440,450,452,486,488,498,502,534,536,544,548)
430      symafm(2:3)=-1
431    case(365,377,413,425,461,473,509,521)
432      symafm(2)=-1
433      symafm(5:6)=-1
434    case(366,378,414,426,462,474,510,522,554,564)
435      symafm(3:5)=-1
436    case(367,379,415,427,463,475,511,523,555,565)
437      symafm(3:4)=-1
438    case(368,380,416,428,464,476,512,524,556,566)
439      symafm(3:4)=-1
440      symafm(6)=-1
441    case(369,381,417,429,465,477,513,525)
442      symafm(2)=-1
443      symafm(5)=-1
444    case(370,382,418,430,466,478,514,526,558,568)
445      symafm(3:6)=-1
446    case(371,383,419,431,467,479,515,527,559,569)
447      symafm(6)=-1
448    case(553,563)
449      symafm(5:6)=-1
450    case(557,567)
451      symafm(5)=-1
452    end select
453  end if
454 
455 !DEBUG
456 !write(std_out,*)' symsgtetra : symrel(:,:,6)=',symrel(:,:,6)
457 !ENDDEBUG
458 
459  if (nogen>1) then
460    call bldgrp(msym,nogen,nsym,symafm,symrel,tnons)
461  end if
462 
463  call spgdata(brvsb,intsb,intsbl,ptintsb,&
464 & ptschsb,schsb,spgaxor,spgroup,sporder,spgorig)
465 
466 !DEBUG
467 !write(std_out,*)' symsgtetra : exit'
468 !do isym=1,nsym
469 !write(std_out,'(i3,2x,9i3,3es13.3,i3)') isym,symrel(:,:,isym),tnons(:,isym),symafm(isym)
470 !end do
471 !ENDDEBUG
472 
473 end subroutine symsgtetra