TABLE OF CONTENTS


ABINIT/symsgcube [ Functions ]

[ Top ] [ Functions ]

NAME

 symsgcube

FUNCTION

 Generate all the symmetry operations starting from the space group symbol
 for the cubic groups (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
 shubnikov= magnetic type of the space group to be generated
 spgaxor = the possible orientation of the axes system
 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

 nsym = the number of symmetry operations
 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,timab

SOURCE

 39 #if defined HAVE_CONFIG_H
 40 #include "config.h"
 41 #endif
 42 
 43 #include "abi_common.h"
 44 
 45 
 46 subroutine symsgcube(msym,nsym,shubnikov,spgaxor,spgorig,spgroup,&
 47 &   spgroupma,symafm,symrel,tnons)
 48 
 49  use defs_basis
 50  use m_profiling_abi
 51 
 52 !This section has been created automatically by the script Abilint (TD).
 53 !Do not modify the following lines by hand.
 54 #undef ABI_FUNC
 55 #define ABI_FUNC 'symsgcube'
 56  use interfaces_18_timing
 57  use interfaces_41_geometry, except_this_one => symsgcube
 58 !End of the abilint section
 59 
 60  implicit none
 61 
 62 !Arguments ------------------------------------
 63 !scalars
 64  integer,intent(in) :: msym,shubnikov,spgaxor,spgorig,spgroup,spgroupma
 65  integer,intent(inout) :: nsym !vz_i
 66 !arrays
 67  integer,intent(inout) :: symafm(msym),symrel(3,3,msym) !vz_i
 68  real(dp),intent(out) :: tnons(3,msym)
 69 
 70 !Local variables -----------------------------
 71 !scalars
 72  integer :: ii,nogen,sporder
 73  character(len=1) :: brvsb
 74  character(len=15) :: intsb,ptintsb,ptschsb,schsb
 75  character(len=35) :: intsbl
 76 !arrays
 77  integer :: gen1(3,3),gen2(3,3),gen3(3,3),gen4(3,3),gen5(3,3),gen6(3,3)
 78  integer :: gen7(3,3),gen8(3,3),gen9(3,3),genmmm(3,3),genmmp(3,3),genmpm(3,3)
 79  integer :: genmpp(3,3),genpmm(3,3),genpmp(3,3),genppm(3,3),genrot(3,3)
 80  integer :: genswm(3,3),genswp(3,3)
 81  real(dp) :: tsec(2)
 82 
 83 !*************************************************************************
 84 
 85 !DEBUG
 86 !write(std_out,*) ' symsgcube : enter with spgroup ',spgroup,' and ',' origin choice ',spgorig
 87 !write(std_out,*) ' msym,nsym=',msym,nsym
 88 !ENDDEBUG
 89 
 90 
 91 !The identity operation belongs to all space groups
 92  symrel(:,:,1)=0 ; symrel(1,1,1)=1 ; symrel(2,2,1)=1 ; symrel(3,3,1)=1
 93 
 94 !Initialize the associated translations matrix to 0
 95  do ii=1,msym
 96    tnons(:,ii)= 0.0d0
 97  end do
 98  nogen=0
 99 
100 !Predefine some generators
101  genmpp(:,:)=0 ; genmpp(1,1)=-1 ; genmpp(2,2)= 1 ; genmpp(3,3)= 1
102  genpmp(:,:)=0 ; genpmp(1,1)= 1 ; genpmp(2,2)=-1 ; genpmp(3,3)= 1
103  genppm(:,:)=0 ; genppm(1,1)= 1 ; genppm(2,2)= 1 ; genppm(3,3)=-1
104  genpmm(:,:)=0 ; genpmm(1,1)= 1 ; genpmm(2,2)=-1 ; genpmm(3,3)=-1
105  genmpm(:,:)=0 ; genmpm(1,1)=-1 ; genmpm(2,2)= 1 ; genmpm(3,3)=-1
106  genmmp(:,:)=0 ; genmmp(1,1)=-1 ; genmmp(2,2)=-1 ; genmmp(3,3)= 1
107  genmmm(:,:)=0 ; genmmm(1,1)=-1 ; genmmm(2,2)=-1 ; genmmm(3,3)=-1
108  genrot(:,:)=0 ; genrot(1,3)=1 ; genrot(3,2)=1 ; genrot(2,1)=1  !reshape((/0,0,1,1,0,0,0,1,0/),(/3,3/),(/0,0/),(/2,1/) )
109  genswm(:,:)=0 ; genswm(2,1)=1 ; genswm(1,2)=1 ; genswm(3,3)=-1 !reshape((/0,1,0,1,0,0,0,0,-1/),(/3,3/),(/0,0/),(/2,1/) )
110  genswp(:,:)=0 ; genswp(2,1)=1 ; genswp(1,2)=1 ; genswp(3,3)=1  !reshape((/0,1,0,1,0,0,0,0,1/),(/3,3/),(/0,0/),(/2,1/) )
111 
112 !Because of problems with the IBM compiler, that does not like reshape
113 !operations, define 9 basic matrices
114  gen1(:,:)=0 ; gen1(1,1)=1
115  gen2(:,:)=0 ; gen2(1,2)=1
116  gen3(:,:)=0 ; gen3(1,3)=1
117  gen4(:,:)=0 ; gen4(2,1)=1
118  gen5(:,:)=0 ; gen5(2,2)=1
119  gen6(:,:)=0 ; gen6(2,3)=1
120  gen7(:,:)=0 ; gen7(3,1)=1
121  gen8(:,:)=0 ; gen8(3,2)=1
122  gen9(:,:)=0 ; gen9(3,3)=1
123 
124 !Default non-magnetic behaviour
125  symafm(1:msym)=1
126 
127 !*************************************************************************
128 
129 !Treat CUBIC groups
130 
131  if(195<=spgroup .and. spgroup<=230)then
132 
133    select case(spgroup)
134    case (195,196,197)        !P23, F23, I23
135      symrel(:,:,2) = genrot(:,:)
136      symrel(:,:,3) = genmmp(:,:)
137      symrel(:,:,4) = genmpm(:,:)
138      symrel(:,:,5) = genpmm(:,:)
139      nogen=5
140    case (200,202,204)                !PmB3, FmB3, ImB3
141      symrel(:,:,2) = genrot(:,:)
142      symrel(:,:,3) = genmmp(:,:)
143      symrel(:,:,4) = genmpm(:,:)
144      nogen=4
145    case (198,199)                !P213, I213
146      symrel(:,:,2) = genrot(:,:)
147      symrel(:,:,3) = genmmp(:,:)
148      tnons(:,3)=(/0.5d0,0.d0,0.5d0/)
149      symrel(:,:,4) = genmpm(:,:)
150      tnons(:,4)=(/0.d0,0.5d0,0.5d0/)
151      symrel(:,:,5) = genpmm(:,:)
152      tnons(:,5)=(/0.5d0,0.5d0,0.d0/)
153      symrel(:,:,6) = gen3(:,:)-gen4(:,:)-gen8(:,:)
154      tnons(:,6)=(/0.5d0,0.5d0,0.d0/)
155      symrel(:,:,7) = -gen3(:,:)-gen4(:,:)+gen8(:,:)
156      tnons(:,7)=(/0.5d0,0.d0,0.5d0/)
157      symrel(:,:,8) = -gen3(:,:)+gen4(:,:)-gen8(:,:)
158      tnons(:,8)=(/0.d0,0.5d0,0.5d0/)
159      symrel(:,:,10) = -gen2(:,:)+gen6(:,:)-gen7(:,:)
160      tnons(:,10)=(/0.d0,0.5d0,0.5d0/)
161      symrel(:,:,11) =  gen2(:,:)-gen6(:,:)-gen7(:,:)
162      tnons(:,11)=(/0.5d0,0.5d0,0.0d0/)
163      symrel(:,:,12) = -gen2(:,:)-gen6(:,:)+gen7(:,:)
164      tnons(:,12)=(/0.5d0,0.d0,0.5d0/)
165      symrel(:,:,9) =  gen2(:,:)+gen6(:,:)+gen7(:,:)
166    case (201)                !Pn-3
167      if (spgorig==1) then
168        symrel(:,:,2) = genrot(:,:)
169        symrel(:,:,3) = genmmp(:,:)
170        symrel(:,:,4) = genmpm(:,:)
171        symrel(:,:,5) = genmmm(:,:)
172        tnons(:,5)=(/0.5d0,0.5d0,0.5d0/)
173        nogen=5
174        if(shubnikov==3)symafm(5)=-1
175        call bldgrp(msym,nogen,nsym,symafm,symrel,tnons)
176      else
177        symrel(:,:,2) = genrot(:,:)
178        symrel(:,:,3) = genmmp(:,:)
179        tnons(:,3)=(/0.5d0,0.5d0,0.d0/)
180        symrel(:,:,4) = genmpm(:,:)
181        tnons(:,4)=(/0.5d0,0.d0,0.5d0/)
182        symrel(:,:,5) = genpmm(:,:)
183        tnons(:,5)=(/0.d0,0.5d0,0.5d0/)
184        symrel(:,:,6) =  gen3(:,:)-gen4(:,:)-gen8(:,:)
185        tnons(:,6)=(/0.d0,0.5d0,0.5d0/)
186        symrel(:,:,7) = -gen3(:,:)-gen4(:,:)+gen8(:,:)
187        tnons(:,7)=(/0.5d0,0.5d0,0.d0/)
188        symrel(:,:,8) = -gen3(:,:)+gen4(:,:)-gen8(:,:)
189        tnons(:,8)=(/0.5d0,0.d0,0.5d0/)
190        symrel(:,:,10) = -gen2(:,:)+gen6(:,:)-gen7(:,:)
191        tnons(:,10)=(/0.5d0,0.d0,0.5d0/)
192        symrel(:,:,11) =  gen2(:,:)-gen6(:,:)-gen7(:,:)
193        tnons(:,11)=(/0.d0,0.5d0,0.5d0/)
194        symrel(:,:,12) = -gen2(:,:)-gen6(:,:)+gen7(:,:)
195        tnons(:,12)=(/0.5d0,0.5d0,0.d0/)
196        symrel(:,:,9) =  gen2(:,:)+gen6(:,:)+gen7(:,:)
197        do ii=1,12
198          symrel(:,:,ii+12)= - symrel(:,:,ii)
199          tnons(:,ii+12)=tnons(:,ii)
200          if(shubnikov==3)symafm(ii+12)=-1
201        end do
202      end if
203      nogen=0
204    case (205,206)                !Pa-3,Ia-3
205      symrel(:,:,2) = genrot(:,:)
206      symrel(:,:,3) = genmmp(:,:)
207      tnons(:,3)=(/0.5d0,0.d0,0.5d0/)
208      symrel(:,:,4) = genmpm(:,:)
209      tnons(:,4)=(/0.d0,0.5d0,0.5d0/)
210      symrel(:,:,5) = genpmm(:,:)
211      tnons(:,5)=(/0.5d0,0.5d0,0.d0/)
212      symrel(:,:,6) =  gen3(:,:)-gen4(:,:)-gen8(:,:)
213      tnons(:,6)=(/0.5d0,0.5d0,0.d0/)
214      symrel(:,:,7) = -gen3(:,:)-gen4(:,:)+gen8(:,:)
215      tnons(:,7)=(/0.5d0,0.d0,0.5d0/)
216      symrel(:,:,8) = -gen3(:,:)+gen4(:,:)-gen8(:,:)
217      tnons(:,8)=(/0.0d0,0.5d0,0.5d0/)
218      symrel(:,:,9) =  gen2(:,:)+gen6(:,:)+gen7(:,:)
219      symrel(:,:,10) = -gen2(:,:)+gen6(:,:)-gen7(:,:)
220      tnons(:,10)=(/0.d0,0.5d0,0.5d0/)
221      symrel(:,:,11) =  gen2(:,:)-gen6(:,:)-gen7(:,:)
222      tnons(:,11)=(/0.5d0,0.5d0,0.d0/)
223      symrel(:,:,12) = -gen2(:,:)-gen6(:,:)+gen7(:,:)
224      tnons(:,12)=(/0.5d0,0.d0,0.5d0/)
225      nogen=0
226    case (203)                !FdB3
227      if (spgorig==1) then
228        symrel(:,:,2) = genrot(:,:)
229        symrel(:,:,3) = genmmp(:,:)
230        symrel(:,:,4) = genmpm(:,:)
231        nogen=4
232        nsym=12
233        call bldgrp(msym,nogen,nsym,symafm,symrel,tnons)
234        do ii=1,12
235          symrel(:,:,ii+12) = - symrel(:,:,ii)
236          tnons(:,ii+12)=(/0.25d0,0.25d0,0.25d0/)
237          if(shubnikov==3)symafm(ii+12)=-1
238        end do
239        nsym=24
240        nogen=0
241      else
242        symrel(:,:,2) = genrot(:,:)
243        symrel(:,:,3) = genmmp(:,:)
244        tnons(:,3)=(/0.25d0,0.25d0,0.d0/)
245        symrel(:,:,4) = genmpm(:,:)
246        tnons(:,4)=(/0.25d0,0.d0,0.25d0/)
247        symrel(:,:,5) = genpmm(:,:)
248        tnons(:,5)=(/0.d0,0.25d0,0.25d0/)
249        symrel(:,:,6) =  gen3(:,:)-gen4(:,:)-gen8(:,:)
250        tnons(:,6)=(/0.d0,0.25d0,0.25d0/)
251        symrel(:,:,7) = -gen3(:,:)-gen4(:,:)+gen8(:,:)
252        tnons(:,7)=(/0.25d0,0.25d0,0.d0/)
253        symrel(:,:,8) = -gen3(:,:)+gen4(:,:)-gen8(:,:)
254        tnons(:,8)=(/0.25d0,0.d0,0.25d0/)
255        symrel(:,:,10) = -gen2(:,:)+gen6(:,:)-gen7(:,:)
256        tnons(:,10)=(/0.25d0,0.d0,0.25d0/)
257        symrel(:,:,11) =  gen2(:,:)-gen6(:,:)-gen7(:,:)
258        tnons(:,11)=(/0.d0,0.25d0,0.25d0/)
259        symrel(:,:,12) = -gen2(:,:)-gen6(:,:)+gen7(:,:)
260        tnons(:,12)=(/0.25d0,0.25d0,0.d0/)
261        symrel(:,:,9) =  gen2(:,:)+gen6(:,:)+gen7(:,:)
262        do ii=1,12
263          symrel(:,:,ii+12) = - symrel(:,:,ii)
264          tnons(:,ii+12)=tnons(:,ii)
265          if(shubnikov==3)symafm(ii+12)=-1
266        end do
267      end if
268      nsym=24
269      nogen=0
270    case (207,209,211)        !P432, F432, I432, PmB3m, FmB3m, ImB3m
271      symrel(:,:,2) = genrot(:,:)
272      symrel(:,:,3) = genmmp(:,:)
273      symrel(:,:,4) = genmpm(:,:)
274      symrel(:,:,5) = genswm(:,:)
275      if(shubnikov==3)symafm(5)=-1
276      nogen=5
277    case (208)                !P4232
278      symrel(:,:,2) = genrot(:,:)
279      symrel(:,:,3) = genmmp(:,:)
280      symrel(:,:,4) = genmpm(:,:)
281      symrel(:,:,5) = genswm(:,:)
282      tnons(:,5)=(/0.5d0,0.5d0,0.5d0/)
283      if(shubnikov==3)symafm(5)=-1
284      nogen=5
285    case (210)                !F4132
286      symrel(:,:,2) = genrot(:,:)
287      symrel(:,:,3) = genmmp(:,:)
288      tnons(:,3)=(/0.d0,0.5d0,0.5d0/)
289      symrel(:,:,4) = genmpm(:,:)
290      tnons(:,4)=(/0.5d0,0.5d0,0.d0/)
291      symrel(:,:,5) = genswm(:,:)
292      tnons(:,5)=(/0.75d0,0.25d0,0.75d0/)
293      if(shubnikov==3)symafm(5)=-1
294      nogen=5
295    case (212)                !P4332
296      symrel(:,:,2) = genrot(:,:)
297      symrel(:,:,3) = genmmp(:,:)
298      tnons(:,3)=(/0.5d0,0.d0,0.5d0/)
299      symrel(:,:,4) = genmpm(:,:)
300      tnons(:,4)=(/0.d0,0.5d0,0.5d0/)
301      symrel(:,:,5) = genswm(:,:)
302      tnons(:,5)=(/0.25d0,0.75d0,0.75d0/)
303      if(shubnikov==3)symafm(5)=-1
304      nogen=5
305    case (213,214)                !P4132, I4132
306      symrel(:,:,2) = genrot(:,:)
307      symrel(:,:,3) = genmmp(:,:)
308      tnons(:,3)=(/0.5d0,0.d0,0.5d0/)
309      symrel(:,:,4) = genmpm(:,:)
310      tnons(:,4)=(/0.d0,0.5d0,0.5d0/)
311      symrel(:,:,5) = genswm(:,:)
312      tnons(:,5)=(/0.75d0,0.25d0,0.25d0/)
313      if(shubnikov==3)symafm(5)=-1
314      nogen=5
315    case (215,216,217)        !PB43m, FB43m, IB43m
316      symrel(:,:,2) = genrot(:,:)
317      symrel(:,:,3) = genmmp(:,:)
318      symrel(:,:,4) = genmpm(:,:)
319      symrel(:,:,5) = genswp(:,:)
320      if(shubnikov==3)symafm(5)=-1
321      nogen=5
322    case (218,219)                !PB43n, FB343c
323      symrel(:,:,2) = genrot(:,:)
324      symrel(:,:,3) = genmmp(:,:)
325      symrel(:,:,4) = genmpm(:,:)
326      symrel(:,:,5) = genswp(:,:)
327      tnons(:,5)=(/0.5d0,0.5d0,0.5d0/)
328      if(shubnikov==3)symafm(5)=-1
329      nogen=5
330    case (220)                !IB43d
331      symrel(:,:,2) = genrot(:,:)
332      symrel(:,:,3) = genmmp(:,:)
333      tnons(:,3)=(/0.5d0,0.d0,0.5d0/)
334      symrel(:,:,4) = genmpm(:,:)
335      tnons(:,4)=(/0.d0,0.5d0,0.5d0/)
336      symrel(:,:,5) = genswp(:,:)
337      tnons(:,5)=(/0.25d0,0.25d0,0.25d0/)
338      if(shubnikov==3)symafm(5)=-1
339      nogen=5
340    case (221,225,229)        !Pm3m,Fm3m,Im3m
341      symrel(:,:,2) = genrot(:,:)
342      symrel(:,:,3) = genmmp(:,:)
343      symrel(:,:,4) = genmpm(:,:)
344      symrel(:,:,5) = genswm(:,:)
345      if(shubnikov==3)then
346        if(spgroupma==94 .or. spgroupma==95 .or. spgroupma==118 .or. &
347 &       spgroupma==119 .or. spgroupma==142 .or. spgroupma==143   )symafm(5)=-1
348      end if
349      nogen=5
350    case (222)                !PnB3n
351      if (spgorig==1) then
352        symrel(:,:,2) = genrot(:,:)
353        symrel(:,:,3) = genmmp(:,:)
354        tnons(:,3)=(/0.5d0,0.5d0,0.d0/)
355        symrel(:,:,4) = genmpm(:,:)
356        tnons(:,4)=(/0.5d0,0.d0,0.5d0/)
357        symrel(:,:,5) = genswm(:,:)
358        tnons(:,5)=(/0.d0,0.d0,0.5d0/)
359        if(shubnikov==3)then
360          if(spgroupma==100 .or. spgroupma==101)symafm(5)=-1
361        end if
362        nogen=5
363        nsym=24
364        call bldgrp(msym,nogen,nsym,symafm,symrel,tnons)
365        do ii=1,24
366          symrel(:,:,ii+24)=-symrel(:,:,ii)
367          tnons(:,ii+24)=tnons(:,ii)
368          if(shubnikov==3)then
369            if(spgroupma==100 .or. spgroupma==102)symafm(ii+24)=-symafm(ii)
370            if(spgroupma==101)symafm(ii+24)=symafm(ii)
371          end if
372        end do
373        nogen=0 ; nsym=48
374      else
375        symrel(:,:,2) = genrot(:,:)
376        symrel(:,:,3) = genmmp(:,:)
377        symrel(:,:,4) = genmpm(:,:)
378        symrel(:,:,5) = genswm(:,:)
379        if(shubnikov==3)then
380          if(spgroupma==100 .or. spgroupma==101)symafm(5)=-1
381        end if
382        nogen=5
383        nsym=24
384        call bldgrp(msym,nogen,nsym,symafm,symrel,tnons)
385        do ii=1,24
386          symrel(:,:,ii+24)=-symrel(:,:,ii)
387          tnons(:,ii+24)=(/0.5d0,0.5d0,0.5d0/)
388          if(shubnikov==3)then
389            if(spgroupma==100 .or. spgroupma==102)symafm(ii+24)=-symafm(ii)
390            if(spgroupma==101)symafm(ii+24)=symafm(ii)
391          end if
392        end do
393        nogen=0 ; nsym=48
394      end if
395    case (223,226)           ! PmB3n, FmB3c
396      symrel(:,:,2) = genrot(:,:)
397      symrel(:,:,3) = genmmp(:,:)
398      symrel(:,:,4) = genmpm(:,:)
399      symrel(:,:,5) = genswm(:,:)
400      tnons(:,5)=(/0.5d0,0.5d0,0.5d0/)
401      if(shubnikov==3)then
402        if(spgroupma==106 .or. spgroupma==124 .or. &
403 &       spgroupma==107 .or. spgroupma==125     )symafm(5)=-1
404      end if
405      nogen=5
406    case (224)                !PnB3m
407      if (spgorig==1) then
408        symrel(:,:,2) = genrot(:,:)
409        symrel(:,:,3) = genmmp(:,:)
410        symrel(:,:,4) = genmpm(:,:)
411        symrel(:,:,5) = genswm(:,:)
412        tnons(:,5)=(/0.5d0,0.5d0,0.5d0/)
413        symrel(:,:,6) = genmmm(:,:)
414        tnons(:,6)=(/0.5d0,0.5d0,0.5d0/)
415      else
416        symrel(:,:,2) = genrot(:,:)
417        symrel(:,:,3) = genmmp(:,:)
418        tnons(:,3)=(/0.5d0,0.5d0,0.0d0/)
419        symrel(:,:,4) = genmpm(:,:)
420        tnons(:,4)=(/0.5d0,0.0d0,0.5d0/)
421        symrel(:,:,5) = genswm(:,:)
422        tnons(:,5)=(/0.5d0,0.5d0,0.d0/)
423        symrel(:,:,6) = genmmm(:,:)
424      end if
425      if(shubnikov==3)then
426        if(spgroupma==112 .or. spgroupma==113)symafm(5)=-1
427        if(spgroupma==112 .or. spgroupma==114)symafm(6)=-1
428      end if
429      nogen=6
430    case (227)                !FdB3m
431      if (spgorig==1) then
432        symrel(:,:,2) = genrot(:,:)
433        symrel(:,:,3) = genmmp(:,:)
434        tnons(:,3)=(/0.d0,0.5d0,0.5d0/)
435        symrel(:,:,4) = genmpm(:,:)
436        tnons(:,4)=(/0.5d0,0.5d0,0.d0/)
437        symrel(:,:,5) = genswm(:,:)
438        tnons(:,5)=(/0.75d0,0.25d0,0.75d0/)
439        symrel(:,:,6) = genmmm(:,:)
440        tnons(:,6)=(/0.25d0,0.25d0,0.25d0/)
441      else
442        symrel(:,:,2) = genrot(:,:)
443        symrel(:,:,3) = genmmp(:,:)
444        tnons(:,3)=(/0.75d0,0.25d0,0.5d0/)
445        symrel(:,:,4) = genmpm(:,:)
446        tnons(:,4)=(/0.25d0,0.5d0,0.75d0/)
447        symrel(:,:,5) = genswm(:,:)
448        tnons(:,5)=(/0.75d0,0.25d0,0.5d0/)
449        symrel(:,:,6) = genmmm(:,:)
450      end if
451      if(shubnikov==3)then
452        if(spgroupma==130 .or. spgroupma==131)symafm(5)=-1
453        if(spgroupma==130 .or. spgroupma==132)symafm(6)=-1
454      end if
455      nogen=6
456    case (228)                !FdB3c
457      if (spgorig==1) then
458        symrel(:,:,2) = genrot(:,:)
459        symrel(:,:,3) = genmmp(:,:)
460        tnons(:,3)=(/0.d0,0.5d0,0.5d0/)
461        symrel(:,:,4) = genmpm(:,:)
462        tnons(:,4)=(/0.5d0,0.5d0,0.d0/)
463        symrel(:,:,5) = genswm(:,:)
464        tnons(:,5)=(/0.75d0,0.25d0,0.75d0/)
465        symrel(:,:,6) = genmmm(:,:)
466        tnons(:,6)=(/0.75d0,0.75d0,0.75d0/)
467      else
468        symrel(:,:,2) = genrot(:,:)
469        symrel(:,:,3) = genmmp(:,:)
470        tnons(:,3)=(/0.25d0,0.75d0,0.5d0/)
471        symrel(:,:,4) = genmpm(:,:)
472        tnons(:,4)=(/0.75d0,0.5d0,0.25d0/)
473        symrel(:,:,5) = genswm(:,:)
474        tnons(:,5)=(/0.75d0,0.25d0,0.d0/)
475        symrel(:,:,6) = genmmm(:,:)
476      end if
477      if(shubnikov==3)then
478        if(spgroupma==136 .or. spgroupma==137)symafm(5)=-1
479        if(spgroupma==136 .or. spgroupma==138)symafm(6)=-1
480      end if
481      nogen=6
482    case (230)                !IaB3d
483      symrel(:,:,2) = genrot(:,:)
484      symrel(:,:,3) = genmmp(:,:)
485      tnons(:,3)=(/0.5d0,0.d0,0.5d0/)
486      symrel(:,:,4) = genmpm(:,:)
487      tnons(:,4)=(/0.d0,0.5d0,0.5d0/)
488      symrel(:,:,5) = genswm(:,:)
489      tnons(:,5)=(/0.75d0,0.25d0,0.25d0/)
490      if(shubnikov==3)then
491        if(spgroupma==147 .or. spgroupma==148)symafm(5)=-1
492      end if
493      nogen=5
494    end select
495 
496 !  End CUBIC space groups
497  end if
498 
499 !***************************************************************************
500 
501  call timab(47,1,tsec)
502 
503  if (nogen>0) then
504    call bldgrp(msym,nogen,nsym,symafm,symrel,tnons)
505  end if
506 
507  call timab(47,2,tsec)
508 
509  call spgdata(brvsb,intsb,intsbl,ptintsb,&
510 & ptschsb,schsb,spgaxor,spgroup,sporder,spgorig)
511 
512 end subroutine symsgcube