SYMLIB (CCP4: Library)
NAME
symlib - Software Library for handling symmetry operations
DESCRIPTION
From CCP4 5.0, the core handling of symmetry information
is done by a set of C functions. Separate
documentation describes the structures and functions used, and the
API for C/C++ programs.
For Fortran programs, the original set of subroutines (held in symlib.f)
has been replaced by an interface to
the C library. From the point of view of an application programmer,
this interface should be identical to the original set of subroutines. This
document originates from the original Fortran library, but should be
applicable to the new library.
The available Fortran calls have been arranged as much as possible into groups by
function. There are often several versions of calls apparently performing the same
or very similar tasks, reflecting the policy of never removing existing
functionality in order to maintain compatibility with programs written using older
versions of the library.
In CCP4, the primary symmetry information is held in a library data
file defined by the logical name SYMINFO. In the standard configuration of CCP4 this
is the file syminfo.lib which normally resides in $CLIBD. This file replaces
the previous data file SYMOP (i.e. $CLIBD/symop.lib).
syminfo.lib holds information for all the standard spacegroups
in the International Tables. For each spacegroup, several alternative settings
are included (e.g. "P 1 2 1", "P 1 1 2" (a.k.a. 1003) and
"P 2 1 1" for spacegroup 3).
Format of the symmetry library file
Each setting of a spacegroup is delimited by begin_spacegroup / end_spacegroup
records, and contains the following items:
number = standard spacegroup number
basisop = change of basis operator
symbol ccp4 = CCP4 spacegroup number e.g. 1003
(0 if not a CCP4 group)
symbol Hall = Hall symbol
symbol xHM = extended Hermann Mauguin symbol
symbol old = CCP4 spacegroup name
(blank if not a CCP4 group)
symbol laue = Laue group symbol
symbol patt = Patterson group symbol
symbol pgrp = Point group symbol
hklasu ccp4 = reciprocal space asymmetric unit
(with respect to standard setting)
mapasu ccp4 = CCP4 real space asymmetric unit
(with respect to standard setting)
(negative ranges if not a CCP4 group)
mapasu zero = origin based real space asymmetric unit
(with respect to current setting)
mapasu nonz = non-origin based real space asymmetric unit
(with respect to current setting)
cheshire = Cheshire cell
(with respect to standard setting)
symop = list of primitive symmetry operators
cenop = list of centering operators
For example:
begin_spacegroup
number 3
basisop z,x,y
symbol ccp4 1003
symbol Hall ' P 2y (z,x,y)'
symbol xHM 'P 1 1 2'
symbol old 'P 1 1 2'
symbol laue '-P 2' '2/m'
symbol patt '-P 2' '2/m'
symbol pgrp ' P 2' '2'
hklasu ccp4 'k>=0 and (l>0 or (l=0 and h>=0))'
mapasu ccp4 0<=x<=1/2; 0<=y<1; 0<=z<1
mapasu zero 0<=x<1; 0<=y<=1/2; 0<=z<1
mapasu nonz 0<=x<1; 0<=y<=1/2; 0<=z<1
cheshire 0<=x<=1/2; 0<=y<=0; 0<=z<=1/2
symop x,y,z
symop -x,-y,z
cenop x,y,z
end_spacegroup
A call to MSYMLB3 should be used to retrieve
information from the symmetry library. Note that not all the data items are
compulsory for MSYMLB3, although older versions of the routine (MSYMLB2, MSYMLB,
MSYGET) still need them.
2. HOW SYMMETRY OPERATIONS ARE STORED AND APPLIED
Storage of symmetry operations
In syminfo.lib the symmetry operations in each spacegroup are listed as
strings, for example X,Y,Z or -Y,X-Y,Z+1/3 etc. To be useful within a
program these string representations have to be converted to some
mathematical representation.
Typically a symmetry operation RSym will consist of a rotation
operation R and a translation operation T (basically a vector). These
are applied to a vector x to obtain x':
x' = Rx + T
It is convenient to represent the rotation by a 3*3 matrix:
( R11 R12 R13 )
[R] = ( R21 R22 R23 )
( R31 R32 R33 )
and the translation by a column vector with 3 elements:
( T1 )
[T] = ( T2 )
( T3 ).
CCP4 uses 4x4 arrays to store these symmetry operations as follows:
RSym = ( R11 R12 R13 T1 )
( R21 R22 R23 T2 )
( R31 R32 R33 T3 )
( 0 0 0 1 )
or
RSym = ( [R] | [T] )
( 0 0 0 | 1 )
Essentially this is a 4x4 matrix holding 3x3 transformation matrix in
the "top-left-hand corner", the 3-element column (translation) vector
in the "top-right-hand corner", and then (0 0 0 1) in the bottom row.
The subroutine MSYMLB3 will obtain the set of symmetry
matrices in this representation for a given spacegroup, whilst
SYMFR2 or SYMFR3 will obtain
individual matrices from the string representation mentioned above.
(SYMTR4 will perform the inverse operation, converting
matrices to string representation.)
Application of symmetry operations
1. Real Space Coordinates
Using the convention outlined above, post-multiplying the 4x4 matrix by a column
vector as follows:
RSym . [xf]
[yf]
[zf]
[1 ]
will apply both the symmetry and the translation operations to real space
coordinates with a single matrix multiplication. The CCP4
MODLIB library provides
matrix-vector routines MATVEC4 and
TRANSFRM which can be used to perform this operation.
2. Reciprocal Space
The inverse of real-space symmetry matrices are applied to reflection
indices by pre-multiplying them by a row vector representing hkl,
ie. h' = h R-1
or,
(h' k' l') = (h k l) R-1
Note that only the operations in the appropriate Laue group are applied to
reflection indices, so there are no translational components (i.e. the
vector part of the operation, [T], is zero).
The subroutine INVSYM
will invert a 4x4 matrix stored in this representation, for the purpose
of applying symmetry operations to reflection indices.
3. Axis Vectors
Real space axis vectors are transformed like reciprocal space
vectors, i.e.
(a' b' c') = (a b c) R-1
while reciprocal space axis vectors are transformed like real space
coordinates, i.e.
(a*' b*' c*') = R (a* b* c*)
(See also the REINDEX documentation.)
3. DESCRIPTION OF THE GROUPS OF ROUTINES
The routines have been broken down into groups according to function.
Group 1: Subroutines for deriving and manipulating symmetry operations
This group contains routines for obtaining the lists of symmetry
operators from the library, and converting between the string (eg Y,X,-Z
etc) and matrix representations of symmetry operators.
Group 1a: Deriving symmetry operator matrices
- MSYMLB3
- Gets the symmetry operators in matrix representation for a
specified spacegroup from the symmetry library. The spacegroup is
specified either by the spacegroup number (which could be a
non-standard CCP4 number such as 4005) or a spacegroup name. It will
match to any valid s/g name but always returns the longest
name.
Replaces: MSYMLB2,
MSYMLB,
MSYGET
- SYMFR2, SYMFR3
- Translates a character string containing symmetry operator(s) into matrix representation, stored in a 4*4 matrix/array.
NB: SYMFR2 will translate real space coordinate operations (e.g.
x+z,z,-y), reciprocal space operations (e.g. h,l-h,-k) and
reciprocal- and real-space axis vector operations (e.g. a*+c*,c*,-b* and
a,c-a,-b respectively). SYMFR3 only translates real space coordinate
operations.
- SYMTR4
- Translates symmetry matrices into character strings with the equivalent symmetry operations.
Replaces: SYMTRN, SYMTR3
- INVSYM
- Invert the 4*4 array holding the symmetry matrices, to get the inverse symmetry operation.
Use MSYMLB3 to obtain the set of symmetry operator matrices given the spacegroup name
or number. SYMFR2/3 will generate individual symmetry operator matrices from their string
representation (useful if the operators are a subset of a spacegroup). SYMTR4
performs the opposite action, and generates string representations of individual
symmetry operations from the matrices.
INVSYM will generate the inverse matrix of a real space symmetry operation, to be
applied to reflection indices as described in section 2.
Internal routines:
- DETERM
- Calculate the determinant of 4*4 matrix.
Group 1b: Deriving information from symmetry operators
- PGDEFN
- Obtain/guess pointgroup and primitive set of symmetry operators from analysis of all symmetry operators.
- PGMDF
- Gronigen subroutine: determine the nature of the rotation and screw axes from the symmetry matrices.
- PGNLAU
- Determine the Laue group from pointgroup name.
- PATSGP
- Determine the Patterson spacegroup from true spacegroup.
- HKLRANGE
- Return HKL ranges chosen in PGNLAU
These routines all derive additional information from the symmetry operators or
the spacegroup name. The subroutine HKLRANGE returns the information stored in
the common block which it shares with PGNLAU
Group 2: Subroutines for testing reflection data
a) Centric reflections
- CENTRIC
- Sets up symmetry elements; must be called first.
- CENTR
- Tests if a reflection is centric
Nb: routines CENTRC and CENPHS both appear to be unused.
Call CENTRIC once to set up the symmetry elements in common blocks shared with CENTR.
This defines the set of centric reflections. Then for each reflection, a call to
CENTR will return whether it is centric.
b) Epsilon zones
- EPSLN
- Sets up tables of epsilon zones; must be called first.
- EPSLON
- Returns the epsilon zone of a given reflection, as well as whether the reflection is systematically absent (using a call to SYSAB).
- SYSAB
- Function: determines if a reflection is systematically absent.
Call EPSLN once to generate the epsilon zones (general sets of reflections eg 00l or 0k0) and determine the multiplicity/fold. For each reflection a call to EPSLON returns the zone and if the reflection is systematically absent. SYSAB is not called directly.
HKLEQ - used in SCALA to test if two reflections have equal indices.
Group 3: Subroutines for choosing asymmetric units
Remember that the choice of asymmetric unit is NOT UNIQUE. These routines define the
set of CCP4 asymmetric units. The limits for these definitions are given in
section 3.
a) Subroutines for choosing asymmetric units for reflection data
- ASUSET
- Set up symmetry for ASUPUT and ASUGET; must be called first. Calls PRTRSM.
- ASUPUT
- Put reflection into asymmetric unit defined in ASUSET (reverse operation of ASUGET). Calls INASU.
- ASUGET
- Recover original indices of a reflection in the asymmetric unit defined in ASUSET (reverse operation of ASUPUT).
- ASUPHP
- Change phase for symmetry related reflection.
Call ASUSET first to set up symmetry operations in common blocks shared with the other
routines. For each reflection calls can then be made to ASUPUT (return the unique hkl
indices in the asymmetric unit and symmetry number) or ASUGET (obtain real space indices
given unique hkl's and symmetry number). INASU will determine whether a given reflection lies in the asymmetric unit and ASUPHP will convert the phase.
Internal routines:
- INASU
- Function: test if reflection is in the asymmetric unit defined by ASUSET.
- PRTRSM
- Print reciprocal space symmetry operations.
Both these routines are called from within other routines, although they can also be
called independently. ASUSET must be called before INASU can be used.
b) Subroutines for choosing asymmetric units in real space consistent with FFT expectations; FFT grids etc.
- SETLIM
- Set the appropriate box (=asymmetric unit) for the true spacegroup (ie not the FFT spacegroup). For cubic symmetry spacegroups, this will be more than one asymmetric unit.
- SETLIM_ZERO
- Set the appropriate box (=asymmetric unit) for the true spacegroup (ie not the FFT spacegroup). For cubic symmetry spacegroups, this will be more than one asymmetric unit.
NB: This s/r differs from SETLIM in using asu limits derived from cctbx.
- SETGRD
- Set up a suitable sampling grid for FFT. Calls FNSMP.
Internal routines:
- FACTRZ
- Function: returns TRUE if N has all prime factors < 19.
- FNDSMP
- Find suitable grid sample.
Group 4: Subroutines for testing coordinate data
- CALC_ORIG_PS
- Creates a list of equivalent origins for a given spacegroup
- XSPECIALS
- Finds which coordinates occupy special positions (i.e. have occupancies less than 1.0) from consideration of the symmetry operations.
- KROT
- Function: apply symmetry operation to coordinate triplet and check if the result lies in the asymmetric unit.
Neither of the routines XSPECIALS or KROT appear to be used in supported CCP4
programs.
Group 5: Subroutines for permuting symmetry operators
Three subroutines for permuting symmetry operations. They do not really belong here
in symlib, but are widely used invisibly in FFT routines using symmetry operations to permute axes for easier fast fourier calculations.
- PRMVCI
- Permutes specified column of an integer input matrix using another matrix.
- PRMVCR
- Equivalent to PRMVCI but operates on a real input matrix.
- ROTFIX
- Permutes inverse symmetry operations.
Group 6: Subroutines for generating and accessing a hash table
A set of routines used in SCALA, POSTREF and REBATCH.
- CCP4_HASH_SETUP
-
Places a value in the internal look-up table.
- CCP4_HASH_LOOKUP
-
Access a value stored in the table.
- CCP4_HASH_ZEROIT
-
Initialise contents of the table to zero.
These routines are not directly related to symmetry operations. Hashing
is a method of storing data value pairs in such a way that they can be
be efficiently retrieved later on; the hash table is the resulting data
structure.
Group 7: Miscellaneous subroutines
- SETRSL
-
Routine to calculate set coefficients for calculation of (sin(theta)/lambda)**2,
from cell dimensions and angles.
- STHLSQ
-
Calculate (sin(theta)/lambda)**2 from h,k,l, using coefficients set by a
call to SETRSL.
- STS3R4
-
Calculate (sin(theta)/lambda)**2 from h,k,l, using coefficients set by a
call to SETRSL. Duplicates STHLSQ exactly.
These three routines share the common block RECPLT. SETRSL and STHLSQ are used only
in CAD, whilst STS3R4 does not appear in any supported program.
This is how the routines are used in CAD. A call to SETRSL with the cell dimensions
and angles sets up coefficients in RECPLT, which are then used by the function STHLSQ
to calculate the quantity "(sin(theta)/lambda)**2" for any given set of
h, k, l indices. From Bragg's Law, this quantity is equal to 1/(4*d**2), that is,
one-quarter of the resolution. Within CAD, multiplication by 4 yields the resolution
1/d**2.
- PSTOPH
-
Phase angle conversion routine.
The exact function of this routine is unclear and it does not appear in any
supported program.
Group 8: Multiple Spacegroups
The set of Fortran calls which mimic the original symlib.f assume
you are working within a single spacegroup. All calls access the
same spacegroup data structure, in analogy with the COMMON blocks
of symlib.f For cases where you wish to work with multiple
spacegroups (e.g. in the program REINDEX,
a different set of calls is provided (the names of which generally
start with "CCP4SPG_F_"). These identify the spacegroup of interest
via an index "sindx" (by analogy with the "mindx" of mtzlib).
- CCP4SPG_F_ASUPUT
- Put reflection in asymmetric unit of spacegroup on index sindx.
- CCP4SPG_F_CENTPHASE
- CCP4SPG_F_EPSLON
- CCP4SPG_F_EQUAL_OPS_ORDER
- Compare two sets of symmetry operators to see if they are
in the same order.
- CCP4SPG_F_GET_LAUE
- Return Laue number and name for a spacegroup onto index "sindx".
- CCP4SPG_F_INASU
- Test whether reflection or it's Friedel mate is in the asymmetric unit of the spacegroup
on index "sindx".
- CCP4SPG_F_IS_CENTRIC
- CCP4SPG_F_IS_SYSABS
- CCP4SPG_F_LOAD_BY_NAME
- Loads a spacegroup onto index "sindx". The spacegroup is identified by the spacegroup name.
- CCP4SPG_F_LOAD_BY_OPS
- Loads a spacegroup onto index "sindx". The spacegroup is identified by the set of symmetry matrices.
4. FORTRAN API
The following calls are available to Fortran programs. They are arranged alphabetically.
Subroutine ASUGET(IHKL,JHKL,ISYM)
Get original indices of reflection from asymmetric unit,
i.e. reverse operation of ASUPUT. Symmetry defined by call to ASUSET.
On input:
- IHKL(3)
-
input unique indices hkl
- ISYM
-
symmetry number for output
odd numbers are for I+
even numbers are for I-
real-space symmetry operation number L = (ISYM-1)/2 + 1
On output:
- JHKL(3)
-
output original indices hkl
The real-space symmetry matrices are applied in ASUPUT by
premultiplying them by a row vector hkl, ie (h'k'l') = (hkl)R.
So here we calculate (hkl) = (h'k'l') R**-1
Subroutine ASUPHP(JHKL,LSYM,ISIGN,PHASIN,PHSOUT)
Generate phase of symmetry equivalent JHKL from that of IHKL.
On input:
- JHKL(3)
-
indices hkl generated in ASUPUT
- LSYM
-
symmetry number for generating JHKL ( found by ASUPUT)
- ISIGN
-
= 1 for I+
= -1 for I-
- PHASIN
-
phase for reflection IHKL(3)
On output:
- PHSOUT
-
phase for reflection JHKL(3)
Subroutine ASUPUT(IHKL,JHKL,ISYM)
Put reflection into asymmetric unit defined by call to ASUSET
On input:
- IHKL(3)
-
input indices hkl
On output:
- JHKL(3)
-
output indices hkl
- ISYM
-
symmetry number for output
odd numbers are for I+
even numbers are for I-
real-space symmetry operation number L = (ISYM-1)/2 + 1
The real-space symmetry matrices are applied by premultiplying them
by a row vector hkl, ie (h'k'l') = (hkl)R
Subroutine ASUSET(SPGNAM, NUMSPG, PGNAME, MSYM,
RSYM, MSYMP, MLAUE, LPRINT)
Set up & store symmetry information for later use in ASUPUT or ASUGET.
Should be used with subroutine LRSYMI to get PGNAME and subroutine LRSYMM
(both in mtzlib) to get RSYM and MSYM.
On input:
- SPGNAM
-
space-group name (not used) ( character)
- NUMSGP
-
space-group number (not used)
- PGNAME
-
point-group name ( character)
- MSYM
-
total number of symmetry operations
- RRSYM(4,4,MSYM)
-
symmetry matrices (real-space)
- LPRINT
-
- printing flag. ( logical)
On output:
- PGNAME
-
point-group name ( character)
- MSYMP
-
number of primitive symmetry operations
- MLAUE
-
Laue group number - See PGNLAU for details
Subroutine CALC_ORIG_PS(NAMSPG_CIF,NSYM,RSYM,NORIG,ORIG,LPAXISX,LPAXISY,LPAXISZ)
Creates a list of equivalent origins for the named spacegroup.
ARGUMENTS
- (I) NAMSPG_CIF
-
spacegroup name (character)
- (I) NSYM
-
number of symmetry operations
- (I) RSYM(4,4,NSYM)
-
symmetry ops stored as 4x4 matrices
- (O) NORIG
-
number of origins.
- (O) ORIG(3,i)
-
vector of alternate origin (for example : 0.5,0.0,0.5)
only positive components. include vector: (0,0,0)
- (O) LPAXISX
-
logical; set true if s/grp is polar along x axis
- (O) LPAXISY
-
logical; set true if s/grp is polar along y axis
- (O) LPAXISZ
-
logical; set true if s/grp is polar along z axis
Taken from Alexei Vagin
Integer Function CCP4_HASH_LOOKUP(NSER)
The function CCP4_HASH_LOOKUP returns the value NFIND (which was input when
setting up the function in the subroutine CCP4_HASH_SETUP) for the large
range variable NSER. Uses hashing. (see comments for
CCP4_HASH_SETUP for description of hashing method).
Subroutine CCP4_HASH_SETUP(NSER,NFIND)
This subroutine sets up a value for the function ccp4_hash_lookup.
When ccp4_hash_lookup(nser) is later evaluated it will return nfind
This function will allow the efficient retrieval of an identifier
for a large range variable (such as a crystal number). The values
of the function ccp4_hash_lookup(nser) are stored in the array
it(2, kpri) where kpri is the prime number used to generate the
function.
The array 'it' lives in the common block which is shared by
ccp4_hash_setup and the function ccp4_hash_lookup
NOTES: A hash table is a way of storing information so that it
easily be retrieved without the need for indexing or long searches.
NSER is referred to as the "key", which is "hashed" (computer-
science speak for "messed up") by the hashing function (in this
case MOD(NSER4,KPRI) + 1) to determine where the value pair will
be stored. The function LOOKUP can then search on the same basis
when supplied with the key, to retrieve the pair in (at most) 3
calculations. Note that KPRI (the table size) MUST BE A PRIME in
order for this method to work.
Subroutine CCP4_HASH_ZEROIT()
Initialises elements of array 'it' used in ccp4_hash_setup and
ccp4_hash_lookup to zero.
Subroutine CCP4SPG_F_ASUPUT(sindx,ihkl[3],jhkl[3],isym)
Put reflection in asymmetric unit of spacegroup on index sindx.
Subroutine CCP4SPG_F_CENTPHASE(sindx,ih[3],cenphs)
Subroutine CCP4SPG_F_EPSLON(sindx,ih[3],epsi,isysab)
Integer Function CCP4SPG_F_EQUAL_OPS_ORDER(msym1,rrsym1,msym2,rrsym2)
Compare two sets of symmetry operators to see if they are
in the same order. This is important for the consistent use
of ISYM which encodes the operator position in the list.
Returns 1 if operator lists are equal and in the same order, 0 otherwise.
Subroutine CCP4SPG_F_GET_LAUE(sindx,nlaue,launam)
Return Laue number and name for a spacegroup onto index "sindx".
Integer Function CCP4SPG_F_INASU(sindx,ihkl[3])
Test whether reflection or it's Friedel mate is in the asymmetric unit of the spacegroup
on index "sindx". Return 1 if in asu, -1 if -h -k -l is in asu, 0 otherwise.
Subroutine CCP4SPG_F_IS_CENTRIC(sindx,ih[3],ic)
Subroutine CCP4SPG_F_IS_SYSABS(sindx,in[3],isysab)
Subroutine CCP4SPG_F_LOAD_BY_NAME(sindx,namspg)
Loads a spacegroup onto index "sindx". The spacegroup is identified by the spacegroup name.
Subroutine CCP4SPG_F_LOAD_BY_OPS(sindx,msym,rrsym)
Loads a spacegroup onto index "sindx". The spacegroup is identified by the set of symmetry matrices.
Subroutine CENTR(IH,IC)
Input IH(3) - reflection indices
Returns IC
Determine whether a reflection is centric (return IC=1)
or not (IC=0). If none of the zone tests is satisfied,
the reflection is non-centric.
Logical Function CENTRC(KHKL,ICENT)
Returns value as true if reflection khkl is centric, false otherwise.
It is general for all point groups - but only for the unique set of
indices which conforms to the criterion of maximising the value of
(khkl(3)*256 + khkl(2))*256 + khkl(1)
as produced by e.g. subroutine turnip in protin and ulysses.
In this case the required tests are controlled by 7 flags in
icent for
0KL H0L HK0 HKK HKH HHL H,-2H,L
(the last is needed in pg312)
Subroutine CENTRIC(NSM,RSMT,IPRINT)
This is Randy Read's method of defining centric reflections.
It uses NSM and the symmetry operators stored in RSMT(4,4,NSM)
It decides how many centric zones there are, and flags them.
set up tests for 0kl h0l hk0 hhl hkh hkk h,-hl hk-h hk-k
-h 2h l 2h -h l hkl
Zones are encoded using an idea from a program by Bricogne.
If h*zone(1) + k*zone(2) + l*zone(3) is equal to 0.0,
that reflection is in that zone. All that is needed is the
most general conditions - a reflection is either centric or
not.
Subroutine DETERM(det,a)
Gets determinant of a matrix
- Input A
-
4*4 matrix (real)
- Output DET
-
determinant of A.
Subroutine EPSLN(NSM,NSMP,RSMT,IPRINT)
It works out the epsilon cards using NSM and the symmetry operators stored in
RSMT(4,4,NSM).
This is Randy's program description:
zones defined as for centric zones, but
fourth number on each line is the multiplicity corresponding
to this zone. last card should always be 0 0 0 n, where n is
appropriate for the lattice (1 for primitive, 2 for face-
centred, etc.), so that general reflections get a value
for epsilon. be very careful of the order of the zones. cards
for reciprocal lattice rows should be given before those for
planes, because the first test that is satisfied defines
the zone.
set up tests for
h00 0k0 00l hh0 h0h 0kk h,-h0 h0-h 0k-k -h2h0 2h-h0 hhh hkl
Subroutine EPSLON(IH,EPSI,ISYSAB)
Input IH(3) - reflection indices
Returns EPSI ( epsilon zone) , and ISYSAB flag.
Systematic absences flagged with ISYSAB = 1
Find the zone a reflection falls into, and return the
appropriate value for the reflection multiplicity factor.
Each reflection must have a zone.
Logical Function FACTRZ(N)
Returns true if N has all prime factors <= 19
Subroutine FNDSMP(MINSMP, NMUL, SAMPLE, NSAMPL)
Find suitable grid sample, approximately = SAMPLE/2 * maximum index,
with required factor, and no prime factor > 19
On entry:
- MINSMP
-
minimum sample, approximately 2 * maximum index
- NMUL
-
required factor
- SAMPLE
-
desired sample factor, ie if = 1.0 (minimum), try to get sample close to MINSMP
On exit:
- nsampl
-
grid sample; if MINSMP<=0, nsampl=nmul
Logical Function HKLEQ(IH,KH)
Checks if indices are equal.
Returns true if indices ih = kh
Subroutine HANDCHANGE(lspgrp,cx,cy,cz)
Used in program phistats.
Subroutine HKLRANGE(IHRNG0,IKRNG1,IKRNG0,IKRNG1,ILRNG0,ILRNG1)
Return HKL ranges chosen in PGNLAUE
Arguments: (INTEGER) HRNG0,HRNG1,KRNG0,KRNG1,LRNG0,LRNG1
Integer Function INASU(IH, NLAUE)
Arguments:
- NLAUE
-
- code number for this pointgroup
- IH(3)
-
- indices
Returns:
INASU = +1 if h k l chosen
INASU = -1 if -h-k-l chosen
INASU = 0 if reflection is out-of-bounds
Subroutine INVSYM(S,ST)
Inverts a 4*4 matrix. Used here to get inverse symmetry operation for
generating equivalent h k l, i.e.
[h'] = [h][St]
h'(j) =Sum(I=1,3)[ h(i)*St(I,J,NS)]
Arguments:
- Input S
-
4*4 matrix to be inverted
- Output ST
-
4*4 matrix (inverse of S)
Integer Function KROT(NS)
Apply NS'th symmetry operation to JP to get LP,
check if lies in asymmetric unit given by NAU.
Returns KROT=0 correct operation, =1 if not.
Subroutine MSYGET(IST,LSPGRP,NSYM,ROT)
Get symmetry operations for space-group LSPGRP
from library file, logical name SYMINFO.
Arguments:
- IST
-
dummy parameter for backwards compatibility
- LSPGRP (input)
-
Name of spacegroup
- IST (input)
-
Stream of library file
- NSYM (output)
-
Number of symmetry operations
- ROT(4,4,NSYM) (output)
-
Rotation/translation matrices
Subroutine MSYMLB(IST,LSPGRP,NAMSPG,NAMPG,NSYMP,NSYM,ROT)
Get symmetry operations from library file, logical name SYMINFO.
Space group defined by LSPGRP - spacegroup number or NAMSPG - spacegroup name.
This routine is obsolete now, and simply wraps MSYMLB3.
Arguments:
- IST
-
dummy parameter for backwards compatibility
- LSPGRP (I/O)
-
spacegroup number
- NAMSPG (I/O)
-
spacegroup name
- NAMPG (O)
-
pointgroup name
- NSYMP (O)
-
number of primitive symmetry operations
- NSYM (O)
-
number of symmetry operations
- ROT(4,4,NSYM)
-
rotation/translation matrices
Subroutine MSYMLB2(IST,LSPGRP,NAMSPG_CIF,NAMPG,NSYMP,NSYM,ROT)
Identical to MSYMLB, except that on output NAMSPG_CIF
has correct CIF format, e.g. 'P 21 21 21'. This routine is obsolete now, and
simply wraps MSYMLB3.
Subroutine MSYMLB3(IST,LSPGRP,NAMSPG_CIF,NAMSPG_CIFS,NAMPG,NSYMP,NSYM,RSYM)
Get symmetry operations from library file, logical name SYMINFO.
- The routine first tries to identify the spacegroup by its number.
This requires that the number is UNIQUE, so alternate settings are
numbered n000 + Int Tab number.
- If the spacegroup number is set to 0 on input, the routine will try to match the
assigned NAMSPG_CIF. This uses the C function ccp4spg_load_by_ccp4_spgname
which will cope with a variety of possible formats for the name, e.g.
presence or absence of spaces.
On entry:
- IST
-
stream number to read file
- LSPGRP
-
spacegroup number, can be set to 0
- NAMSPG_CIF
-
any acceptable spacegroup name: this will be used to
identify the spacegroup if possible
Returns:
- LSPGRP
-
spacegroup number
- NAMSPG_CIF
-
full spacegroup name
- NAMSPG_CIFS
-
name without any spaces
- NAMPG
-
pointgroup name
- NSYMP
-
number of primitive symmetry operations,
only different from NSYM in non-primitive spacegroups
- NSYM
-
total number of symmetry operations
- RSYM(4,4,NSYM)
-
Symmetry rotation/translation matrices
Subroutine PATSGP(SPGNAM, PGNAME, PATNAM, LPATSG)
Determine Patterson spacegroup from true space-group
On entry:
- SPGNAM
-
space-group name. Only used to determine lattice centering
- PGNAME
-
point-group name
On exit:
- PATNAM
-
name of Patterson spacegroup
- LPATSG
-
number of Patterson spacegroup
Subroutine PGDEFN(NAMPG,NSYMP,NSYM,RSYMT,LPRINT)
Arguments:
- Input NSYM
-
- number of symmetry operators. ( integer)
- Input RSYMT
-
- 4*4 symmetry matrices. ( real)
- Input LPRINT
-
- printing flag. ( logical)
- Returns NAMPG
-
- name of point group. ( character)
- Returns NSYMP
-
- number of primitive symmetry operators. ( integer)
- Returns RSYMT
-
- possibly reordered.
This subroutine chooses the primitive set of symmetry operators.
If necessary it re-orders the symmetry operators to give the
primitive ones first.
This subroutine works out the point group name NAMPG. That is ; it checks
rotation axes, etc etc and recognises these point groups. (It DOES NOT cope
with mirror planes etc)
Gronigen MDF stuff: It now sets up the common block MDFPAR for
MDF file mods and fills in the symmetry info. See subroutine for
details.
Subroutine PGMDF(JLASS,JCENTR,JSCREW)
Gronigen subroutine.
Use this subroutine to transfer information to and from MDFPAR.
If JLASS eq 0 then fill JLASS JCENTR JSCREW from common block.
If JLASS gt 0 then fill KLASS ICENTR ISCREW in common block.
Subroutine PGNLAU(NAMPG,NLAUE,LAUNAM)
Choose Laue group from PG name.
On entry:
- NAMPG
-
point-group name ( character)
On exit:
- NLAUE
-
Laue group number ( integer)
- LAUNAM
-
Laue group name ( character)
This subroutine returns a laue code number used to choose
the unique region of reciprocal space for each point group.
The number nlaue is the same as the one set in CAD for this purpose.
Pointgroup Laue group Limits
3 pg1 1bar hkl:l>=0 hk0:h>=0 0k0:k>=0 1,2
pg1bar
4 pg2 (b) 2/m hkl:k>=0, l>=0 hk0:h>=0 3/b,4/b....
pgm pg2/m
5 pg2 (c) 2/m hkl:k>=0, l>=0 h0l:h>=0 1003,1004
6 pg222 mmm hkl:h>=0, k>=0, l>=0 16 ...
pgmm2 pgmmm
7 pg4 4/m hkl:h>=0, l>=0 with k>=0 if h=0 and
pg4bar pg4/m k>0 if h>0
8 pg422 4/mmm hkl:h>=0, k>=0, l>=0 89..
pg4mm pg4bar2m pg4barm2 pg4/mmm
9 pg3 3bar hkl:h>=0, k>0 00l:l>0 143..
pg3bar
10 pg312 3/m hkl:h>=0, k>=0 with k<=h for all l.
pg32 pg3m pg3m1 pg3barm1 if k = 0 l>=0
Space group numbers : 149-151-153 157 159 162 163
11 pg321 3bar1m hkl:h>=0, k>=0 with k<=h for all l.
pg31m pg3bar1m if h = k l>=0
Space group numbers : 150-152-154
12 pg6 6/m hkl:h>=0, k>=0, l>=0 with k>=0 if h=0
pg6bar 6/m and k> 0 if h>0
13 pg622 6/mmm hkl:h>=0, k>=0, l>=0 with h>=k 177..
pg6mm pg6barm2 pg6bar2m pg 6/mmm
14 pg23 m3 hkl:h>=0, k>=0, l>=0 with l>=h, k>=h
pgm3bar
15 pg432 m3m hkl:h>=0, k>=0, l>=0 with k>=l
pg4bar3m pgm3barm
Subroutine PRMVCI(PERM,JV,N,N1)
- Input PERM
-
- 4*4 matrix (real)
- Input JV
-
- N1*3 matrix (integer)
- Output JV
-
- N1*3 matrix (integer)
This has been modified by permuting the Nth column by matrix PERM.
Here is the code for clarity:
C---- Permute
C
C DO 10 I = 1,3
C BV(I) = PERM(I,1)*JV(N,1) + PERM(I,2)*JV(N,2) +
C + PERM(I,3)*JV(N,3)
C 10 CONTINUE
C
C---- Copy back
C
C DO 20 I = 1,3
C JV(N,I) = NINT(BV(I))
C 20 CONTINUE
Subroutine PRMVCR(PERM,AV,N,N1)
- Input PERM
-
- 4*4 matrix (real)
- Input AV
-
- N1*3 matrix (real)
- Output AV
-
- N1*3 matrix (real)
This has been modified by permuting the Nth column by matrix PERM.
See PRMVCI - this routine is its real equivalent.
Subroutine PRTRSM(PGNAME, NSYMP, RSYMIV)
Print reciprocal space symmetry operations
The real-space symmetry matrices are applied by premultiplying them
by a row vector hkl, ie (h'k'l') = (hkl)R
Subroutine PSTOPH (PSIX,PSIY,PSIZ,PHIX,PHIY,PHIZ,AVPHI)
Convert PSIX,PSIY,PSIZ (= epsx,epsy,epsz) to PHIX,PHIY,PHIZ, using AVPHI
All angles in radians
Subroutine ROTFIX
Permutes inverse symmetry operations
Matrices passed in Common block ATSYM
Subroutine SETGRD(NLAUE,SAMPLE,NXMIN,NYMIN,NZMIN,NX,NY,NZ)
Set up a suitable sampling grid for FFT
Input:
- NLAUE
-
Laue-group for FFT/SF calculation
- SAMPLE
-
default fineness of sample, ie if = 1.0 (minimum),
try to get sampling as close to minimum as possible.
Typically = 1.5 to get sample at traditional 3 * maximum index
- NXMIN NYMIN NZMIN
-
minimum sampling (true XYZ)
Output:
- NX,NY,NZ
-
sampling intervals along X,Y,Z
The sampling intervals must satisfy the following conditions:
- approximately SAMPLE * minimum sampling
- no prime factor > 19
- special restrictions for particular space-groups
This is ALL the point groups.
PG1 PG1bar PG2 PGm PG2/m PG222 PGmm2 PGmmm
PG4 PG4bar PG4/m PG422 PG4mm PG4bar2m PG4/mmm
PG3 PG3bar PG32 PG3m PG3barm
PG6 PG6bar PG6/m PG622 PG6mm PG6bar2m PG6/mmm
PG23 PGm/3bar PG432 PG4bar3m PGm3bar m
We use:
PG1 PG1bar PG2 PG2/m PG222 PGmmm
PG4 PG4/m PG422 PG4/mmm
PG3 PG3bar PG32 PG3bar/m
PG6 PG6/m PG622 PG6/mmm
PG23 PGm/3bar PG432 PGm3barm
For grid restrictions we only need to know the laue number.
Here is the table:
3 pg1 1bar hkl:l>=0 hk0:h>=0 0k0:k>=0 1,2
4 pg2 2/m hkl:k>=0, l>=0 hk0:h>=0 3/b,4/b....
5 pg2(c) 2/m hkl:k>=0, l>=0 h0l:h>=0 1003,1004
6 pg222 mmm hkl:h>=0, k>=0, l>=0 16 ...
7 pg4 4/m hkl:h>=0, l>=0 with k>=0 if h=0 and
8 pg422 4/mmm hkl:h>=0, k>=0, l>=0 89..
9 pg3 3bar hkl:h>=0, k>0 00l:l>0 143..
10 pg312 3/m hkl:h>=0, k>=0 with k<=h for all l.
if k = 0 l>=0
Space group numbers : 149-151-153
11 pg321 3/m hkl:h>=0, k>=0 with k<=h for all l.
if h = k l>=0
Space group numbers : 150-152-154
12 pg6 6/m hkl:h>=0, k>=0, l>=0 with k=0 if h=0
13 pg622 6/mmm
14 pg23 m3
15 pg432 m3m
Subroutine SETLIM(LSPGRP,XYZLIM)
Set appropriate box (asymmetric unit) for spacegroup (true spacegroup rather than
FFT spacegroup) LSPGRP. For cubic symmetry spacegroups, this will be more than
one asymmetric unit.
On entry:
- lspgrp
-
true spacegroup (not FFT spacegroup)
On exit:
- xyzlim(2,3)
-
minimum, maximum limits on x,y,z (fractions of cell); if spacegroup not
recognized, returns xzylim(1,1) = -1.0
Note that the minimum limits (xyzlim(1,)) will always = 0.0
Subroutine SETLIM_ZERO(LSPGRP,XYZLIM)
Set appropriate box (asymmetric unit) for spacegroup (true spacegroup rather than
FFT spacegroup) LSPGRP. For cubic symmetry spacegroups, this will be more than
one asymmetric unit.
NB This s/r differs from SETLIM in that the limits are taken from cctbx
via CCP4's syminfo.lib file.
On entry:
- lspgrp
-
true spacegroup (not FFT spacegroup)
On exit:
- xyzlim(2,3)
-
minimum, maximum limits on x,y,z (fractions of cell); if spacegroup not
recognized, returns xzylim(1,1) = -1.0
Note that the minimum limits (xyzlim(1,)) will always = 0.0
Subroutine SETRSL(A,B,C,ALPHA,BETA,GAMMA)
Routine to calculate coefficients for (sin(theta)/lambda)**2 from
h,k,l for general axes.
First calculates the components of input axes in an orthonormal
basis, then calculate components of reciprocal axes in same basis.
The input angles are in degrees.
Real Function STHLSQ(IH,IK,IL)
Calculate (sin(theta)/lambda)**2 from h,k,l. The coefficients are set by
a previous call to SETRSL. Works for any kind of axes.
Real Function STS3R4(IH,IK,IL)
Calculate (sin(theta)/lambda)**2 from h,k,l. The coefficients are set by a
call to SETRSL. Works for any kind of axes.
Subroutine SYMFR2(ICOL,I1,NS,ROT)
Read and interpret symmetry operations
SYMFR2 recognises the following types of input:
real space symmetry operations, e.g. X+1/2,Y-X,Z
reciprocal space operations, e.g. h,l-h,-k
reciprocal axis vectors, e.g. a*+c*,c*,-b*
real space axis vectors, e.g. a,c-a,-b
The subroutine returns the appropriate 4x4 transformation
matrix for each operation. The calling program must
interpret the resulting matrix(ces) correctly.
On entry I1 is the first character of ICOL to look at (say after
keyword 'SYMM')
NS is the number of the first symmetry operation to be read, & returns
with the number of the last one read.
On exit, ROT(4,4,NS) contains the real-space symmetry matrices, in standard
convention, i.e.
[x'] = [s][x]
x'(I)=Sum(J=1,3)ROT(I,J,NS)*x(J) + ROT(I,4,NS)
Input:
- ICOL
-
character string containing symmetry operations
- I1
-
first character in ICOL to interpret from.
Output:
- ROT(I,4,NS)
-
contains the fractional translations.
Subroutine SYMFR3(ICOL,I1,NS,ROT,EFLAG)
Read and interpret symmetry operations.
Arguments:
- ICOL (I) CHARACTER*80
-
Line containing the symmetry operations
- I1 (I) INTEGER
-
First character to look at in ICOL (say after keyword 'SYM')
- NS (I/O) INTEGER
-
is the number of the first symmetry operation to be read, & returns with the
number of the last one read (ie you can have more than one on a line!)
- ROT (O) REAL
-
Array (4,4,at_least_NS), on exit contains the real-space
symmetry matrices, in standard convention, i.e.
[x'] = [s][x]
x'(I)=Sum(J=1,3)ROT(I,J,NS)*x(J) + ROT(I,4,NS)
ROT(I,4,NS) contains the fractional translations
- EFLAG (O) INTEGER
-
Error flag - on exit, if 0 then OK, > 0, an error occurred.
Subroutine SYMTRN(NSM,RSM)
Symmetry translation from matrix back to characters
This translates the symmetry matrices RSM(4,4,NSM) into INT TAB
character strings
It gives the real and reciprocal space operations.
eg X,Y,Z H , K, L
eg -Y,X-Y, Z -H-K, H, L etc
That is more complicated than you might think!!
Subroutine SYMTR3(NSM,RSM)
Symmetry translation from matrix back to characters
This translates the symmetry matrices RSM(4,4,NSM) into INT TAB
character strings
It gives the real and reciprocal space operations.
eg X,Y,Z H , K, L
eg -Y,X-Y, Z -H-K, H, L etc
That is more complicated than you might think!!
Arguments :
- NSM (I) INTEGER
-
Number of Symmetry operations
- RSM (I) REAL
-
Array of dimension (4,4,at least NSM) containing symmetry operations on input
- SYMCHS (O) CHARACTER*(*)
-
Array of dimension at least NSM containing int tab char strings on output
- IPRINT (I) INTEGER
-
Print flag
=0 No printing
=1 Print the int tab strings
Subroutine SYMTR4(NSYM,RSM,SYMCHS)
Symmetry translation from matrix back to characters
This translates the symmetry matrices RSM(4,4,NSM) into INT TAB
character strings
It gives the real and reciprocal space operations.
eg X,Y,Z H , K, L
eg -Y,X-Y, Z -H-K, H, L etc
That is more complicated than you might think!!
Arguments :
- Nsym (I) INTEGER
-
Number of Symmetry operations
- Rsm (I) REAL
-
Array of dimension (4,4,at least Nsym) containing symmetry operations on input
- Symchs (O) CHARACTER*(*)
-
Array of dimension at least Nsym containing int tab char strings on output
Subroutine SYSAB(IN,ISYSAB)
Input IN(3) - reflection indices
Returns ISYSAB flag.
Systematic absences flagged with ISYSAB = 1
Only reflns with EPSI > 1 need be considered.
Subroutine XSPECIALS(NSYM,RSYM,XF,YF,ZF,NSPEC)
- Input NSYM
-
- number of symmetry operators. ( integer)
- Input RSYM
-
- 4*4*NSYM symmetry matrices. ( real)
- Input XF YF ZF
-
- a coordinate in fractional coordinates.
- Output NSPEC
-
- the multiplicity of the coordinate.
eg: NSPEC = 3 for an atom on a 3fod axis.
This subroutine finds what coordinates occupy special positions
i.e. have occupancies less than 1.0
from consideration of the symmetry operations.
5. DEFINITION OF THE CCP4 ASYMMETRIC UNIT
There is no standard defined asymmetric unit so the definitions are arbitrary and
may differ between differ packages. The subroutines in
group 3.a are used to define the CCP4 asymmetric unit, and
to determine whether a reflection falls within this definition.
Below are the definitions of the reciprocal space
and the real space asymmetric units under the CCP4
convention.
a. Reciprocal Space Asymmetric Unit Definitions
The reciprocal space asymmetric unit is defined in the subroutine ASUSET from the
Laue group using calls to the s/r's PGDEFN and PGNLAU. The limits of the CCP4
asymmetric unit are (from PGNLAU):
Pointgroup
| Laue group
| Limits
| Spacegroup Nos
|
3
| pg1 pg1bar
| 1bar
| hkl:l>=0 hk0:h>=0 0k0:k>=0
| 1,2
|
4
| pg2 (b) pgm pg2/m
| 2/m
| hkl:k>=0, l>=0 hk0:h>=0
| 3,4....
|
5
| pg2 (c)
| 2/m
| hkl:k>=0, l>=0 h0l:h>=0
| 1003, 1004
|
6
| pg222 pgmm2 pgmmm
| mmm
| hkl:h>=0, k>=0, l>=0
| 16 ...
|
7
| pg4 pg4bar pg4/m
| 4/m
| hkl:h>=0, l>=0 with k>=0 if h=0
and k>0 if h>0
| 75,..
|
8
| pg422 pg4mm pg4bar2m pg4barm2 pg4/mm
| 4/mmm
| hkl:h>=0, k>=0, l>=0
| 89,..
|
9
| pg3 pg3bar
| 3bar
| hkl:h>=0, k>0 00l:l>0
| 143,..
|
10
| pg312 pg32 pg3m pg3m1 pg3barm1
| 3/m
| hkl:h>=0, k>=0 with k<=h for all l.
if k=0 l>=0
| 149 151 153 157 159 162 163
|
11
| pg321 pg31m pg3bar1m
| 3bar1m
| hkl:h>=0, k>=0 with k<=h for all l.
if k=h l>=0
| 150 152 154
|
12
| pg6 pg6bar
| 6/m
| hkl:h>=0, k>=0, l>=0 with k>=0 if h=0
and k>0 if h>0
| 168..
|
13
| pg622 pg6mm pg6barm2 pg6bar2m pg6/mmm
| 6/mmm
| hkl:h>=0, k>=0, l>=0 with h>=k
| 177..
|
14
| pg23 pgm3bar
| m3
| hkl:h>=0, k>=0, l>=0 with l>=h, k>=h
| 195..
|
15
| pg432 pg4bar3m pgm3barm
| m3m
| hkl:h>=0, k>=0, l>=0 with k>=1
| 209..
|
b. Real Space Asymmetric Unit Definitions
The subroutine SETLIM contains the definitions of the
real space asymmetric unit. Note that not all of the spacegroups have a
definition within SETLIM.
No.
| Spacegroup
| Upper limits on x, y, z (*)
|
1
| P 1
| x < 1, y < 1, z < 1,
|
2
| P -1
| x < 1, y <= 1/2, z < 1,
|
3
| P 1 2 1
| x <= 1/2, y < 1, z < 1,
|
4
| P 1 21 1
| x < 1, y < 1/2, z < 1,
|
5
| C 1 2 1
| x <= 1/2, y < 1/2, z < 1,
|
10
| P 1 2/M 1
| x <= 1/2, y <= 1/2, z < 1,
|
16
| P 2 2 2
| x <= 1/2, y <= 1/2, z < 1,
|
17
| P 2 2 21
| x <= 1/2, y <= 1/2, z < 1,
|
18
| P 21 21 2
| x < 1, y <= 1/4, z < 1,
|
19
| P 21 21 21
| x < 1, y < 1, z <= 1/4,
|
20
| C 2 2 21
| x <= 1/2, y <= 1/4, z < 1,
|
21
| C 2 2 2
| x <= 1/2, y <= 1/4, z < 1,
|
22
| F 2 2 2
| x <= 1/4, y <= 1/4, z < 1,
|
23
| I 2 2 2
| x <= 1/2, y <= 1/4, z <= 1,
|
24
| I 21 21 21
| x <= 1/2, y <= 1/4, z < 1,
|
47
| P 2/M 2/M 2/M
| x <= 1/2, y <= 1/2, z <= 1/2,
|
65
| C 2/M 2/M 2/M
| x <= 1/2, y <= 1/4, z <= 1/2,
|
69
| F 2/M 2/M 2/M
| x <= 1/4, y <= 1/4, z <= 1/2,
|
71
| I 2/M 2/M 2/M
| x <= 1/2, y <= 1/4, z <= 1/2,
|
75
| P 4
| x <= 1/2, y <= 1/2, z < 1,
|
76
| P 41
| x < 1, y < 1, z < 1/4,
|
77
| P 42
| x <= 1/2, y < 1, z < 1/2,
|
78
| P 43
| x < 1, y < 1, z < 1/4,
|
79
| I 4
| x <= 1/2, y <= 1/2, z <= 1/2,
|
80
| I 41
| x <= 1/2, y < 1, z < 1/4,
|
83
| P 4/M
| x <= 1/2, y <= 1/2, z <= 1/2,
|
87
| I 4/M
| x <= 1/2, y <= 1/2, z <= 1/4,
|
89
| P 4 2 2
| x <= 1/2, y <= 1/2, z <= 1/2,
|
90
| P 4 21 2
| x <= 1/2, y <= 1/2, z <= 1/2,
|
91
| P 41 2 2
| x < 1, y < 1, z <= 1/8,
|
92
| P 41 21 2
| x < 1, y < 1, z <= 1/8,
|
93
| P 42 2 2
| x <= 1/2, y < 1, z <= 1/4,
|
94
| P 42 21 2
| x <= 1/2, y <= 1/2, z <= 1/2,
|
95
| P 43 2 2
| x < 1, y < 1, z <= 1/8,
|
96
| P 43 21 2
| x < 1, y < 1, z <= 1/8,
|
97
| I 4 2 2
| x <= 1/2, y <= 1/2, z <= 1/4,
|
98
| I 41 2 2
| x <= 1/2, y < 1, z <= 1/8,
|
123
| P 4/M 2/M 2/M
| x <= 1/2, y <= 1/2, z <= 1/2,
|
139
| I 4/M 2/M 2/M
| x <= 1/2, y <= 1/2, z <= 1/4,
|
143
| P 3
| x <= 2/3, y <= 2/3, z < 1,
|
144
| P 31
| x < 1, y < 1, z < 1/3,
|
145
| P 32
| x < 1, y < 1, z < 1/3,
|
146
| H 3
| x <= 2/3, y <= 2/3, z < 1/3,
|
147
| P -3
| x <= 2/3, y <= 2/3, z <= 1/2,
|
148
| R -3
| x <= 2/3, y <= 2/3, z <= 1/6,
|
149
| P 3 1 2
| x <= 2/3, y <= 2/3, z <= 1/2,
|
150
| P 3 2 1
| x <= 2/3, y <= 2/3, z <= 1/2,
|
151
| P 31 1 2
| x < 1, y < 1, z <= 1/6,
|
152
| P 31 2 1
| x < 1, y < 1, z <= 1/6,
|
153
| P 32 1 2
| x < 1, y < 1, z <= 1/6,
|
154
| P 32 2 1
| x < 1, y < 1, z <= 1/6,
|
155
| H 3 2
| x <= 2/3, y <= 2/3, z <= 1/6,
|
162
| P -31 2/M
| x <= 2/3, y <= 1/2, z <= 1/2,
|
164
| P -3 2/M 1
| x <= 2/3, y <= 1/3, z <= 1,
|
166
| R -3 2/M
| x <= 2/3, y <= 2/3, z <= 1/6,
|
168
| P 6
| x <= 2/3, y <= 1/2, z < 1,
|
169
| P 61
| x < 1, y < 1, z < 1/6,
|
170
| P 65
| x < 1, y < 1, z < 1/6,
|
171
| P 62
| x < 1, y < 1, z < 1/3,
|
172
| P 64
| x < 1, y < 1, z < 1/3,
|
173
| P 63
| x <= 2/3, y <= 2/3, z < 1/2,
|
175
| P 6/M
| x <= 2/3, y <= 2/3, z <= 1/2,
|
177
| P 6 2 2
| x <= 2/3, y <= 1/2, z <= 1/2,
|
178
| P 61 2 2
| x < 1, y < 1, z <= 1/12,
|
179
| P 65 2 2
| x < 1, y < 1, z <= 1/12,
|
180
| P 62 2 2
| x < 1, y < 1, z <= 1/6,
|
181
| P 64 2 2
| x < 1, y < 1, z <= 1/6,
|
182
| P 63 2 2
| x <= 2/3, y <= 2/3, z <= 1/4,
|
191
| P 6/M 2/M 2/M
| x <= 2/3, y <= 1/3, z <= 1/2,
|
195
| P 2 3
| x < 1, y < 1, z <= 1/2,
|
196
| F 2 3
| x <= 1/4, y <= 1/4, z < 1,
|
197
| I 2 3
| x < 1, y < 1, z <= 1/2,
|
198
| P 21 3
| x <= 1/2, y <= 1/2, z < 1,
|
199
| I 21 3
| x <= 1/2, y <= 1/2, z <= 1/2,
|
200
| P 2/M -3
| x <= 1/2, y <= 1/2, z <= 1/2,
|
202
| F 2/M -3
| x <= 1/2, y <= 1/2, z <= 1/4,
|
204
| I 2/M -3
| x <= 1/2, y <= 1/2, z <= 1/2,
|
207
| P 4 3 2
| x < 1, y <= 1/2, z <= 1/2,
|
208
| P 42 3 2
| x <= 1/2, y < 1, z <= 1/4,
|
209
| F 4 3 2
| x <= 1/2, y <= 1/2, z <= 1/2,
|
210
| F 41 3 2
| x <= 1/2, y < 1, z <= 1/8,
|
211
| I 4 3 2
| x <= 1/2, y <= 1/2, z <= 1/4,
|
212
| P 43 3 2
| x < 1, y < 1, z <= 1/8,
|
213
| P 41 3 2
| x < 1, y < 1, z <= 1/8,
|
214
| I 41 3 2
| x <= 1/2, y < 1, z <= 1/8,
|
221
| P 4/M -3 2/M
| x <= 1/2, y <= 1/2, z <= 1/2,
|
225
| F 4/M -3 2/M
| x <= 1/2, y <= 1/4, z <= 1/4,
|
229
| I 4/M -3 2/M
| x <= 1/2, y <= 1/2, z <= 1/4,
|
(*) The limits are in fractional coordinates,
and the lower limits are always x=0, y=0, z=0.