1!
2! Copyright (C) 2006-2014 Quantum ESPRESSO group
3! This file is distributed under the terms of the
4! GNU General Public License. See the file `License'
5! in the root directory of the present distribution,
6! or http://www.gnu.org/copyleft/gpl.txt .
7!
8!
9!-----------------------------------------------------------------------------
10SUBROUTINE divide_class_so( code_group, nrot, smat,d_spin, has_e, nclass, &
11                           nelem,elem, which_irr )
12!-----------------------------------------------------------------------------
13!! This subroutine receives as input a set of nrot 3x3 matrices smat,
14!! and nrot complex 2x2 matrices d_spin, which  are assumed to be the
15!! operations of the point group given by code_group. Only the operations
16!! that do not contain the 2\pi rotation (-E) are given in input.
17!! smat are in cartesian coordinates.
18!! This routine divides the double group in classes and find:
19
20!! * nclass:         the number of classes of the double group;
21!! * nelem(iclass):  for each class, the number of elements of the class;
22!! * elem(i,iclass): 1<i<nelem(iclass) for each class tells which matrices
23!!                   smat belong to that class;
24!! * has_e(i,iclass):   equal to -1 if the operation is multiplied by -E, 1 otherwise;
25!! * which_irr(iclass): for each class gives the position of that class in the
26!!                      character table associated with the group and provided
27!!                      by the routine set_irr_rap_so.
28
29!! NB: changing the order of
30!! the elements in the character table must reflect in
31!! a change to which_irr. Presently the character tables
32!! are those given by G.F. Koster, Space Group and their
33!! representations.
34!! Several equivalent names for the irreducible representation
35!! are given. D, G, L, S are used for Delta, Gamma, Lambda
36!! and Sigma.
37!
38USE kinds, ONLY : DP
39!
40IMPLICIT NONE
41!
42INTEGER :: code_group
43!! The code of the point group
44INTEGER :: nrot
45!! The number of symmetry operation
46INTEGER :: nclass
47!! The number of classes
48INTEGER :: nelem(24)
49!! The elements of each class
50INTEGER :: elem(12,24)
51!! Which elements in the smat list for each class
52INTEGER :: has_e(12,24)
53!! If -1 the element is multiplied by -E
54INTEGER :: which_irr(24)
55!! See main description
56!
57! ... local variables
58!
59REAL(DP) :: smat(3,3,nrot), cmat(3,3), ax(3), bx(3), cx(3)
60REAL(DP) :: smate(3,3,2*nrot)
61COMPLEX(DP) :: d_spin(2,2,48), d_spine(2,2,96), c_spin(2,2)
62!
63INTEGER :: done(96), irot, jrot, krot, iclass, i, other, other1
64INTEGER :: tipo_sym, set_e, ipol, axis, axis1, axis2, ts, nused, iaxis(4), &
65           iax, ibx, icx, aclass, bclass, cclass,  &
66           imax, imbx, imcx, amclass, bmclass, cmclass, ind2(3)
67REAL(DP), PARAMETER :: eps = 1.d-7
68REAL(DP) :: angle_rot, angle_rot_s, ars, ax_save(3,3:5), angle_vectors
69LOGICAL :: compare_mat_so, is_axis, isok, isok1, is_parallel
70LOGICAL :: done_ax(6)
71!
72! Divide the group in classes.
73!
74DO irot=1,nrot
75   smate(:,:,irot)=smat(:,:,irot)
76   smate(:,:,irot+nrot)=smat(:,:,irot)
77   d_spine(:,:,irot)=d_spin(:,:,irot)
78   d_spine(:,:,irot+nrot)=-d_spin(:,:,irot)
79END DO
80!
81!  If there are doubts that the input matrices are not a point group uncomment
82!  the call to this routine.
83!
84!CALL check_tgroup(2*nrot,d_spine,smate)
85!
86!
87nclass=0
88nelem=0
89done=0
90DO irot=1,2*nrot
91   IF (done(irot)==0) THEN
92      nclass=nclass+1
93      DO jrot=1,2*nrot
94         CALL coniug_mat_so(smate(1,1,jrot),d_spine(1,1,jrot), &
95                            smate(1,1,irot),d_spine(1,1,irot), &
96                            cmat,c_spin)
97         DO krot=1,2*nrot
98            IF (compare_mat_so(cmat,c_spin,smate(1,1,krot),d_spine(1,1,krot)) &
99                          .AND.done(krot)==0) THEN
100               nelem(nclass)=nelem(nclass)+1
101               IF (krot.le.nrot) THEN
102                  elem(nelem(nclass),nclass)=krot
103                  has_e(nelem(nclass),nclass)=1
104               ELSE
105                  elem(nelem(nclass),nclass)=krot-nrot
106                  has_e(nelem(nclass),nclass)=-1
107               ENDIF
108               done(krot)=1
109            ENDIF
110         ENDDO
111      ENDDO
112   ENDIF
113ENDDO
114!
115!  For each class we should now decide which_irr. This depends on the group
116!  and on the tables of characters of the irreducible representation.
117!
118which_irr(1)=1
119IF (code_group==1) THEN
120!
121!  C_1
122!
123   IF (nclass /= 2) CALL errore('divide_class_so','Wrong classes for C_1',1)
124   which_irr(2)=2
125ELSEIF (code_group==2.OR.code_group==3.OR.code_group==4) THEN
126!
127!  C_i, C_s, C_2
128!
129   IF (nclass /= 4) &
130            CALL errore('divide_class_so','Wrong classes for C_i, C_s or C_2',1)
131
132   DO iclass=2,nclass
133      IF (tipo_sym(smat(1,1,elem(1,iclass)))==1) THEN
134         which_irr(iclass)=2
135      ELSE
136         which_irr(iclass)=set_e(has_e(1,iclass),3)
137      END IF
138   END DO
139ELSEIF (code_group==5) THEN
140!
141!  C_3
142!
143! The function angle_rot(smat) provides the rotation angle of the matrix smat
144!
145   IF (nclass /= 6) CALL errore('divide_class_so','Wrong classes for C_3',1)
146   DO iclass=2,nclass
147      IF (tipo_sym(smat(1,1,elem(1,iclass)))==1) THEN
148         which_irr(iclass)=2
149      ELSE
150         IF (ABS(angle_rot(smat(1,1,elem(1,iclass)))-120.d0)<eps) THEN
151            which_irr(iclass)=set_e(has_e(1,iclass),3)
152         ELSE
153            which_irr(iclass)=set_e(has_e(1,iclass),5)
154         ENDIF
155      ENDIF
156   ENDDO
157ELSEIF (code_group==6) THEN
158!
159!  C_4
160!
161   IF (nclass /= 8) CALL errore('divide_class_so','Wrong classes for C_4',1)
162   DO iclass=2,nclass
163      ts=tipo_sym(smat(1,1,elem(1,iclass)))
164      IF (ts==1) THEN
165         which_irr(iclass)=2
166      ELSEIF (ts==4) THEN
167         which_irr(iclass)=set_e(has_e(1,iclass),5)
168      ELSEIF (ts==3) THEN
169         ars=angle_rot(smat(1,1,elem(1,iclass)))
170         IF (ABS(ars-90.d0)<eps) THEN
171            which_irr(iclass)=set_e(has_e(1,iclass),3)
172         ELSEIF (ABS(ars-270.d0)<eps) THEN
173            which_irr(iclass)=set_e(has_e(1,iclass),7)
174         ELSE
175            CALL errore('divide_class_so','wrong angle',1)
176         ENDIF
177      ELSE
178         CALL errore('divide_class_so','wrong sym_type',1)
179      ENDIF
180   ENDDO
181ELSEIF (code_group==7) THEN
182!
183!  C_6
184!
185   IF (nclass /= 12) CALL errore('divide_class_so','Wrong classes for C_6',1)
186   DO iclass=2,nclass
187      ts=tipo_sym(smat(1,1,elem(1,iclass)))
188      IF (ts==1) THEN
189         which_irr(iclass)=2
190      ELSE IF (ts==4) THEN
191         which_irr(iclass)=set_e(has_e(1,iclass),7)
192      ELSEIF (ts==3) THEN
193         ars=angle_rot(smat(1,1,elem(1,iclass)))
194         IF (ABS(ars-60.d0)<eps) THEN
195            which_irr(iclass)=set_e(has_e(1,iclass),3)
196         ELSEIF (ABS(ars-120.d0)<eps) THEN
197            which_irr(iclass)=set_e(has_e(1,iclass),5)
198         ELSEIF (ABS(ars-240.d0)<eps) THEN
199            which_irr(iclass)=set_e(has_e(1,iclass),9)
200         ELSEIF (ABS(ars-300.d0)<eps) THEN
201            which_irr(iclass)=set_e(has_e(1,iclass),11)
202         ELSE
203            CALL errore('divide_class_so','wrong angle',1)
204         ENDIF
205      ELSE
206         CALL errore('divide_class_so','wrong sym_type',1)
207      ENDIF
208   ENDDO
209ELSEIF (code_group==8) THEN
210!
211!  D_2
212!
213   IF (nclass /= 5) CALL errore('divide_class_so','Wrong classes for D_2',1)
214!
215!  first search -E
216!
217   nused=1
218   DO iclass=2, nclass
219      ts=tipo_sym(smat(1,1,elem(1,iclass)))
220      IF (ts==1) THEN
221         which_irr(iclass)=2
222      ELSE
223         iaxis(nused)=iclass
224         nused=nused+1
225      ENDIF
226   ENDDO
227
228   CALL versor(smat(1,1,elem(1,iaxis(1))),ax)
229   CALL which_c2(ax,iax)
230   CALL versor(smat(1,1,elem(1,iaxis(2))),bx)
231   CALL which_c2(bx,ibx)
232   CALL versor(smat(1,1,elem(1,iaxis(3))),cx)
233   CALL which_c2(cx,icx)
234
235   CALL is_d2(iax, ibx, icx, ind2)
236
237   which_irr(iaxis(1))=ind2(1)+2
238   which_irr(iaxis(2))=ind2(2)+2
239   which_irr(iaxis(3))=ind2(3)+2
240
241ELSEIF (code_group==9) THEN
242!
243!  D_3
244!
245   IF (nclass /= 6) CALL errore('divide_class_so','Wrong classes for D_3',1)
246   DO iclass=2,nclass
247      ts=tipo_sym(smat(1,1,elem(1,iclass)))
248      IF (ts==1) THEN
249         which_irr(iclass)=2
250      ELSE IF (ts==4) THEN
251         which_irr(iclass)=set_e(has_e(1,iclass),5)
252      ELSEIF (ts==3) THEN
253         which_irr(iclass)=set_e(has_e(1,iclass),3)
254      ELSE
255         CALL errore('divide_class_so','wrong sym_type',1)
256      ENDIF
257   ENDDO
258ELSEIF (code_group==10) THEN
259!
260!  D_4
261!
262   IF (nclass /= 7) CALL errore('divide_class_so','Wrong classes for D_4',1)
263   DO iclass=2,nclass
264      ts=tipo_sym(smat(1,1,elem(1,iclass)))
265      IF (ts==3) THEN
266         which_irr(iclass)=set_e(has_e(1,iclass),3)
267         CALL versor(smat(1,1,elem(1,iclass)),ax)
268         axis=0
269         DO ipol=1,3
270            IF (is_axis(ax,ipol)) axis=ipol
271         ENDDO
272         axis1=MOD(ipol,3)+1
273         axis2=MOD(ipol+1,3)+1
274         IF (axis==0) call errore('divide_class_so','unknown D_4 axis ',1)
275      ENDIF
276   END DO
277   DO iclass=2,nclass
278      ts=tipo_sym(smat(1,1,elem(1,iclass)))
279      IF (ts==1) THEN
280         which_irr(iclass)=2
281      ELSEIF (ts==4) THEN
282         CALL versor(smat(1,1,elem(1,iclass)),ax)
283         IF (is_axis(ax,axis)) THEN
284            which_irr(iclass)=5
285         ELSEIF (is_axis(ax,axis1).or.is_axis(ax,axis2)) THEN
286            which_irr(iclass)=6
287         ELSE
288            which_irr(iclass)=7
289         END IF
290      ELSEIF (ts.ne.3) THEN
291         CALL errore('divide_class_so','wrong sym_type for D_4',1)
292      END IF
293   END DO
294ELSEIF (code_group==11) THEN
295!
296!  D_6
297!
298   IF (nclass /= 9) CALL errore('divide_class_so','Wrong classes for D_6',1)
299   DO iclass=2,nclass
300      ts=tipo_sym(smat(1,1,elem(1,iclass)))
301      IF (ts==1) THEN
302         which_irr(iclass)=2
303      ELSE IF (ts==3) THEN
304         ars=angle_rot(smat(1,1,elem(1,iclass)))
305         IF ((ABS(ars-60.d0)<eps).OR.(ABS(ars-300.d0)<eps) ) THEN
306            which_irr(iclass)=set_e(has_e(1,iclass),3)
307         ELSE
308            which_irr(iclass)=set_e(has_e(1,iclass),5)
309         ENDIF
310      ELSEIF (ts==4) THEN
311         CALL versor(smat(1,1,elem(1,iclass)),ax)
312         IF (is_axis(ax,3)) THEN
313            which_irr(iclass)=7
314         ELSE
315            CALL which_c2(ax, iax)
316            IF (iax==1 .OR. iax==10 .OR. iax==11) THEN
317               which_irr(iclass)=8
318            ELSEIF (iax==2 .OR. iax==12 .OR. iax==13) THEN
319               which_irr(iclass)=9
320            ELSE
321               CALL errore('divide_sym_so','D_6, C_2 axis not recognized',1)
322            END IF
323         END IF
324      ELSE
325         CALL errore('divide_class_so','wrong sym_type',1)
326      END IF
327   END DO
328ELSEIF (code_group==12) THEN
329!
330!  C_2v
331!
332   IF (nclass /= 5) CALL errore('divide_class_so','Wrong classes for C_2v',1)
333   iax=0
334   ibx=0
335   icx=0
336   DO iclass=2,nclass
337      ts=tipo_sym(smat(1,1,elem(1,iclass)))
338      IF (ts==1) THEN
339         which_irr(iclass)=2
340      ELSEIF (ts==4) THEN
341         CALL versor(smat(1,1,elem(1,iclass)), ax)
342         CALL which_c2( ax, iax)
343         which_irr(iclass)=3
344      ELSEIF (ts==5) THEN
345         IF (ibx==0) THEN
346            CALL mirror_axis(smat(1,1,elem(1,iclass)), bx)
347            CALL which_c2( bx, ibx)
348            bclass=iclass
349         ELSE
350            CALL mirror_axis(smat(1,1,elem(1,iclass)), bx)
351            CALL which_c2( bx, icx)
352            cclass=iclass
353         ENDIF
354      ENDIF
355   ENDDO
356   CALL is_c2v(iax, ibx, icx, isok)
357   IF (isok) THEN
358      which_irr(bclass)=4
359      which_irr(cclass)=5
360   ELSE
361      CALL is_c2v(iax, icx, ibx, isok1)
362      IF (.NOT.isok1) CALL errore('divide_class_so','problem with C_2v',1)
363      which_irr(bclass)=5
364      which_irr(cclass)=4
365   ENDIF
366
367ELSEIF (code_group==13) THEN
368!
369!  C_3v
370!
371   IF (nclass /= 6) CALL errore('divide_class_so','Wrong classes for C_3v',1)
372   DO iclass=2,nclass
373      ts=tipo_sym(smat(1,1,elem(1,iclass)))
374      IF (ts==1) THEN
375         which_irr(iclass)=2
376      ELSEIF (ts==3) THEN
377         which_irr(iclass)=set_e(has_e(1,iclass),3)
378      ELSEIF (ts==5) THEN
379         which_irr(iclass)=set_e(has_e(1,iclass),5)
380      ELSE
381         CALL errore('divide_class_so','wrong operation',1)
382      ENDIF
383   ENDDO
384ELSEIF (code_group==14) THEN
385!
386!  C_4v
387!
388   IF (nclass /= 7) CALL errore('divide_class_so','Wrong classes for C_4v',1)
389
390   DO iclass=2,nclass
391      ts=tipo_sym(smat(1,1,elem(1,iclass)))
392      IF (ts==1) THEN
393         which_irr(iclass)=2
394      ELSEIF (ts==3) THEN
395         which_irr(iclass)=set_e(has_e(1,iclass),3)
396      ELSEIF (ts==4) THEN
397         which_irr(iclass)=5
398      ELSEIF (ts==5) THEN
399         CALL mirror_axis(smat(1,1,elem(1,iclass)), ax)
400         CALL which_c2(ax, iax)
401         IF (iax < 4) THEN
402            which_irr(iclass)=6
403         ELSE
404            which_irr(iclass)=7
405         ENDIF
406      ENDIF
407   ENDDO
408
409ELSEIF (code_group==15) THEN
410!
411!  C_6v
412!
413   IF (nclass /= 9) CALL errore('divide_class_so','Wrong classes for C_6v',1)
414   DO iclass=2,nclass
415      ts=tipo_sym(smat(1,1,elem(1,iclass)))
416      IF (ts==1) THEN
417         which_irr(iclass)=2
418      ELSE IF (ts==3) THEN
419         ars=angle_rot(smat(1,1,elem(1,iclass)))
420         IF ((ABS(ars-60.d0)<eps).OR.(ABS(ars-300.d0)<eps)) THEN
421            which_irr(iclass)=set_e(has_e(1,iclass),3)
422         ELSE
423            which_irr(iclass)=set_e(has_e(1,iclass),5)
424         ENDIF
425      ELSEIF (ts==4) THEN
426         which_irr(iclass)=7
427      ELSEIF (ts==5) THEN
428         CALL mirror_axis(smat(1,1,elem(1,iclass)), ax)
429         CALL which_c2(ax,iax)
430         IF (iax==2 .OR. iax==12 .OR. iax==13) THEN
431            which_irr(iclass)=9
432         ELSEIF (iax==1 .OR. iax==10 .OR. iax==11) THEN
433            which_irr(iclass)=8
434         ELSE
435            CALL errore('divide_class_so','C_6v mirror symmetry not recognized',1)
436         ENDIF
437      ENDIF
438   ENDDO
439ELSEIF (code_group==16) THEN
440!
441!  C_2h
442!
443   IF (nclass /= 8) CALL errore('divide_class_so','Wrong classes for C_2h',1)
444   DO iclass=2,nclass
445      ts=tipo_sym(smat(1,1,elem(1,iclass)))
446      IF (ts==1) THEN
447         which_irr(iclass)=2
448      ELSEIF (ts==4) THEN
449         which_irr(iclass)=set_e(has_e(1,iclass),3)
450      ELSEIF (ts==2) THEN
451         which_irr(iclass)=set_e(has_e(1,iclass),5)
452      ELSEIF (ts==5) THEN
453         which_irr(iclass)=set_e(has_e(1,iclass),7)
454      ELSE
455         CALL errore('divide_class_so','wrong sym_type',1)
456      ENDIF
457   ENDDO
458ELSEIF (code_group==17) THEN
459!
460!  C_3h
461!
462   IF (nclass /= 12) CALL errore('divide_class_so','Wrong classes for C_3h',1)
463   DO iclass=2,nclass
464      ts=tipo_sym(smat(1,1,elem(1,iclass)))
465      IF (ts==1) THEN
466         which_irr(iclass)=2
467      ELSE IF (ts==3) THEN
468         ars=angle_rot(smat(1,1,elem(1,iclass)))
469         IF (ABS(ars-120.d0)<eps) THEN
470            which_irr(iclass)=set_e(has_e(1,iclass),3)
471         ELSE
472            which_irr(iclass)=set_e(has_e(1,iclass),5)
473         END IF
474      ELSEIF (ts==5) THEN
475         which_irr(iclass)=set_e(has_e(1,iclass),7)
476      ELSEIF (ts==6) THEN
477         IF (ABS(angle_rot_s(smat(1,1,elem(1,iclass)))-120.d0)<eps) THEN
478            which_irr(iclass)=set_e(has_e(1,iclass),9)
479         ELSE
480            which_irr(iclass)=set_e(has_e(1,iclass),11)
481         END IF
482      ELSE
483         CALL errore('divide_class_so','wrong sym_type',1)
484      ENDIF
485   ENDDO
486ELSEIF (code_group==18) THEN
487!
488!  C_4h
489!
490   IF (nclass /= 16) CALL errore('divide_class_so','Wrong classes for C_4h',1)
491   DO iclass=2,nclass
492      ts=tipo_sym(smat(1,1,elem(1,iclass)))
493      IF (ts==1) THEN
494         which_irr(iclass)=2
495      ELSE IF (ts==3) THEN
496         IF (angle_rot(smat(1,1,elem(1,iclass)))-90.d0<eps) THEN
497            which_irr(iclass)=set_e(has_e(1,iclass),3)
498         ELSE
499            which_irr(iclass)=set_e(has_e(1,iclass),7)
500         END IF
501      ELSEIF (ts==4) THEN
502         which_irr(iclass)=set_e(has_e(1,iclass),5)
503      ELSEIF (ts==2) THEN
504         which_irr(iclass)=set_e(has_e(1,iclass),9)
505      ELSEIF (ts==5) THEN
506         which_irr(iclass)=set_e(has_e(1,iclass),13)
507      ELSEIF (ts==6) THEN
508         IF (angle_rot_s(smat(1,1,elem(1,iclass)))-90.d0<eps) THEN
509            which_irr(iclass)=set_e(has_e(1,iclass),15)
510         ELSE
511            which_irr(iclass)=set_e(has_e(1,iclass),11)
512         END IF
513      ELSE
514         CALL errore('divide_class_so','wrong operation',1)
515      END IF
516   END DO
517ELSEIF (code_group==19) THEN
518!
519!  C_6h
520!
521   IF (nclass /= 24) CALL errore('divide_class_so','Wrong classes for C_6h',1)
522   DO iclass=2,nclass
523      ts=tipo_sym(smat(1,1,elem(1,iclass)))
524      IF (ts==1) THEN
525         which_irr(iclass)=2
526      ELSE IF (ts==3) THEN
527         ars=angle_rot(smat(1,1,elem(1,iclass)))
528         IF (ABS(ars-60.d0)<eps) THEN
529            which_irr(iclass)=set_e(has_e(1,iclass),3)
530         ELSEIF (ABS(ars-120.d0)<eps) THEN
531            which_irr(iclass)=set_e(has_e(1,iclass),5)
532         ELSEIF (ABS(ars-240.d0)<eps) THEN
533            which_irr(iclass)=set_e(has_e(1,iclass),9)
534         ELSEIF (ABS(ars-300.d0)<eps) THEN
535            which_irr(iclass)=set_e(has_e(1,iclass),11)
536         END IF
537      ELSEIF (ts==4) THEN
538         which_irr(iclass)=set_e(has_e(1,iclass),7)
539      ELSEIF (ts==2) THEN
540         which_irr(iclass)=set_e(has_e(1,iclass),13)
541      ELSEIF (ts==5) THEN
542         which_irr(iclass)=set_e(has_e(1,iclass),19)
543      ELSEIF (ts==6) THEN
544         ars=angle_rot_s(smat(1,1,elem(1,iclass)))
545         IF (ABS(ars-60.d0)<eps) THEN
546            which_irr(iclass)=set_e(has_e(1,iclass),21)
547         ELSEIF (ABS(ars-120.d0)<eps) THEN
548            which_irr(iclass)=set_e(has_e(1,iclass),23)
549         ELSEIF (ABS(ars-240.d0)<eps) THEN
550            which_irr(iclass)=set_e(has_e(1,iclass),15)
551         ELSEIF (ABS(ars-300.d0)<eps) THEN
552            which_irr(iclass)=set_e(has_e(1,iclass),17)
553         END IF
554      ELSE
555         CALL errore('divide_class_so','wrong operation',1)
556      ENDIF
557   ENDDO
558ELSEIF (code_group==20) THEN
559!
560!  D_2h
561!
562!  mirror_axis gives the normal to the mirror plane
563!
564   IF (nclass /= 10) CALL errore('divide_class_so','Wrong classes for D_2h',1)
565!
566!  First check if the axis are parallel to x, y or z
567!
568   iax=0
569   ibx=0
570   icx=0
571   imax=0
572   imbx=0
573   imcx=0
574   DO iclass=2,nclass
575      ts=tipo_sym(smat(1,1,elem(1,iclass)))
576      IF (ts==1) THEN
577         which_irr(iclass)=2
578      ELSEIF (ts==4) THEN
579         CALL versor(smat(1,1,elem(1,iclass)),ax)
580         IF (iax==0) THEN
581            CALL which_c2(ax, iax)
582            aclass=iclass
583         ELSEIF (ibx==0) THEN
584            CALL which_c2(ax, ibx)
585            bclass=iclass
586         ELSEIF (icx==0) THEN
587            CALL which_c2(ax, icx)
588            cclass=iclass
589         ELSE
590            CALL errore('divide_class_so','D_2h too many C_2 axis',1)
591         ENDIF
592      ELSEIF (ts==2) THEN
593         which_irr(iclass)=set_e(has_e(1,iclass),6)
594      ELSEIF (ts==5) THEN
595         CALL mirror_axis(smat(1,1,elem(1,iclass)),ax)
596         IF (imax==0) THEN
597            CALL which_c2(ax, imax)
598            amclass=iclass
599         ELSEIF (imbx==0) THEN
600            CALL which_c2(ax, imbx)
601            bmclass=iclass
602         ELSEIF (imcx==0) THEN
603            CALL which_c2(ax, imcx)
604            cmclass=iclass
605         ELSE
606            CALL errore('divide_class_so','D_2h too many mirrors',1)
607         ENDIF
608      ELSE
609         CALL errore('divide_class_so','D_2h operation not recognized',1)
610      ENDIF
611   ENDDO
612
613   CALL is_d2( iax, ibx, icx, ind2)
614
615   which_irr(aclass)=ind2(1)+2
616   which_irr(bclass)=ind2(2)+2
617   which_irr(cclass)=ind2(3)+2
618
619   IF (imax==iax) which_irr(amclass) = which_irr(aclass) + 5
620   IF (imax==ibx) which_irr(amclass) = which_irr(bclass) + 5
621   IF (imax==icx) which_irr(amclass) = which_irr(cclass) + 5
622   IF (imbx==iax) which_irr(bmclass) = which_irr(aclass) + 5
623   IF (imbx==ibx) which_irr(bmclass) = which_irr(bclass) + 5
624   IF (imbx==icx) which_irr(bmclass) = which_irr(cclass) + 5
625   IF (imcx==iax) which_irr(cmclass) = which_irr(aclass) + 5
626   IF (imcx==ibx) which_irr(cmclass) = which_irr(bclass) + 5
627   IF (imcx==icx) which_irr(cmclass) = which_irr(cclass) + 5
628
629ELSEIF (code_group==21) THEN
630!
631!  D_3h
632!
633   IF (nclass /= 9) CALL errore('divide_class_so','Wrong classes for D_3h',1)
634   DO iclass=2,nclass
635      ts=tipo_sym(smat(1,1,elem(1,iclass)))
636      IF (ts==1) THEN
637         which_irr(iclass)=2
638      ELSE IF (ts==3) THEN
639         which_irr(iclass)=set_e(has_e(1,iclass),3)
640      ELSE IF (ts==4) THEN
641         which_irr(iclass)=5
642      ELSE IF (ts==5) THEN
643         IF (nelem(iclass)>1) THEN
644            which_irr(iclass)=9
645         ELSE
646            which_irr(iclass)=6
647         END IF
648      ELSE IF (ts==6) THEN
649         which_irr(iclass)=set_e(has_e(1,iclass),7)
650      END IF
651   END DO
652ELSEIF (code_group==22) THEN
653!
654!  D_4h
655!
656!  First search the order 4 axis
657!
658   IF (nclass /= 14) CALL errore('divide_class_so','Wrong classes for D_4h',1)
659   DO iclass=2,nclass
660      ts=tipo_sym(smat(1,1,elem(1,iclass)))
661      IF (ts==3) THEN
662         which_irr(iclass)=set_e(has_e(1,iclass),3)
663         CALL versor(smat(1,1,elem(1,iclass)),ax)
664         axis=0
665         DO ipol=1,3
666            IF (is_axis(ax,ipol)) axis=ipol
667         ENDDO
668         IF (axis==0) call errore('divide_class_so','unknown D_4h axis ',1)
669      ENDIF
670   END DO
671   DO iclass=2,nclass
672      ts=tipo_sym(smat(1,1,elem(1,iclass)))
673      IF (ts==1) THEN
674         which_irr(iclass)=2
675      ELSE IF (ts==4) THEN
676         which_irr(iclass)=0
677         CALL versor(smat(1,1,elem(1,iclass)),ax)
678         IF (is_axis(ax,axis)) THEN
679            which_irr(iclass)=5
680         ELSE
681           DO ipol=1,3
682              IF (is_axis(ax,ipol)) which_irr(iclass)=6
683           ENDDO
684           IF (which_irr(iclass)==0) which_irr(iclass)=7
685         END IF
686      ELSEIF (ts==2) THEN
687         which_irr(iclass)=set_e(has_e(1,iclass),8)
688      ELSEIF (ts==5) THEN
689         which_irr(iclass)=0
690         CALL mirror_axis(smat(1,1,elem(1,iclass)),ax)
691         IF (is_axis(ax,axis)) THEN
692            which_irr(iclass)=12
693         ELSE
694            DO ipol=1,3
695               IF (is_axis(ax,ipol)) which_irr(iclass)=13
696            ENDDO
697            IF (which_irr(iclass)==0) which_irr(iclass)=14
698         END IF
699      ELSEIF (ts==6) THEN
700         which_irr(iclass)=set_e(has_e(1,iclass),10)
701      END IF
702   END DO
703ELSEIF (code_group==23) THEN
704!
705!  D_6h
706!
707   IF (nclass /= 18) CALL errore('divide_class_so','Wrong classes for D_6h',1)
708   DO iclass=2,nclass
709      ts=tipo_sym(smat(1,1,elem(1,iclass)))
710      IF (ts==1) THEN
711         which_irr(iclass)=2
712      ELSE IF (ts==3) THEN
713         ars=angle_rot(smat(1,1,elem(1,iclass)))
714         IF ((ABS(ars-60.d0)<eps).OR.(ABS(ars-300.d0)<eps)) THEN
715            which_irr(iclass)=set_e(has_e(1,iclass),3)
716         ELSE
717            which_irr(iclass)=set_e(has_e(1,iclass),5)
718         END IF
719      ELSE IF (ts==4) THEN
720         IF (nelem(iclass)==2) THEN
721            which_irr(iclass)=7
722         ELSE
723            CALL versor(smat(1,1,elem(1,iclass)),ax)
724            CALL which_c2(ax,iax)
725            IF (iax==1 .OR. iax==10 .OR. iax==11) THEN
726               which_irr(iclass)=8
727            ELSEIF (iax==2 .OR. iax==12 .OR. iax==13) THEN
728               which_irr(iclass)=9
729            ELSE
730               CALL errore('divide_class_so','Problem with C_2 of D_6h',1)
731            ENDIF
732         END IF
733      ELSE IF (ts==2) THEN
734         which_irr(iclass)=set_e(has_e(1,iclass),10)
735      ELSE IF (ts==5) THEN
736          IF (nelem(iclass)==2) THEN
737             which_irr(iclass)=16
738          ELSE
739             CALL mirror_axis(smat(1,1,elem(1,iclass)),ax)
740             CALL which_c2(ax, iax)
741             IF (iax==1 .OR. iax==10 .OR. iax==11) THEN
742!
743!   sigma_d
744!
745                which_irr(iclass)=17
746             ELSEIF (iax==2 .OR. iax==12 .OR. iax==13) THEN
747!
748!   sigma_v
749!
750                which_irr(iclass)=18
751             ELSE
752               CALL errore('divide_class_so','Problem with mirror of D_6h',1)
753             ENDIF
754          END IF
755      ELSE IF (ts==6) THEN
756         ars=angle_rot_s(smat(1,1,elem(1,iclass)))
757         IF ((ABS(ars-60.d0)<eps).OR.(ABS(ars-300.d0)<eps)) THEN
758            which_irr(iclass)=set_e(has_e(1,iclass),14)
759         ELSE
760            which_irr(iclass)=set_e(has_e(1,iclass),12)
761         END IF
762      END IF
763   END DO
764ELSEIF (code_group==24) THEN
765!
766!  D_2d
767!
768   IF (nclass /= 7) CALL errore('divide_class_so','Wrong classes for D_2d',1)
769   DO iclass=2,nclass
770      ts=tipo_sym(smat(1,1,elem(1,iclass)))
771      IF (ts==1) THEN
772          which_irr(iclass)=2
773      ELSE IF (ts==6) THEN
774         which_irr(iclass)=set_e(has_e(1,iclass),3)
775      ELSE IF (ts==4) THEN
776         IF (nelem(iclass)==2) THEN
777            which_irr(iclass)=5
778         ELSE
779            which_irr(iclass)=6
780         END IF
781      ELSE IF (ts==5) THEN
782         which_irr(iclass)=7
783      ELSE
784         CALL errore('divide_class_so','wrong operation',1)
785      END IF
786   END DO
787ELSEIF (code_group==25) THEN
788!
789!  D_3d
790!
791   IF (nclass /= 12) CALL errore('divide_class_so','Wrong classes for D_3d',1)
792   DO iclass=2,nclass
793      ts=tipo_sym(smat(1,1,elem(1,iclass)))
794      IF (ts==1) THEN
795         which_irr(iclass)=2
796      ELSEIF (ts==3) THEN
797         which_irr(iclass)=set_e(has_e(1,iclass),3)
798      ELSE IF (ts==4) THEN
799         which_irr(iclass)=set_e(has_e(1,iclass),5)
800      ELSE IF (ts==2) THEN
801         which_irr(iclass)=set_e(has_e(1,iclass),7)
802      ELSE IF (ts==6) THEN
803         which_irr(iclass)=set_e(has_e(1,iclass),9)
804      ELSE IF (ts==5) THEN
805         which_irr(iclass)=set_e(has_e(1,iclass),11)
806      ELSE
807         CALL errore('divide_class_so','wrong operation',1)
808      END IF
809   END DO
810ELSEIF (code_group==26) THEN
811!
812!  S_4
813!
814   IF (nclass /= 8) CALL errore('divide_class_so','Wrong classes for S_4',1)
815   DO iclass=2,nclass
816      ts=tipo_sym(smat(1,1,elem(1,iclass)))
817      IF (ts==1) THEN
818         which_irr(iclass)=2
819      ELSE IF (ts==4) THEN
820         which_irr(iclass)=set_e(has_e(1,iclass),5)
821      ELSE IF (ts==6) THEN
822         IF (ABS(angle_rot_s(smat(1,1,elem(1,iclass)))-90.d0)<eps) THEN
823            which_irr(iclass)=set_e(has_e(1,iclass),7)
824         ELSE
825            which_irr(iclass)=set_e(has_e(1,iclass),3)
826         END IF
827      ELSE
828         CALL errore('divide_class_so','wrong operation',1)
829      END IF
830   END DO
831ELSE IF (code_group==27) THEN
832!
833!  S_6
834!
835   IF (nclass /= 12) CALL errore('divide_class_so','Wrong classes for S_6',1)
836   DO iclass=2,nclass
837      ts=tipo_sym(smat(1,1,elem(1,iclass)))
838      IF (ts==1) THEN
839         which_irr(iclass)=2
840      ELSE IF (ts==3) THEN
841         IF (ABS(angle_rot(smat(1,1,elem(1,iclass)))-120.d0)<eps) THEN
842            which_irr(iclass)=set_e(has_e(1,iclass),3)
843         ELSE
844            which_irr(iclass)=set_e(has_e(1,iclass),5)
845         END IF
846      ELSE IF (ts==2) THEN
847         which_irr(iclass)=set_e(has_e(1,iclass),7)
848      ELSE IF (ts==6) THEN
849         IF (ABS(angle_rot_s(smat(1,1,elem(1,iclass)))-60.d0)<eps) THEN
850            which_irr(iclass)=set_e(has_e(1,iclass),11)
851         ELSE
852            which_irr(iclass)=set_e(has_e(1,iclass),9)
853         END IF
854      ELSE
855         CALL errore('divide_class_so','wrong operation',1)
856      END IF
857   END DO
858ELSEIF (code_group==28) THEN
859!
860!  T
861!
862   IF (nclass /= 7) CALL errore('divide_class_so','Wrong classes for T',1)
863   DO iclass=2,nclass
864      ts=tipo_sym(smat(1,1,elem(1,iclass)))
865      IF (ts==1) THEN
866         which_irr(iclass)=2
867      ELSE IF (ts==4) THEN
868         which_irr(iclass)=3
869      ELSE IF (ts==3) THEN
870         IF (ABS(angle_rot(smat(1,1,elem(1,iclass)))-120.d0)<eps) THEN
871            which_irr(iclass)=set_e(has_e(1,iclass),4)
872         ELSE IF (ABS(angle_rot(smat(1,1,elem(1,iclass)))-240.d0)<eps) THEN
873            which_irr(iclass)=set_e(has_e(1,iclass),6)
874         ENDIF
875      ELSE
876         CALL errore('divide_class_so','wrong sym type',1)
877      END IF
878   END DO
879ELSE IF (code_group==29) THEN
880!
881!  T_h
882!
883   IF (nclass /= 14) CALL errore('divide_class_so','Wrong classes for T_h',1)
884   DO iclass=2,nclass
885      ts=tipo_sym(smat(1,1,elem(1,iclass)))
886      IF (ts==1) THEN
887         which_irr(iclass)=2
888      ELSE IF (ts==3) THEN
889         IF (ABS(angle_rot(smat(1,1,elem(1,iclass)))-120.d0)<eps) THEN
890            which_irr(iclass)=set_e(has_e(1,iclass),4)
891         ELSE
892            which_irr(iclass)=set_e(has_e(1,iclass),6)
893         END IF
894      ELSE IF (ts==4) THEN
895         which_irr(iclass)=3
896      ELSE IF (ts==2) THEN
897         which_irr(iclass)=set_e(has_e(1,iclass),8)
898      ELSE IF (ts==6) THEN
899         IF (ABS(angle_rot_s(smat(1,1,elem(1,iclass)))-60.d0)<eps) THEN
900            which_irr(iclass)=set_e(has_e(1,iclass),13)
901         ELSE
902            which_irr(iclass)=set_e(has_e(1,iclass),11)
903         END IF
904      ELSE IF (ts==5) THEN
905         which_irr(iclass)=10
906      ELSE
907         CALL errore('divide_class_so','wrong operation',1)
908      END IF
909   END DO
910ELSEIF (code_group==30) THEN
911!
912!  T_d
913!
914   IF (nclass /= 8) CALL errore('divide_class_so','Wrong classes for T_d',1)
915   DO iclass=2,nclass
916      ts=tipo_sym(smat(1,1,elem(1,iclass)))
917      IF (ts==1) THEN
918         which_irr(iclass)=2
919      ELSE IF (ts==3) THEN
920         which_irr(iclass)=set_e(has_e(1,iclass),3)
921      ELSE IF (ts==4) THEN
922         which_irr(iclass)=5
923      ELSE IF (ts==6) THEN
924         which_irr(iclass)=set_e(has_e(1,iclass),6)
925      ELSE IF (ts==5) THEN
926         which_irr(iclass)=8
927      ELSE
928         CALL errore('divide_class_so','wrong operation',1)
929      END IF
930   END DO
931ELSEIF (code_group==31) THEN
932!
933!  O
934!
935   IF (nclass /= 8) CALL errore('divide_class_so','Wrong classes for O',1)
936   DO iclass=2,nclass
937      ts=tipo_sym(smat(1,1,elem(1,iclass)))
938      IF (ts==1) THEN
939         which_irr(iclass)=2
940      ELSE IF (ts==4) THEN
941         IF (nelem(iclass)==3) THEN
942            which_irr(iclass)=5
943         ELSE
944            which_irr(iclass)=8
945         ENDIF
946      ELSE IF (ts==3) THEN
947         IF (nelem(iclass)==8) THEN
948            which_irr(iclass)=set_e(has_e(1,iclass),3)
949         ELSE
950            which_irr(iclass)=set_e(has_e(1,iclass),6)
951         END IF
952      ELSE
953         CALL errore('divide_class_so','wrong operation',1)
954      END IF
955   END DO
956ELSEIF (code_group==32) THEN
957!
958!  O_h
959!
960   IF (nclass /= 16) CALL errore('divide_class_so','Wrong classes for O_h',1)
961   DO iclass=2,nclass
962      ts=tipo_sym(smat(1,1,elem(1,iclass)))
963      IF (ts==4) THEN
964         IF (nelem(iclass)==6) THEN
965            which_irr(iclass)=5
966         ELSE
967            which_irr(iclass)=8
968         END IF
969      ELSE IF (ts==3) THEN
970         IF (nelem(iclass)==8) THEN
971            which_irr(iclass)=set_e(has_e(1,iclass),3)
972         ELSE
973            which_irr(iclass)=set_e(has_e(1,iclass),6)
974         END IF
975      ELSE IF (ts==2) THEN
976         which_irr(iclass)=set_e(has_e(1,iclass),9)
977      ELSE IF (ts==5) THEN
978         IF (nelem(iclass)==12) THEN
979            which_irr(iclass)=16
980         ELSE
981            which_irr(iclass)=13
982         END IF
983      ELSE IF (ts==6) THEN
984         IF (nelem(iclass)==8) THEN
985            which_irr(iclass)=set_e(has_e(1,iclass),11)
986         ELSE
987            which_irr(iclass)=set_e(has_e(1,iclass),14)
988         END IF
989      ELSE IF (ts==1) THEN
990           which_irr(iclass)=2
991      ELSE
992         CALL errore('divide_class','wrong operation',1)
993      END IF
994   ENDDO
995ELSE
996 CALL errore('divide_class_so','code_group not correct',1)
997ENDIF
998!
999RETURN
1000!
1001END SUBROUTINE divide_class_so
1002
1003
1004!-----------------------------------------------------------------------------
1005SUBROUTINE coniug_mat_so( a, a_spin, b, b_spin, c, c_spin )
1006!-----------------------------------------------------------------------------
1007USE kinds, ONLY : DP
1008!
1009IMPLICIT NONE
1010!
1011REAL(DP) :: a(3,3), b(3,3), c(3,3)
1012COMPLEX(DP) :: a_spin(2,2), b_spin(2,2), c_spin(2,2)
1013!
1014 c=MATMUL(a,MATMUL(b,TRANSPOSE(a)))
1015 c_spin=MATMUL(a_spin,MATMUL(b_spin,TRANSPOSE(CONJG(a_spin))))
1016!
1017RETURN
1018!
1019END SUBROUTINE coniug_mat_so
1020!
1021!
1022!-----------------------------------------------------------------------------
1023FUNCTION compare_mat_so( a, a_spin, b, b_spin )
1024!-----------------------------------------------------------------------------
1025!!  This function compare two 3x3 matrices and two 2x2 matrices
1026!!  and returns .true. if they coincide.
1027!
1028USE kinds, ONLY : DP
1029!
1030IMPLICIT NONE
1031!
1032REAL(DP) :: a(3,3), b(3,3), csum
1033REAL(DP), PARAMETER :: eps=1.d-7
1034COMPLEX(DP) :: a_spin(2,2), b_spin(2,2)
1035LOGICAL :: compare_mat_so
1036INTEGER :: i, j
1037!
1038 csum=0.d0
1039DO i=1,2
1040   DO j=1,2
1041      csum=csum+ABS(a_spin(i,j)-b_spin(i,j))
1042   END DO
1043END DO
1044 compare_mat_so=((ABS(MAXVAL(a-b))<eps).AND.             &
1045                 (ABS(MINVAL(a-b))<eps).AND. (csum < eps ))
1046!
1047RETURN
1048END FUNCTION compare_mat_so
1049!
1050!-----------------------------------------------------------------------------
1051SUBROUTINE set_irr_rap_so( code_group, nclass_ref, nrap_ref, char_mat, &
1052                          name_rap, name_class, name_class1 )
1053!-----------------------------------------------------------------------------
1054!! This subroutine collects the character tables of the 32 crystallographic
1055!! double point groups. Various names have been used in the litterature
1056!! to identify D, G, L, S are used for Delta, Gamma, Lambda and Sigma.
1057!
1058USE kinds, ONLY : DP
1059!
1060IMPLICIT NONE
1061!
1062INTEGER :: nclass_ref
1063!! Output: number of classes
1064INTEGER :: nrap_ref
1065!! Output: number of irreducible representation
1066INTEGER :: code_group
1067!! Input: code of the group
1068CHARACTER(LEN=15) :: name_rap(12)
1069!! Output: name of the representations
1070CHARACTER(LEN=5) :: name_class(24)
1071!! Output: name of the classes
1072CHARACTER(LEN=5) :: name_class1(24)
1073!! Output: name of the classes
1074COMPLEX(DP) :: char_mat(12,24)
1075!! Output: character matrix
1076!
1077! ... local variables
1078!
1079REAL(DP) :: sqr3d2, sqrt2, sqrt3, dsq2
1080!
1081sqrt2 =SQRT(2.d0)
1082sqrt3 =SQRT(3.d0)
1083sqr3d2=sqrt3*0.5d0
1084dsq2  =sqrt2*0.5d0
1085
1086 char_mat=(1.d0,0.d0)
1087 char_mat(:,2)=(-1.d0,0.d0)
1088
1089name_class1="     "
1090name_class(1)="E   "
1091name_class(2)="-E  "
1092
1093
1094IF (code_group==1) THEN
1095!
1096! C_1
1097!
1098   nclass_ref=2
1099   nrap_ref=1
1100
1101   name_rap(1)="G_2  "
1102
1103ELSEIF (code_group==2) THEN
1104!
1105! C_i
1106!
1107   nclass_ref=4
1108
1109   name_class(3)="i    "
1110   name_class(4)="-i   "
1111
1112   nrap_ref=2
1113
1114   name_rap(1)="G_2+"
1115   char_mat(1,4)=(-1.d0,0.d0)
1116
1117   name_rap(2)="G_2-"
1118   char_mat(2,3)=(-1.d0,0.d0)
1119
1120ELSEIF (code_group==3) THEN
1121!
1122! C_s
1123!
1124   nclass_ref=4
1125   name_class(3)="s    "
1126   name_class(4)="-s   "
1127
1128   nrap_ref=2
1129
1130   name_rap(1)="G_3   "
1131   char_mat(1,3)=(0.d0, 1.d0)
1132   char_mat(1,4)=(0.d0,-1.d0)
1133
1134   name_rap(2)="G_4   "
1135   char_mat(2,3)=(0.d0,-1.d0)
1136   char_mat(2,4)=(0.d0, 1.d0)
1137
1138ELSEIF (code_group==4) THEN
1139!
1140! C_2
1141!
1142   nclass_ref=4
1143   name_class(3)="C2   "
1144   name_class(4)="-C2  "
1145
1146   nrap_ref=2
1147
1148   name_rap(1)="G_3   "
1149   char_mat(1,3)=(0.d0, 1.d0)
1150   char_mat(1,4)=(0.d0,-1.d0)
1151
1152   name_rap(2)="G_4   "
1153   char_mat(2,3)=(0.d0,-1.d0)
1154   char_mat(2,4)=(0.d0, 1.d0)
1155
1156
1157ELSEIF (code_group==5) THEN
1158!
1159! C_3    NB: The signs of the characters of the classes C3^2 -C3^2
1160!            are changed with respect to Koster, Space groups and
1161!            their representations. They match the table in Koster,
1162!            Dimmock, Wheeler, Statz, Properties of the 32 point groups.
1163!
1164   nclass_ref=6
1165   name_class(3)="C3   "
1166   name_class(4)="-C3  "
1167   name_class(5)="C3^2 "
1168   name_class(6)="-C3^2"
1169
1170   nrap_ref=3
1171
1172   name_rap(1)="G_4  "
1173   char_mat(1,3)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1174   char_mat(1,4)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1175!   char_mat(1,5)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1176!   char_mat(1,6)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1177   char_mat(1,5)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1178   char_mat(1,6)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1179
1180   name_rap(2)="G_5  "
1181   char_mat(2,3)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1182   char_mat(2,4)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1183!   char_mat(2,5)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1184!   char_mat(2,6)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1185   char_mat(2,5)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1186   char_mat(2,6)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1187
1188   name_rap(3)="G_6  "
1189   char_mat(3,3)=(-1.0d0,0.d0)
1190!   char_mat(3,6)=(-1.0d0,0.d0)
1191   char_mat(3,5)=(-1.0d0,0.d0)
1192
1193
1194ELSEIF (code_group==6) THEN
1195!
1196! C_4   NB: The signs of the characters of the class C4^3 -C4^3
1197!            are changed with respect to Koster, Space groups and
1198!            their representation.
1199!
1200   nclass_ref=8
1201   name_class(3)="C4   "
1202   name_class(4)="-C4  "
1203   name_class(5)="C2   "
1204   name_class(6)="-C2  "
1205   name_class(7)="C4^3 "
1206   name_class(8)="-C4^3"
1207
1208   nrap_ref=4
1209
1210   name_rap(1)="G_5  "
1211   char_mat(1,3)=CMPLX( dsq2, dsq2,kind=DP)
1212   char_mat(1,4)=CMPLX(-dsq2,-dsq2,kind=DP)
1213   char_mat(1,5)=( 0.d0,1.d0)
1214   char_mat(1,6)=( 0.d0,-1.d0)
1215!   char_mat(1,7)=CMPLX(-dsq2, dsq2,kind=DP)
1216!   char_mat(1,8)=CMPLX( dsq2,-dsq2,kind=DP)
1217   char_mat(1,7)=CMPLX( dsq2,-dsq2,kind=DP)
1218   char_mat(1,8)=CMPLX(-dsq2, dsq2,kind=DP)
1219
1220   name_rap(2)="G_6  "
1221   char_mat(2,3)=CMPLX( dsq2,-dsq2,kind=DP)
1222   char_mat(2,4)=CMPLX(-dsq2, dsq2,kind=DP)
1223   char_mat(2,5)=( 0.d0,-1.d0)
1224   char_mat(2,6)=( 0.d0, 1.d0)
1225!   char_mat(2,7)=CMPLX(-dsq2,-dsq2,kind=DP)
1226!   char_mat(2,8)=CMPLX( dsq2, dsq2,kind=DP)
1227   char_mat(2,7)=CMPLX( dsq2, dsq2,kind=DP)
1228   char_mat(2,8)=CMPLX(-dsq2,-dsq2,kind=DP)
1229
1230   name_rap(3)="G_7  "
1231   char_mat(3,3)=CMPLX(-dsq2,-dsq2,kind=DP)
1232   char_mat(3,4)=CMPLX( dsq2, dsq2,kind=DP)
1233   char_mat(3,5)=( 0.d0,1.d0)
1234   char_mat(3,6)=( 0.d0,-1.d0)
1235!   char_mat(3,7)=CMPLX( dsq2,-dsq2,kind=DP)
1236!   char_mat(3,8)=CMPLX(-dsq2, dsq2,kind=DP)
1237   char_mat(3,7)=CMPLX(-dsq2, dsq2,kind=DP)
1238   char_mat(3,8)=CMPLX( dsq2,-dsq2,kind=DP)
1239
1240   name_rap(4)="G_8  "
1241   char_mat(4,3)=CMPLX(-dsq2, dsq2,kind=DP)
1242   char_mat(4,4)=CMPLX( dsq2,-dsq2,kind=DP)
1243   char_mat(4,5)=( 0.d0,-1.d0)
1244   char_mat(4,6)=( 0.d0, 1.d0)
1245!   char_mat(4,7)=CMPLX( dsq2, dsq2,kind=DP)
1246!   char_mat(4,8)=CMPLX(-dsq2,-dsq2,kind=DP)
1247   char_mat(4,7)=CMPLX(-dsq2,-dsq2,kind=DP)
1248   char_mat(4,8)=CMPLX( dsq2, dsq2,kind=DP)
1249
1250
1251ELSEIF (code_group==7) THEN
1252!
1253! C_6   NB: The signs of characters of the C3^2 -C3^2 C6^5 -C6^5 classes
1254!           are changed with respect to Koster, Space groups and their
1255!           representation. They match the table in Koster, Dimmock,
1256!           Wheeler, Statz, Properties of the 32 point groups.
1257!
1258   nclass_ref=12
1259   name_class(3)="C6  "
1260   name_class(4)="-C6 "
1261   name_class(5)="C3 "
1262   name_class(6)="-C3 "
1263   name_class(7)="C2 "
1264   name_class(8)="-C2 "
1265   name_class(9)="C3^2"
1266   name_class(10)="-C3^2"
1267   name_class(11)="C6^5"
1268   name_class(12)="-C6^5"
1269
1270   nrap_ref=6
1271
1272   name_rap(1)="G_7  "
1273   char_mat(1,2)=(-1.d0, 0.d0)
1274   char_mat(1,3)=CMPLX( sqr3d2, 0.5d0,kind=DP)
1275   char_mat(1,4)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
1276   char_mat(1,5)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1277   char_mat(1,6)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1278   char_mat(1,7)=( 0.d0, 1.d0)
1279   char_mat(1,8)=( 0.d0,-1.d0)
1280!   char_mat(1,9)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1281!   char_mat(1,10)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1282!   char_mat(1,11)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1283!   char_mat(1,12)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1284   char_mat(1,9)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1285   char_mat(1,10)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1286   char_mat(1,11)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1287   char_mat(1,12)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1288
1289   name_rap(2)="G_8  "
1290   char_mat(2,2)=(-1.d0, 0.d0)
1291   char_mat(2,3)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1292   char_mat(2,4)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1293   char_mat(2,5)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1294   char_mat(2,6)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1295   char_mat(2,7)=( 0.d0,-1.d0)
1296   char_mat(2,8)=( 0.d0, 1.d0)
1297!   char_mat(2,9)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1298!   char_mat(2,10)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1299!   char_mat(2,11)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
1300!   char_mat(2,12)=CMPLX( sqr3d2, 0.5d0,kind=DP)
1301   char_mat(2,9)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1302   char_mat(2,10)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1303   char_mat(2,11)=CMPLX( sqr3d2, 0.5d0,kind=DP)
1304   char_mat(2,12)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
1305
1306   name_rap(3)="G_9   "
1307   char_mat(3,2)=(-1.d0, 0.d0)
1308   char_mat(3,3)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
1309   char_mat(3,4)=CMPLX( sqr3d2, 0.5d0,kind=DP)
1310   char_mat(3,5)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1311   char_mat(3,6)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1312   char_mat(3,7)=( 0.d0,-1.d0)
1313   char_mat(3,8)=( 0.d0, 1.d0)
1314!   char_mat(5,9)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1315!   char_mat(5,10)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1316!   char_mat(5,11)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1317!   char_mat(5,12)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1318   char_mat(3,9)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1319   char_mat(3,10)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1320   char_mat(3,11)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1321   char_mat(3,12)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1322
1323   name_rap(4)="G_10  "
1324   char_mat(4,2)=(-1.d0, 0.d0)
1325   char_mat(4,3)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1326   char_mat(4,4)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1327   char_mat(4,5)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1328   char_mat(4,6)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1329   char_mat(4,7)=( 0.d0, 1.d0)
1330   char_mat(4,8)=( 0.d0,-1.d0)
1331!   char_mat(4,9)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1332!   char_mat(4,10)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1333!   char_mat(4,11)=CMPLX( sqr3d2, 0.5d0,kind=DP)
1334!   char_mat(4,12)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
1335   char_mat(4,9)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1336   char_mat(4,10)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1337   char_mat(4,11)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
1338   char_mat(4,12)=CMPLX( sqr3d2, 0.5d0,kind=DP)
1339
1340   name_rap(5)="G_11 "
1341   char_mat(5,2)=(-1.d0, 0.d0)
1342   char_mat(5,3)=( 0.d0, 1.d0)
1343   char_mat(5,4)=( 0.d0,-1.d0)
1344   char_mat(5,5)=(-1.d0, 0.d0)
1345   char_mat(5,7)=( 0.d0,-1.d0)
1346   char_mat(5,8)=( 0.d0, 1.d0)
1347!   char_mat(5,10)=(-1.d0, 0.d0)
1348!   char_mat(5,11)=( 0.d0, 1.d0)
1349!   char_mat(5,12)=( 0.d0,-1.d0)
1350   char_mat(5,9)=(-1.d0, 0.d0)
1351   char_mat(5,11)=( 0.d0,-1.d0)
1352   char_mat(5,12)=( 0.d0, 1.d0)
1353
1354   name_rap(6)="G_12  "
1355   char_mat(6,2)=(-1.d0, 0.d0)
1356   char_mat(6,3)=( 0.d0,-1.d0)
1357   char_mat(6,4)=( 0.d0, 1.d0)
1358   char_mat(6,5)=(-1.d0, 0.d0)
1359   char_mat(6,7)=( 0.d0, 1.d0)
1360   char_mat(6,8)=( 0.d0,-1.d0)
1361!   char_mat(6,10)=(-1.d0, 0.d0)
1362!   char_mat(6,11)=( 0.d0,-1.d0)
1363!   char_mat(6,12)=( 0.d0, 1.d0)
1364   char_mat(6,9)=(-1.d0, 0.d0)
1365   char_mat(6,11)=( 0.d0, 1.d0)
1366   char_mat(6,12)=( 0.d0,-1.d0)
1367
1368ELSEIF (code_group==8) THEN
1369!
1370! D_2
1371!
1372   nclass_ref=5
1373   name_class(3)="C2 "
1374   name_class1(3)="-C2 "
1375   name_class(4)="C2' "
1376   name_class1(4)="-C2' "
1377   name_class(5)="C2''"
1378   name_class1(5)="-C2''"
1379
1380   nrap_ref=1
1381
1382   name_rap(1)="G_5  "
1383   char_mat(1,1)=( 2.d0, 0.d0)
1384   char_mat(1,2)=(-2.d0, 0.d0)
1385   char_mat(1,3)=( 0.d0, 0.d0)
1386   char_mat(1,4)=( 0.d0, 0.d0)
1387   char_mat(1,5)=( 0.d0, 0.d0)
1388
1389ELSEIF (code_group==9) THEN
1390!
1391! D_3
1392!
1393   nclass_ref=6
1394   name_class(3)="2C3  "
1395   name_class(4)="-2C3  "
1396   name_class(5)=" 3C2'"
1397   name_class(6)="-3C2'"
1398
1399   nrap_ref=3
1400
1401   name_rap(1)="G_4  "
1402   char_mat(1,1)=( 2.d0, 0.d0)
1403   char_mat(1,2)=(-2.d0, 0.d0)
1404   char_mat(1,4)=(-1.d0, 0.d0)
1405   char_mat(1,5)=( 0.d0, 0.d0)
1406   char_mat(1,6)=( 0.d0, 0.d0)
1407
1408   name_rap(2)="G_5  "
1409   char_mat(2,2)=(-1.d0, 0.d0)
1410   char_mat(2,3)=(-1.d0, 0.d0)
1411   char_mat(2,5)=( 0.d0, 1.d0)
1412   char_mat(2,6)=( 0.d0,-1.d0)
1413
1414   name_rap(3)="G_6  "
1415   char_mat(3,2)=(-1.d0, 0.d0)
1416   char_mat(3,3)=(-1.d0, 0.d0)
1417   char_mat(3,5)=( 0.d0,-1.d0)
1418   char_mat(3,6)=( 0.d0, 1.d0)
1419
1420
1421ELSEIF (code_group==10) THEN
1422!
1423! D_4
1424!
1425   nclass_ref=7
1426   name_class(3)="2C4  "
1427   name_class(4)="-2C4  "
1428   name_class(5)=" C2 "
1429   name_class1(5)="-C2 "
1430   name_class(6)="2C2'"
1431   name_class1(6)="-2C2'"
1432   name_class(7)="2C2'' "
1433   name_class1(7)="-2C2''"
1434
1435   nrap_ref=2
1436
1437   name_rap(1)="G_6  "
1438   char_mat(1,1)=( 2.d0, 0.d0)
1439   char_mat(1,2)=(-2.d0, 0.d0)
1440   char_mat(1,3)=CMPLX( sqrt2, 0.d0,kind=DP)
1441   char_mat(1,4)=CMPLX(-sqrt2, 0.d0,kind=DP)
1442   char_mat(1,5)=( 0.d0, 0.d0)
1443   char_mat(1,6)=( 0.d0, 0.d0)
1444   char_mat(1,7)=( 0.d0, 0.d0)
1445
1446   name_rap(2)="G_7  "
1447   char_mat(2,1)=( 2.d0, 0.d0)
1448   char_mat(2,2)=(-2.d0, 0.d0)
1449   char_mat(2,3)=CMPLX(-sqrt2, 0.d0,kind=DP)
1450   char_mat(2,4)=CMPLX( sqrt2, 0.d0,kind=DP)
1451   char_mat(2,5)=( 0.d0, 0.d0)
1452   char_mat(2,6)=( 0.d0, 0.d0)
1453   char_mat(2,7)=( 0.d0, 0.d0)
1454
1455
1456ELSEIF (code_group==12) THEN
1457!
1458! C_2v
1459!
1460   nclass_ref=5
1461   name_class(3)=" C2"
1462   name_class1(3)="-C2"
1463   name_class(4)=" s_v"
1464   name_class1(4)="-s_v"
1465   name_class(5)=" s_v'"
1466   name_class1(5)="-s_v'"
1467
1468   nrap_ref=1
1469
1470   name_rap(1)="G_5  D_5"
1471   char_mat(1,1)=( 2.d0, 0.d0)
1472   char_mat(1,2)=(-2.d0, 0.d0)
1473   char_mat(1,3)=( 0.d0, 0.d0)
1474   char_mat(1,4)=( 0.d0, 0.d0)
1475   char_mat(1,5)=( 0.d0, 0.d0)
1476
1477ELSEIF (code_group==13) THEN
1478!
1479! C_3v
1480!
1481   nclass_ref=6
1482   name_class(3)="2C3  "
1483   name_class(4)="-2C3 "
1484   name_class(5)="3s_v"
1485   name_class(6)="-3s_v"
1486
1487   nrap_ref=3
1488
1489   name_rap(1)="G_4   L_6 "
1490   char_mat(1,1)=( 2.d0, 0.d0)
1491   char_mat(1,2)=(-2.d0, 0.d0)
1492   char_mat(1,4)=(-1.d0, 0.d0)
1493   char_mat(1,5)=( 0.d0, 0.d0)
1494   char_mat(1,6)=( 0.d0, 0.d0)
1495
1496   name_rap(2)="G_5   L_4 "
1497   char_mat(2,3)=(-1.d0, 0.d0)
1498   char_mat(2,5)=( 0.d0, 1.d0)
1499   char_mat(2,6)=( 0.d0,-1.d0)
1500
1501   name_rap(3)="G_6   L_5 "
1502   char_mat(3,3)=(-1.d0, 0.d0)
1503   char_mat(3,5)=( 0.d0,-1.d0)
1504   char_mat(3,6)=( 0.d0, 1.d0)
1505
1506ELSEIF (code_group==14) THEN
1507!
1508! C_4v
1509!
1510   nclass_ref=7
1511   name_class(3)="2C4  "
1512   name_class(4)="-2C4 "
1513   name_class(5)=" C2"
1514   name_class1(5)="-C2"
1515   name_class(6)=" 2s_v"
1516   name_class1(6)="-2s_v"
1517   name_class(7)=" 2s_d"
1518   name_class1(7)="-2s_d"
1519
1520   nrap_ref=2
1521
1522   name_rap(1)="G_6   D_6"
1523   char_mat(1,1)=( 2.d0, 0.d0)
1524   char_mat(1,2)=(-2.d0, 0.d0)
1525   char_mat(1,3)=CMPLX( sqrt2, 0.d0,kind=DP)
1526   char_mat(1,4)=CMPLX(-sqrt2, 0.d0,kind=DP)
1527   char_mat(1,5)=( 0.d0, 0.d0)
1528   char_mat(1,6)=( 0.d0, 0.d0)
1529   char_mat(1,7)=( 0.d0, 0.d0)
1530
1531   name_rap(2)="G_7   D_7"
1532   char_mat(2,1)=( 2.d0, 0.d0)
1533   char_mat(2,2)=(-2.d0, 0.d0)
1534   char_mat(2,3)=CMPLX(-sqrt2, 0.d0,kind=DP)
1535   char_mat(2,4)=CMPLX( sqrt2, 0.d0,kind=DP)
1536   char_mat(2,5)=( 0.d0, 0.d0)
1537   char_mat(2,6)=( 0.d0, 0.d0)
1538   char_mat(2,7)=( 0.d0, 0.d0)
1539
1540ELSEIF (code_group==15.OR.code_group==11) THEN
1541!
1542! C_6v, D_6
1543!
1544   nclass_ref=9
1545   name_class(3)=" 2C6"
1546   name_class(4)="-2C6"
1547   name_class(5)=" 2C3 "
1548   name_class(6)="-2C3"
1549   name_class(7)=" C2 "
1550   name_class1(7)="-C2 "
1551   IF (code_group==15) THEN
1552      name_class(8)=" 3s_v"
1553      name_class1(8)="-3s_v"
1554      name_class(9)=" 3s_d"
1555      name_class1(9)="-3s_d"
1556   ELSE
1557      name_class(8)=" 3C2'"
1558      name_class1(8)="-3C2'"
1559      name_class(9)="3C2''"
1560      name_class1(9)="-3C2''"
1561   END IF
1562
1563   nrap_ref=3
1564
1565   name_rap(1)="G_7  "
1566   char_mat(1,1)=( 2.d0, 0.d0)
1567   char_mat(1,2)=(-2.d0, 0.d0)
1568   char_mat(1,3)=CMPLX( sqrt3, 0.d0,kind=DP)
1569   char_mat(1,4)=CMPLX(-sqrt3, 0.d0,kind=DP)
1570   char_mat(1,6)=(-1.d0, 0.d0)
1571   char_mat(1,7)=( 0.d0, 0.d0)
1572   char_mat(1,8)=( 0.d0, 0.d0)
1573   char_mat(1,9)=( 0.d0, 0.d0)
1574
1575   name_rap(2)="G_8  "
1576   char_mat(2,1)=( 2.d0, 0.d0)
1577   char_mat(2,2)=(-2.d0, 0.d0)
1578   char_mat(2,3)=CMPLX(-sqrt3, 0.d0,kind=DP)
1579   char_mat(2,4)=CMPLX( sqrt3, 0.d0,kind=DP)
1580   char_mat(2,6)=(-1.d0, 0.d0)
1581   char_mat(2,7)=( 0.d0, 0.d0)
1582   char_mat(2,8)=( 0.d0, 0.d0)
1583   char_mat(2,9)=( 0.d0, 0.d0)
1584
1585   name_rap(3)="G_9  "
1586   char_mat(3,1)=( 2.d0, 0.d0)
1587   char_mat(3,2)=(-2.d0, 0.d0)
1588   char_mat(3,3)=( 0.d0, 0.d0)
1589   char_mat(3,4)=( 0.d0, 0.d0)
1590   char_mat(3,5)=(-2.d0, 0.d0)
1591   char_mat(3,6)=( 2.d0, 0.d0)
1592   char_mat(3,7)=( 0.d0, 0.d0)
1593   char_mat(3,8)=( 0.d0, 0.d0)
1594   char_mat(3,9)=( 0.d0, 0.d0)
1595
1596ELSEIF (code_group==16) THEN
1597!
1598! C_2h
1599!
1600   nclass_ref=8
1601   name_class(3)="C2  "
1602   name_class(4)="-C2 "
1603   name_class(5)="i   "
1604   name_class(6)="-i  "
1605   name_class(7)="s_h "
1606   name_class(8)="-s_h"
1607
1608   nrap_ref=4
1609
1610   name_rap(1)="G_3+"
1611   char_mat(1,3)=( 0.d0, 1.d0)
1612   char_mat(1,4)=( 0.d0,-1.d0)
1613   char_mat(1,6)=(-1.d0, 0.d0)
1614   char_mat(1,7)=( 0.d0, 1.d0)
1615   char_mat(1,8)=( 0.d0,-1.d0)
1616
1617   name_rap(2)="G_4+"
1618   char_mat(2,3)=( 0.d0,-1.d0)
1619   char_mat(2,4)=( 0.d0, 1.d0)
1620   char_mat(2,6)=(-1.d0, 0.d0)
1621   char_mat(2,7)=( 0.d0,-1.d0)
1622   char_mat(2,8)=( 0.d0, 1.d0)
1623
1624   name_rap(3)="G_3-"
1625   char_mat(3,3)=( 0.d0, 1.d0)
1626   char_mat(3,4)=( 0.d0,-1.d0)
1627   char_mat(3,5)=(-1.d0, 0.d0)
1628   char_mat(3,7)=( 0.d0,-1.d0)
1629   char_mat(3,8)=( 0.d0, 1.d0)
1630
1631   name_rap(4)="G_4- "
1632   char_mat(4,3)=( 0.d0,-1.d0)
1633   char_mat(4,4)=( 0.d0, 1.d0)
1634   char_mat(4,5)=(-1.d0, 0.d0)
1635   char_mat(4,7)=( 0.d0, 1.d0)
1636   char_mat(4,8)=( 0.d0,-1.d0)
1637
1638ELSEIF (code_group==17) THEN
1639!
1640!     Changed several signs. Now the character table is equal to the table
1641!     of C_6 and the signs match the table in Koster, Dimmock,
1642!     Wheeler, Statz, Properties of the 32 point groups.
1643!
1644! C_3h
1645!
1646   nclass_ref=12
1647   name_class(3)="C3  "
1648   name_class(4)="-C3 "
1649   name_class(5)="C3^2"
1650   name_class(6)="-C3^2"
1651   name_class(7)=" s_h "
1652   name_class(8)="-s_h "
1653   name_class(9)=" S3 "
1654   name_class(10)="-S3 "
1655   name_class(11)=" S3^5"
1656   name_class(12)="-S3^5"
1657
1658   nrap_ref=6
1659
1660   name_rap(1)="G_7  "
1661   char_mat(1,2)=(-1.d0, 0.d0)
1662   char_mat(1,3)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1663   char_mat(1,4)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1664!   char_mat(1,5)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1665!   char_mat(1,6)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1666   char_mat(1,5)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1667   char_mat(1,6)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1668   char_mat(1,7)=( 0.d0, 1.d0)
1669   char_mat(1,8)=( 0.d0,-1.d0)
1670!   char_mat(1,9)= CMPLX(-sqr3d2, 0.5d0,kind=DP)
1671!   char_mat(1,10)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1672   char_mat(1,9)= CMPLX( sqr3d2,-0.5d0,kind=DP)
1673   char_mat(1,10)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1674   char_mat(1,11)=CMPLX( sqr3d2, 0.5d0,kind=DP)
1675   char_mat(1,12)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
1676
1677   name_rap(2)="G_8  "
1678   char_mat(2,2)=(-1.d0, 0.d0)
1679   char_mat(2,3)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1680   char_mat(2,4)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1681!   char_mat(2,5)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1682!   char_mat(2,6)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1683   char_mat(2,5)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1684   char_mat(2,6)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1685   char_mat(2,7)=( 0.d0,-1.d0)
1686   char_mat(2,8)=( 0.d0, 1.d0)
1687!   char_mat(2,9)= CMPLX(-sqr3d2,-0.5d0,kind=DP)
1688!   char_mat(2,10)=CMPLX( sqr3d2, 0.5d0,kind=DP)
1689   char_mat(2,9)= CMPLX( sqr3d2, 0.5d0,kind=DP)
1690   char_mat(2,10)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
1691   char_mat(2,11)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1692   char_mat(2,12)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1693
1694   name_rap(3)="G_9  "
1695   char_mat(3,2)=(-1.d0, 0.d0)
1696   char_mat(3,3)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1697   char_mat(3,4)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1698!   char_mat(3,5)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1699!   char_mat(3,6)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1700   char_mat(3,5)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1701   char_mat(3,6)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1702   char_mat(3,7)=( 0.d0,-1.d0)
1703   char_mat(3,8)=( 0.d0, 1.d0)
1704!   char_mat(3,9)= CMPLX( sqr3d2,-0.5d0,kind=DP)
1705!   char_mat(3,10)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1706   char_mat(3,9)= CMPLX(-sqr3d2, 0.5d0,kind=DP)
1707   char_mat(3,10)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1708   char_mat(3,11)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
1709   char_mat(3,12)=CMPLX( sqr3d2, 0.5d0,kind=DP)
1710
1711   name_rap(4)="G_10  "
1712   char_mat(4,2)=(-1.d0, 0.d0)
1713   char_mat(4,3)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1714   char_mat(4,4)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1715!   char_mat(4,5)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1716!   char_mat(4,6)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1717   char_mat(4,5)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1718   char_mat(4,6)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1719   char_mat(4,7)=( 0.d0, 1.d0)
1720   char_mat(4,8)=( 0.d0,-1.d0)
1721!   char_mat(4,9)= CMPLX( sqr3d2, 0.5d0,kind=DP)
1722!   char_mat(4,10)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
1723   char_mat(4,9)= CMPLX(-sqr3d2,-0.5d0,kind=DP)
1724   char_mat(4,10)=CMPLX( sqr3d2, 0.5d0,kind=DP)
1725   char_mat(4,11)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1726   char_mat(4,12)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1727
1728   name_rap(5)="G_11 "
1729   char_mat(5,2)=(-1.d0, 0.d0)
1730   char_mat(5,3)=(-1.d0, 0.d0)
1731!   char_mat(5,5)=( 1.d0, 0.d0)
1732!   char_mat(5,6)=(-1.d0, 0.d0)
1733!   char_mat(5,7)=( 0.d0, 1.d0)
1734!   char_mat(5,8)=( 0.d0,-1.d0)
1735   char_mat(5,5)=(-1.d0, 0.d0)
1736   char_mat(5,6)=( 1.d0, 0.d0)
1737   char_mat(5,7)=( 0.d0,-1.d0)
1738   char_mat(5,8)=( 0.d0, 1.d0)
1739   char_mat(5,9)=( 0.d0,-1.d0)
1740   char_mat(5,10)=( 0.d0, 1.d0)
1741!   char_mat(5,11)=( 0.d0,-1.d0)
1742!   char_mat(5,12)=( 0.d0, 1.d0)
1743   char_mat(5,11)=( 0.d0, 1.d0)
1744   char_mat(5,12)=( 0.d0,-1.d0)
1745
1746   name_rap(6)="G_12 "
1747   char_mat(6,2)=(-1.d0, 0.d0)
1748   char_mat(6,3)=(-1.d0, 0.d0)
1749!   char_mat(6,5)=( 1.d0, 0.d0)
1750!   char_mat(6,6)=(-1.d0, 0.d0)
1751!   char_mat(6,7)=( 0.d0,-1.d0)
1752!   char_mat(6,8)=( 0.d0, 1.d0)
1753   char_mat(6,5)=(-1.d0, 0.d0)
1754   char_mat(6,6)=( 1.d0, 0.d0)
1755   char_mat(6,7)=( 0.d0, 1.d0)
1756   char_mat(6,8)=( 0.d0,-1.d0)
1757   char_mat(6,9)=( 0.d0, 1.d0)
1758   char_mat(6,10)=( 0.d0,-1.d0)
1759!   char_mat(6,11)=( 0.d0, 1.d0)
1760!   char_mat(6,12)=( 0.d0,-1.d0)
1761   char_mat(6,11)=( 0.d0,-1.d0)
1762   char_mat(6,12)=( 0.d0, 1.d0)
1763
1764ELSEIF (code_group==18) THEN
1765!
1766! C_4h   NB: The signs of the characters of the classes C4^3 -C4^2 S4 -S4
1767!            are changed with respect to Koster, Space groups and
1768!            their representations. They match the table in Koster,
1769!            Dimmock, Wheeler, Statz, Properties of the 32 point groups.
1770!
1771   nclass_ref=16
1772   name_class(3)="C4 "
1773   name_class(4)="-C4 "
1774   name_class(5)="C2  "
1775   name_class(6)="-C2 "
1776   name_class(7)="C4^3"
1777   name_class(8)="-C4^3"
1778   name_class(9)="i "
1779   name_class(10)="-i "
1780   name_class(11)="S4^3"
1781   name_class(12)="-S4^3"
1782   name_class(13)="s_h "
1783   name_class(14)="-s_h"
1784   name_class(15)="S4 "
1785   name_class(16)="-S4"
1786
1787   nrap_ref=8
1788
1789   name_rap(1)="G_5+ "
1790   char_mat(1,3)=CMPLX( dsq2, dsq2,kind=DP)
1791   char_mat(1,4)=CMPLX(-dsq2,-dsq2,kind=DP)
1792   char_mat(1,5)=( 0.d0, 1.d0)
1793   char_mat(1,6)=( 0.d0,-1.d0)
1794   char_mat(1,7)=CMPLX( dsq2,-dsq2,kind=DP)
1795   char_mat(1,8)=CMPLX(-dsq2, dsq2,kind=DP)
1796   char_mat(1,10)=(-1.d0, 0.d0)
1797   char_mat(1,11)=CMPLX( dsq2, dsq2,kind=DP)
1798   char_mat(1,12)=CMPLX(-dsq2,-dsq2,kind=DP)
1799   char_mat(1,13)=( 0.d0, 1.d0)
1800   char_mat(1,14)=( 0.d0,-1.d0)
1801   char_mat(1,15)=CMPLX( dsq2,-dsq2,kind=DP)
1802   char_mat(1,16)=CMPLX(-dsq2, dsq2,kind=DP)
1803
1804   name_rap(2)="G_6+ "
1805   char_mat(2,3)=CMPLX( dsq2,-dsq2,kind=DP)
1806   char_mat(2,4)=CMPLX(-dsq2, dsq2,kind=DP)
1807   char_mat(2,5)=( 0.d0,-1.d0)
1808   char_mat(2,6)=( 0.d0, 1.d0)
1809   char_mat(2,7)=CMPLX( dsq2, dsq2,kind=DP)
1810   char_mat(2,8)=CMPLX(-dsq2,-dsq2,kind=DP)
1811   char_mat(2,10)=(-1.d0, 0.d0)
1812   char_mat(2,11)=CMPLX( dsq2,-dsq2,kind=DP)
1813   char_mat(2,12)=CMPLX(-dsq2, dsq2,kind=DP)
1814   char_mat(2,13)=( 0.d0,-1.d0)
1815   char_mat(2,14)=( 0.d0, 1.d0)
1816   char_mat(2,15)=CMPLX( dsq2, dsq2,kind=DP)
1817   char_mat(2,16)=CMPLX(-dsq2,-dsq2,kind=DP)
1818
1819   name_rap(3)="G_7+ "
1820   char_mat(3,3)=CMPLX(-dsq2,-dsq2,kind=DP)
1821   char_mat(3,4)=CMPLX( dsq2, dsq2,kind=DP)
1822   char_mat(3,5)=( 0.d0, 1.d0)
1823   char_mat(3,6)=( 0.d0,-1.d0)
1824   char_mat(3,7)=CMPLX(-dsq2, dsq2,kind=DP)
1825   char_mat(3,8)=CMPLX( dsq2,-dsq2,kind=DP)
1826   char_mat(3,10)=(-1.d0, 0.d0)
1827   char_mat(3,11)=CMPLX(-dsq2,-dsq2,kind=DP)
1828   char_mat(3,12)=CMPLX( dsq2, dsq2,kind=DP)
1829   char_mat(3,13)=( 0.d0, 1.d0)
1830   char_mat(3,14)=( 0.d0,-1.d0)
1831   char_mat(3,15)=CMPLX(-dsq2, dsq2,kind=DP)
1832   char_mat(3,16)=CMPLX( dsq2,-dsq2,kind=DP)
1833
1834   name_rap(4)="G_8+ "
1835   char_mat(4,3)=CMPLX(-dsq2, dsq2,kind=DP)
1836   char_mat(4,4)=CMPLX( dsq2,-dsq2,kind=DP)
1837   char_mat(4,5)=( 0.d0,-1.d0)
1838   char_mat(4,6)=( 0.d0, 1.d0)
1839   char_mat(4,7)=CMPLX(-dsq2,-dsq2,kind=DP)
1840   char_mat(4,8)=CMPLX( dsq2, dsq2,kind=DP)
1841   char_mat(4,10)=(-1.d0, 0.d0)
1842   char_mat(4,11)=CMPLX(-dsq2, dsq2,kind=DP)
1843   char_mat(4,12)=CMPLX( dsq2,-dsq2,kind=DP)
1844   char_mat(4,13)=( 0.d0,-1.d0)
1845   char_mat(4,14)=( 0.d0, 1.d0)
1846   char_mat(4,15)=CMPLX(-dsq2,-dsq2,kind=DP)
1847   char_mat(4,16)=CMPLX( dsq2, dsq2,kind=DP)
1848
1849   name_rap(5)="G_5- "
1850   char_mat(5,3)=CMPLX( dsq2, dsq2,kind=DP)
1851   char_mat(5,4)=CMPLX(-dsq2,-dsq2,kind=DP)
1852   char_mat(5,5)=( 0.d0, 1.d0)
1853   char_mat(5,6)=( 0.d0,-1.d0)
1854   char_mat(5,7)=CMPLX( dsq2,-dsq2,kind=DP)
1855   char_mat(5,8)=CMPLX(-dsq2, dsq2,kind=DP)
1856   char_mat(5,9)=(-1.d0, 0.d0)
1857   char_mat(5,11)=CMPLX(-dsq2,-dsq2,kind=DP)
1858   char_mat(5,12)=CMPLX( dsq2, dsq2,kind=DP)
1859   char_mat(5,13)=( 0.d0,-1.d0)
1860   char_mat(5,14)=( 0.d0, 1.d0)
1861   char_mat(5,15)=CMPLX(-dsq2, dsq2,kind=DP)
1862   char_mat(5,16)=CMPLX( dsq2,-dsq2,kind=DP)
1863
1864   name_rap(6)="G_6- "
1865   char_mat(6,3)=CMPLX( dsq2,-dsq2,kind=DP)
1866   char_mat(6,4)=CMPLX(-dsq2, dsq2,kind=DP)
1867   char_mat(6,5)=( 0.d0,-1.d0)
1868   char_mat(6,6)=( 0.d0, 1.d0)
1869   char_mat(6,7)=CMPLX( dsq2, dsq2,kind=DP)
1870   char_mat(6,8)=CMPLX(-dsq2,-dsq2,kind=DP)
1871   char_mat(6,9)=(-1.d0, 0.d0)
1872   char_mat(6,11)=CMPLX(-dsq2, dsq2,kind=DP)
1873   char_mat(6,12)=CMPLX( dsq2,-dsq2,kind=DP)
1874   char_mat(6,13)=( 0.d0, 1.d0)
1875   char_mat(6,14)=( 0.d0,-1.d0)
1876   char_mat(6,15)=CMPLX(-dsq2,-dsq2,kind=DP)
1877   char_mat(6,16)=CMPLX( dsq2, dsq2,kind=DP)
1878
1879   name_rap(7)="G_7- "
1880   char_mat(7,3)=CMPLX(-dsq2,-dsq2,kind=DP)
1881   char_mat(7,4)=CMPLX( dsq2, dsq2,kind=DP)
1882   char_mat(7,5)=( 0.d0, 1.d0)
1883   char_mat(7,6)=( 0.d0,-1.d0)
1884   char_mat(7,7)=CMPLX(-dsq2, dsq2,kind=DP)
1885   char_mat(7,8)=CMPLX( dsq2,-dsq2,kind=DP)
1886   char_mat(7,9)=(-1.d0, 0.d0)
1887   char_mat(7,11)=CMPLX( dsq2, dsq2,kind=DP)
1888   char_mat(7,12)=CMPLX(-dsq2,-dsq2,kind=DP)
1889   char_mat(7,13)=( 0.d0,-1.d0)
1890   char_mat(7,14)=( 0.d0, 1.d0)
1891   char_mat(7,15)=CMPLX( dsq2,-dsq2,kind=DP)
1892   char_mat(7,16)=CMPLX(-dsq2, dsq2,kind=DP)
1893
1894   name_rap(8)="G_8- "
1895   char_mat(8,3)=CMPLX(-dsq2, dsq2,kind=DP)
1896   char_mat(8,4)=CMPLX( dsq2,-dsq2,kind=DP)
1897   char_mat(8,5)=( 0.d0,-1.d0)
1898   char_mat(8,6)=( 0.d0, 1.d0)
1899   char_mat(8,7)=CMPLX(-dsq2,-dsq2,kind=DP)
1900   char_mat(8,8)=CMPLX( dsq2, dsq2,kind=DP)
1901   char_mat(8,9)=(-1.d0, 0.d0)
1902   char_mat(8,11)=CMPLX( dsq2,-dsq2,kind=DP)
1903   char_mat(8,12)=CMPLX(-dsq2, dsq2,kind=DP)
1904   char_mat(8,13)=( 0.d0, 1.d0)
1905   char_mat(8,14)=( 0.d0,-1.d0)
1906   char_mat(8,15)=CMPLX( dsq2, dsq2,kind=DP)
1907   char_mat(8,16)=CMPLX(-dsq2,-dsq2,kind=DP)
1908
1909ELSEIF (code_group==19) THEN
1910!
1911! C_6h  NB: The signs of the characters of C3^2 -C3^2 C6^5 -C6^5 S6 -S6 S3 -S3
1912!           are changed with respect to Koster, Space groups and their
1913!           representation. They match the table in Koster, Dimmock,
1914!           Wheeler, Statz, Properties of the 32 point groups.
1915!           They are also consistent with the fact that C_6h=C_6xC_i
1916!
1917   nclass_ref=24
1918   name_class(3)="C6  "
1919   name_class(4)="-C6 "
1920   name_class(5)="C3  "
1921   name_class(6)="-C3 "
1922   name_class(7)="C2  "
1923   name_class(8)="-C2 "
1924   name_class(9)="C3^2"
1925   name_class(10)="-C3^2"
1926   name_class(11)="C6^5"
1927   name_class(12)="-C6^5"
1928   name_class(13)="i  "
1929   name_class(14)="-i "
1930   name_class(15)="S3^5"
1931   name_class(16)="-S3^5"
1932   name_class(17)="S6^5"
1933   name_class(18)="-S6^5"
1934   name_class(19)="s_h "
1935   name_class(20)="-s_h"
1936   name_class(21)="S6 "
1937   name_class(22)="-S6 "
1938   name_class(23)="S3 "
1939   name_class(24)="-S3 "
1940
1941   nrap_ref=12
1942
1943   name_rap(1)="G_7+ "
1944   char_mat(1,3)=CMPLX( sqr3d2, 0.5d0,kind=DP)
1945   char_mat(1,4)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
1946   char_mat(1,5)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1947   char_mat(1,6)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1948   char_mat(1,7)=( 0.d0, 1.0d0)
1949   char_mat(1,8)=( 0.d0,-1.0d0)
1950!   char_mat(1,9)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1951!   char_mat(1,10)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1952!   char_mat(1,11)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1953!   char_mat(1,12)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1954   char_mat(1,9)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1955   char_mat(1,10)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1956   char_mat(1,11)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1957   char_mat(1,12)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1958   char_mat(1,14)=(-1.0d0, 0.d0)
1959   char_mat(1,15)=CMPLX( sqr3d2, 0.5d0,kind=DP)
1960   char_mat(1,16)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
1961   char_mat(1,17)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1962   char_mat(1,18)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1963   char_mat(1,19)=( 0.d0, 1.0d0)
1964   char_mat(1,20)=( 0.d0,-1.0d0)
1965!   char_mat(1,21)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1966!   char_mat(1,22)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1967!   char_mat(1,23)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1968!   char_mat(1,24)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1969   char_mat(1,21)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1970   char_mat(1,22)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1971   char_mat(1,23)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1972   char_mat(1,24)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1973
1974   name_rap(2)="G_8+ "
1975   char_mat(2,3)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1976   char_mat(2,4)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1977   char_mat(2,5)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1978   char_mat(2,6)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1979   char_mat(2,7)=( 0.d0,-1.0d0)
1980   char_mat(2,8)=( 0.d0, 1.0d0)
1981!   char_mat(2,9)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1982!   char_mat(2,10)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1983!   char_mat(2,11)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
1984!   char_mat(2,12)=CMPLX( sqr3d2, 0.5d0,kind=DP)
1985   char_mat(2,9)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1986   char_mat(2,10)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1987   char_mat(2,11)=CMPLX( sqr3d2, 0.5d0,kind=DP)
1988   char_mat(2,12)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
1989   char_mat(2,14)=(-1.0d0, 0.d0)
1990   char_mat(2,15)=CMPLX( sqr3d2,-0.5d0,kind=DP)
1991   char_mat(2,16)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
1992   char_mat(2,17)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
1993   char_mat(2,18)=CMPLX(-0.5d0, sqr3d2,kind=DP)
1994   char_mat(2,19)=( 0.d0,-1.0d0)
1995   char_mat(2,20)=( 0.d0, 1.0d0)
1996!   char_mat(2,21)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
1997!   char_mat(2,22)=CMPLX( 0.5d0, sqr3d2,kind=DP)
1998!   char_mat(2,23)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
1999!   char_mat(2,24)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2000   char_mat(2,21)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2001   char_mat(2,22)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2002   char_mat(2,23)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2003   char_mat(2,24)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2004
2005
2006   name_rap(3)="G_9+ "
2007   char_mat(3,3)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2008   char_mat(3,4)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2009   char_mat(3,5)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2010   char_mat(3,6)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2011   char_mat(3,7)=( 0.d0,-1.0d0)
2012   char_mat(3,8)=( 0.d0, 1.0d0)
2013!   char_mat(3,9)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2014!   char_mat(3,10)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2015!   char_mat(3,11)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2016!   char_mat(3,12)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2017   char_mat(3,9)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2018   char_mat(3,10)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2019   char_mat(3,11)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2020   char_mat(3,12)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2021   char_mat(3,14)=(-1.0d0, 0.d0)
2022   char_mat(3,15)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2023   char_mat(3,16)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2024   char_mat(3,17)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2025   char_mat(3,18)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2026   char_mat(3,19)=( 0.d0,-1.0d0)
2027   char_mat(3,20)=( 0.d0, 1.0d0)
2028!   char_mat(3,21)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2029!   char_mat(3,22)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2030!   char_mat(3,23)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2031!   char_mat(3,24)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2032   char_mat(3,21)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2033   char_mat(3,22)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2034   char_mat(3,23)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2035   char_mat(3,24)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2036
2037   name_rap(4)="G_10+ "
2038   char_mat(4,3)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2039   char_mat(4,4)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2040   char_mat(4,5)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2041   char_mat(4,6)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2042   char_mat(4,7)=( 0.d0, 1.0d0)
2043   char_mat(4,8)=( 0.d0,-1.0d0)
2044!   char_mat(4,9)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2045!   char_mat(4,10)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2046!   char_mat(4,11)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2047!   char_mat(4,12)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2048   char_mat(4,9)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2049   char_mat(4,10)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2050   char_mat(4,11)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2051   char_mat(4,12)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2052   char_mat(4,14)=(-1.0d0, 0.d0)
2053   char_mat(4,15)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2054   char_mat(4,16)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2055   char_mat(4,17)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2056   char_mat(4,18)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2057   char_mat(4,19)=( 0.d0, 1.0d0)
2058   char_mat(4,20)=( 0.d0,-1.0d0)
2059!   char_mat(4,21)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2060!   char_mat(4,22)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2061!   char_mat(4,23)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2062!   char_mat(4,24)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2063   char_mat(4,21)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2064   char_mat(4,22)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2065   char_mat(4,23)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2066   char_mat(4,24)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2067
2068   name_rap(5)="G_11+ "
2069   char_mat(5,3)=( 0.d0, 1.0d0)
2070   char_mat(5,4)=( 0.d0,-1.0d0)
2071   char_mat(5,5)=(-1.d0, 0.d0 )
2072   char_mat(5,7)=( 0.d0,-1.0d0)
2073   char_mat(5,8)=( 0.d0, 1.0d0)
2074!   char_mat(5,10)=(-1.d0, 0.d0 )
2075!   char_mat(5,11)=( 0.d0, 1.0d0)
2076!   char_mat(5,12)=( 0.d0,-1.0d0)
2077   char_mat(5,9)=(-1.d0, 0.d0 )
2078   char_mat(5,11)=( 0.d0,-1.0d0)
2079   char_mat(5,12)=( 0.d0, 1.0d0)
2080   char_mat(5,14)=(-1.0d0, 0.d0)
2081   char_mat(5,15)=( 0.d0, 1.0d0)
2082   char_mat(5,16)=( 0.d0,-1.0d0)
2083   char_mat(5,17)=(-1.0d0, 0.d0)
2084   char_mat(5,19)=( 0.d0,-1.0d0)
2085   char_mat(5,20)=( 0.d0, 1.0d0)
2086!   char_mat(5,22)=(-1.0d0, 0.d0)
2087!   char_mat(5,23)=( 0.d0, 1.0d0)
2088!   char_mat(5,24)=( 0.d0,-1.0d0)
2089   char_mat(5,21)=(-1.0d0, 0.d0)
2090   char_mat(5,23)=( 0.d0,-1.0d0)
2091   char_mat(5,24)=( 0.d0, 1.0d0)
2092
2093   name_rap(6)="G_12+ "
2094   char_mat(6,3)=( 0.d0,-1.0d0)
2095   char_mat(6,4)=( 0.d0, 1.0d0)
2096   char_mat(6,5)=(-1.d0, 0.d0 )
2097   char_mat(6,7)=( 0.d0, 1.0d0)
2098   char_mat(6,8)=( 0.d0,-1.0d0)
2099!   char_mat(6,10)=(-1.d0, 0.d0 )
2100!   char_mat(6,11)=( 0.d0,-1.0d0)
2101!   char_mat(6,12)=( 0.d0, 1.0d0)
2102   char_mat(6,9)=(-1.d0, 0.d0 )
2103   char_mat(6,11)=( 0.d0, 1.0d0)
2104   char_mat(6,12)=( 0.d0,-1.0d0)
2105   char_mat(6,14)=(-1.0d0, 0.d0)
2106   char_mat(6,15)=( 0.d0,-1.0d0)
2107   char_mat(6,16)=( 0.d0, 1.0d0)
2108   char_mat(6,17)=(-1.0d0, 0.d0)
2109   char_mat(6,19)=( 0.d0, 1.0d0)
2110   char_mat(6,20)=( 0.d0,-1.0d0)
2111!   char_mat(6,22)=(-1.0d0, 0.d0)
2112!   char_mat(6,23)=( 0.d0,-1.0d0)
2113!   char_mat(6,24)=( 0.d0, 1.0d0)
2114   char_mat(6,21)=(-1.0d0, 0.d0)
2115   char_mat(6,23)=( 0.d0, 1.0d0)
2116   char_mat(6,24)=( 0.d0,-1.0d0)
2117
2118
2119   name_rap(7)="G_7- "
2120   char_mat(7,3)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2121   char_mat(7,4)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2122   char_mat(7,5)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2123   char_mat(7,6)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2124   char_mat(7,7)=( 0.d0, 1.0d0)
2125   char_mat(7,8)=( 0.d0,-1.0d0)
2126!   char_mat(7,9)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2127!   char_mat(7,10)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2128!   char_mat(7,11)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2129!   char_mat(7,12)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2130   char_mat(7,9)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2131   char_mat(7,10)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2132   char_mat(7,11)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2133   char_mat(7,12)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2134   char_mat(7,13)=(-1.0d0, 0.d0)
2135   char_mat(7,15)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2136   char_mat(7,16)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2137   char_mat(7,17)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2138   char_mat(7,18)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2139   char_mat(7,19)=( 0.d0,-1.0d0)
2140   char_mat(7,20)=( 0.d0, 1.0d0)
2141!   char_mat(7,21)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2142!   char_mat(7,22)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2143!   char_mat(7,23)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2144!   char_mat(7,24)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2145   char_mat(7,21)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2146   char_mat(7,22)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2147   char_mat(7,23)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2148   char_mat(7,24)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2149
2150   name_rap(8)="G_8- "
2151   char_mat(8,3)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2152   char_mat(8,4)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2153   char_mat(8,5)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2154   char_mat(8,6)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2155   char_mat(8,7)=( 0.d0,-1.0d0)
2156   char_mat(8,8)=( 0.d0, 1.0d0)
2157!   char_mat(8,9)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2158!   char_mat(8,10)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2159!   char_mat(8,11)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2160!   char_mat(8,12)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2161   char_mat(8,9)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2162   char_mat(8,10)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2163   char_mat(8,11)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2164   char_mat(8,12)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2165   char_mat(8,13)=(-1.0d0, 0.d0)
2166   char_mat(8,15)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2167   char_mat(8,16)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2168   char_mat(8,17)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2169   char_mat(8,18)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2170   char_mat(8,19)=( 0.d0, 1.0d0)
2171   char_mat(8,20)=( 0.d0,-1.0d0)
2172!   char_mat(8,21)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2173!   char_mat(8,22)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2174!   char_mat(8,23)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2175!   char_mat(8,24)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2176   char_mat(8,21)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2177   char_mat(8,22)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2178   char_mat(8,23)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2179   char_mat(8,24)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2180
2181
2182   name_rap(9)="G_9- "
2183   char_mat(9,3)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2184   char_mat(9,4)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2185   char_mat(9,5)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2186   char_mat(9,6)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2187   char_mat(9,7)=( 0.d0,-1.0d0)
2188   char_mat(9,8)=( 0.d0, 1.0d0)
2189!   char_mat(9,9)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2190!   char_mat(9,10)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2191!   char_mat(9,11)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2192!   char_mat(9,12)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2193   char_mat(9,9)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2194   char_mat(9,10)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2195   char_mat(9,11)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2196   char_mat(9,12)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2197   char_mat(9,13)=(-1.0d0, 0.d0)
2198   char_mat(9,15)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2199   char_mat(9,16)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2200   char_mat(9,17)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2201   char_mat(9,18)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2202   char_mat(9,19)=( 0.d0, 1.0d0)
2203   char_mat(9,20)=( 0.d0,-1.0d0)
2204!   char_mat(9,21)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2205!   char_mat(9,22)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2206!   char_mat(9,23)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2207!   char_mat(9,24)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2208   char_mat(9,21)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2209   char_mat(9,22)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2210   char_mat(9,23)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2211   char_mat(9,24)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2212
2213   name_rap(10)="G_10- "
2214   char_mat(10,3)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2215   char_mat(10,4)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2216   char_mat(10,5)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2217   char_mat(10,6)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2218   char_mat(10,7)=( 0.d0, 1.0d0)
2219   char_mat(10,8)=( 0.d0,-1.0d0)
2220!   char_mat(10,9)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2221!   char_mat(10,10)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2222!   char_mat(10,11)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2223!   char_mat(10,12)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2224   char_mat(10,9)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2225   char_mat(10,10)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2226   char_mat(10,11)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2227   char_mat(10,12)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2228   char_mat(10,13)=(-1.0d0, 0.d0)
2229   char_mat(10,15)=CMPLX( sqr3d2,-0.5d0,kind=DP)
2230   char_mat(10,16)=CMPLX(-sqr3d2, 0.5d0,kind=DP)
2231   char_mat(10,17)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2232   char_mat(10,18)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2233   char_mat(10,19)=( 0.d0,-1.0d0)
2234   char_mat(10,20)=( 0.d0, 1.0d0)
2235!   char_mat(10,21)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2236!   char_mat(10,22)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2237!   char_mat(10,23)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2238!   char_mat(10,24)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2239   char_mat(10,21)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2240   char_mat(10,22)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2241   char_mat(10,23)=CMPLX( sqr3d2, 0.5d0,kind=DP)
2242   char_mat(10,24)=CMPLX(-sqr3d2,-0.5d0,kind=DP)
2243
2244   name_rap(11)="G_11-"
2245   char_mat(11,3)=( 0.d0, 1.0d0)
2246   char_mat(11,4)=( 0.d0,-1.0d0)
2247   char_mat(11,5)=(-1.d0, 0.d0 )
2248   char_mat(11,7)=( 0.d0,-1.0d0)
2249   char_mat(11,8)=( 0.d0, 1.0d0)
2250!   char_mat(11,10)=(-1.d0, 0.d0 )
2251!   char_mat(11,11)=( 0.d0, 1.0d0)
2252!   char_mat(11,12)=( 0.d0,-1.0d0)
2253   char_mat(11,9)=(-1.d0, 0.d0 )
2254   char_mat(11,11)=( 0.d0,-1.0d0)
2255   char_mat(11,12)=( 0.d0, 1.0d0)
2256   char_mat(11,13)=(-1.0d0, 0.d0)
2257   char_mat(11,15)=( 0.d0,-1.0d0)
2258   char_mat(11,16)=( 0.d0, 1.0d0)
2259   char_mat(11,18)=(-1.0d0, 0.d0)
2260   char_mat(11,19)=( 0.d0, 1.0d0)
2261   char_mat(11,20)=( 0.d0,-1.0d0)
2262!   char_mat(11,21)=(-1.0d0, 0.d0)
2263!   char_mat(11,23)=( 0.d0,-1.0d0)
2264!   char_mat(11,24)=( 0.d0, 1.0d0)
2265   char_mat(11,22)=(-1.0d0, 0.d0)
2266   char_mat(11,23)=( 0.d0, 1.0d0)
2267   char_mat(11,24)=( 0.d0,-1.0d0)
2268
2269   name_rap(12)="G_12-"
2270   char_mat(12,3) =( 0.d0,-1.0d0)
2271   char_mat(12,4) =( 0.d0, 1.0d0)
2272   char_mat(12,5) =(-1.d0, 0.d0 )
2273   char_mat(12,7) =( 0.d0, 1.0d0)
2274   char_mat(12,8) =( 0.d0,-1.0d0)
2275!   char_mat(12,10)=(-1.d0, 0.d0 )
2276!   char_mat(12,11)=( 0.d0,-1.0d0)
2277!   char_mat(12,12)=( 0.d0, 1.0d0)
2278   char_mat(12,9)=(-1.d0, 0.d0 )
2279   char_mat(12,11)=( 0.d0, 1.0d0)
2280   char_mat(12,12)=( 0.d0,-1.0d0)
2281   char_mat(12,13)=(-1.0d0, 0.d0)
2282   char_mat(12,15)=( 0.d0, 1.0d0)
2283   char_mat(12,16)=( 0.d0,-1.0d0)
2284   char_mat(12,18)=(-1.0d0, 0.d0)
2285   char_mat(12,19)=( 0.d0,-1.0d0)
2286   char_mat(12,20)=( 0.d0, 1.0d0)
2287!   char_mat(12,21)=(-1.0d0, 0.d0)
2288!   char_mat(12,23)=( 0.d0, 1.0d0)
2289!   char_mat(12,24)=( 0.d0,-1.0d0)
2290   char_mat(12,22)=(-1.0d0, 0.d0)
2291   char_mat(12,23)=( 0.d0,-1.0d0)
2292   char_mat(12,24)=( 0.d0, 1.0d0)
2293
2294
2295
2296ELSEIF (code_group==20) THEN
2297!
2298! D_2h
2299!
2300   nclass_ref=10
2301   name_class(3)=" C2 "
2302   name_class1(3)="-C2 "
2303   name_class(4)=" C2'  "
2304   name_class1(4)="-C2'"
2305   name_class(5)="C2''  "
2306   name_class1(5)="-C2''"
2307   name_class(6)="i   "
2308   name_class(7)="-i  "
2309   name_class(8)=" s_v"
2310   name_class1(8)="-s_v"
2311   name_class(9)=" s_v'"
2312   name_class1(9)="-s_v'"
2313   name_class(10)="s_v''"
2314   name_class1(10)="-s_v''"
2315
2316   nrap_ref=2
2317
2318   name_rap(1)="G_5+ "
2319   char_mat(1,1)=( 2.d0, 0.d0)
2320   char_mat(1,2)=(-2.d0, 0.d0)
2321   char_mat(1,3)=( 0.d0, 0.d0)
2322   char_mat(1,4)=( 0.d0, 0.d0)
2323   char_mat(1,5)=( 0.d0, 0.d0)
2324   char_mat(1,6)=( 2.d0, 0.d0)
2325   char_mat(1,7)=(-2.d0, 0.d0)
2326   char_mat(1,8)=( 0.d0, 0.d0)
2327   char_mat(1,9)=( 0.d0, 0.d0)
2328   char_mat(1,10)=( 0.d0, 0.d0)
2329
2330   name_rap(2)="G_5- "
2331   char_mat(2,1)=( 2.d0, 0.d0)
2332   char_mat(2,2)=(-2.d0, 0.d0)
2333   char_mat(2,3)=( 0.d0, 0.d0)
2334   char_mat(2,4)=( 0.d0, 0.d0)
2335   char_mat(2,5)=( 0.d0, 0.d0)
2336   char_mat(2,6)=(-2.d0, 0.d0)
2337   char_mat(2,7)=( 2.d0, 0.d0)
2338   char_mat(2,8)=( 0.d0, 0.d0)
2339   char_mat(2,9)=( 0.d0, 0.d0)
2340   char_mat(2,10)=( 0.d0, 0.d0)
2341
2342ELSEIF (code_group==21) THEN
2343!
2344! D_3h
2345!
2346   nclass_ref=9
2347   name_class(3)="2C3 "
2348   name_class(4)="-2C3"
2349   name_class(5)=" 3C2'"
2350   name_class1(5)="-3C2'"
2351   name_class(6)=" s_h"
2352   name_class1(6)="-s_h'"
2353   name_class(7)="2S3 "
2354   name_class(8)="-2S3"
2355   name_class(9)=" 3s_v"
2356   name_class1(9)="-3s_v"
2357
2358   nrap_ref=3
2359
2360   name_rap(1)="G_7  "
2361   char_mat(1,1)=( 2.d0, 0.d0)
2362   char_mat(1,2)=(-2.d0, 0.d0)
2363   char_mat(1,4)=(-1.d0, 0.d0)
2364   char_mat(1,5)=( 0.d0, 0.d0)
2365   char_mat(1,6)=( 0.d0, 0.d0)
2366   char_mat(1,7)=CMPLX( sqrt3, 0.d0,kind=DP)
2367   char_mat(1,8)=CMPLX(-sqrt3, 0.d0,kind=DP)
2368   char_mat(1,9)=( 0.d0, 0.d0)
2369
2370   name_rap(2)="G_8  "
2371   char_mat(2,1)=( 2.d0, 0.d0)
2372   char_mat(2,2)=(-2.d0, 0.d0)
2373   char_mat(2,4)=(-1.d0, 0.d0)
2374   char_mat(2,5)=( 0.d0, 0.d0)
2375   char_mat(2,6)=( 0.d0, 0.d0)
2376   char_mat(2,7)=CMPLX(-sqrt3, 0.d0,kind=DP)
2377   char_mat(2,8)=CMPLX( sqrt3, 0.d0,kind=DP)
2378   char_mat(2,9)=( 0.d0, 0.d0)
2379
2380   name_rap(3)="G_9  "
2381   char_mat(3,1)=( 2.d0, 0.d0)
2382   char_mat(3,2)=(-2.d0, 0.d0)
2383   char_mat(3,3)=(-2.d0, 0.d0)
2384   char_mat(3,4)=( 2.d0, 0.d0)
2385   char_mat(3,5)=( 0.d0, 0.d0)
2386   char_mat(3,6)=( 0.d0, 0.d0)
2387   char_mat(3,7)=( 0.d0, 0.d0)
2388   char_mat(3,8)=( 0.d0, 0.d0)
2389   char_mat(3,9)=( 0.d0, 0.d0)
2390
2391ELSEIF (code_group==22) THEN
2392!
2393! D_4h
2394!
2395   nclass_ref=14
2396   name_class(3)="2C4 "
2397   name_class(4)="-2C4"
2398   name_class(5)=" C2 "
2399   name_class1(5)="-C2 "
2400   name_class(6)=" 2C2'"
2401   name_class1(6)="-2C2' "
2402   name_class(7)="2C2''"
2403   name_class1(7)="-2C2''"
2404   name_class(8)="i   "
2405   name_class(9)="-i  "
2406   name_class(10)="2S4 "
2407   name_class(11)="-2S4"
2408   name_class(12)=" s_h"
2409   name_class1(12)="-s_h"
2410   name_class(13)=" 2s_v"
2411   name_class1(13)="-2s_v"
2412   name_class(14)=" 2s_d"
2413   name_class1(14)="-2s_d"
2414
2415   nrap_ref=4
2416
2417   name_rap(1)="G_6+  M_6+"
2418   char_mat(1,1)=( 2.d0, 0.d0)
2419   char_mat(1,2)=(-2.d0, 0.d0)
2420   char_mat(1,3)=CMPLX( sqrt2, 0.d0,kind=DP)
2421   char_mat(1,4)=CMPLX(-sqrt2, 0.d0,kind=DP)
2422   char_mat(1,5)=( 0.d0, 0.d0)
2423   char_mat(1,6)=( 0.d0, 0.d0)
2424   char_mat(1,7)=( 0.d0, 0.d0)
2425   char_mat(1,8)=( 2.d0, 0.d0)
2426   char_mat(1,9)=(-2.d0, 0.d0)
2427   char_mat(1,10)=CMPLX( sqrt2, 0.d0,kind=DP)
2428   char_mat(1,11)=CMPLX(-sqrt2, 0.d0,kind=DP)
2429   char_mat(1,12)=( 0.d0, 0.d0)
2430   char_mat(1,13)=( 0.d0, 0.d0)
2431   char_mat(1,14)=( 0.d0, 0.d0)
2432
2433   name_rap(2)="G_7+  M_7+"
2434   char_mat(2,1)=( 2.d0, 0.d0)
2435   char_mat(2,2)=(-2.d0, 0.d0)
2436   char_mat(2,3)=CMPLX(-sqrt2, 0.d0,kind=DP)
2437   char_mat(2,4)=CMPLX( sqrt2, 0.d0,kind=DP)
2438   char_mat(2,5)=( 0.d0, 0.d0)
2439   char_mat(2,6)=( 0.d0, 0.d0)
2440   char_mat(2,7)=( 0.d0, 0.d0)
2441   char_mat(2,8)=( 2.d0, 0.d0)
2442   char_mat(2,9)=(-2.d0, 0.d0)
2443   char_mat(2,10)=CMPLX(-sqrt2, 0.d0,kind=DP)
2444   char_mat(2,11)=CMPLX( sqrt2, 0.d0,kind=DP)
2445   char_mat(2,12)=( 0.d0, 0.d0)
2446   char_mat(2,13)=( 0.d0, 0.d0)
2447   char_mat(2,14)=( 0.d0, 0.d0)
2448
2449   name_rap(3)="G_6-  M_6- "
2450   char_mat(3,1)=( 2.d0, 0.d0)
2451   char_mat(3,2)=(-2.d0, 0.d0)
2452   char_mat(3,3)=CMPLX( sqrt2, 0.d0,kind=DP)
2453   char_mat(3,4)=CMPLX(-sqrt2, 0.d0,kind=DP)
2454   char_mat(3,5)=( 0.d0, 0.d0)
2455   char_mat(3,6)=( 0.d0, 0.d0)
2456   char_mat(3,7)=( 0.d0, 0.d0)
2457   char_mat(3,8)=(-2.d0, 0.d0)
2458   char_mat(3,9)=( 2.d0, 0.d0)
2459   char_mat(3,10)=CMPLX(-sqrt2, 0.d0,kind=DP)
2460   char_mat(3,11)=CMPLX( sqrt2, 0.d0,kind=DP)
2461   char_mat(3,12)=( 0.d0, 0.d0)
2462   char_mat(3,13)=( 0.d0, 0.d0)
2463   char_mat(3,14)=( 0.d0, 0.d0)
2464
2465   name_rap(4)="G_7-  M_7- "
2466   char_mat(4,1)=( 2.d0, 0.d0)
2467   char_mat(4,2)=(-2.d0, 0.d0)
2468   char_mat(4,3)=CMPLX(-sqrt2, 0.d0,kind=DP)
2469   char_mat(4,4)=CMPLX( sqrt2, 0.d0,kind=DP)
2470   char_mat(4,5)=( 0.d0, 0.d0)
2471   char_mat(4,6)=( 0.d0, 0.d0)
2472   char_mat(4,7)=( 0.d0, 0.d0)
2473   char_mat(4,8)=(-2.d0, 0.d0)
2474   char_mat(4,9)=( 2.d0, 0.d0)
2475   char_mat(4,10)=CMPLX( sqrt2, 0.d0,kind=DP)
2476   char_mat(4,11)=CMPLX(-sqrt2, 0.d0,kind=DP)
2477   char_mat(4,12)=( 0.d0, 0.d0)
2478   char_mat(4,13)=( 0.d0, 0.d0)
2479   char_mat(4,14)=( 0.d0, 0.d0)
2480
2481ELSEIF (code_group==23) THEN
2482!
2483! D_6h
2484!
2485   nclass_ref=18
2486   name_class(3)="2C6 "
2487   name_class(4)="-2C6 "
2488   name_class(5)="2C3  "
2489   name_class(6)="-2C3 "
2490   name_class(7)="C2  "
2491   name_class1(7)="-C2 "
2492   name_class(8)=" 3C2'"
2493   name_class1(8)="-3C2'"
2494   name_class(9)="3C2''"
2495   name_class1(9)="-3C2''"
2496   name_class(10)=" i "
2497   name_class(11)="-i "
2498   name_class(12)="2S3"
2499   name_class(13)="-2S3"
2500   name_class(14)="2S6"
2501   name_class(15)="-2S6"
2502   name_class(16)=" s_h"
2503   name_class1(16)="-s_h"
2504   name_class(17)=" 3s_v"
2505   name_class1(17)="-3s_v"
2506   name_class(18)=" 3s_d"
2507   name_class1(18)="-3s_d"
2508
2509   nrap_ref=6
2510
2511   name_rap(1)="G_7+ "
2512   char_mat(1,1)=( 2.d0, 0.d0)
2513   char_mat(1,2)=(-2.d0, 0.d0)
2514   char_mat(1,3)=CMPLX( sqrt3, 0.d0,kind=DP)
2515   char_mat(1,4)=CMPLX(-sqrt3, 0.d0,kind=DP)
2516   char_mat(1,5)=( 1.d0, 0.d0)
2517   char_mat(1,6)=(-1.d0, 0.d0)
2518   char_mat(1,7)=( 0.d0, 0.d0)
2519   char_mat(1,8)=( 0.d0, 0.d0)
2520   char_mat(1,9)=( 0.d0, 0.d0)
2521   char_mat(1,10)=( 2.d0, 0.d0)
2522   char_mat(1,11)=(-2.d0, 0.d0)
2523   char_mat(1,12)=CMPLX( sqrt3, 0.d0,kind=DP)
2524   char_mat(1,13)=CMPLX(-sqrt3, 0.d0,kind=DP)
2525   char_mat(1,14)=( 1.d0, 0.d0)
2526   char_mat(1,15)=(-1.d0, 0.d0)
2527   char_mat(1,16)=( 0.d0, 0.d0)
2528   char_mat(1,17)=( 0.d0, 0.d0)
2529   char_mat(1,18)=( 0.d0, 0.d0)
2530
2531   name_rap(2)="G_8+   "
2532   char_mat(2,1)=( 2.d0, 0.d0)
2533   char_mat(2,2)=(-2.d0, 0.d0)
2534   char_mat(2,3)=CMPLX(-sqrt3, 0.d0,kind=DP)
2535   char_mat(2,4)=CMPLX( sqrt3, 0.d0,kind=DP)
2536   char_mat(2,5)=( 1.d0, 0.d0)
2537   char_mat(2,6)=(-1.d0, 0.d0)
2538   char_mat(2,7)=( 0.d0, 0.d0)
2539   char_mat(2,8)=( 0.d0, 0.d0)
2540   char_mat(2,9)=( 0.d0, 0.d0)
2541   char_mat(2,10)=( 2.d0, 0.d0)
2542   char_mat(2,11)=(-2.d0, 0.d0)
2543   char_mat(2,12)=CMPLX(-sqrt3, 0.d0,kind=DP)
2544   char_mat(2,13)=CMPLX( sqrt3, 0.d0,kind=DP)
2545   char_mat(2,14)=( 1.d0, 0.d0)
2546   char_mat(2,15)=(-1.d0, 0.d0)
2547   char_mat(2,16)=( 0.d0, 0.d0)
2548   char_mat(2,17)=( 0.d0, 0.d0)
2549   char_mat(2,18)=( 0.d0, 0.d0)
2550
2551   name_rap(3)="G_9+   "
2552   char_mat(3,1)=( 2.d0, 0.d0)
2553   char_mat(3,2)=(-2.d0, 0.d0)
2554   char_mat(3,3)=( 0.d0, 0.d0)
2555   char_mat(3,4)=( 0.d0, 0.d0)
2556   char_mat(3,5)=(-2.d0, 0.d0)
2557   char_mat(3,6)=( 2.d0, 0.d0)
2558   char_mat(3,7)=( 0.d0, 0.d0)
2559   char_mat(3,8)=( 0.d0, 0.d0)
2560   char_mat(3,9)=( 0.d0, 0.d0)
2561   char_mat(3,10)=( 2.d0, 0.d0)
2562   char_mat(3,11)=(-2.d0, 0.d0)
2563   char_mat(3,12)=( 0.d0, 0.d0)
2564   char_mat(3,13)=( 0.d0, 0.d0)
2565   char_mat(3,14)=(-2.d0, 0.d0)
2566   char_mat(3,15)=( 2.d0, 0.d0)
2567   char_mat(3,16)=( 0.d0, 0.d0)
2568   char_mat(3,17)=( 0.d0, 0.d0)
2569   char_mat(3,18)=( 0.d0, 0.d0)
2570
2571   name_rap(4)="G_7- "
2572   char_mat(4,1)=( 2.d0, 0.d0)
2573   char_mat(4,2)=(-2.d0, 0.d0)
2574   char_mat(4,3)=CMPLX( sqrt3, 0.d0,kind=DP)
2575   char_mat(4,4)=CMPLX(-sqrt3, 0.d0,kind=DP)
2576   char_mat(4,5)=( 1.d0, 0.d0)
2577   char_mat(4,6)=(-1.d0, 0.d0)
2578   char_mat(4,7)=( 0.d0, 0.d0)
2579   char_mat(4,8)=( 0.d0, 0.d0)
2580   char_mat(4,9)=( 0.d0, 0.d0)
2581   char_mat(4,10)=(-2.d0, 0.d0)
2582   char_mat(4,11)=( 2.d0, 0.d0)
2583   char_mat(4,12)=CMPLX(-sqrt3, 0.d0,kind=DP)
2584   char_mat(4,13)=CMPLX( sqrt3, 0.d0,kind=DP)
2585   char_mat(4,14)=(-1.d0, 0.d0)
2586   char_mat(4,15)=( 1.d0, 0.d0)
2587   char_mat(4,16)=( 0.d0, 0.d0)
2588   char_mat(4,17)=( 0.d0, 0.d0)
2589   char_mat(4,18)=( 0.d0, 0.d0)
2590
2591   name_rap(5)="G_8-   "
2592   char_mat(5,1)=( 2.d0, 0.d0)
2593   char_mat(5,2)=(-2.d0, 0.d0)
2594   char_mat(5,3)=CMPLX(-sqrt3, 0.d0,kind=DP)
2595   char_mat(5,4)=CMPLX( sqrt3, 0.d0,kind=DP)
2596   char_mat(5,5)=( 1.d0, 0.d0)
2597   char_mat(5,6)=(-1.d0, 0.d0)
2598   char_mat(5,7)=( 0.d0, 0.d0)
2599   char_mat(5,8)=( 0.d0, 0.d0)
2600   char_mat(5,9)=( 0.d0, 0.d0)
2601   char_mat(5,10)=(-2.d0, 0.d0)
2602   char_mat(5,11)=( 2.d0, 0.d0)
2603   char_mat(5,12)=CMPLX( sqrt3, 0.d0,kind=DP)
2604   char_mat(5,13)=CMPLX(-sqrt3, 0.d0,kind=DP)
2605   char_mat(5,14)=(-1.d0, 0.d0)
2606   char_mat(5,15)=( 1.d0, 0.d0)
2607   char_mat(5,16)=( 0.d0, 0.d0)
2608   char_mat(5,17)=( 0.d0, 0.d0)
2609   char_mat(5,18)=( 0.d0, 0.d0)
2610
2611   name_rap(6)="G_9-   "
2612   char_mat(6,1)=( 2.d0, 0.d0)
2613   char_mat(6,2)=(-2.d0, 0.d0)
2614   char_mat(6,3)=( 0.d0, 0.d0)
2615   char_mat(6,4)=( 0.d0, 0.d0)
2616   char_mat(6,5)=(-2.d0, 0.d0)
2617   char_mat(6,6)=( 2.d0, 0.d0)
2618   char_mat(6,7)=( 0.d0, 0.d0)
2619   char_mat(6,8)=( 0.d0, 0.d0)
2620   char_mat(6,9)=( 0.d0, 0.d0)
2621   char_mat(6,10)=(-2.d0, 0.d0)
2622   char_mat(6,11)=( 2.d0, 0.d0)
2623   char_mat(6,12)=( 0.d0, 0.d0)
2624   char_mat(6,13)=( 0.d0, 0.d0)
2625   char_mat(6,14)=( 2.d0, 0.d0)
2626   char_mat(6,15)=(-2.d0, 0.d0)
2627   char_mat(6,16)=( 0.d0, 0.d0)
2628   char_mat(6,17)=( 0.d0, 0.d0)
2629   char_mat(6,18)=( 0.d0, 0.d0)
2630
2631ELSEIF (code_group==24) THEN
2632!
2633! D_2d
2634!
2635   nclass_ref=7
2636   name_class(3)="2S4 "
2637   name_class(4)="-2S4"
2638   name_class(5)="C2  "
2639   name_class1(5)="-C2 "
2640   name_class(6)="2C2' "
2641   name_class1(6)="-2C2'"
2642   name_class(7)="2s_d "
2643   name_class1(7)="-2s_d"
2644
2645   nrap_ref=2
2646
2647   name_rap(1)="G_6  X_6  W_6"
2648   char_mat(1,1)=( 2.d0, 0.d0)
2649   char_mat(1,2)=(-2.d0, 0.d0)
2650   char_mat(1,3)=CMPLX( sqrt2, 0.d0,kind=DP)
2651   char_mat(1,4)=CMPLX(-sqrt2, 0.d0,kind=DP)
2652   char_mat(1,5)=( 0.d0, 0.d0)
2653   char_mat(1,6)=( 0.d0, 0.d0)
2654   char_mat(1,7)=( 0.d0, 0.d0)
2655
2656   name_rap(2)="G_7  X_7  W_7"
2657   char_mat(2,1)=( 2.d0, 0.d0)
2658   char_mat(2,2)=(-2.d0, 0.d0)
2659   char_mat(2,3)=CMPLX(-sqrt2, 0.d0,kind=DP)
2660   char_mat(2,4)=CMPLX( sqrt2, 0.d0,kind=DP)
2661   char_mat(2,5)=( 0.d0, 0.d0)
2662   char_mat(2,6)=( 0.d0, 0.d0)
2663   char_mat(2,7)=( 0.d0, 0.d0)
2664
2665ELSEIF (code_group==25) THEN
2666!
2667! D_3d
2668!
2669   nclass_ref=12
2670   name_class(3)="2C3 "
2671   name_class(4)="-2C3"
2672   name_class(5)="3C2'"
2673   name_class(6)="-3C2'"
2674   name_class(7)="i  "
2675   name_class(8)="-i  "
2676   name_class(9)="2S6 "
2677   name_class(10)="-2S6"
2678   name_class(11)="3s_v "
2679   name_class(12)="-3s_v"
2680
2681   nrap_ref=6
2682
2683   name_rap(1)="G_4+  L_6+"
2684   char_mat(1,1)=( 2.d0, 0.d0)
2685   char_mat(1,2)=(-2.d0, 0.d0)
2686   char_mat(1,4)=(-1.d0, 0.d0)
2687   char_mat(1,5)=( 0.d0, 0.d0)
2688   char_mat(1,6)=( 0.d0, 0.d0)
2689   char_mat(1,7)=( 2.d0, 0.d0)
2690   char_mat(1,8)=(-2.d0, 0.d0)
2691   char_mat(1,10)=(-1.d0, 0.d0)
2692   char_mat(1,11)=( 0.d0, 0.d0)
2693   char_mat(1,12)=( 0.d0, 0.d0)
2694
2695   name_rap(2)="G_5+  L_4+"
2696   char_mat(2,2)=(-1.d0, 0.d0)
2697   char_mat(2,3)=(-1.d0, 0.d0)
2698   char_mat(2,5)=( 0.d0, 1.d0)
2699   char_mat(2,6)=( 0.d0,-1.d0)
2700   char_mat(2,8)=(-1.d0, 0.d0)
2701   char_mat(2,9)=(-1.d0, 0.d0)
2702   char_mat(2,11)=( 0.d0, 1.d0)
2703   char_mat(2,12)=( 0.d0,-1.d0)
2704
2705   name_rap(3)="G_6+  L_5+"
2706   char_mat(3,2)=(-1.d0, 0.d0)
2707   char_mat(3,3)=(-1.d0, 0.d0)
2708   char_mat(3,5)=( 0.d0,-1.d0)
2709   char_mat(3,6)=( 0.d0, 1.d0)
2710   char_mat(3,8)=(-1.d0, 0.d0)
2711   char_mat(3,9)=(-1.d0, 0.d0)
2712   char_mat(3,11)=( 0.d0,-1.d0)
2713   char_mat(3,12)=( 0.d0, 1.d0)
2714
2715   name_rap(4)="G_4-  L_6-"
2716   char_mat(4,1)=( 2.d0, 0.d0)
2717   char_mat(4,2)=(-2.d0, 0.d0)
2718   char_mat(4,4)=(-1.d0, 0.d0)
2719   char_mat(4,5)=( 0.d0, 0.d0)
2720   char_mat(4,6)=( 0.d0, 0.d0)
2721   char_mat(4,7)=(-2.d0, 0.d0)
2722   char_mat(4,8)=( 2.d0, 0.d0)
2723   char_mat(4,9)=(-1.d0, 0.d0)
2724   char_mat(4,11)=( 0.d0, 0.d0)
2725   char_mat(4,12)=( 0.d0, 0.d0)
2726
2727   name_rap(5)="G_5-  L_4-"
2728   char_mat(5,2)=(-1.d0, 0.d0)
2729   char_mat(5,3)=(-1.d0, 0.d0)
2730   char_mat(5,5)=( 0.d0, 1.d0)
2731   char_mat(5,6)=( 0.d0,-1.d0)
2732   char_mat(5,7)=(-1.d0, 0.d0)
2733   char_mat(5,10)=(-1.d0, 0.d0)
2734   char_mat(5,11)=( 0.d0,-1.d0)
2735   char_mat(5,12)=( 0.d0, 1.d0)
2736
2737   name_rap(6)="G_6-  L_5-"
2738   char_mat(6,2)=(-1.d0, 0.d0)
2739   char_mat(6,3)=(-1.d0, 0.d0)
2740   char_mat(6,5)=( 0.d0,-1.d0)
2741   char_mat(6,6)=( 0.d0, 1.d0)
2742   char_mat(6,7)=(-1.d0, 0.d0)
2743   char_mat(6,10)=(-1.d0, 0.d0)
2744   char_mat(6,11)=( 0.d0, 1.d0)
2745   char_mat(6,12)=( 0.d0,-1.d0)
2746
2747ELSEIF (code_group==26) THEN
2748!
2749! S_4 ! This character table has been found to be working in at least one case.
2750!       NB: The signs of the characters reported by Koster, Space groups and
2751!            their representations are in the comment. They do not work at
2752!            least in one case. The characters in the tables in Koster,
2753!            Dimmock, Wheeler, Statz, Properties of the 32 point groups
2754!            do not match neither those used here nor those of Koster.
2755!            They do not work at least in one case.
2756!            Please report any problem that you might find with S_4.
2757!
2758   nclass_ref=8
2759   name_class(3)=" S4^3"
2760   name_class(4)="-S4^3"
2761   name_class(5)=" C2 "
2762   name_class(6)="-C2 "
2763   name_class(7)=" S4 "
2764   name_class(8)="-S4 "
2765
2766   nrap_ref=4
2767
2768   name_rap(1)="G_5 "
2769   char_mat(1,3)=CMPLX( dsq2, dsq2,kind=DP)
2770   char_mat(1,4)=CMPLX(-dsq2,-dsq2,kind=DP)
2771   char_mat(1,5)=( 0.d0, 1.d0)
2772   char_mat(1,6)=( 0.d0,-1.d0)
2773!   char_mat(1,7)=CMPLX(-dsq2, dsq2,kind=DP)
2774!   char_mat(1,8)=CMPLX( dsq2,-dsq2,kind=DP)
2775   char_mat(1,7)=CMPLX( dsq2,-dsq2,kind=DP)
2776   char_mat(1,8)=CMPLX(-dsq2, dsq2,kind=DP)
2777
2778   name_rap(2)="G_6 "
2779   char_mat(2,3)=CMPLX( dsq2,-dsq2,kind=DP)
2780   char_mat(2,4)=CMPLX(-dsq2, dsq2,kind=DP)
2781   char_mat(2,5)=( 0.d0,-1.d0)
2782   char_mat(2,6)=( 0.d0, 1.d0)
2783!   char_mat(2,7)=CMPLX(-dsq2,-dsq2,kind=DP)
2784!   char_mat(2,8)=CMPLX( dsq2, dsq2,kind=DP)
2785   char_mat(2,7)=CMPLX( dsq2, dsq2,kind=DP)
2786   char_mat(2,8)=CMPLX(-dsq2,-dsq2,kind=DP)
2787
2788   name_rap(3)="G_7 "
2789   char_mat(3,3)=CMPLX(-dsq2,-dsq2,kind=DP)
2790   char_mat(3,4)=CMPLX( dsq2, dsq2,kind=DP)
2791   char_mat(3,5)=( 0.d0, 1.d0)
2792   char_mat(3,6)=( 0.d0,-1.d0)
2793!   char_mat(3,7)=CMPLX( dsq2,-dsq2,kind=DP)
2794!   char_mat(3,8)=CMPLX(-dsq2, dsq2,kind=DP)
2795   char_mat(3,7)=CMPLX(-dsq2, dsq2,kind=DP)
2796   char_mat(3,8)=CMPLX( dsq2,-dsq2,kind=DP)
2797
2798   name_rap(4)="G_8 "
2799   char_mat(4,3)=CMPLX(-dsq2, dsq2,kind=DP)
2800   char_mat(4,4)=CMPLX( dsq2,-dsq2,kind=DP)
2801   char_mat(4,5)=( 0.d0,-1.d0)
2802   char_mat(4,6)=( 0.d0, 1.d0)
2803!   char_mat(4,7)=CMPLX( dsq2, dsq2,kind=DP)
2804!   char_mat(4,8)=CMPLX(-dsq2,-dsq2,kind=DP)
2805   char_mat(4,7)=CMPLX(-dsq2,-dsq2,kind=DP)
2806   char_mat(4,8)=CMPLX( dsq2, dsq2,kind=DP)
2807
2808ELSEIF (code_group==27) THEN
2809!
2810! S_6    NB: The signs of the characters of the classes C3^2 -C3^2 S6 -S6
2811!            are changed with respect to Koster, Space groups and
2812!            their representations. They match the table in Koster,
2813!            Dimmock, Wheeler, Statz, Properties of the 32 point groups.
2814!
2815!
2816   nclass_ref=12
2817   name_class(3)="C3  "
2818   name_class(4)="-C3 "
2819   name_class(5)="C3^2"
2820   name_class(6)="-C3^2"
2821   name_class(7)=" i  "
2822   name_class(8)="-i  "
2823   name_class(9)="S6^5"
2824   name_class(10)="-S6^5"
2825   name_class(11)="S6"
2826   name_class(12)="-S6"
2827
2828   nrap_ref=6
2829
2830   name_rap(1)="G_4+"
2831   char_mat(1,3)=CMPLX( 0.5d0, sqr3d2 ,kind=DP)
2832   char_mat(1,4)=CMPLX(-0.5d0,-sqr3d2 ,kind=DP)
2833
2834!   char_mat(1,5)=CMPLX(-0.5d0, sqr3d2 ,kind=DP)
2835!   char_mat(1,6)=CMPLX( 0.5d0,-sqr3d2 ,kind=DP)
2836   char_mat(1,5)=CMPLX( 0.5d0,-sqr3d2 ,kind=DP)
2837   char_mat(1,6)=CMPLX(-0.5d0, sqr3d2 ,kind=DP)
2838
2839   char_mat(1,8)=(-1.0d0, 0.d0 )
2840   char_mat(1,9)=CMPLX( 0.5d0, sqr3d2 ,kind=DP)
2841   char_mat(1,10)=CMPLX(-0.5d0,-sqr3d2 ,kind=DP)
2842!   char_mat(1,11)=CMPLX(-0.5d0, sqr3d2 ,kind=DP)
2843!   char_mat(1,12)=CMPLX( 0.5d0,-sqr3d2 ,kind=DP)
2844   char_mat(1,11)=CMPLX( 0.5d0,-sqr3d2 ,kind=DP)
2845   char_mat(1,12)=CMPLX(-0.5d0, sqr3d2 ,kind=DP)
2846
2847   name_rap(2)="G_5+"
2848   char_mat(2,3)=CMPLX( 0.5d0,-sqr3d2 ,kind=DP)
2849   char_mat(2,4)=CMPLX(-0.5d0, sqr3d2 ,kind=DP)
2850
2851!   char_mat(2,5)=CMPLX(-0.5d0,-sqr3d2 ,kind=DP)
2852!   char_mat(2,6)=CMPLX( 0.5d0, sqr3d2 ,kind=DP)
2853   char_mat(2,5)=CMPLX( 0.5d0, sqr3d2 ,kind=DP)
2854   char_mat(2,6)=CMPLX(-0.5d0,-sqr3d2 ,kind=DP)
2855
2856   char_mat(2,8)=(-1.0d0, 0.d0 )
2857   char_mat(2,9)=CMPLX( 0.5d0,-sqr3d2 ,kind=DP)
2858   char_mat(2,10)=CMPLX(-0.5d0, sqr3d2 ,kind=DP)
2859!   char_mat(2,11)=CMPLX(-0.5d0,-sqr3d2 ,kind=DP)
2860!   char_mat(2,12)=CMPLX( 0.5d0, sqr3d2 ,kind=DP)
2861   char_mat(2,11)=CMPLX( 0.5d0, sqr3d2 ,kind=DP)
2862   char_mat(2,12)=CMPLX(-0.5d0,-sqr3d2 ,kind=DP)
2863
2864   name_rap(3)="G_6+"
2865   char_mat(3,3)=(-1.0d0, 0.d0 )
2866!   char_mat(3,6)=(-1.0d0, 0.d0 )
2867   char_mat(3,5)=(-1.0d0, 0.d0 )
2868   char_mat(3,8)=(-1.0d0, 0.d0 )
2869   char_mat(3,9)=(-1.0d0, 0.d0 )
2870!   char_mat(3,12)=(-1.0d0, 0.d0 )
2871   char_mat(3,11)=(-1.0d0, 0.d0 )
2872
2873   name_rap(4)="G_4-"
2874   char_mat(4,3)=CMPLX( 0.5d0, sqr3d2 ,kind=DP)
2875   char_mat(4,4)=CMPLX(-0.5d0,-sqr3d2 ,kind=DP)
2876!   char_mat(4,5)=CMPLX(-0.5d0, sqr3d2 ,kind=DP)
2877!   char_mat(4,6)=CMPLX( 0.5d0,-sqr3d2 ,kind=DP)
2878   char_mat(4,5)=CMPLX( 0.5d0,-sqr3d2 ,kind=DP)
2879   char_mat(4,6)=CMPLX(-0.5d0, sqr3d2 ,kind=DP)
2880   char_mat(4,7)=(-1.0d0, 0.d0 )
2881   char_mat(4,9)=CMPLX(-0.5d0,-sqr3d2 ,kind=DP)
2882   char_mat(4,10)=CMPLX( 0.5d0, sqr3d2 ,kind=DP)
2883!   char_mat(4,11)=CMPLX( 0.5d0,-sqr3d2 ,kind=DP)
2884!   char_mat(4,12)=CMPLX(-0.5d0, sqr3d2 ,kind=DP)
2885   char_mat(4,11)=CMPLX(-0.5d0, sqr3d2 ,kind=DP)
2886   char_mat(4,12)=CMPLX( 0.5d0,-sqr3d2 ,kind=DP)
2887
2888   name_rap(5)="G_5-"
2889   char_mat(5,3)=CMPLX( 0.5d0,-sqr3d2 ,kind=DP)
2890   char_mat(5,4)=CMPLX(-0.5d0, sqr3d2 ,kind=DP)
2891!   char_mat(5,5)=CMPLX(-0.5d0,-sqr3d2 ,kind=DP)
2892!   char_mat(5,6)=CMPLX( 0.5d0, sqr3d2 ,kind=DP)
2893   char_mat(5,5)=CMPLX( 0.5d0, sqr3d2 ,kind=DP)
2894   char_mat(5,6)=CMPLX(-0.5d0,-sqr3d2 ,kind=DP)
2895   char_mat(5,7)=(-1.0d0, 0.d0 )
2896   char_mat(5,9)=CMPLX(-0.5d0, sqr3d2 ,kind=DP)
2897   char_mat(5,10)=CMPLX( 0.5d0,-sqr3d2 ,kind=DP)
2898!   char_mat(5,11)=CMPLX( 0.5d0, sqr3d2 ,kind=DP)
2899!   char_mat(5,12)=CMPLX(-0.5d0,-sqr3d2 ,kind=DP)
2900   char_mat(5,11)=CMPLX(-0.5d0,-sqr3d2 ,kind=DP)
2901   char_mat(5,12)=CMPLX( 0.5d0, sqr3d2 ,kind=DP)
2902
2903   name_rap(6)="G_6-"
2904   char_mat(6,3)=(-1.0d0, 0.d0 )
2905!   char_mat(6,6)=(-1.0d0, 0.d0 )
2906   char_mat(6,5)=(-1.0d0, 0.d0 )
2907   char_mat(6,7)=(-1.0d0, 0.d0 )
2908   char_mat(6,10)=(-1.0d0, 0.d0 )
2909!   char_mat(6,11)=(-1.0d0, 0.d0 )
2910   char_mat(6,12)=(-1.0d0, 0.d0 )
2911
2912
2913ELSEIF (code_group==28) THEN
2914!
2915! NB: The signs of the characters of the classes C3 -C3 C3^2 -C3^2
2916!            are changed with respect to Koster, Space groups and
2917!            their representations. They match the table in Koster,
2918!            Dimmock, Wheeler, Statz, Properties of the 32 point groups.
2919
2920!
2921! T
2922!
2923   nclass_ref=7
2924   name_class(3)="3C2  "
2925   name_class1(3)="-3C2 "
2926   name_class(4)=" 4C3 "
2927   name_class(5)="-4C3 "
2928   name_class(6)=" 4C3'"
2929   name_class(7)="-4C3'"
2930
2931   nrap_ref=3
2932
2933   name_rap(1)="G_5 "
2934   char_mat(1,1)=( 2.d0, 0.d0)
2935   char_mat(1,2)=(-2.d0, 0.d0)
2936   char_mat(1,3)=( 0.d0, 0.d0)
2937   char_mat(1,5)=(-1.d0, 0.d0)
2938   char_mat(1,7)=(-1.d0, 0.d0)
2939
2940   name_rap(2)="G_6 "
2941   char_mat(2,1)=( 2.d0, 0.d0)
2942   char_mat(2,2)=(-2.d0, 0.d0)
2943   char_mat(2,3)=( 0.d0, 0.d0)
2944!   char_mat(2,4)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2945!   char_mat(2,5)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2946!   char_mat(2,6)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2947!   char_mat(2,7)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2948   char_mat(2,4)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2949   char_mat(2,5)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2950   char_mat(2,6)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2951   char_mat(2,7)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2952
2953   name_rap(3)="G_7 "
2954   char_mat(3,1)=( 2.d0, 0.d0)
2955   char_mat(3,2)=(-2.d0, 0.d0)
2956   char_mat(3,3)=( 0.d0, 0.d0)
2957!   char_mat(3,4)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2958!   char_mat(3,5)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2959!   char_mat(3,6)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2960!   char_mat(3,7)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2961   char_mat(3,4)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
2962   char_mat(3,5)=CMPLX( 0.5d0, sqr3d2,kind=DP)
2963   char_mat(3,6)=CMPLX(-0.5d0, sqr3d2,kind=DP)
2964   char_mat(3,7)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
2965
2966ELSEIF (code_group==29) THEN
2967! NB: The signs of the characters of the some classes
2968!            are changed with respect to Koster, Space groups and
2969!            their representations. They match the table in Koster,
2970!            Dimmock, Wheeler, Statz, Properties of the 32 point groups.
2971!
2972! T_h
2973!
2974   nclass_ref=14
2975   name_class(3)=" 3C2 "
2976   name_class1(3)="-3C2 "
2977   name_class(4)=" 4C3 "
2978   name_class(5)="-4C3 "
2979   name_class(6)=" 4C3'"
2980   name_class(7)="-4C3'"
2981   name_class(8)="i   "
2982   name_class(9)="-i   "
2983   name_class(10)=" 3s_h"
2984   name_class1(10)="-3s_h"
2985   name_class(11)="4S6'"
2986   name_class(12)="-4S6'"
2987   name_class(13)=" 4S6 "
2988   name_class(14)="-4S6 "
2989
2990   nrap_ref=6
2991
2992   name_rap(1)="G_5+"
2993   char_mat(1,1)=( 2.d0, 0.d0)
2994   char_mat(1,2)=(-2.d0, 0.d0)
2995   char_mat(1,3)=( 0.d0, 0.d0)
2996   char_mat(1,5)=(-1.d0, 0.d0)
2997   char_mat(1,7)=(-1.d0, 0.d0)
2998   char_mat(1,8)=( 2.d0, 0.d0)
2999   char_mat(1,9)=(-2.d0, 0.d0)
3000   char_mat(1,10)=( 0.d0, 0.d0)
3001   char_mat(1,12)=(-1.d0, 0.d0)
3002   char_mat(1,14)=(-1.d0, 0.d0)
3003
3004   name_rap(2)="G_6+"
3005   char_mat(2,1)=( 2.d0, 0.d0)
3006   char_mat(2,2)=(-2.d0, 0.d0)
3007   char_mat(2,3)=( 0.d0, 0.d0)
3008!   char_mat(2,4)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3009!   char_mat(2,5)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3010!   char_mat(2,6)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3011!   char_mat(2,7)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3012   char_mat(2,4)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3013   char_mat(2,5)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3014   char_mat(2,6)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3015   char_mat(2,7)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3016   char_mat(2,8)=( 2.d0, 0.d0)
3017   char_mat(2,9)=(-2.d0, 0.d0)
3018   char_mat(2,10)=( 0.d0, 0.d0)
3019!   char_mat(2,11)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3020!   char_mat(2,12)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3021!   char_mat(2,13)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3022!   char_mat(2,14)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3023   char_mat(2,11)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3024   char_mat(2,12)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3025   char_mat(2,13)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3026   char_mat(2,14)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3027
3028   name_rap(3)="G_7+"
3029   char_mat(3,1)=( 2.d0, 0.d0)
3030   char_mat(3,2)=(-2.d0, 0.d0)
3031   char_mat(3,3)=( 0.d0, 0.d0)
3032!   char_mat(3,4)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3033!   char_mat(3,5)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3034!   char_mat(3,6)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3035!   char_mat(3,7)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3036   char_mat(3,4)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3037   char_mat(3,5)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3038   char_mat(3,6)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3039   char_mat(3,7)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3040   char_mat(3,8)=( 2.d0, 0.d0)
3041   char_mat(3,9)=(-2.d0, 0.d0)
3042   char_mat(3,10)=( 0.d0, 0.d0)
3043!   char_mat(3,11)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3044!   char_mat(3,12)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3045!   char_mat(3,13)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3046!   char_mat(3,14)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3047   char_mat(3,11)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3048   char_mat(3,12)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3049   char_mat(3,13)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3050   char_mat(3,14)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3051
3052   name_rap(4)="G_5-"
3053   char_mat(4,1)=( 2.d0, 0.d0)
3054   char_mat(4,2)=(-2.d0, 0.d0)
3055   char_mat(4,3)=( 0.d0, 0.d0)
3056   char_mat(4,5)=(-1.d0, 0.d0)
3057   char_mat(4,7)=(-1.d0, 0.d0)
3058   char_mat(4,8)=(-2.d0, 0.d0)
3059   char_mat(4,9)=( 2.d0, 0.d0)
3060   char_mat(4,10)=( 0.d0, 0.d0)
3061   char_mat(4,11)=(-1.d0, 0.d0)
3062   char_mat(4,13)=(-1.d0, 0.d0)
3063
3064
3065   name_rap(5)="G_6-"
3066   char_mat(5,1)=( 2.d0, 0.d0)
3067   char_mat(5,2)=(-2.d0, 0.d0)
3068   char_mat(5,3)=( 0.d0, 0.d0)
3069!   char_mat(5,4)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3070!   char_mat(5,5)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3071!   char_mat(5,6)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3072!   char_mat(5,7)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3073   char_mat(5,4)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3074   char_mat(5,5)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3075   char_mat(5,6)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3076   char_mat(5,7)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3077   char_mat(5,8)=(-2.d0, 0.d0)
3078   char_mat(5,9)=( 2.d0, 0.d0)
3079   char_mat(5,10)=( 0.d0, 0.d0)
3080!   char_mat(5,11)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3081!   char_mat(5,12)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3082!   char_mat(5,13)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3083!   char_mat(5,14)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3084   char_mat(5,11)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3085   char_mat(5,12)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3086   char_mat(5,13)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3087   char_mat(5,14)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3088
3089   name_rap(6)="G_7-"
3090   char_mat(6,1)=( 2.d0, 0.d0)
3091   char_mat(6,2)=(-2.d0, 0.d0)
3092   char_mat(6,3)=( 0.d0, 0.d0)
3093!   char_mat(6,4)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3094!   char_mat(6,5)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3095!   char_mat(6,6)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3096!   char_mat(6,7)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3097   char_mat(6,4)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3098   char_mat(6,5)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3099   char_mat(6,6)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3100   char_mat(6,7)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3101   char_mat(6,8)=(-2.d0, 0.d0)
3102   char_mat(6,9)=( 2.d0, 0.d0)
3103   char_mat(6,10)=( 0.d0, 0.d0)
3104!   char_mat(6,11)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3105!   char_mat(6,12)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3106!   char_mat(6,13)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3107!   char_mat(6,14)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3108   char_mat(6,11)=CMPLX( 0.5d0, sqr3d2,kind=DP)
3109   char_mat(6,12)=CMPLX(-0.5d0,-sqr3d2,kind=DP)
3110   char_mat(6,13)=CMPLX( 0.5d0,-sqr3d2,kind=DP)
3111   char_mat(6,14)=CMPLX(-0.5d0, sqr3d2,kind=DP)
3112
3113ELSEIF (code_group==30) THEN
3114!
3115! T_d
3116!
3117   nclass_ref=8
3118   name_class(3)="8C3  "
3119   name_class(4)="-8C3 "
3120   name_class(5)=" 3C2 "
3121   name_class1(5)="-3C2 "
3122   name_class(6)="6S4 "
3123   name_class(7)="-6S4 "
3124   name_class(8)="6s_d"
3125   name_class1(8)="-6s_d"
3126
3127   nrap_ref=3
3128
3129   name_rap(1)="G_6  P_6"
3130   char_mat(1,1)=( 2.d0, 0.d0)
3131   char_mat(1,2)=(-2.d0, 0.d0)
3132   char_mat(1,4)=(-1.d0, 0.d0)
3133   char_mat(1,5)=( 0.d0, 0.d0)
3134   char_mat(1,6)=CMPLX( sqrt2, 0.d0,kind=DP)
3135   char_mat(1,7)=CMPLX(-sqrt2, 0.d0,kind=DP)
3136   char_mat(1,8)=( 0.d0, 0.d0)
3137
3138   name_rap(2)="G_7  P_7"
3139   char_mat(2,1)=( 2.d0, 0.d0)
3140   char_mat(2,2)=(-2.d0, 0.d0)
3141   char_mat(2,4)=(-1.d0, 0.d0)
3142   char_mat(2,5)=( 0.d0, 0.d0)
3143   char_mat(2,6)=CMPLX(-sqrt2, 0.d0,kind=DP)
3144   char_mat(2,7)=CMPLX( sqrt2, 0.d0,kind=DP)
3145   char_mat(2,8)=( 0.d0, 0.d0)
3146
3147   name_rap(3)="G_8  P_8"
3148   char_mat(3,1)=( 4.d0, 0.d0)
3149   char_mat(3,2)=(-4.d0, 0.d0)
3150   char_mat(3,3)=(-1.d0, 0.d0)
3151   char_mat(3,5)=( 0.d0, 0.d0)
3152   char_mat(3,6)=( 0.d0, 0.d0)
3153   char_mat(3,7)=( 0.d0, 0.d0)
3154   char_mat(3,8)=( 0.d0, 0.d0)
3155
3156
3157ELSEIF (code_group==31) THEN
3158!
3159! O
3160!
3161   nclass_ref=8
3162   name_class(3)="8C3  "
3163   name_class(4)="-8C3 "
3164   name_class(5)=" 3C2"
3165   name_class1(5)="-3C2"
3166   name_class(6)="6C4  "
3167   name_class(7)="-6C4 "
3168   name_class(8)=" 6C2'"
3169   name_class1(8)="-6C2'"
3170
3171   nrap_ref=3
3172
3173   name_rap(1)="G_6  "
3174   char_mat(1,1)=(  2.d0, 0.d0)
3175   char_mat(1,2)=( -2.d0, 0.d0)
3176   char_mat(1,4)=( -1.d0, 0.d0)
3177   char_mat(1,5)=(  0.d0, 0.d0)
3178   char_mat(1,6)=CMPLX( sqrt2, 0.d0,kind=DP)
3179   char_mat(1,7)=CMPLX(-sqrt2, 0.d0,kind=DP)
3180   char_mat(1,8)=(  0.d0, 0.d0)
3181
3182   name_rap(2)="G_7  "
3183   char_mat(2,1)=(  2.d0, 0.d0)
3184   char_mat(2,2)=( -2.d0, 0.d0)
3185   char_mat(2,4)=( -1.d0, 0.d0)
3186   char_mat(2,5)=(  0.d0, 0.d0)
3187   char_mat(2,6)=CMPLX(-sqrt2, 0.d0,kind=DP)
3188   char_mat(2,7)=CMPLX( sqrt2, 0.d0,kind=DP)
3189   char_mat(2,8)=(  0.d0, 0.d0)
3190
3191   name_rap(3)="G_8  "
3192   char_mat(3,1)=(  4.d0, 0.d0)
3193   char_mat(3,2)=( -4.d0, 0.d0)
3194   char_mat(3,3)=( -1.d0, 0.d0)
3195   char_mat(3,5)=(  0.d0, 0.d0)
3196   char_mat(3,6)=(  0.d0, 0.d0)
3197   char_mat(3,7)=(  0.d0, 0.d0)
3198   char_mat(3,8)=(  0.d0, 0.d0)
3199
3200
3201ELSEIF (code_group==32) THEN
3202!
3203! O_h
3204!
3205   nclass_ref=16
3206   name_class(3)="8C3  "
3207   name_class(4)="-8C3 "
3208   name_class(5)=" 3C2"
3209   name_class1(5)="-3C2"
3210   name_class(6)="6C4  "
3211   name_class(7)="-6C4 "
3212   name_class(8)=" 6C2'"
3213   name_class1(8)="-6C2'"
3214   name_class(9)="i "
3215   name_class(10)="-i  "
3216   name_class(11)="8S6  "
3217   name_class(12)="-8S6 "
3218   name_class(13)=" 3s_h"
3219   name_class1(13)="-3s_h"
3220   name_class(14)="6S4  "
3221   name_class(15)="-6S4 "
3222   name_class(16)=" 6s_d"
3223   name_class1(16)="-6s_d"
3224
3225   nrap_ref=6
3226
3227   name_rap(1)="G_6+  "
3228   char_mat(1,1)=( 2.d0, 0.d0)
3229   char_mat(1,2)=(-2.d0, 0.d0)
3230   char_mat(1,4)=(-1.d0, 0.d0)
3231   char_mat(1,5)=( 0.d0, 0.d0)
3232   char_mat(1,6)=CMPLX( sqrt2, 0.d0,kind=DP)
3233   char_mat(1,7)=CMPLX(-sqrt2, 0.d0,kind=DP)
3234   char_mat(1,8)=( 0.d0, 0.d0)
3235   char_mat(1,9)=( 2.d0, 0.d0)
3236   char_mat(1,10)=(-2.d0, 0.d0)
3237   char_mat(1,12)=(-1.d0, 0.d0)
3238   char_mat(1,13)=( 0.d0, 0.d0)
3239   char_mat(1,14)=CMPLX( sqrt2, 0.d0,kind=DP)
3240   char_mat(1,15)=CMPLX(-sqrt2, 0.d0,kind=DP)
3241   char_mat(1,16)=( 0.d0, 0.d0)
3242
3243   name_rap(2)="G_7+  "
3244   char_mat(2,1)=( 2.d0, 0.d0)
3245   char_mat(2,2)=(-2.d0, 0.d0)
3246   char_mat(2,4)=(-1.d0, 0.d0)
3247   char_mat(2,5)=( 0.d0, 0.d0)
3248   char_mat(2,6)=CMPLX(-sqrt2, 0.d0,kind=DP)
3249   char_mat(2,7)=CMPLX( sqrt2, 0.d0,kind=DP)
3250   char_mat(2,8)=( 0.d0, 0.d0)
3251   char_mat(2,9)=( 2.d0, 0.d0)
3252   char_mat(2,10)=(-2.d0, 0.d0)
3253   char_mat(2,12)=(-1.d0, 0.d0)
3254   char_mat(2,13)=( 0.d0, 0.d0)
3255   char_mat(2,14)=CMPLX(-sqrt2, 0.d0,kind=DP)
3256   char_mat(2,15)=CMPLX( sqrt2, 0.d0,kind=DP)
3257   char_mat(2,16)=( 0.d0, 0.d0)
3258
3259   name_rap(3)="G_8+  "
3260   char_mat(3,1)=( 4.d0, 0.d0)
3261   char_mat(3,2)=(-4.d0, 0.d0)
3262   char_mat(3,3)=(-1.d0, 0.d0)
3263   char_mat(3,5)=( 0.d0, 0.d0)
3264   char_mat(3,6)=( 0.d0, 0.d0)
3265   char_mat(3,7)=( 0.d0, 0.d0)
3266   char_mat(3,8)=( 0.d0, 0.d0)
3267   char_mat(3,9)=( 4.d0, 0.d0)
3268   char_mat(3,10)=(-4.d0, 0.d0)
3269   char_mat(3,11)=(-1.d0, 0.d0)
3270   char_mat(3,13)=( 0.d0, 0.d0)
3271   char_mat(3,14)=( 0.d0, 0.d0)
3272   char_mat(3,15)=( 0.d0, 0.d0)
3273   char_mat(3,16)=( 0.d0, 0.d0)
3274
3275   name_rap(4)="G_6-  "
3276   char_mat(4,1)=( 2.d0, 0.d0)
3277   char_mat(4,2)=(-2.d0, 0.d0)
3278   char_mat(4,4)=(-1.d0, 0.d0)
3279   char_mat(4,5)=( 0.d0, 0.d0)
3280   char_mat(4,6)=CMPLX(sqrt2, 0.d0,kind=DP)
3281   char_mat(4,7)=CMPLX(-sqrt2, 0.d0,kind=DP)
3282   char_mat(4,8)=( 0.d0, 0.d0)
3283   char_mat(4,9)=(-2.d0, 0.d0)
3284   char_mat(4,10)=( 2.d0, 0.d0)
3285   char_mat(4,11)=(-1.d0, 0.d0)
3286   char_mat(4,13)=( 0.d0, 0.d0)
3287   char_mat(4,14)=CMPLX(-sqrt2, 0.d0,kind=DP)
3288   char_mat(4,15)=CMPLX( sqrt2, 0.d0,kind=DP)
3289   char_mat(4,16)=( 0.d0, 0.d0)
3290
3291   name_rap(5)="G_7-  "
3292   char_mat(5,1)=( 2.d0, 0.d0)
3293   char_mat(5,2)=(-2.d0, 0.d0)
3294   char_mat(5,4)=(-1.d0, 0.d0)
3295   char_mat(5,5)=( 0.d0, 0.d0)
3296   char_mat(5,6)=CMPLX(-sqrt2, 0.d0,kind=DP)
3297   char_mat(5,7)=CMPLX( sqrt2, 0.d0,kind=DP)
3298   char_mat(5,8)=( 0.d0, 0.d0)
3299   char_mat(5,9)=(-2.d0, 0.d0)
3300   char_mat(5,10)=( 2.d0, 0.d0)
3301   char_mat(5,11)=(-1.d0, 0.d0)
3302   char_mat(5,13)=( 0.d0, 0.d0)
3303   char_mat(5,14)=CMPLX( sqrt2, 0.d0,kind=DP)
3304   char_mat(5,15)=CMPLX(-sqrt2, 0.d0,kind=DP)
3305   char_mat(5,16)=( 0.d0, 0.d0)
3306
3307
3308   name_rap(6)="G_8-  "
3309   char_mat(6,1)=( 4.d0, 0.d0)
3310   char_mat(6,2)=(-4.d0, 0.d0)
3311   char_mat(6,3)=(-1.d0, 0.d0)
3312   char_mat(6,5)=( 0.d0, 0.d0)
3313   char_mat(6,6)=( 0.d0, 0.d0)
3314   char_mat(6,7)=( 0.d0, 0.d0)
3315   char_mat(6,8)=( 0.d0, 0.d0)
3316   char_mat(6,9)=(-4.d0, 0.d0)
3317   char_mat(6,10)=( 4.d0, 0.d0)
3318   char_mat(6,12)=(-1.d0, 0.d0)
3319   char_mat(6,13)=( 0.d0, 0.d0)
3320   char_mat(6,14)=( 0.d0, 0.d0)
3321   char_mat(6,15)=( 0.d0, 0.d0)
3322   char_mat(6,16)=( 0.d0, 0.d0)
3323ELSE
3324   CALL errore('set_irr_rap_so','code number not allowed',1)
3325END IF
3326!
3327RETURN
3328!
3329END SUBROUTINE set_irr_rap_so
3330
3331
3332!--------------------------------------------------------------------------
3333FUNCTION is_complex_so( code )
3334!--------------------------------------------------------------------------
3335!! This function receives a code of the group and provide .true. or
3336!! .false. if the double group HAS or HAS NOT complex irreducible
3337!! representations. The order is the following:
3338!
3339!!   1  "C_1 " F    11 "D_6 " F    21 "D_3h" F    31 "O   " F
3340
3341!!   2  "C_i " F    12 "C_2v" F    22 "D_4h" F    32 "O_h " F
3342
3343!!   3  "C_s " T    13 "C_3v" T    23 "D_6h" F
3344
3345!!   4  "C_2 " T    14 "C_4v" F    24 "D_2d" F
3346
3347!!   5  "C_3 " T    15 "C_6v" F    25 "D_3d" T
3348
3349!!   6  "C_4 " T    16 "C_2h" T    26 "S_4 " T
3350
3351!!   7  "C_6 " T    17 "C_3h" T    27 "S_6 " T
3352
3353!!   8  "D_2 " F    18 "C_4h" T    28 "T   " T
3354
3355!!   9  "D_3 " T    19 "C_6h" T    29 "T_h " T
3356
3357!!   10 "D_4 " F    20 "D_2h" F    30 "T_d " F
3358!
3359IMPLICIT NONE
3360!
3361INTEGER :: code
3362LOGICAL :: is_complex_so
3363
3364LOGICAL :: complex_aux(32)
3365
3366DATA complex_aux  / .FALSE., .FALSE., .TRUE.,  .TRUE.,  .TRUE. , &
3367                    .TRUE. , .TRUE. , .FALSE., .TRUE.,  .FALSE., &
3368                    .FALSE., .FALSE., .TRUE.,  .FALSE., .FALSE., &
3369                    .TRUE.,  .TRUE. , .TRUE.,  .TRUE.,  .FALSE., &
3370                    .FALSE., .FALSE., .FALSE., .FALSE., .TRUE.,  &
3371                    .TRUE. , .TRUE. , .TRUE. , .TRUE. , .FALSE., &
3372                    .FALSE., .FALSE.  /
3373
3374IF (code < 1 .OR. code > 32 ) CALL errore('is_complex', &
3375                                          'code is out of range',1)
3376!
3377is_complex_so= complex_aux(code)
3378!
3379RETURN
3380!
3381END FUNCTION is_complex_so
3382
3383
3384!----------------------------------------------------------------------------
3385SUBROUTINE write_group_info( flag )
3386!----------------------------------------------------------------------------
3387!! This routine writes on output the main information on the point group
3388!! If flag is .false. writes only the character table. If flag is .true.
3389!! writes also the elements of each class.
3390!
3391USE rap_point_group,      ONLY : code_group, nclass, nelem, elem, which_irr, &
3392                                 char_mat, name_rap, name_class, gname,      &
3393                                 elem_name
3394USE rap_point_group_so,   ONLY : nrap, nelem_so, elem_so, has_e, &
3395                                 which_irr_so, char_mat_so, name_rap_so, &
3396                                 name_class_so, d_spin, name_class_so1,  &
3397                                 elem_name_so
3398USE rap_point_group_is,   ONLY : code_group_is, gname_is
3399USE spin_orb,             ONLY : domag
3400USE noncollin_module,     ONLY : noncolin
3401USE io_global,            ONLY : stdout
3402!
3403IMPLICIT NONE
3404!
3405INTEGER :: iclass, irot, i, idx, irap
3406LOGICAL :: is_complex, is_complex_so, flag
3407!
3408IF (noncolin) THEN
3409   IF (domag) THEN
3410      WRITE(stdout,'(/,5x,"the magnetic double point group is ", &
3411                     & a11," [",a11,"]")') &
3412                      gname, gname_is
3413      WRITE(stdout,'(5x,"using the double point group ",a11)') &
3414                      gname_is
3415   ELSE
3416      WRITE(stdout,'(/,5x,"double point group ",a11)') gname
3417   END IF
3418   WRITE(stdout,'(5x, "there are", i3," classes and",i3, &
3419                     &   " irreducible representations")') nclass, nrap
3420ELSE
3421   WRITE(stdout,'(/,5x,"point group ",a11)') gname
3422   WRITE(stdout,'(5x, "there are", i3," classes")') nclass
3423ENDIF
3424WRITE(stdout,'(5x, "the character table:")')
3425IF (noncolin) THEN
3426   WRITE(stdout,'(/,7x,12(a5,1x))') (name_class_so(irot), &
3427                                     irot=1,MIN(12,nclass))
3428   WRITE(stdout,'(7x,12(a5,1x))') (name_class_so1(irot), &
3429                                     irot=1,MIN(12,nclass))
3430
3431   DO iclass=1,nrap
3432      WRITE(stdout,'(a5,12f6.2)') name_rap_so(iclass), &
3433               (REAL(char_mat_so(iclass,irot)), irot=1,MIN(nclass,12))
3434   END DO
3435   IF (nclass > 12 ) THEN
3436      WRITE(stdout,'(/,7x,12(a5,1x))') (name_class_so(irot), &
3437                                        irot=13,nclass)
3438      WRITE(stdout,'(7x,12(a5,1x))') (name_class_so1(irot), &
3439                                        irot=13,nclass)
3440      DO iclass=1,nrap
3441         WRITE(stdout,'(a5,12f6.2)') name_rap_so(iclass), &
3442                   (REAL(char_mat_so(iclass,irot)), irot=13,nclass)
3443      END DO
3444   END IF
3445   idx=code_group
3446   IF (noncolin.and.domag) idx=code_group_is
3447   IF (is_complex_so(idx)) THEN
3448      WRITE(stdout,'(/,5x,"imaginary part")')
3449      WRITE(stdout,'(/,7x,12(a5,1x))') (name_class_so(irot), &
3450                                        irot=1,MIN(12,nclass))
3451      WRITE(stdout,'(7x,12(a5,1x))') (name_class_so1(irot), &
3452                                        irot=1,MIN(12,nclass))
3453      DO iclass=1,nrap
3454         WRITE(stdout,'(a5,12f6.2)') name_rap_so(iclass), &
3455              (AIMAG(char_mat_so(iclass,irot)),irot=1, MIN(nclass,12))
3456      END DO
3457      IF (nclass > 12 ) THEN
3458         WRITE(stdout,'(/,7x,12(a5,1x))') (name_class_so(irot), &
3459                                        irot=13,nclass)
3460         WRITE(stdout,'(7x,12(a5,1x))') (name_class_so1(irot), &
3461                                        irot=13,nclass)
3462         DO iclass=1,nrap
3463            WRITE(stdout,'(a5,12f6.2)') name_rap_so(iclass), &
3464                 (AIMAG(char_mat_so(iclass,irot)),irot=13, nclass)
3465         END DO
3466      END IF
3467   END IF
3468   IF (flag) THEN
3469      WRITE(stdout,'(/5x, "the symmetry operations in each class and &
3470                          &the name of the first element:",/)')
3471      DO irap = 1, nclass
3472         DO iclass=1,nclass
3473            IF ( which_irr_so(iclass) /= irap ) CYCLE
3474            WRITE(stdout,'(5x,2a5,12i5)') &
3475                              name_class_so(which_irr_so(iclass)), &
3476                              name_class_so1(which_irr_so(iclass)), &
3477               (elem_so(i,iclass)*has_e(i,iclass), i=1,nelem_so(iclass))
3478            WRITE(stdout,'(10x,a)') elem_name_so(1,iclass)
3479         END DO
3480      ENDDO
3481   ENDIF
3482ELSE
3483   WRITE(stdout,'(/,7x,12(a5,1x))') (name_class(irot),irot=1,nclass)
3484   DO iclass=1,nclass
3485      WRITE(stdout,'(a5,12f6.2)') name_rap(iclass), &
3486         (REAL(char_mat(iclass,irot)),irot=1,nclass)
3487   ENDDO
3488   idx=code_group
3489   IF (noncolin.and.domag) idx=code_group_is
3490   IF (is_complex(idx)) THEN
3491      WRITE(stdout,'(5x,"imaginary part")')
3492      DO iclass=1,nclass
3493         WRITE(stdout,'(a5,12f6.2)') name_rap(iclass), &
3494              (AIMAG(char_mat(iclass,irot)),irot=1,nclass)
3495      ENDDO
3496   ENDIF
3497   IF (flag) THEN
3498      WRITE(stdout,'(/5x, "the symmetry operations in each class and &
3499                           &the name of the first element:",/)')
3500      DO irap = 1, nclass
3501         DO iclass=1,nclass
3502            IF (which_irr(iclass)/=irap) CYCLE
3503            WRITE(stdout,'(5x,a5,12i5)') name_class(which_irr(iclass)), &
3504               (elem(i,iclass), i=1,nelem(iclass))
3505!
3506!    The name of the first element of each class is written on output
3507!
3508            WRITE(stdout,'(10x,a)') elem_name(1,iclass)
3509         ENDDO
3510      ENDDO
3511   END IF
3512END IF
3513RETURN
3514END SUBROUTINE write_group_info
3515
3516
3517!---------------------------------------------------------------------------
3518SUBROUTINE find_u( s, u )
3519!---------------------------------------------------------------------------
3520!!  This subroutine receives as input a 3x3 rotation matrix s, and gives
3521!!  as output the matrix u which represents the same rotation in the spin
3522!!  space. Only one of the two u matrices is given. See below for the
3523!!  definition of the sign.
3524!
3525USE kinds,     ONLY : DP
3526USE constants, ONLY : pi
3527!
3528IMPLICIT NONE
3529!
3530REAL(DP) :: s(3,3)
3531!
3532COMPLEX(DP) :: u(2,2)
3533!
3534REAL(DP), PARAMETER :: eps=1.d-8
3535REAL(DP)  :: det, saux(3,3), ax(3), angle, cosa, sina, angle_rot
3536!
3537!  For consistency check uncomment here
3538!
3539!COMPLEX(DP) :: a, as, b, bs
3540!REAL(DP) :: r(3,3), r1(3,3), diff
3541
3542det = s(1,1) * ( s(2,2) * s(3,3) - s(3,2) * s(2,3) )-   &
3543      s(1,2) * ( s(2,1) * s(3,3) - s(3,1) * s(2,3) )+   &
3544      s(1,3) * ( s(2,1) * s(3,2) - s(3,1) * s(2,2) )
3545!
3546!  inversion has no effect in spin space, so improper rotations are
3547!  multiplied by inversion
3548!
3549IF (ABS(det+1.d0) < eps) THEN
3550   saux=-s
3551ELSE
3552   saux=s
3553ENDIF
3554!
3555! Check for identity or inversion
3556!
3557IF ((ABS(saux(1,1)-1.d0) < eps).AND. &
3558    (ABS(saux(2,2)-1.d0) < eps).AND. &
3559    (ABS(saux(3,3)-1.d0) < eps).AND. &
3560    (ABS(saux(1,2)) < eps).AND.(ABS(saux(2,1)) < eps) &
3561.AND.(ABS(saux(2,3)) < eps).AND. &
3562    (ABS(saux(3,2)) < eps).AND.(ABS(saux(1,3)) < eps) &
3563.AND.(ABS(saux(3,1)) < eps)) THEN
3564   u(1,1)=(1.d0,0.d0)
3565   u(1,2)=(0.d0,0.d0)
3566   u(2,1)=(0.d0,0.d0)
3567   u(2,2)=(1.d0,0.d0)
3568   RETURN
3569ENDIF
3570!
3571!   Find the rotation axis and the rotation angle
3572!
3573CALL versor(saux,ax)
3574angle=angle_rot(saux)
3575!write(6,'(5x,"find u",3f12.5,5x,f12.5)') ax(1), ax(2), ax(3), angle
3576angle=0.5d0*angle*pi/180.d0
3577cosa=COS(angle)
3578sina=SIN(angle)
3579!write(6,'(2f12.5)') cosa, sina
3580!
3581!  set the spin space rotation matrix elements
3582!
3583u(1,1)=CMPLX(cosa,-ax(3)*sina,kind=DP)
3584u(1,2)=CMPLX(-ax(2)*sina, -ax(1)*sina,kind=DP)
3585u(2,1)=-CONJG(u(1,2))
3586u(2,2)=CONJG(u(1,1))
3587!
3588!  To each 3x3 rotation one can associate two 2x2 rotation matrices in spin
3589!  space. This function returns the U matrix with positive cosa term
3590!
3591IF (cosa < -eps ) u=-u
3592
3593!IF (ABS(cosa) < eps) THEN
3594!
3595!  Special case when cosa=0. For this rotation we must take the negative sign.
3596!
3597!   IF (ax(1)*ax(3) < -eps) u=-u
3598!ENDIF
3599!
3600!   Here compute the 3x3 rotation matrix starting form the axis, angle
3601!   or from the rotation in spin space for consistency check.
3602!
3603!angle=angle*2.d0
3604!cosa=COS(angle)
3605!sina=SIN(angle)
3606!r1(1,1)=1.d0+(1.d0-cosa)*(ax(1)**2-1)
3607!r1(1,2)=-ax(3)*sina+(1.d0-cosa)*ax(1)*ax(2)
3608!r1(1,3)=ax(2)*sina+(1.d0-cosa)*ax(1)*ax(3)
3609!r1(2,1)=ax(3)*sina+(1.d0-cosa)*ax(1)*ax(2)
3610!r1(2,2)=1.d0+(1.d0-cosa)*(ax(2)**2-1)
3611!r1(2,3)=-ax(1)*sina+(1.d0-cosa)*ax(2)*ax(3)
3612!r1(3,1)=-ax(2)*sina+(1.d0-cosa)*ax(1)*ax(3)
3613!r1(3,2)=ax(1)*sina+(1.d0-cosa)*ax(2)*ax(3)
3614!r1(3,3)=1.d0+(1.d0-cosa)*(ax(3)**2-1)
3615
3616!a=u(1,1)
3617!as=u(2,2)
3618!b=u(1,2)
3619!bs=-u(2,1)
3620
3621!r(1,1)=0.5d0*(a**2+as**2-b**2-bs**2)
3622!r(1,2)=0.5d0*(0.d0,1.d0)*(as**2+bs**2-a**2-b**2)
3623!r(1,3)=-(a*b+as*bs)
3624
3625!r(2,1)=-0.5d0*(0.d0,1.d0)*(as**2-a**2+b**2-bs**2)
3626!r(2,2)=0.5d0*(a**2+b**2+as**2+bs**2)
3627!r(2,3)=(0.d0,1.d0)*(as*bs-a*b)
3628
3629!r(3,1)=(bs*a+as*b)
3630!r(3,2)=(0.d0,1.d0)*(as*b-bs*a)
3631!r(3,3)=(a*as-b*bs)
3632
3633!diff=ABS(r(1,1)-saux(1,1))+ &
3634!     ABS(r(1,2)-saux(1,2))+ &
3635!     ABS(r(1,3)-saux(1,3))+ &
3636!     ABS(r(2,1)-saux(2,1))+ &
3637!     ABS(r(2,2)-saux(2,2))+ &
3638!     ABS(r(2,3)-saux(2,3))+ &
3639!     ABS(r(3,1)-saux(3,1))+ &
3640!     ABS(r(3,2)-saux(3,2))+ &
3641!     ABS(r(3,3)-saux(3,3))
3642
3643
3644!write(6,*) diff
3645!write(6,'(3f15.5)') r1(1,1),r1(1,2),r1(1,3)
3646!write(6,'(3f15.5)') r1(2,1),r1(2,2),r1(2,3)
3647!write(6,'(3f15.5)') r1(3,1),r1(3,2),r1(3,3)
3648!write(6,*)
3649!write(6,'(3f15.5)') r(1,1),r(1,2),r(1,3)
3650!write(6,'(3f15.5)') r(2,1),r(2,2),r(2,3)
3651!write(6,'(3f15.5)') r(3,1),r(3,2),r(3,3)
3652!write(6,*)
3653!write(6,'(4f15.5)') u(1,1),u(1,2)
3654!write(6,'(4f15.5)') u(2,1),u(2,2)
3655!
3656RETURN
3657!
3658END SUBROUTINE find_u
3659
3660
3661!-----------------------------------------------------------------------------
3662FUNCTION set_e( hase, ind )
3663!-----------------------------------------------------------------------------
3664IMPLICIT NONE
3665!
3666INTEGER :: set_e, hase, ind
3667!
3668IF (hase==-1) THEN
3669   set_e=ind+1
3670ELSE
3671   set_e=ind
3672ENDIF
3673!
3674RETURN
3675!
3676END FUNCTION set_e
3677
3678
3679!-----------------------------------------------------------------------------
3680SUBROUTINE check_tgroup( nsym, a, b )
3681!-----------------------------------------------------------------------------
3682!! This subroutine receives a set of 2x2 and 3x3 rotation matrices and
3683!! checks if they are a group.
3684!
3685USE kinds, ONLY : DP
3686!
3687IMPLICIT NONE
3688!
3689INTEGER,     INTENT(IN) :: nsym
3690COMPLEX(DP), INTENT(IN) :: a(2,2,96)
3691REAL(DP),    INTENT(IN) :: b(3,3,nsym)
3692REAL(DP) ::   d(3,3), b1(3,3), b2(3,3), b3(3,3)
3693COMPLEX(dp):: c(2,2), a1(2,2), a2(2,2), a3(2,2)
3694INTEGER :: done
3695LOGICAL :: compare_mat_so
3696
3697INTEGER :: i, j, k
3698
3699DO i=1,nsym
3700   a1(:,:)=a(:,:,i)
3701   b1(:,:)=b(:,:,i)
3702   DO j=1,nsym
3703      a2(:,:)=a(:,:,j)
3704      b2(:,:)=b(:,:,j)
3705      c=MATMUL(a1,a2)
3706      d=MATMUL(b1,b2)
3707      done=0
3708      do k=1,nsym
3709         a3(:,:)=a(:,:,k)
3710         b3(:,:)=b(:,:,k)
3711         IF (compare_mat_so(d,c,b3,a3)) THEN
3712            done=done+1
3713         ENDIF
3714      ENDDO
3715      IF (done.ne.1) write(6,*) 'problem, i,j',i,j
3716   END DO
3717END DO
3718RETURN
3719END SUBROUTINE check_tgroup
3720
3721
3722SUBROUTINE set_class_el_name_so( nsym, sname, has_e, nclass, &
3723                                     nelem_so, elem_so, elem_name_so )
3724!
3725IMPLICIT NONE
3726!
3727INTEGER :: nsym
3728CHARACTER(LEN=45) :: sname(nsym)
3729CHARACTER(LEN=55) :: elem_name_so(12,24)
3730INTEGER :: nclass, nelem_so(24), elem_so(12,24), has_e(12,24)
3731!
3732INTEGER :: iclass, ielem
3733!
3734DO iclass=1,nclass
3735   DO ielem=1,nelem_so(iclass)
3736      elem_name_so(ielem,iclass)=sname(elem_so(ielem,iclass))
3737      IF (has_e(ielem,iclass)==-1) elem_name_so(ielem,iclass)=&
3738                          TRIM(elem_name_so(ielem,iclass)) // ' E'
3739   ENDDO
3740ENDDO
3741!
3742RETURN
3743!
3744END SUBROUTINE set_class_el_name_so
3745