SgInfo Reference Section


How to #include "sginfo.h"

There are two ways of including sginfo.h. The first way is to simply have #include "sginfo.h", the second is to precede this statement by #define SGCOREDEF__.

Without #define SGCOREDEF__ various declarations, globals, and #define's are not "visible". Please check sginfo.h if you need to know the details.


Handling errors: SgError

sginfo.h defines "const char *SgError;". Whenever a call to one of the SgInfo library routines generates an error condiditon, SgError is set via SetSgError(). An application can now take appropriate actions, e.g. print the message in a pop-up window.

Important: if an error occurred and SgError was set, it has to be reset (SgError = NULL;) prior to more calls to the library!

If SGCOREDEF__ is defined before including sginfo.h, the static character array SgErrorBuffer is visible, too. This buffer is intended to hold "dynamic" error messages and could be useful when expanding the library.


Rotation & Translation Base Factors

Except for two routines in sgio.c (ParseSymXYZ() and FormatFraction()) floating point arithmetic is not used in the SgInfo library routines. By choosing proper Rotation & Translation Base Factors all computations can be done with integer arithmetic. Base factors for two types of matrices are defined in sginfo.h.

For Seitz matrices, always labeled SeitzMx or sometimes SMx, the Seitz Matrix Rotation Base Factor is always 1, there is no explicit definition.
The Seitz Matrix Translation Base Factor STBF is defined to be 12 in sginfo.h.

For Change-of-Basis matrices, always labeled CBMx or InvCBMx, the Change-of-Basis Matrix Rotation Base Factor CRBF is defined to be 12.
The Change-of-Basis Matrix Translation Base Factor CTBF is defined to be 72.


SgInfo "legal" Seitz matrices

A "legal" SgInfo Seitz matrix has a 3x3 rotation part for which GetRotMxInfo() returns not 0, and 3 integer translation components in the range [0...(STBF-1)].

Another possible definition for the 3x3 rotation part: all rotation operations which appear in ITVA Table 7. (The 230 Space groups) are legal SgInfo rotation matrices.

Rotation matrices which are the result of a transformation of a centred space group to a primitive setting are in general not legal SgInfo rotation matrices.

E.g, when transforming space group "I 4" to a primitive setting, the resulting rotation matrices are not accepted by GetRotMxInfo(). However, it is of course possible to transform the space group after its generation. But then it is no longer allowed to use any of the SgInfo library routines.

Remark: the set of legal rotation matrices is derived by applying matrix inversion and appropriate basis transformations to the set of 12 "prototype" matrices defined in TabXtalRotMx in sginfo.h.


SgInfo structures


T_SgInfo

typedef struct
  {
    int                  GenOption;
    int                  Centric;
    int                  InversionOffOrigin;
    const T_LatticeInfo  *LatticeInfo;
    int                  StatusLatticeTr;
    int                  OriginShift[3];
    int                  nList;
    int                  MaxList;
    T_RTMx               *ListSeitzMx;
    T_RotMxInfo          *ListRotMxInfo;
    int                  OrderL;
    int                  OrderP;
    int                  XtalSystem;
    int                  UniqueRefAxis;
    int                  UniqueDirCode;
    int                  ExtraInfo;
    int                  PointGroup;
    int                  nGenerator;
    int                   Generator_iList[4];
    char                 HallSymbol[MaxLenHallSymbol + 1];
    const T_TabSgName    *TabSgName;
    const int            *CCMx_LP;
    int                  n_si_Vector;
    int                  si_Vector[9];
    int                  si_Modulus[3];
  }
  T_SgInfo;
GenOption
takes three possible values:
 0 = full group generation
 1 = trusted: set Centric/InversionOffOrigin/LatticeInfo only
-1 = no group generation
GenOption should always be set to 0. The alternative values are only for SgInfo internal usage (see the description of Add2ListSeitzMx()).
A call to InitSgInfo() will set GenOption to 0.
Centric
takes three possible values:
 0 = acentric
 1 = inversion operation in ListSeitzMx
-1 = inversion operation removed from ListSeitzMx
After a call to CompleteSgInfo() Centric will be either 0 or -1.
InversionOffOrigin
takes two possible values:
 0 = ListSeitzMx does not ...
 1 = ListSeitzMx does contain an inversion operation off origin
LatticeInfo
is a pointer to a constant T_LatticeInfo structure.
StatusLatticeTr
takes three possible values:
 0 = translation matrices removed from ListSeitzMx
 1 = all translation matrices in ListSeitzMx
-1 = some translation matrices could be missing in ListSeitzMx
After a call to CompleteSgInfo() StatusLatticeTr will always be 0.
OriginShift
stores the Hall symbol origin-shift translation vector.
nList
holds the number of valid matrices stored in ListSeitzMx (and corresponding entries in ListRotMxInfo).
MaxList
holds the maximum number of matrices which can be stored in ListSeitzMx (and entries which can be stored in ListRotMxInfo). MaxList should always be set to 192, the largest possible number of symmetry operations.
ListSeitzMx
is a pointer to an array of T_RTMx unions. These unions hold the rotation parts of the Seitz matrices (with factor 1) and the translation parts, multiplied by the Seitz-matrix-Translation-Base-Factor STBF.
ListRotMxInfo
is a pointer to an array of T_RotMxInfo structures. For each entry in ListSeitzMx, ListRotMxInfo stores the corresponding set of EigenVector, Order, and two codes for the reference axis.
OrderL
holds the order of the space group, i.e. the maximum number of symmetry equivalent positions. For primitive space groups OrderL = OrderP.
Only valid after a call to CompleteSgInfo()!
OrderP
holds the order of the primitive subgroup of the space group.
Only valid after a call to CompleteSgInfo()!
XtalSystem
is one of the constants:
XS_Unknown
XS_Triclinic
XS_Monoclinic
XS_Orthorhombic
XS_Tetragonal
XS_Trigonal
XS_Hexagonal
defined in sginfo.h.
Only valid after a call to CompleteSgInfo()!
The string, e.g., "Orthorhombic" can be obtained by XS_Name[XS_Orthorhombic]. XS_Name is also defined in sginfo.h. Printing the crystal system is therefore done like:
printf("This space group belongs to the %s crystal system\n",
XS_Name[SgInfo.XtalSystem]);
UniqueRefAxis
holds a code for the unique axis (see T_RotMxInfo). For triclinic and orthorhombic space groups UniqueRefAxis is 0.
Only valid after a call to CompleteSgInfo()!
UniqueDirCode
holds a code for the unique axis (see T_RotMxInfo). For triclinic and orthorhombic space groups UniqueDirCode is 0.
Only valid after a call to CompleteSgInfo()!
ExtraInfo
is one of the constants:
EI_Unknown
EI_Enantiomorphic
EI_Obverse
EI_Reverse
defined in sginfo.h.
Only valid after a call to CompleteSgInfo()!
The string, e.g., "Obverse" can be obtained by EI_Name[EI_Obverse]. EI_Name is also defined in sginfo.h. Printing ExtraInfo is therefore done like:
printf("This space group is %s\n",
EI_Name[SgInfo.ExtraInfo]);
PointGroup
is one of the constants:
PG_Unknown
PG_1
PG_1b
PG_2
PG_m
PG_2_m
PG_222
PG_mm2
PG_mmm
PG_4
PG_4b
PG_4_m
PG_422
PG_4mm
PG_4b2m
PG_4bm2
PG_4_mmm
PG_3
PG_3b
PG_321
PG_312
PG_32
PG_3m1
PG_31m
PG_3m
PG_3bm1
PG_3b1m
PG_3bm
PG_6
PG_6b
PG_6_m
PG_622
PG_6mm
PG_6bm2
PG_6b2m
PG_6_mmm
PG_23
PG_m3b
PG_432
PG_4b3m
PG_m3bm
defined in sginfo.h.
Only valid after a call to CompleteSgInfo()!
To evaluate PointGroup, sginfo.h defines several macros and arrays:
PG_Index(PG_Code)
returns an index which can be used to access PG_Names and LG_Code_of_PG_Index.
PG_Names[PG_Index(SgInfo.PointGroup)] is a pointer to a string, e.g. "m-3m" if SgInfo.PointGroup = PG_m3bm.
LG_Code_of_PG_Index[PG_Index(SgInfo.PointGroup)] is the point group code of the Laue group of SgInfo.PointGroup, e.g. PG_mmm if SgInfo.PointGroup is PG_mm2 or PG_222.
PG_Number(PG_Code)
returns the point group number (range = 1...32) of PG_Code, in the order of Table 10.6.1 of ITVA 1983.
LG_Number(PG_Code)
returns the Laue group number (range = 1...11) of PG_Code, in the order of Table 10.6.1 of ITVA 1983.
nGenerator
holds the number of generator indices stored in Generator_iList.
Only valid after a call to CompleteSgInfo()!
Generator_iList
holds the indices of the Seitz matrices in ListSeitzMx which have been selected for the construction of HallSymbol.
Only valid after a call to CompleteSgInfo()!
HallSymbol
stores the Hall symbol built by SgInfo.
Only valid after a call to CompleteSgInfo()!
TabSgName
is a pointer to a constant T_TabSgName structure.
After building Sginfo.HallSymbol, the internal table of space group symbols is scanned for the same Hall symbol. If there is a match, TabSgName will point to the corresponding entry, otherwise TabSgName will hold the NULL pointer.
Only valid after a call to CompleteSgInfo()!
CCMx_LP
is a pointer to a constant integer array holding the "change-of-centering" (rotation) matrix from the centred to the primitive setting. In contrast to the rotation part of a "CBMx" labeled change-of-basis matrix the factor for the elements of CCMx_LP is 1.
Only valid after a call to CompleteSgInfo()!
n_si_Vector
holds the number (range 0...3) of semi-invariant vectors and moduli stored in si_Vector and si_Modulus.
Only valid after a call to Set_si()!
si_Vector
stores the semi-invariant vectors. The first vector is si_Vector[0..2], the second si_Vector[3..5], and so on.
Only valid after a call to Set_si()!
si_Modulus
stores the semi-invariant moduli.
Only valid after a call to Set_si()!

T_LatticeInfo

typedef struct
  {
    const int  Code;
    const int  nTrVector;
    const int  *TrVector;
  }
  T_LatticeInfo;

T_LatticeInfo is exclusively defined to hold the data of Table 1 in the Hall symbol definition.

Code is one of the letters in column 1 of the table, e.g., 'P' or 'F'. nTrVector is the corresponding entry in column 2. TrVector points to an array of 3 *  nTrVector integers. The first three values are always {0,0,0}, the next three values denote the second lattice translation vector, mulitplied by STBF, and so on. Since the translation parts of the Seitz matrices in SgInfo.ListSeitzMx are also multiplied by STBF, Seitz matrices and lattice translation vectors are compatible.


T_RTMx

typedef union
  {
    struct { int R[9], T[3]; } s;
    int                        a[12];
  }
  T_RTMx;

This union is used for two types of matrices throughout SgInfo, namely Seitz matrices and change-of-basis matrices. For a discussion on the base factors see above.

The matrix elements can be accessed in two ways:

Sometimes it is advantageous to access the matrix elements in an continuous fashion, e.g. to multiply all elements by -1, sometimes it is more convenient to have separate pointers for the rotation and translation parts. The definition of T_RTMx allows for both.

Using the ITVA "augmented 4*4 matrix" notation, the mapping of a position (x,y,z) to the symmetry equivalent position (xs,ys,zs) by a SgInfo Seitz matrix can be summerized by:

        ( xs )     ( R[0]  R[1]  R[2]  T[0]/STBF )  ( x )
        ( ys )  =  ( R[3]  R[4]  R[5]  T[1]/STBF )  ( y )
        ( zs )     ( R[6]  R[7]  R[8]  T[2]/STBF )  ( z )
        ( 1  )     (  0     0     0     1        )  ( 1 )

The definition for SgInfo change-of-basis matrices which transform the coordinates (xa,ya,za) of a position in basis system A to coordinates (xb,yb,zb) in system B, is:

        ( xb )     ( R[0]/CRBF  R[1]/CRBF  R[2]/CRBF  T[0]/CTBF )  ( xa )
        ( yb )  =  ( R[3]/CRBF  R[4]/CRBF  R[5]/CRBF  T[1]/CTBF )  ( ya )
        ( zb )     ( R[6]/CRBF  R[7]/CRBF  R[8]/CRBF  T[2]/CTBF )  ( za )
        ( 1  )     (  0          0          0          1        )  ( 1  )

T_RotMxInfo

typedef struct
  {
    int  EigenVector[3];
    int  Order;
    int  Inverse;
    int  RefAxis;
    int  DirCode;
  }
  T_RotMxInfo;

T_RotMxInfo is a structure which holds the characteristics of an affiliated rotation matrix. (Note: throughout SgInfo there is no explicit affiliation of rotation matrices and T_RotMxInfo structures. This is mainly because all routines of SgInfo are intended to also work with SgInfo.ListRotMxInfo = NULL.)
The only recommended way to set a T_RotMxInfo structure is a call to GetRotMxInfo().

EigenVector
holds the rotation axis of the affiliated rotation matrix, e.g. [0,0,1] for the z-axis.
Order
is one of { 1,2,3,4,6,-1,-2,-3,-4,-6 } (the turn angle of the rotation matrix is 360/abs(Order)).
Inverse
If Inverse is 0, the rotation matrix has a positive turn angle (counter-clock-wise).
If Inverse is 1, the rotation matrix has a negative turn angle (clock-wise).
RefAxis
is one of { o,x,y,z }.
Note: In general RefAxis does not coincide with the direction given by EigenVector.
DirCode
is one of { .,=,",',|,\,* }.

The values of RefAxis and DirCode form a superset of the Hall symbol axis symbols. These symbols are meant to be only for SgInfo internal purposes. (The definition of TabXtalRotMx in sginfo.h and a look at GetRotMxInfo() could help to get an idea of the meanings of RefAxis and DirCode.)


T_TabSgName

typedef struct
  {
    const char  *HallSymbol;
    const int   SgNumber;
    const char  *Extension;
    const char  *SgLabels;
  }
  T_TabSgName;

T_TabSgName is exclusively defined to hold the data of Table 6 in the Hall symbol definition. However, the data are arranged in a different way. A sample entry of TabSgName as defined in sginfo.h is:

     " A 2 2 -1ac",        68, "1cab",  "A_b_a_a",

HallSymbol and SgNumber need no comment. Extension is a separate string (remember, in the lists of space group symbols the extension appears like 68:1cab). SgLabels points to the Herman-Mauguin symbol(s). The symbols for the different directions are separated by underscores. Another TabSgName entry is:

     " B 2",                5, "c2",    "C_2 = B_1_1_2 = B_2",

Here we have three alternatives for the Herman-Mauguin symbol, separated by " = ". The underscores (instead of blanks) help the SgInfo routines to decide what belongs a to a specific Herman-Mauguin symbol and where it ends.

For most purposes there should be no need to access the table entries directly. See the descriptions for PrintTabSgNameEntry() and PrintFullHM_SgName().


T_Eq_hkl

typedef struct
  {
    int  M;      /* Multiplicity */
    int  N;      /* Number of equivalent hkl to follow */
    int  h[24];  /* If hkl == 000 M = N = 1 */
    int  k[24];  /* If hkl != 000 M = 2 * N */
    int  l[24];  /* List of hkl does not contain friedel mates */
    int  TH[24]; /* Phase shift relative to h[0], k[0], l[0] */
  }
  T_Eq_hkl;

For more information see the description of BuildEq_hkl().


SgInfo functions

Functions marked with "@" are less important to the "end user". Unless you don't want to add functions to the library you're not going to miss anything if you skip the marked entries.


void SetSgError(const char *msg);

@

This function is used to set the global SgError variable. SgError will only be set if it is NULL. Successive calls to SetSgError() without clearing SgError (setting SgError to NULL) will preserve the first error message.

Example:

        SetSgError("Internal Error: This should never happen");

int iModPositive(int ix, int iy);

If iy > 0, this function returns ix % iy if ix is positive, and ix % iy + iy if ix is negative. Therefore the result will always be a positive integer in the range [0...(iy-1)].

Example:

        SMx.s.T[i] = iModPositive(f * lsmx->s.T[i] + TrV[i], STBF);

int traceRotMx(const int *RotMx);

@

traceRotMx() returns the sum of the diagonal elements of RotMx, i.e., RotMx[0] + RotMx[4] + RotMx[8].

Example:

        trace = traceRotMx(SeitzMx.s.R);

int deterRotMx(const int *RotMx);

@

deterRotMx() returns the determinant of RotMx.

Example:

        deter = deterRotMx(SeitzMx.s.R);

void RotMx_t_Vector(int *R_t_V, const int *RotMx, const int *Vector,
                    int FacTr);

"Rotation matrix times vector" computes the matrix product:

        ( R_t_V[0] )     ( RotMx[0]  RotMx[1]  RotMx[2] )  ( Vector[0] )
        ( R_t_V[1] )  =  ( RotMx[3]  RotMx[4]  RotMx[5] )  ( Vector[1] )
        ( R_t_V[2] )     ( RotMx[6]  RotMx[7]  RotMx[8] )  ( Vector[2] )

In addition, if FacTr > 0, iModPositive(R_t_V[i], FacTr) will be called for i=0...2.

Caveat: R_t_V and Vector may not point to the same memory area.

Example:

        int  RV[3], RotMx[9], V[3];
        ...
        RotMx_t_Vector(RV, RotMx, V, STBF);

void RotMxMultiply(int *rmxab, const int *rmxa, const int *rmxb);

RotMxMultiply() computes the matrix product rmxab = rmxa * rmxb.

Caveat: rmxab and rmxa or rmxb may not point to the same memory area.

Example:

        int  RAB[9], RA[9], RB[9];
        ...
        RotMxMultiply(RAB, RA, RB);

void RotateRotMx(int *RotMx, const int *RMx, const int *InvRMx);

@

RotateRotMx computes the double matrix product RotMx = RMx * RotMx * InvRMx. This function is used to, e.g., transform a two-fold symmetry operation parallel z to a two-fold symmetry operation parallel x. E.g. GetRotMxInfo() makes use of RotateRotMx() to generate all allowed crystallographic rotation operations by applying appropriate transformation matrices to the rotation matrices tabulated in TabXtalRotMx, defined in sginfo.h.

Example:

        int  RotMx[9], T[9], InvT[9];
        ...
        RotRotMx(RotMx, T, InvT);

void SeitzMxMultiply(T_RTMx *smxab,
                     const T_RTMx *smxa, const T_RTMx *smxb);

SeitzMxMultiply() computes the matrix product smxab = smxa * smxb.

Caveat: smxab and smxa or smxb may not point to the same memory area.

In addition, iModPositive(smxab->s.T[i], STBF) will be called for i=0...2. I.e., the translation components of smxab are guaranteed to be in the range [0...(STBF-1)].

Example:

        T_RTMx  SAB, SA, SB;
        ...
        SeitzMxMultiply(&SAB, &SA, &SB);

void RTMxMultiply(T_RTMx *rtmxab, const T_RTMx *rtmxa, const T_RTMx *rtmxb,
                  int FacAug, int FacTr);

RTMxMultiply() computes the matrix product rtmxab = rtmxa * rtmxb. The translation components of rtmxa are multiplied by FacAug before they are used in the calculation. This allows for easy multiplication of matrices with different rotation and/or translation base factors. (Using the 4x4 augmented matrix notation, the bottom line of rtmxb would read (0 0 0 FacAug). Hence the name "FacAug".)

Caveat: rtmxab and rtmxa or rtmxb may not point to the same memory area.

In addition, if FacTr > 0, iModPositive(rtmxab->s.T[i], FacTr) will be called for i=0...2. I.e., the translation components of rtmxab are then guaranteed to be in the range [0...(FacTr-1)].

Example:

        T_RTMx  RTAB, RTA, RTB;
        ...
        RTMxMultiply(&RTAB, &RTA, &RTB);

void InverseRotMx(const int *RotMx, int *InvRotMx);

InverseRotMx() computes the inverse of the RotMx 3x3 rotation matrix, multiplied by the determinant of RotMx, and stores the result in InvRotMx.
(Speaking more correctly: InverseRotMx() computes the inverse of RotMx, not divided by the determinant of RotMx.)

Caveat: RotMx and InvRotMx may not point to the same memory area.

Example:

        int  R[9], InvR[9];
        ...
        InverseRotMx(R, InvR);

void InverseRTMx(const T_RTMx *RTMx, T_RTMx *InvRTMx);

InverseRTMx() computes the inverse of RTMx, multiplied by the determinant of RTMx->s.R, and stores the result in InvRTMx.
(Speaking more correctly: InverseRTMx() computes the inverse of RTMx, not divided by the determinant of RTMx->s.R.)

Caveat: RTMx and InvRTMx may not point to the same memory area.

Example:

        T_RTMx  RT, InvRT;
        ...
        InverseRTMx(&RT, &InvRT);

int IsSMxTransl0(const T_LatticeInfo *LatticeInfo, const int *SeitzMxT);

@

If SeitzMxT = 0 mod STBF after addition of one of the lattice translation vectors stored in LatticeInfo->TrVector, IsSMxTransl0() returns 1; otherwise IsSMxTransl0() returns 0.

Example:

        const T_LatticeInfo  *LatticeInfo;
        T_RTMx               SeitzMx;
        int                  Tr0;
        ...
        Tr0 = IsSMxTransl0(LatticeInfo, SeitzMx.s.T);

int CompareSeitzMx(const T_LatticeInfo *LatticeInfo,
                   const T_RTMx *SeitzMxA, const T_RTMx *SeitzMxB);
@

If the rotation matrices SeitzMxA->s.R and SeitzMxB->s.R are equal, and if SeitzMxA->s.T - SeitzMxB->s.T = 0 mod STBF after addition of one of the lattice translation vectors stored in LatticeInfo->TrVector, CompareSeitzMx() returns 0; otherwise CompareSeitzMx() returns -1.

Example:

        const T_LatticeInfo  *LatticeInfo;
        T_RTMx               SA, SB;
        int                  CmpResult;
        ...
        CmpResult = CompareSeitzMx(LatticeInfo, &SA, &SB);

int GetRotMxOrder(const int *RotMx);

@

GetRotMxOrder() returns the rotational order of the RotMx 3x3 rotation matrix, i.e. one of { 1,2,3,4,6,-1,-2,-3,-4,-6 } (the turn angle of the rotation matrix is 360/abs(Order)).
If GetRotMxOrder() returns 0, RotMx is not a legal SgInfo rotation matrix, or not a rotation matrix at all.
The return value is obtained by evaluation of the trace and the determinant of RotMx (see e.g. Boisen & Gibbs [2]).

Example:

        int  R[9], Order;
        ...
        Order = GetRotMxOrder(R);

int GetRotMxInfo(const int *RotMx, T_RotMxInfo *RotMxInfo);

@

GetRotMxInfo() sets the RotMxInfo structure (see the description of T_RotMxInfo) for the RotMx 3x3 rotation matrix. GetRotMxInfo() returns the rotational order of RotMx (see GetRotMxOrder()).

Example:

        int          R[9], Order;
        T_RotMxInfo  RMxI;
        ...
        Order = GetRotMxInfo(R, &RMxI);

const T_RotMxInfo *ListOrBufRotMxInfo(const T_SgInfo *SgInfo, int iList,
                                      T_RotMxInfo *BufRotMxInfo);
@

If SgInfo->ListRotMxInfo != NULL, ListOrBufRotMxInfo() returns a pointer to the iList'th element of Sginfo->ListRotMxInfo.
If SgInfo->ListRotMxInfo == NULL, ListOrBufRotMxInfo() calls
GetRotMxInfo(SgInfo->ListSeitzMx[iList].s.R, BufRotMxInfo)
and returns BufRotMxInfo.
If ListOrBufRotMxInfo() returns NULL, the iList'th element of SgInfo->ListSeitzMx is not a legal Sginfo Seitz matrix, and SgError is set.

Example:

        T_SgInfo     SgInfo;
        int          iList;
        T_RotMxInfo  BufRMxI, *RMxI;
        ...
        RMxI = ListOrBufRotMxInfo(&SgInfo, iList, &BufRMxI);

int Add2ListSeitzMx(T_SgInfo *SgInfo, const T_RTMx *NewSMx);

Add2ListSeitzMx() adds the Seitz matrix NewSMx to SgInfo->ListSeitzMx if NewSMx is not already in the list and NewSMx is a legal SgInfo Seitz matrix.
If SgInfo->ListSeitzMx is empty, Add2ListSeitzMx() will insert the identity operation before inserting NewSMx; i.e., Add2ListSeitzMx() guarantees that SgInfo->ListSeitzMx[0] will always be the identity operation.
If SgInfo->ListRotMxInfo != NULL, Add2ListSeitzMx() will also set the appropriate ListRotMxInfo entry.

The behaviour of Add2ListSeitzMx() depends on the value of SgInfo->GenOption:

GenOption ==  0
Full group generation, this is the recommended setting.
Add2ListSeitzMx() will multiply the new Seitz matrices with any Seitz matrix already in the list. If the result of a multiplication is again a new Seitz matrix, Add2ListSeitzMx() will add this matrix, too, and also start the multiplication loop. This process continues until no further Seitz matrices are obtained.
Finally SgInfo->ListSeitzMx will contain the whole set of symmetry operations, including lattice translation operations and the inversion operation at the origin, if present.
Using GenOption == 0 is the only way to ensure, that the Seitz matrices passed to Add2ListSeitzMx() actually form a closed group.
GenOption ==  1
Trusted: set Centric/InversionOffOrigin/LatticeInfo only.
This mode is used to translate the "well-known" Hall symbols of the internal list of space group symbols. Any other use of this mode is discouraged.
Add2ListSeitzMx() will not include lattice translation operations and the inversion operation at the origin into SgInfo->ListSeitzMx, but will only set SgInfo->Centric or SgInfo->InversionOffOrigin and SgInfo->LatticeInfo, thereby reducing the time needed for group generation.
GenOption == -1
No group generation.
This mode is used by the space group finding functions to translate the Hall symbols of the internal list to generator matrices.
Add2ListSeitzMx() will only add NewSMx, but not include lattice translation operations and the inversion operation at the origin. Like for GenOption == 1, only SgInfo->Centric or SgInfo->InversionOffOrigin and SgInfo->LatticeInfo will be set.
Add2ListSeitzMx() returns:
 1, if NewSMx has actually been added to Sginfo->ListSeitzMx.
 0, if NewSMx was already in the list.
-1, if an error occurred. In this case SgError points to the error message.

Example:

        int       Result;
        T_SgInfo  SgInfo;
        T_RTMx    SMx;
        ...
        Result = Add2ListSeitzMx(&SgInfo, &SMx);

int AddInversion2ListSeitzMx(T_SgInfo *SgInfo);

AddInversion2ListSeitzMx() sets up a buffer Seitz matrix with the inversion operation at the origin and adds this matrix to SgInfo->ListSeitzMx via a call to Add2ListSeitzMx().
This is the recommended way to add the inversion operation. Setting Sginfo->Centric directly is discouraged.

AddInversion2ListSeitzMx() returns:
 1, if the inversion matrix has actually been added to Sginfo->ListSeitzMx.
 0, if the inversion matrix was already in the list.
-1, if an error occurred. In this case SgError points to the error message.

Example:

        int       Result;
        T_SgInfo  SgInfo;
        ...
        Result = AddInversion2ListSeitzMx(&SgInfo);

int AddLatticeTr2ListSeitzMx(T_SgInfo *SgInfo,
                             const T_LatticeInfo *LatticeInfo);

For each lattice translation vector of LatticeInfo, AddLatticeTr2ListSeitzMx() sets up a buffer Seitz matrix with the rotation part being the identity matrix, and adds this Seitz matrix to SgInfo->ListSeitzMx via a call to Add2ListSeitzMx().
This is the recommended way to introduce lattice translations. Setting Sginfo->LatticeInfo directly is discouraged.

AddLatticeTr2ListSeitzMx() returns:
 0, if no error occurred.
-1, if an error occurred. In this case SgError points to the error message.

Example:

        int                  Result;
        T_SgInfo             SgInfo;
        const T_LatticeInfo  *LatticeInfo;
        ...
        Result = AddLatticeTr2ListSeitzMx(&SgInfo, LatticeInfo);

int ApplyOriginShift(T_SgInfo *SgInfo);

@

ApplyOriginShift()

Note that - if the origin shift vector is different from (0,0,0) - an inversion operation at the origin becomes an inversion off origin. The other way round, an inversion off origin can become an inversion at the origin. ApplyOriginShift() partly keeps track of this; only CompleteSgInfo() ensures that SgInfo->Centric and SgInfo->InversionOffOrigin are correct.

ApplyOriginShift() returns:
 1, if a non (0,0,0) shift was successfully applied.
 0, if the shift is (0,0,0).
-1, if an error occurred. In this case SgError points to the error message.

Example:

        int       Result;
        T_SgInfo  SgInfo;
        ...
        Result = ApplyOriginShift(&SgInfo);

int FindSeitzMx(const T_SgInfo *SgInfo,
                int Order, int HonorSign, int RefAxis, int DirCode);
@

FindSeitzMx() scans SgInfo->ListSeitzMx for a matrix with rotational order Order if HonorSign is 1, or, if HonorSign is 0, Order or -Order. Furthermore, if one or both of RefAxis and DirCode (see T_RotMxInfo) is not 0, these have to match, too.

FindSeitzMx() returns the index of the first matching Seitz matrix. If there is no match, -1 is returned.

Example:

        int       iList;
        T_SgInfo  SgInfo;
        ...
        iList = FindSeitzMx(&SgInfo, 4, 0, 0, '=');

void InitSgInfo(T_SgInfo *SgInfo);

InitSgInfo() initializes (resets) all members of SgInfo, but not SgInfo->MaxList, SgInfo->ListSeitzMx, and SgInfo->ListRotMxInfo.

Prior to the first call to InitSgInfo(), MaxList has to be set, followed by memory allocation for ListSeitzMx and ListRotMxInfo.
It is also possible to set SgInfo->ListRotMxInfo = NULL;. All functions of SgInfo will still work, though less fast.

The recommended value for MaxList is 192, the maximum number of symmetry operations for all space groups.

Example:

        T_SgInfo  SgInfo;
        ...
        SgInfo.MaxList       = 192;
        SgInfo.ListSeitzMx   = malloc(...);
        SgInfo.ListRotMxInfo = malloc(...);
        InitSgInfo(&SgInfo);

int CompleteSgInfo(T_SgInfo *SgInfo);

CompleteSgInfo() applies the origin shift vector stored in Sginfo->OriginShift, reduces SgInfo->ListSeitzMx to the acentric, primitive subgroup, sorts the list, and updates all members of SgInfo, but not SgInfo->n_si_Vector, SgInfo->si_Vector, and SgInfo->si_Modulus (these are set by calling Set_si()).

CompleteSgInfo() returns 0 if no error occured, and -1 otherwise. In the latter case SgError points to the error message.

Example:

        int       Result;
        T_SgInfo  SgInfo;
        ...
        Result = CompleteSgInfo(&SgInfo);

int CB_SMx(T_RTMx *CSiC,
           const T_RTMx *CBMx, const T_RTMx *SMx, const T_RTMx *InvCBMx);
@

CB_SMx() ("Change-of-Basis for Seitz Matrix") computes the double matrix product:

        CSiC = CBMx * SMx * InvCBMx

For a discussion on basis transformations of symmetry operations see, e.g., Boisen & Gibbs [2].

Note: CSiC and SMX may point to the same memory location.

CB_SMx() returns 0 if no error occured, and -1 otherwise. In the latter case SgError points to the error message.

Example:

        int     Result;
        T_RTMx  SMx, TfSMx;
        R_RTMx  CBMx, InvCBMx;
        ...
        Result = CB_SMx(&TfSMx, &CBMx, &SMx, &InvCBMx);

int TransformSgInfo(const T_SgInfo *SgInfo,
                    const T_RTMx *CBMx, const T_RTMx *InvCBMx,
                    T_SgInfo *BC_SgInfo);

TransformSgInfo() generates the full group for SgInfo, transforms each Seitz matrix by calling CB_SMx(), followed by a call to Add2ListSeitzMx() for BC_SgInfo.
Before calling TransformSgInfo(), InitSgInfo() has to be called for BC_SgInfo. Afterwards CompleteSgInfo() has to be called.

TransformSgInfo() returns 0 if no error occured, and -1 otherwise. In the latter case SgError points to the error message.

Example:

        T_SgInfo  SgInfo, BC_SgInfo;
        T_RTMx    CBMx, InvCBMx;
        ...
        InitSgInfo(&BC_SgInfo);

        if (TransformSgInfo(&SgInfo, &CBMx, &InvCBMx, &BC_SgInfo) == 0)
          CompleteSgInfo(&BC_SgInfo);

const T_TabSgName *FindTabSgNameEntry(const char *UserSgName,
                                      int VolLetter);

FindTabSgNameEntry() is able to find the proper entry in the internal table of space group symbols, given either a space group number, a Schönflies symbol, or a Hermann-Mauguin symbol, passed by UserSgName.
If VolLetter equals 'I', the conventions of ITVI apply (e.g., with UserSgName = "5" setting "B112" is obtained), otherwise the conventions of ITVA apply (e.g., with UserSgName = "5" setting "C121" is obtained).

If there is no proper entry in the internal table, FindTabSgNameEntry() returns NULL.

Note: FindTabSgNameEntry() tries very hard to find a proper entry. In an attempt to cover the most common historical computer notations, several "extras" have been introduced; and, of course, the search is not case-sensitive.
A few examples:

        Fm3m    =>  F m -3 m
        FM3-M   =>  F m -3 m
        Fd3m S  =>  F d -3 m :1
        Fd3m Z  =>  F d -3 m :2
        P2b     =>  P 1 2 1
        P2c     =>  P 1 1 2
        R32R    =>  R 3 2 :R

Example:

        const T_TabSgName  *TSgN;
        char               *UserSgName;
        ...
        TSgN = FindTabSgNameEntry(UserSgName, 'A');

unsigned int SgID_Number(const T_TabSgName *tsgn);

SgID_Number() returns an unique integer number in the range [1...36015] for the entry of the internal table pointed to by tsgn.
Of a five digit number, the right 3 digits give the space group number, while the first and the second digit on the left encode cell choice, unique axis, origin choice, or basis system (rhombohedral or hexagonal), respectively.

Upon error SgID_Number() returns 0.

Example:

        unsigned int       SgID;
        const T_TabSgName  *TSgN;
        ...
        SgID = SgID_Number(TSgN);

int ParseSymXYZ(const char *SymXYZ, T_RTMx *SeitzMx, int FacTr);

ParseSymXYZ() translates (ignoring letter case) the string SymXYZ, e.g. "-y, x-y, z+1/3", to a T_RTMx structure SeitzMx.
Both, fractional and decimal notation for the translation components is possible. Translational components are converted to integer through muliplication by FacTr, followed by rounding. If the deviation introduced by this manipulation is greater than 1% of FacTr, an error condition is generated.

Caveat: Factors for x,y,z, e.g. "2*x" are currently not supported.

ParseSymXYZ() returns 0 if no error occured, and -1 otherwise.

Example:

        int     Result;
        char    *SymXYZ;
        T_RTMx  SeitzMx;
        ...
        Result = ParseSymXYZ(SymXYZ, &SeitzMx, STBF);

int ParseHallSymbol(const char *hsym, T_SgInfo *SgInfo);

ParseHallSymbol() translates (ignoring letter case) the hsym Hall symbol to Seitz matrices and adds them to SgInfo->ListSeitzMx by calling the lower level function Add2ListSeitzMx().

In order to prevent problems with exotic operating systems, several "character-aliases" have been introduced:
_ => white-space
. => white-space
q => '
+ => "

ParseHallSymbol() always returns the index of the position of the last parsed character of hsym. If an error occurred, SgError points the the error message.

Example:

        int       i;
        int       PosHallSymbol;
        T_SgInfo  SgInfo;
        char      *HallSymbol;
        ...
        PosHallSymbol = ParseHallSymbol(HallSymbol, &SgInfo);

        if (SgError != NULL)
        {
          fprintf(stderr, "    %s\n", HallSymbol);
          for (i = 0; i < PosHallSymbol; i++) putc('-', stderr);
          fprintf(stderr, "---^\n");
          fprintf(stderr, "%s\n", SgError);
          exit(1);
        }

int PrintFullHM_SgName(const T_TabSgName *tsgn, int space, FILE *fpout);

PrintFullHM_SgName() prints a Hermann-Mauguin symbol of the entry of the internal table of space group symbols, pointed to by tsgn, on the file pointed to by fpout. With space = 0 output is, e.g., "Fm-3m", with space = '_' output is "F_m_-3_m"; any character is permitted.

PrintFullHM_SgName() returns the number of characters printed. (Caveat: output error diagnostics are currently not supported.)

Example:

        int                nout;
        const T_TabSgName  *TSgN;
        FILE               *fpout;
        ...
        nout = PrintFullHM_SgName(TSgN, ' ', fpout);

void PrintTabSgNameEntry(const T_TabSgName *tsgn, int Style, int space,
                         FILE *fpout);

PrintTabSgNameEntry() prints an entry of the internal table of space group symbols, pointed to by tsgn, on the file pointed to by fpout.
With space = 0 output is, e.g.:

        "225  Oh^5  Fm-3m  -F 4 2 3"
With space = '_' output is:
        "225  Oh^5  F_m_-3_m  -F 4 2 3"
Any character is permitted.
With Style = 0 the four output fields will be separated by exactly two white-space characters. With Style = 1 the fields will be "tabulated" at positions 1,12,26, and 54, respectively.

Caveat: output error diagnostics are currently not supported.

Example:

        const T_TabSgName  *TSgN;
        FILE               *fpout;
        ...
        PrintTabSgNameEntry(TSgN, 0, 0, fpout);

const char *FormatFraction(int nume, int deno, int Decimal,
                           char *Buffer, int SizeBuffer);

FormatFraction() takes the numerator nume and the denominator deno. For Decimal = 0 the numerator and denominator are made relatively prime, e.g. 6/12 becomes 1/2. Then the result is printed to a string (see below). If the simplified denominator is 1, only the numerator will be printed and no "/" appears.
For Decimal = 1 FormatFraction() simply prints the decimal number, e.g. ".5".

If Buffer is not the NULL pointer, SizeBuffer indicates the number of characters allocted for Buffer, output goes to Buffer, and Buffer will also be the pointer returned.
If Buffer is the NULL pointer, output goes to a static internal character array, and a pointer to this array is returned. (Be careful with Buffer = NULL. Any new call to FormatFraction() with Buffer = NULL destroys the previous result, of course!)

FormatFraction() returns NULL if Buffer (or the internal character array) is too small. In addition SgError is set.
However, with SizeBuffer >= 40 this should never happen.

Example:

        const char  *ff;
        T_RTMx      SeitzMx;
        ...
        ff = FormatFraction(SeitzMx.s.T[0], STBF, NULL, 0);

const char *RTMx2XYZ(const T_RTMx *RTMx, int FacRo, int FacTr,
                     int Decimal, int TrFirst, int Low,
                     const char *Seperator,
                     char *BufferXYZ, int SizeBufferXYZ);

RTMx2XYZ() prints the T_RTMx structure RTMx with the rotation base factor FacRo (e.g. 1 or CRBF) and the translation base factor FacTr (e.g. STBF or CTBF) in "XYZ" format (e.g. "-y,x-y,z+1/3" or "1/2*x+1/2*y,z,1/2*x-1/2*y") to a string (see below).

For Decimal = 0 output is, e.g., "x,y+1/2,-z", for Decimal = 1 output is "x,y+.5,-z".

For TrFirst = 0 ("translation first") output is, e.g., "x,y+1/2,-z", for TrFirst = 1 output is "x,1/2+y,-z".

For Low = 0 output is, e.g., "X,Y+1/2,-Z", for Low = 1 output is "x,y+1/2,-z".

For Seperator = NULL output is, e.g., "x,y+1/2,-z", for, e.g., Seperator = " # " output is "x # y+1/2 # -z". Any string is permitted.

If BufferXYZ is not the NULL pointer, SizeBufferXYZ indicates the number of characters allocted for BufferXYZ, output goes to BufferXYZ, and BufferXYZ will also be the pointer returned.
If BufferXYZ is the NULL pointer, output goes to a static internal character array, and a pointer to this array is returned. (Be careful with BufferXYZ = NULL. Any new call to RTMx2XYZ() with BufferXYZ = NULL destroys the previous result, of course!)

RTMx2XYZ() returns NULL if BufferXYZ (or the internal character array) is too small. In addition SgError is set.
However, with SizeBufferXYZ >= 80 this should never happen.

Example:

        const char  *xyz;
        T_RTMx      CBMx;
        ...
        xyz = RTMx2XYZ(&CBMx, CRBF, CTBF, 0, 0, 1, ", ", NULL, 0);

void PrintMapleRTMx(const T_RTMx *RTMx, int FacRo, int FacTr,
                    const char *Label, FILE *fpout);

PrintMapleRTMx() prints the T_RTMx structure RTMx with the rotation base factor FacRo (e.g. 1 or CRBF) and the translation base factor FacTr (e.g. STBF or CTBF) in Maple V matrix format
on the file pointed to by fpout.

For Label = NULL output is, e.g., " := matrix(4,4, [ -1,0,0,0, 0,1,0,1/2, 0,0,-1,0, 0,0,0,1]);",
for, e.g., Label = "m1" output is "m1 := matrix(4,4, [ -1,0,0,0, 0,1,0,1/2, 0,0,-1,0, 0,0,0,1]);".
This is the resulting Maple V "Pretty Print" output:

             [ -1  0   0   0  ]
             [                ]
             [  0  1   0  1/2 ]
       m1 := [                ]
             [  0  0  -1   0  ]
             [                ]
             [  0  0   0   1  ]

Note: the bottom line of the matrix is always "[ 0 0 0 1 ]".

Caveat: output error diagnostics are currently not supported.

Example:

        T_RTMx  CBMx;
        FILE    *fpout;
        ...
        PrintMapleRTMx(&CBMx, CRBF, CTBF, "CBMx", fpout);

void ListSgInfo(const T_SgInfo *SgInfo, int F_XYZ, int F_Verbose,
                FILE *fpout);

ListSgInfo() prints several items of SgInfo on the file pointed to by fpout, e.g.:

        Point Group  3
        Laue  Group  -3
        Trigonal
        Unique Axis  z
        Obverse

        Order     9
        Order P   3

        s.i.Vector  Modulus
          0  0  1   0

If F_XYZ != 0 the entries of Sginfo->ListSeitzMx will be appended in "XYZ" format, e.g.:

        #List     3

        x, y, z
        -y, x-y, z
        -x+y, -x, z

If F_Verbose != 0 the entries of Sginfo->ListSeitzMx and Sginfo->ListRotMxInfo will be appended in the following form:

        #List     3

        (1)    1    [ 0  0  0] 'o' '.'    x, y, z
          1  0  0      0
          0  1  0      0
          0  0  1      0

        (2)    3    [ 0  0  1] 'z' '='    -y, x-y, z
          0 -1  0      0
          1 -1  0      0
          0  0  1      0

        (3)    3^-1 [ 0  0  1] 'z' '='    -x+y, -x, z
         -1  1  0      0
         -1  0  0      0
          0  0  1      0

ListSgInfo() is meant to be for debugging purposes. Furthermore, the source code of ListSgInfo() can serve as an example on how to access and interpret the T_SgInfo structure members.

If an error occurrs SgError points to the error message.

Caveat: output error diagnostics are currently not supported.

Example:

        T_SgInfo  SgInfo;
        FILE      *fpout;
        ...
        ListSgInfo(&SgInfo, 0, 1, fpout);

const T_TabSgName *FindReferenceSpaceGroup(T_SgInfo *SgInfo,
                                           T_RTMx *CBMx, T_RTMx *InvCBMx);

FindReferenceSpaceGroup() scans the internal table of space group symbols and tries to match the Hall generators to symmetry operations of SgInfo - thereby building the change-of-basis matrix CBMx and its inverse InvCBMx - until a correct entry of the internal table has been found.

CBMx transforms coordinates of the "SgInfo-setting" to coordinates of the reference setting.

Upon successful completition, a pointer to an entry of TabSgName is returned. Otherwise NULL is returned and SgError points to the error message.

Example:

        const T_TabSgName  *ReferenceTSgN;
        T_SgInfo           SgInfo;
        T_RTMx             CBMX, InvCBMx;
        ...
        ReferenceTSgN = FindReferenceSpaceGroup(&SgInfo, &CBMx, &InvCBMx);

The result of FindReferenceSpaceGroup() is not necessarily a "standard" setting, but a more or less arbitrary "reference"setting. To obtain the standard setting further action has to be taken. A nice trick to get the standard setting is to use the list of Schönflies symbols in a call to FindTabSgNameEntry():

        #define SGCOREDEF__
        #include "sginfo.h"
        ...
        StandardTSgN = FindTabSgNameEntry(
                         SchoenfliesSymbols[ReferenceTSgN->SgNumber], 'A');

int IsSysAbsent_hkl(const T_SgInfo *SgInfo,
                    int h, int k, int l, int *TH_Restriction);

IsSysAbsent_hkl() tests if the reflection with indices hkl is systematic absent for the space group as described by SgInfo. If so, a code for the symmetry operation which causes the systematic absence is returned. The code is +-(iList + iTrV * SgInfo->nList + 1), where iList is the index of ListSeitzMx, and iTrv is the index of the lattice translation vector as stored in SgInfo->LatticeInfo. If the symmetry operation was multiplied by the inversion operation at the origin (centro-symmetric space groups only), the sign is negative, positive otherwise.
For non-absent, i.e. permitted reflections 0 is returned.
In addition, if TH_Restriction is not the NULL pointer, and the reflection hkl has no phase restriction, *TH_Restriction is set to -1. For values >= 0 the restricted phase angles in degrees are given by P = (*TH_Restriction) * (180 / STBF), and P + 180, respectively.

If an error occurrs SgError points to the error message.

Example:

        int       h, k, l, is_absent, restriction;
        T_SgInfo  SgInfo;
        ...
        is_absent = IsSysAbsent_hkl(&SgInfo, h, k, l, &restriction);

int BuildEq_hkl(const T_SgInfo *SgInfo,
                T_Eq_hkl *Eq_hkl, int h, int k, int l);

BuildEq_hkl() generates the symmetry equivalent indices (for the space group as described by SgInfo) of reflection hkl and stores them in the Eq_hkl structure. Of the Friedel mates only one is included in the list.

Eq_hkl->M
holds the multiplicity of hkl.
Eq_hkl->N
holds the number of symmetry equivalent indices stored.
For hkl == 000 M = N = 1.
For hkl != 000 M = 2 * N.
Eq_hkl->h, Eq_hkl->k, Eq_hkl->l
store the symmetry equivalent indices.
Eq_hkl->TH
stores the list of phase shifts relative to Eq_hkl->h|k|l[0], which are equal to hkl.
The phase shift Eq_hkl->TH[i] of Eq_hkl->h|k|l[i] adds to the phase angle of hkl.
The phase shift in degrees is Eq_hkl->TH[i] * (360 / STBF).

Upon success BuildEq_hkl() returns Eq_hkl->M, otherwise 0 is returned and SgError points to the error message.

Example:

        int       M, h, k, l;
        T_SgInfo  SgInfo;
        T_Eq_hkl  Eq_hkl;
        ...
        M = BuildEq_hkl(&SgInfo, &Eq_hkl, h, k, l);

int AreSymEquivalent_hkl(const T_SgInfo *SgInfo, int h1, int k1, int l1,
                                                 int h2, int k2, int l2);

AreSymEquivalent_hkl() tests whether the reflection with indices h1 k1 l1 is symmetry equivalent to h2 k2 l2 for the space group as described by SgInfo. If so, a code for the symmetry operation which maps the reflections is returned. The code is +-(iList + 1), where iList is the index of ListSeitzMx. If the symmetry operation was multiplied by the inversion operation at the origin (centro-symmetric space groups only), the sign is negative, positive otherwise.
For symmetry independent reflections 0 is returned.

Example:

        int       h1, k1, l1, h2, k2, l2, are_equiv;
        T_SgInfo  SgInfo;
        ...
        are_equiv = AreSymEquivalent_hkl(&SgInfo, h1, k1, l1, h2, k2, l2);

void SetListMin_hkl(const T_SgInfo *SgInfo,            int  Maxk, int  Maxl,
                                            int *Minh, int *Mink, int *Minl);
int IsSuppressed_hkl(const T_SgInfo *SgInfo, int Minh, int Mink, int Minl,
                                                       int Maxk, int Maxl,
                                             int    h, int    k, int    l);

These two functions are designed to help in constructing a hkl list for the space group as described by SgInfo.
Depending on the crystal system and on the positive values of Maxk and Maxl, SetListMin_hkl() sets Minh, Mink, and Minl (Minh is always 0). Inside a loop over Minh|k|l...Maxh|k|l IsSuppressed_hkl() decides, whether a symmetry equivalent reflex of hklhas already appeared before. If so, a code for the corresponding symmetry operation is returned; otherwise 0 is returned.
The code is +-(iList + 1), where iList is the index of ListSeitzMx. If the symmetry operation was multiplied by the inversion operation at the origin (centro-symmetric space groups only), the sign is negative, positive otherwise.

The Simple_hklList() function of sginfo.c gives an example for the usage of SetListMin_hkl() and IsSuppressed_hkl().


int Is_si(const T_SgInfo *SgInfo, int h, int k, int l);

Is_si() returns 1 if the phase angle of the reflex with index hkl is a structure semi-invariant, 0 otherwise.
Set_si() has to be called before the first call to Is_si().

Example:

        int       h, k, l, si;
        T_SgInfo  SgInfo;
        ...
        si = Is_si(&SgInfo, h, k, l);

int Set_si(T_SgInfo *SgInfo);

Set_si() sets the SgInfo structure members n_si_Vector, si_Vector, and si_Modulus.
CompleteSgInfo() has to be called prior to Set_si().

Upon success Set_si() returns 0; otherwise -1 is returned and SgError points to the error message.

Example:

        int       Result;
        T_SgInfo  SgInfo;
        ...
        Result = Set_si(&SgInfo);

Literature: Carmelo Giacovazzo, Direct Methods in Crystallography, Academic Press Inc. (London) 1980.

Caveat: the results of the Set_si() algorithm are extensively tested and (believed to be) correct for any space group and any setting. However, the algorithm itself is anything but sophisticated and consumes much more time and memory than necessary; a better implemantation is desirable. (On the other hand, even this brute-force implementation is very unlikely to ever cause a noteworthy delay.)


void Set_uvw(const T_SgInfo *SgInfo, int h, int k, int l, int *uvw);

For each of the SgInfo->n_si_Vector semi-invariant vectors Set_uvw() computes the scalar product "s.i.vector * hkl modulo s.i.moduls". The n_si_Vector results go to uvw and are intended to be used to establish the primitivity condition (see Giacovazzo).
Set_si() has to be called before the first call to Set_uvw().

Example:

        int       h, k, l, uvw[3];
        T_SgInfo  SgInfo;
        ...
        Set_uvw(&SgInfo, h, k, l, uvw);

SgInfo macros

sginfo.h defines three macros for simple but frequent tasks.


Sg_nLoopInv(SgInfo_)

Sg_nLoopInv returns the constant 2 if SgInfo->Centric is -1, and the constant 1 otherwise.
Purpose: using this macro makes an application independent of the conventions for SgInfo->Centric.

Example:

        int       nLoopInv, iLoopInv;
        T_SgInfo  SgInfo;
        ...
        nLoopInv = Sg_nLoopInv(&SgInfo);
        ...
        for (iLoopInv = 0; iLoopInv < nLoopInv; iLoopInv++)
          ...

See also: LoopSymOps()


InitRotMx(RotMx, diagonal)

InitRotMx sets all off-diagonal elements of the rotation matrix RotMx to 0. The diagonal elements are set to diagonal.

Example:

        int  RotMx[9];
        ...
        InitRotMx(RotMx, 1);

InitSeitzMx(SeitzMx_, diagonal)

InitSeitzMx sets all off-diagonal elements of the Seitz matrix SeitzMx_ to 0. The diagonal elements are set to diagonal.

Example:

        T_RTMx  SeitzMx;
        ...
        InitSeitzMx(&SeitzMx, 1);

Ralf W. Grosse-Kunstleve <rwgk@cci.lbl.gov>