TABLE OF CONTENTS


ABINIT/m_ptgroups [ Modules ]

[ Top ] [ Modules ]

NAME

 m_ptgroups

FUNCTION

  This module contains the irreducible representations and the
  character tables of the 32 point groups.

COPYRIGHT

 Copyright (C) 2010-2024 ABINIT group (MG)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .


m_ptgroups/copy_irrep [ Functions ]

[ Top ] [ m_ptgroups ] [ Functions ]

NAME

  copy_irrep

FUNCTION

  Perform a copy of a set of irrep_t datatypes. Optionally one can multiply
  by a phase factor.

SOURCE

892 subroutine copy_irrep(In_irreps,Out_irreps,phase_fact)
893 
894  implicit none
895 
896 !Arguments ------------------------------------
897  type(irrep_t),intent(in) :: In_irreps(:)
898  type(irrep_t),intent(inout) :: Out_irreps(:)
899  complex(dpc),optional,intent(in) :: phase_fact(:)
900 
901 !Local variables-------------------------------
902 !scalars
903  integer :: irp,dim1,dim2,in_nsym,in_dim,isym
904  character(len=500) :: msg
905 !arrays
906  complex(dpc) :: my_phase_fact(In_irreps(1)%nsym)
907 ! *********************************************************************
908 
909  !@irrep_t
910  dim1 = SIZE( In_irreps)
911  dim2 = SIZE(Out_irreps)
912  if (dim1/=dim2) then
913    msg = " irreps to be copied have different dimension"
914    ABI_ERROR(msg)
915  end if
916 
917  my_phase_fact=cone
918  if (PRESENT(phase_fact)) then
919    my_phase_fact=phase_fact
920    if (SIZE(phase_fact) /= In_irreps(1)%nsym) then
921      msg = " irreps to be copied have different dimension"
922      ABI_ERROR(msg)
923    end if
924  end if
925 
926  do irp=1,dim1
927    in_dim  = In_irreps(irp)%dim
928    in_nsym = In_irreps(irp)%nsym
929    call init_irrep(Out_irreps(irp),in_nsym,in_dim)
930    Out_irreps(irp)%name = In_irreps(irp)%name
931    do isym=1,in_nsym
932      Out_irreps(irp)%mat(:,:,isym) = In_irreps(irp)%mat(:,:,isym) * my_phase_fact(isym)
933      Out_irreps(irp)%trace(isym) = get_trace(Out_irreps(irp)%mat(:,:,isym))
934    end do
935  end do
936 
937 end subroutine copy_irrep

m_ptgroups/get_classes [ Functions ]

[ Top ] [ m_ptgroups ] [ Functions ]

NAME

 get_classes

FUNCTION

  Given a set of nsym 3x3 operations in reciprocal or real space,
  which are supposed to form a group, this routine divides the group into classes.

INPUTS

 nsym=number of symmetry operation.
 sym(3,3,nsym)=The symmetry operations.

OUTPUT

 nclass=The number of classes
 nelements(1:nclass)=For each class, the number of elements
 elements_idx(ii,1:nclass)=For each class, this table gives the index
   of its elements (ii=1,..,nelements(iclass))

NOTES

  * A class is defined as the set of distinct elements obtained by
    considering for each element, S, of the group all its conjugate
    elements X^-1 S X where X range over all the elements of the group.

  * It does not work in case of non-collinear magnetism.

  * The routine assumes that anti-ferromagnetic symmetries (if any) have been removed by the caller.

SOURCE

245 subroutine get_classes(nsym,sym,nclass,nelements,elements_idx)
246 
247  implicit none
248 
249 !Arguments ------------------------------------
250 !scalars
251  integer,intent(in) :: nsym
252  integer,intent(out) :: nclass
253 !arrays
254  integer,intent(in) :: sym(3,3,nsym)
255  integer,intent(out) :: nelements(nsym),elements_idx(nsym,nsym)
256 
257 !Local variables-------------------------------
258 !scalars
259  integer :: isym,jsym,ksym,identity_idx,ierr
260  character(len=500) :: msg
261 !arrays
262  integer :: cjg(3,3),ss(3,3),xx(3,3),xxm1(3,3),test(3,3)
263  integer :: identity(3,3)
264  integer :: dummy_symafm(nsym)
265  logical :: found(nsym),found_identity
266 
267 !************************************************************************
268 
269  ! === Check if identity is present in the first position ===
270  identity=RESHAPE((/1,0,0,0,1,0,0,0,1/),(/3,3/)); found_identity=.FALSE.
271 
272  do isym=1,nsym
273    if (ALL(sym(:,:,isym)==identity)) then
274      found_identity=.TRUE.; identity_idx=isym; EXIT
275    end if
276  end do
277  if (.not.found_identity.or.identity_idx/=1) then
278   write(msg,'(3a)')&
279 &  'Either identity is not present or it is not the first operation ',ch10,&
280 &  'check set of symmetry operations '
281   ABI_ERROR(msg)
282  end if
283   
284  dummy_symafm=1
285  call sg_multable(nsym,dummy_symafm,sym,ierr)
286  ABI_CHECK(ierr==0,"Error in group closure")
287 
288  nclass=0; nelements(:)=0; elements_idx(:,:)=0; found(:)=.FALSE.
289  do isym=1,nsym
290    if (.not.found(isym)) then
291      nclass=nclass+1
292      ss(:,:)=sym(:,:,isym)
293 
294      do jsym=1,nsym ! Form conjugate.
295        xx(:,:)=sym(:,:,jsym)
296        call mati3inv(xx,xxm1) ; xxm1=TRANSPOSE(xxm1)
297        cjg(:,:)=MATMUL(xxm1,MATMUL(ss,xx))
298        do ksym=1,nsym ! Is it already found?
299          test(:,:)=sym(:,:,ksym)
300          if (.not.found(ksym).and.(ALL((test-cjg)==0))) then
301            found(ksym)=.TRUE.
302            nelements(nclass)=nelements(nclass)+1
303            elements_idx(nelements(nclass),nclass)=ksym
304          end if
305        end do
306      end do
307 
308    end if
309  end do
310 
311 end subroutine get_classes

m_ptgroups/get_point_group [ Functions ]

[ Top ] [ m_ptgroups ] [ Functions ]

NAME

  get_point_group

FUNCTION

INPUTS

 ptg_name=point group name as returned by symptgroup.

OUTPUT

 nsym=Number of symmetries in the point group.
 nclass=Number of classes.
 sym(3,3,nsym)=Elements of the point group ordered by classe.
 class_ids(2,nclass)=Initial and final index in sym, for each
 Irreps(nclass)=Datatype gathering data on the different irreducible representations.

SOURCE

114 subroutine get_point_group(ptg_name,nsym,nclass,sym,class_ids,class_names,Irreps)
115 
116  implicit none
117 
118 !Arguments ------------------------------------
119 !scalars
120  integer,intent(out) :: nclass,nsym
121  character(len=*),intent(in) :: ptg_name
122 !arrays
123  integer,allocatable,intent(out) :: sym(:,:,:),class_ids(:,:)
124  character(len=5),allocatable,intent(out) :: class_names(:)
125  type(irrep_t),allocatable,intent(out) :: Irreps(:)
126 
127 !Local variables-------------------------------
128 !scalars
129  integer :: irp,isym
130  !character(len=500) :: msg
131 
132 ! *************************************************************************
133 
134  SELECT CASE (TRIM(ADJUSTL(ptg_name)))
135  CASE ('1')
136    call ptg_C1  (nsym,nclass,sym,class_ids,class_names,Irreps)
137  CASE ('-1')
138    call ptg_Ci  (nsym,nclass,sym,class_ids,class_names,Irreps)
139  CASE ('2')
140    call ptg_C2  (nsym,nclass,sym,class_ids,class_names,Irreps)
141  CASE ('m',"-2") ! Abinit uses "-2"
142    call ptg_Cs  (nsym,nclass,sym,class_ids,class_names,Irreps)
143  CASE ('2/m')
144    call ptg_C2h (nsym,nclass,sym,class_ids,class_names,Irreps)
145  CASE ('222')
146    call ptg_D2  (nsym,nclass,sym,class_ids,class_names,Irreps)
147  CASE ('mm2')
148    call ptg_C2v (nsym,nclass,sym,class_ids,class_names,Irreps)
149  CASE ('mmm')
150    call ptg_D2h (nsym,nclass,sym,class_ids,class_names,Irreps)
151  CASE ('4')
152    call ptg_C4  (nsym,nclass,sym,class_ids,class_names,Irreps)
153  CASE ('-4')
154    call ptg_S4  (nsym,nclass,sym,class_ids,class_names,Irreps)
155  CASE ('4/m')
156    call ptg_C4h (nsym,nclass,sym,class_ids,class_names,Irreps)
157  CASE ('422')
158    call ptg_D4  (nsym,nclass,sym,class_ids,class_names,Irreps)
159  CASE ('4mm')
160    call ptg_C4v (nsym,nclass,sym,class_ids,class_names,Irreps)
161  CASE ('-42m')
162    call ptg_D2d (nsym,nclass,sym,class_ids,class_names,Irreps)
163  CASE ('4/mmm')
164    call ptg_D4h (nsym,nclass,sym,class_ids,class_names,Irreps)
165  CASE ('3')
166    call ptg_C3  (nsym,nclass,sym,class_ids,class_names,Irreps)
167  CASE ('-3')
168    call ptg_C3i (nsym,nclass,sym,class_ids,class_names,Irreps)
169  CASE ('32')
170    call ptg_D3  (nsym,nclass,sym,class_ids,class_names,Irreps)
171  CASE ('3m')
172    call ptg_C3v (nsym,nclass,sym,class_ids,class_names,Irreps)
173  CASE ('-3m')
174    call ptg_D3d (nsym,nclass,sym,class_ids,class_names,Irreps)
175  CASE ('6')
176    call ptg_C6  (nsym,nclass,sym,class_ids,class_names,Irreps)
177  CASE ('-6')
178    call ptg_C3h (nsym,nclass,sym,class_ids,class_names,Irreps)
179  CASE ('6/m')
180    call ptg_C6h (nsym,nclass,sym,class_ids,class_names,Irreps)
181  CASE ('622')
182    call ptg_D6  (nsym,nclass,sym,class_ids,class_names,Irreps)
183  CASE ('6mm')
184    call ptg_C6v (nsym,nclass,sym,class_ids,class_names,Irreps)
185  CASE ('-62m')
186    call ptg_D3h (nsym,nclass,sym,class_ids,class_names,Irreps)
187  CASE ('6/mmm')
188    call ptg_D6h (nsym,nclass,sym,class_ids,class_names,Irreps)
189  CASE ('23')
190    call ptg_T   (nsym,nclass,sym,class_ids,class_names,Irreps)
191  CASE ('m-3')
192    call ptg_Th  (nsym,nclass,sym,class_ids,class_names,Irreps)
193  CASE ('432')
194    call ptg_O   (nsym,nclass,sym,class_ids,class_names,Irreps)
195  CASE ('-43m')
196    call ptg_Td  (nsym,nclass,sym,class_ids,class_names,Irreps)
197  CASE ('m-3m')
198    call ptg_Oh  (nsym,nclass,sym,class_ids,class_names,Irreps)
199  CASE DEFAULT
200    ABI_BUG(sjoin("Unknown value for ptg_name:", ptg_name))
201  END SELECT
202 
203  ! Calculate the trace of each irreducible representation in order to have the character at hand.
204  do irp=1,SIZE(Irreps)
205    ABI_MALLOC(Irreps(irp)%trace, (nsym))
206    do isym=1,nsym
207      Irreps(irp)%trace(isym) = get_trace(Irreps(irp)%mat(:,:,isym))
208    end do
209  end do
210 
211 end subroutine get_point_group

m_ptgroups/groupk_free [ Functions ]

[ Top ] [ m_ptgroups ] [ Functions ]

NAME

  groupk_free

FUNCTION

  Deallocate all memory allocate in the group_k_t.

SOURCE

1061 subroutine groupk_free(Gk)
1062 
1063  implicit none
1064 
1065 !Arguments ------------------------------------
1066  type(group_k_t),intent(inout) :: Gk
1067 
1068 ! *************************************************************************
1069 
1070 ! integer
1071  ABI_SFREE(Gk%class_ids)
1072  ABI_SFREE(Gk%sym)
1073 
1074 !real
1075  ABI_SFREE(Gk%tnons)
1076 
1077 !character
1078  ABI_SFREE(Gk%class_names)
1079 
1080 !type
1081  if (allocated(Gk%Irreps)) then
1082    call irrep_free(Gk%Irreps)
1083    ABI_FREE(Gk%Irreps)
1084  end if
1085 
1086 end subroutine groupk_free

m_ptgroups/groupk_from_file [ Functions ]

[ Top ] [ m_ptgroups ] [ Functions ]

NAME

  groupk_from_file

FUNCTION

  Initialize the group_k_t datatype from an external database retrieved from
  the Bilbao server via the ptg.py script.

INPUTS

  fname(len=*)=file name

OUTPUT

  ierr=Status error
  Lgrps<group_k_t>=The structure completely initialized.

TODO

   This is a stub. I still have to complete the fileformat for the Bilbao database.

SOURCE

672 subroutine groupk_from_file(Lgrps,spgroup,fname,nkpt,klist,ierr)
673 
674  implicit none
675 
676 !Arguments ------------------------------------
677 !scalars
678  integer,intent(in) :: spgroup
679  integer,intent(out) :: ierr,nkpt
680  character(len=fnlen),intent(in) :: fname
681 !arrays
682  type(group_k_t),target,allocatable :: Lgrps(:)
683  real(dp),pointer :: klist(:,:)
684 
685 !Local variables-------------------------------
686 !scalars
687  integer,parameter :: last_file_version=1
688  integer :: unt,fvers,ik,nsym_ltgk,isym,nirreps_k,irp,icls
689  integer :: irrep_idx,irrep_dim,sym_idx,ita_spgroup,nels,prev,now
690  character(len=IRREPNAME_LEN) :: irrep_name
691  character(len=1) :: basis
692  character(len=500) :: msg
693  type(group_k_t),pointer :: Gk
694 !arrays
695  integer,allocatable :: nelements(:),elements_idx(:,:)
696  real(dp) :: kpt(3)
697  character(len=10),allocatable :: kname(:)
698  type(irrep_t),pointer :: OneIrr
699 
700 ! *************************************************************************
701 
702  ierr=0
703  if (open_file(fname,msg,newunit=unt,form="formatted") /=0) then
704    ABI_ERROR(msg)
705  end if
706 
707  read(unt,*,ERR=10)          ! Skip the header.
708  read(unt,*,ERR=10) fvers    ! File version.
709  if (fvers > last_file_version) then
710    write(msg,"(2(a,i0))")" Found file format= ",fvers," but the latest supported version is: ",last_file_version
711    ABI_ERROR(msg)
712  end if
713  read(unt,*,ERR=10) ita_spgroup
714  read(unt,*,ERR=10) basis
715 
716  if (spgroup/=ita_spgroup) then
717    write(msg,'(a,2i0)')&
718 &   " Input space group does not match with the value reported on file: ",spgroup,ita_spgroup
719    ABI_ERROR(msg)
720  end if
721 
722  if (basis /= "b") then
723    msg=" Wrong value for basis: "//TRIM(basis)
724    ABI_ERROR(msg)
725  end if
726 
727  ! * Read the list of the k-points.
728  read(unt,*,ERR=10) nkpt
729 
730  ABI_MALLOC(Lgrps,(nkpt))
731 
732  ABI_MALLOC(klist,(3,nkpt))
733  ABI_MALLOC(kname,(nkpt))
734  do ik=1,nkpt
735    read(unt,*,ERR=10) klist(:,ik), kname(ik)
736  end do
737 
738  ! * Read tables for each k-point
739  do ik=1,nkpt
740 
741    read(unt,*,ERR=10) kpt
742    read(unt,*,ERR=10) nsym_ltgk
743    Gk  => Lgrps(ik)
744 
745    Gk%spgroup = ita_spgroup
746    Gk%nsym    = nsym_ltgk
747    Gk%point   = kpt
748    ABI_MALLOC(Gk%sym,(3,3,nsym_ltgk))
749    ABI_MALLOC(Gk%tnons,(3,nsym_ltgk))
750 
751    do isym=1,nsym_ltgk ! Read symmetries of the little group.
752      read(unt,*,ERR=10) Gk%sym(:,:,isym)
753      read(unt,*,ERR=10) Gk%tnons(:,isym)
754    end do
755 
756    ABI_MALLOC(nelements,(nsym_ltgk))
757    ABI_MALLOC(elements_idx,(nsym_ltgk,nsym_ltgk))
758 
759    call get_classes(nsym_ltgk,Gk%sym,Gk%nclass,nelements,elements_idx)
760 
761    ! The operations reported on the file are supposed to be packed in classes
762    ! otherwise one should perform a rearrangement of the indices.
763    prev = 0
764    do icls=1,Gk%nclass
765      do isym=1,nelements(icls)
766        now = elements_idx(isym,icls)
767        if ( (now-prev) /= 1 ) then
768          write(msg,"(2(a,i0))")&
769 &          " Symmetries on file are not ordered in classes. icls= ",icls,", isym= ",isym
770          ABI_ERROR(msg)
771        else
772          prev = now
773        end if
774      end do
775    end do
776 
777    ABI_MALLOC(Gk%class_ids,(2,Gk%nclass))
778    do icls=1,Gk%nclass
779      nels = nelements(icls)
780      Gk%class_ids(1,icls) = elements_idx(1,   icls)
781      Gk%class_ids(2,icls) = elements_idx(nels,icls)
782    end do
783 
784    ABI_FREE(nelements)
785    ABI_FREE(elements_idx)
786 
787    ! Read the irreducible representations.
788    read(unt,*,ERR=10) nirreps_k
789    ABI_CHECK(Gk%nclass == nirreps_k,"Gk%nclass /= nirreps_k")
790 
791    !$$ allocate(Gk%class_names(Gk%nclass))
792    ABI_MALLOC(Gk%Irreps,(nirreps_k))
793 
794    do irp=1,nirreps_k
795      OneIrr =>  Gk%Irreps(irp)
796      read(unt,*,ERR=10) irrep_idx, irrep_dim, irrep_name
797      call init_irrep(OneIrr,nsym_ltgk,irrep_dim,irrep_name)
798      do isym=1,nsym_ltgk
799        read(unt,*,ERR=10) sym_idx, OneIrr%mat(:,:,isym)
800        ABI_CHECK(sym_idx==irp,"sym_idx/=irp!")
801        ! Matrix elements on file are in the form (rho, theta) with theta given in degrees.
802        call cmplx_sphcart(OneIrr%mat(:,:,isym),from="Sphere",units="Degrees")
803        OneIrr%trace(isym) = get_trace(OneIrr%mat(:,:,isym))
804      end do
805    end do
806 
807  end do
808 
809  close(unt)
810  RETURN
811  !
812  ! Handle IO-error.
813 10 ierr=1
814  close(unt)
815  RETURN
816 
817 end subroutine groupk_from_file

m_ptgroups/init_irrep [ Functions ]

[ Top ] [ m_ptgroups ] [ Functions ]

NAME

  alloc_irrep

FUNCTION

  Initialize an instance of the irrep_t datatype.

INPUTS

  nsym=The number of symmetries.
  irr_dim=The dimension of the irrep.
  [irr_name]=The name of theirrep. "???" is used if not given

OUTPUT

  Irrep<irrep_t>=

SOURCE

959 subroutine init_irrep(Irrep,nsym,irr_dim,irr_name)
960 
961  implicit none
962 
963 !Arguments ------------------------------------
964 !scalars
965  integer,intent(in) :: nsym
966 !arrays
967  integer,intent(in) :: irr_dim
968  character(len=*),optional,intent(in) :: irr_name
969  type(irrep_t),intent(inout) :: Irrep
970 
971 !Local variables-------------------------------
972  !character(len=500) :: msg
973 ! *********************************************************************
974 
975  !@irrep_t
976  Irrep%dim      = irr_dim
977  Irrep%nsym     = nsym
978  Irrep%name     = "???"
979  if (PRESENT(irr_name)) Irrep%name = irr_name
980 
981  ABI_MALLOC(Irrep%mat,(irr_dim,irr_dim,nsym))
982  Irrep%mat = czero
983 
984  ABI_MALLOC(Irrep%trace,(nsym))
985  Irrep%trace = czero
986 
987 end subroutine init_irrep

m_ptgroups/irrep_free_0d [ Functions ]

[ Top ] [ m_ptgroups ] [ Functions ]

NAME

 irrep_free_0d

FUNCTION

  Deallocate all memory allocated in the irrep_t datatype.

SOURCE

831 subroutine irrep_free_0d(Irrep)
832 
833  implicit none
834 
835 !Arguments ------------------------------------
836  type(irrep_t),intent(inout) :: Irrep
837 
838 ! *********************************************************************
839 
840  !@irrep_t
841  if (allocated(Irrep%trace))  then
842    ABI_FREE(Irrep%trace)
843  end if
844  if (allocated(Irrep%mat))  then
845    ABI_FREE(Irrep%mat)
846  end if
847 
848 end subroutine irrep_free_0d

m_ptgroups/irrep_free_1d [ Functions ]

[ Top ] [ m_ptgroups ] [ Functions ]

NAME

 irrep_free_1d

FUNCTION

  Deallocate all memory allocated in the irrep_t datatype.

SOURCE

862 subroutine irrep_free_1d(Irrep)
863 
864  implicit none
865 
866 !Arguments ------------------------------------
867  type(irrep_t),intent(inout) :: Irrep(:)
868 
869 !Local variables-------------------------------
870  integer :: irp
871 ! *********************************************************************
872 
873  do irp=1,SIZE(Irrep)
874    call irrep_free_0d(Irrep(irp))
875  end do
876 
877 end subroutine irrep_free_1d

m_ptgroups/locate_sym [ Functions ]

[ Top ] [ m_ptgroups ] [ Functions ]

NAME

  locate_sym

FUNCTION

  Given a symmetry operation asym, this routine returns its index in the Ptg%sym
  array and the index of the class it belongs to.

INPUTS

OUTPUT

SOURCE

537 subroutine locate_sym(Ptg,asym,sym_idx,cls_idx,ierr)
538 
539  implicit none
540 
541 !Arguments ------------------------------------
542 !scalars
543  integer,intent(out) :: sym_idx,cls_idx
544  integer,optional,intent(out) :: ierr
545  type(point_group_t),intent(in) :: Ptg
546 !arrays
547  integer,intent(in) :: asym(3,3)
548 
549 !Local variables-------------------------------
550 !scalars
551  integer :: isym,icls
552  character(len=500) :: msg
553 
554 ! *********************************************************************
555 
556  sym_idx = 0
557  do isym=1,Ptg%nsym
558    if (ALL(asym == Ptg%sym(:,:,isym) )) then
559      sym_idx = isym
560      EXIT
561    end if
562  end do
563 
564  cls_idx = 0
565  do icls=1,Ptg%nclass
566    if (sym_idx >= Ptg%class_ids(1,icls) .and. &
567 &      sym_idx <= Ptg%class_ids(2,icls) ) then
568      cls_idx = icls
569      EXIT
570    end if
571  end do
572 
573  if (PRESENT(ierr)) ierr=0
574  if (sym_idx==0 .or. cls_idx==0) then
575    write(msg,'(a,9(i0,1x),3a,i1,a,i1)')&
576 &    " Symmetry: ",asym," not found in point group table ",ch10,&
577 &    " sym_idx= ",sym_idx, " and cls_idx= ",cls_idx
578    if (PRESENT(ierr)) then
579      ierr=1
580      ABI_WARNING(msg)
581    else
582      ABI_ERROR(msg)
583    end if
584  end if
585 
586 end subroutine locate_sym

m_ptgroups/mult_table [ Functions ]

[ Top ] [ m_ptgroups ] [ Functions ]

NAME

 mult_table

FUNCTION

  Given a set of nsym 3x3 operations which are supposed to form a group,
  this routine constructs the multiplication table of the group.

INPUTS

 nsym=number of symmetry operation
 sym(3,3,nsym)=the operations

OUTPUT

  mtab(nsym,nsym)=The index of the product S_i * S_j in the input set sym.

SOURCE

608 subroutine mult_table(nsym,sym,mtab)
609 
610  implicit none
611 
612 !Arguments ------------------------------------
613 !scalars
614  integer,intent(in) :: nsym
615 !arrays
616  integer,intent(in) :: sym(3,3,nsym)
617  integer,intent(out) :: mtab(nsym,nsym)
618 
619 !Local variables-------------------------------
620 !scalars
621  integer :: isym,jsym,ksym
622  !character(len=500) :: msg
623 !arrays
624  integer :: prod_ij(3,3),found(nsym)
625 
626 !************************************************************************
627 
628  do jsym=1,nsym
629    found(:)=0 ! Each symmetry should compare only once in a given (row|col).
630 
631    do isym=1,nsym
632      prod_ij = MATMUL(sym(:,:,isym),sym(:,:,jsym))
633      do ksym=1,nsym
634        if ( ALL(prod_ij == sym(:,:,ksym)) ) then
635          found(ksym)=found(ksym)+1
636          mtab(isym,jsym) = ksym
637        end if
638      end do
639    end do ! jsym
640 
641    if ( ANY(found /= 1)) then
642      write(std_out,*)"found = ",found
643      ABI_ERROR("Input elements do not form a group")
644    end if
645  end do ! isym
646 
647 end subroutine mult_table

m_ptgroups/point_group_free [ Functions ]

[ Top ] [ m_ptgroups ] [ Functions ]

NAME

 point_group_free

FUNCTION

  Deallocate all memory allocated in the point_group_t datatype.

SOURCE

374 subroutine point_group_free(Ptg)
375 
376  implicit none
377 
378 !Arguments ------------------------------------
379  type(point_group_t),intent(inout) :: Ptg
380 
381 !Local variables-------------------------------
382 ! *********************************************************************
383 
384  !@point_group_t
385  if (allocated(Ptg%class_ids))   then
386    ABI_FREE(Ptg%class_ids)
387  end if
388  if (allocated(Ptg%sym)) then
389    ABI_FREE(Ptg%sym)
390  end if
391  if (allocated(Ptg%class_names)) then
392    ABI_FREE(Ptg%class_names)
393  end if
394 
395  if (allocated(Ptg%Irreps)) then
396    call irrep_free(Ptg%Irreps)
397    ABI_FREE(Ptg%Irreps)
398  end if
399 
400 end subroutine point_group_free

m_ptgroups/point_group_init [ Functions ]

[ Top ] [ m_ptgroups ] [ Functions ]

NAME

 point_group_init

FUNCTION

  Creation method for the point_group_t datatype.

INPUTS

  ptg_name=The name of the point group (International conventions).

OUTPUT

  The datatype completely initialized.

SOURCE

421 subroutine point_group_init(Ptg,ptg_name)
422 
423  implicit none
424 
425 !Arguments ------------------------------------
426 !scalars
427  character(len=5),intent(in) :: ptg_name
428  type(point_group_t),intent(inout) :: Ptg
429 ! *********************************************************************
430 
431  !@point_group_t
432  !call wrtout(std_out," Retrieving point group data for: "//TRIM(ptg_name),"COLL")
433 
434  Ptg%gname = ptg_name
435  call get_point_group(Ptg%gname,Ptg%nsym,Ptg%nclass,Ptg%sym,Ptg%class_ids,Ptg%class_names,Ptg%Irreps)
436 
437 end subroutine point_group_init

m_ptgroups/point_group_print [ Functions ]

[ Top ] [ m_ptgroups ] [ Functions ]

NAME

 point_group_print

FUNCTION

INPUTS

OUTPUT

SOURCE

454 subroutine point_group_print(Ptg,header,unit,mode_paral,prtvol)
455 
456  implicit none
457 
458 !Arguments ------------------------------------
459 !scalars
460  integer,optional,intent(in) :: unit,prtvol
461  character(len=4),optional,intent(in) :: mode_paral
462  character(len=*),optional,intent(in) :: header
463  type(point_group_t),target,intent(in) :: Ptg
464 
465 !Local variables-------------------------------
466  integer :: my_unt,my_prtvol,irp,icls,sidx
467  complex(dpc) :: trace
468  character(len=4) :: my_mode
469  character(len=500) :: msg
470  type(irrep_t),pointer :: Row
471 ! *********************************************************************
472 
473  my_unt   =std_out; if (PRESENT(unit      )) my_unt   =unit
474  my_prtvol=0      ; if (PRESENT(prtvol    )) my_prtvol=prtvol
475  my_mode  ='COLL' ; if (PRESENT(mode_paral)) my_mode  =mode_paral
476 
477  msg=' ==== Point Group Table ==== '
478  if (PRESENT(header)) msg=' ==== '//TRIM(ADJUSTL(header))//' ==== '
479  call wrtout(my_unt,msg,my_mode)
480 
481  write(std_out,*)REPEAT("=",80)
482  write(std_out,*)" Point group : ",TRIM(Ptg%gname)," Number of symmetries ",Ptg%nsym," Number of classes    ",Ptg%nclass
483 
484  write(std_out,"(a6)",advance="no")"Class "
485  do icls=1,Ptg%nclass
486    write(std_out,"('|',a10)",advance="no")Ptg%class_names(icls)
487  end do
488  write(std_out,"('|')",advance="no")
489  write(std_out,*)" "
490 
491  write(std_out,"(a6)",advance="no")"Mult  "
492  do icls=1,Ptg%nclass
493    write(std_out,"('|',i10)",advance="no")Ptg%class_ids(2,icls)-Ptg%class_ids(1,icls) + 1
494  end do
495  write(std_out,"('|')",advance="no")
496  write(std_out,*)" "
497 
498  do irp=1,SIZE(Ptg%Irreps)
499    Row =>  Ptg%Irreps(irp)
500    write(std_out,'(a6)',advance="no")TRIM(Row%name)
501 
502    do icls=1,Ptg%nclass
503      sidx = Ptg%class_ids(1,icls)
504      trace = Row%trace(sidx)
505      if (ABS(AIMAG(trace)) > tol6) then
506         write(std_out,"('|',(2f5.2))",advance="no")trace
507       else
508         write(std_out,"('|',(f10.2))",advance="no")REAL(trace)
509       end if
510    end do
511 
512    write(std_out,"('|')",advance="no")
513    write(std_out,*)" "
514  end do
515 
516  write(std_out,*)REPEAT("=",80)
517 
518 end subroutine point_group_print

m_ptgroups/show_character_tables [ Functions ]

[ Top ] [ m_ptgroups ] [ Functions ]

NAME

  show_character_tables

FUNCTION

   Printout of the caracter tables of the 32 point groups.

INPUTS

  [unit]=Unit number of output file. Defaults to std_out

OUTPUT

  Only writing.

SOURCE

331 subroutine show_character_tables(unit)
332 
333  implicit none
334 
335 !Arguments ------------------------------------
336 !scalars
337  integer,optional,intent(in) :: unit
338 
339 !Local variables-------------------------------
340  integer :: igrp,my_unt
341  character(len=5) :: ptg_name
342  type(point_group_t) :: Ptg
343 !arrays
344  !integer,allocatable :: elements_idx(:,:),nelements(:)
345 
346 ! *********************************************************************
347 
348  my_unt = std_out; if (PRESENT(unit)) my_unt=unit
349 
350  do igrp=1,SIZE(ptgroup_names)
351    ptg_name = ptgroup_names(igrp)
352    call point_group_init(Ptg,ptg_name)
353    call point_group_print(Ptg,unit=my_unt)
354    !allocate(nelements(Ptg%nsym),elements_idx(Ptg%nsym,Ptg%nsym))
355    !call get_classes(Ptg%nsym,Ptg%sym,nclass,nelements,elements_idx)
356    !deallocate(nelements,elements_idx)
357    call point_group_free(Ptg)
358  end do
359 
360 end subroutine show_character_tables

m_ptgroups/sum_irreps [ Functions ]

[ Top ] [ m_ptgroups ] [ Functions ]

NAME

  sum_irreps

FUNCTION

INPUTS

OUTPUT

SOURCE

1004 function sum_irreps(Irrep1,Irrep2,ii,jj,kk,ll) result(res)
1005 
1006  implicit none
1007 
1008 !Arguments ------------------------------------
1009 !scalars
1010  integer,intent(in) :: ii,jj,kk,ll
1011 !arrays
1012  type(irrep_t),intent(in) :: Irrep1,Irrep2
1013  complex(dpc) :: res
1014 
1015 !Local variables-------------------------------
1016  integer :: isym,nsym,ierr
1017  !character(len=500) :: msg
1018 ! *********************************************************************
1019 
1020  ierr=0; res = czero
1021 
1022  nsym = Irrep1%nsym
1023  if (nsym /= Irrep2%nsym) then
1024    ABI_WARNING("Irreps have different nsym")
1025    ierr=ierr+1
1026  end if
1027 
1028  if (Irrep1%dim /= Irrep2%dim) then
1029    ABI_WARNING("Irreps have different dimensions")
1030    write(std_out,*)Irrep1%dim,Irrep2%dim
1031    ierr=ierr+1
1032  end if
1033 
1034  if (ii>Irrep2%dim .or. jj>Irrep2%dim .or. &
1035 &    kk>Irrep1%dim .or. ll>Irrep1%dim) then
1036    ABI_WARNING("Wrong indices")
1037    write(std_out,*)ii,Irrep2%dim,jj,Irrep2%dim,kk>Irrep1%dim,ll,Irrep1%dim
1038    ierr=ierr+1
1039  end if
1040 
1041  if (ierr/=0) RETURN
1042 
1043  do isym=1,nsym
1044    res = res + DCONJG(Irrep1%mat(ii,jj,isym)) * Irrep2%mat(kk,ll,isym)
1045  end do
1046 
1047 end function sum_irreps