ProSHADE  0.7.6.6 (JUL 2022)
Protein Shape Detection
ProSHADE_symmetry.hpp
Go to the documentation of this file.
1 
22 //==================================================== ProSHADE
23 #include "ProSHADE_distances.hpp"
24 
25 //==================================================== Overinclusion protection
26 #ifndef PROSHADE_SYMMETRY
27 #define PROSHADE_SYMMETRY
28 
29 //==================================================== ProSHADE_internal_symmetry Namespace
36 {
37  bool isSymmetrySame ( std::vector< proshade_double* >* ret, proshade_double* sym, proshade_double simThres, proshade_signed* matchedPos );
38  bool isSymmetrySame ( std::vector< proshade_double* >* ret, proshade_double* sym, proshade_double simThres, proshade_signed* matchedPos, proshade_double fscVal );
39  bool detectTetrahedralSymmetry ( std::vector< proshade_double* >* CSymList, proshade_double axErr, proshade_double minPeakHeight );
40  void findTetra4C3s ( std::vector< proshade_double* >* CSymList, std::vector< proshade_double* >* ret,
41  proshade_double axErr, ProSHADE_internal_data::ProSHADE_data* dataObj,
42  proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight );
43  bool testGroupAgainstSymmetry ( std::vector< proshade_double* >* CSymList, std::vector< proshade_unsign >* grp,
44  proshade_double* sym, proshade_double axErr, proshade_double angle, bool improve,
45  proshade_unsign pos = 0 );
46  bool findMissingAxes ( std::vector< std::vector< proshade_unsign > >* possibilities,
47  std::vector< proshade_double* >* CSymList, proshade_unsign requiredNoAxes,
48  proshade_double axErr, proshade_double angle, proshade_unsign fold,
49  ProSHADE_internal_data::ProSHADE_data* dataObj, proshade_double minPeakHeight );
50  proshade_double missingAxisHeight ( proshade_double xVal, proshade_double yVal, proshade_double zVal,
51  ProSHADE_internal_data::ProSHADE_data* dataObj, proshade_unsign fold, proshade_double axErr );
52  std::vector < proshade_double* > findMissingAxisPoints ( proshade_double xVal, proshade_double yVal, proshade_double zVal,
53  ProSHADE_internal_data::ProSHADE_data* dataObj, proshade_double axErr );
54  bool sortArrVecHlp ( const proshade_double* a, const proshade_double* b );
55  void saveMissingAxisNewOnly ( std::vector< proshade_double* >* axVec, proshade_double axX, proshade_double axY,
56  proshade_double axZ, proshade_double height, proshade_unsign fold, proshade_double axErr );
57  void searchMissingSymmetrySpace ( ProSHADE_internal_data::ProSHADE_data* dataObj, std::vector< proshade_double* >* CSymList,
58  std::vector< proshade_unsign >* grp, std::vector< proshade_double* >* hlpVec,
59  proshade_double axErr, proshade_double angle, proshade_unsign fold,
60  proshade_double minPeakHeight );
61  void findTetra3C2s ( std::vector< proshade_double* >* CSymList, std::vector< proshade_double* >* ret,
62  proshade_double axErr, ProSHADE_internal_data::ProSHADE_data* dataObj,
63  proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight );
64  bool testGroupAgainstGroup ( std::vector< proshade_double* >* CSymList, std::vector< proshade_unsign >* grp1,
65  std::vector< proshade_double* >* RetList, std::vector< proshade_unsign >* grp2,
66  proshade_double angle, proshade_double axErr );
67  bool detectOctahedralSymmetry ( std::vector< proshade_double* >* CSymList, proshade_double axErr, proshade_double minPeakHeight );
68  void findOcta3C4s ( std::vector< proshade_double* >* CSymList, std::vector< proshade_double* >* ret,
69  proshade_double axErr, ProSHADE_internal_data::ProSHADE_data* dataObj,
70  proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight );
71  void findOcta4C3s ( std::vector< proshade_double* >* CSymList, std::vector< proshade_double* >* ret,
72  proshade_double axErr, ProSHADE_internal_data::ProSHADE_data* dataObj,
73  proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight );
74  void findOcta6C2s ( std::vector< proshade_double* >* CSymList, std::vector< proshade_double* >* ret,
75  proshade_double axErr, ProSHADE_internal_data::ProSHADE_data* dataObj,
76  proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight );
77  bool findMissingAxesDual ( std::vector< proshade_unsign >* possibilities,
78  std::vector< proshade_double* >* CSymList, std::vector< proshade_double* >* ret, std::vector<
79  proshade_unsign >* retGroup, proshade_unsign requiredNoAxes, proshade_double axErr,
80  proshade_unsign noMatchesG1, proshade_double angle1, proshade_unsign noMatchesG2,
81  proshade_double angle2, proshade_unsign fold, ProSHADE_internal_data::ProSHADE_data* dataObj );
82  proshade_signed addAxisUnlessSame ( proshade_unsign fold, proshade_double axX, proshade_double axY, proshade_double axZ,
83  proshade_double axHeight, std::vector< proshade_double* >* prosp, proshade_double axErr );
84  proshade_signed addAxisUnlessSame ( proshade_unsign fold, proshade_double axX, proshade_double axY, proshade_double axZ,
85  proshade_double axHeight, proshade_double averageFSC, std::vector< proshade_double* >* prosp, proshade_double axErr );
86  bool checkFittingAxisDualAndSave ( std::vector< proshade_unsign >* retGroup, std::vector< proshade_double* >* ret,
87  proshade_unsign fold, proshade_double axX, proshade_double axY, proshade_double axZ,
88  std::vector< proshade_double* >* prosp, proshade_double axErr, proshade_unsign noMatchesG1,
89  proshade_double angle1, proshade_unsign noMatchesG2, proshade_double angle2,
91  bool detectIcosahedralSymmetry ( std::vector< proshade_double* >* CSymList, proshade_double axErr, proshade_double minPeakHeight );
92  void findIcos6C5s ( std::vector< proshade_double* >* CSymList, std::vector< proshade_double* >* ret,
93  proshade_double axErr, ProSHADE_internal_data::ProSHADE_data* dataObj,
94  proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight );
95  void findIcos10C3s ( std::vector< proshade_double* >* CSymList, std::vector< proshade_double* >* ret,
96  proshade_double axErr, ProSHADE_internal_data::ProSHADE_data* dataObj,
97  proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight );
98  void findIcos15C2s ( std::vector< proshade_double* >* CSymList, std::vector< proshade_double* >* ret,
99  proshade_double axErr, ProSHADE_internal_data::ProSHADE_data* dataObj,
100  proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight );
101  bool findMissingAxesTriple ( std::vector< proshade_unsign >* possibilities,
102  std::vector< proshade_double* >* CSymList, std::vector< proshade_double* >* ret, std::vector<
103  proshade_unsign >* retGroup, proshade_unsign requiredNoAxes, proshade_double axErr,
104  proshade_unsign noMatchesG1, proshade_double angle1, proshade_unsign noMatchesG2,
105  proshade_double angle2, proshade_unsign noMatchesG3, proshade_double angle3,
106  proshade_unsign fold, ProSHADE_internal_data::ProSHADE_data* dataObj );
107  void checkFittingAxisTripleAndSave ( std::vector< proshade_unsign >* retGroup, std::vector< proshade_double* >* ret,
108  proshade_unsign fold, proshade_double axX, proshade_double axY, proshade_double axZ,
109  std::vector< proshade_double* >* prosp, proshade_double axErr,
110  proshade_unsign noMatchesG1, proshade_double angle1, proshade_unsign noMatchesG2,
111  proshade_double angle2, proshade_unsign noMatchesG3, proshade_double angle3,
113  proshade_double findPredictedSingleAxisHeight ( proshade_double* axis, proshade_double fold, ProSHADE_internal_data::ProSHADE_data* dataObj, ProSHADE_settings* settings );
114  void findPredictedAxesHeights ( std::vector< proshade_double* >* ret, ProSHADE_internal_data::ProSHADE_data* dataObj, ProSHADE_settings* settings );
115  void optimiseDGroupAngleFromAxesHeights ( std::vector < std::vector< proshade_double > >* ret, ProSHADE_internal_data::ProSHADE_data* dataObj, ProSHADE_settings* settings );
116  void optimiseDGroupAngleFromAxesHeights ( std::vector < std::vector< proshade_double > >* allCs, std::vector< proshade_unsign > selection, ProSHADE_internal_data::ProSHADE_data* dataObj,
117  ProSHADE_settings* settings );
118  void predictIcosAxes ( std::vector< proshade_double* >* CSymList, std::vector< std::vector< proshade_double* > >* ret, proshade_double axErr, proshade_double minPeakHeight );
119  void predictOctaAxes ( std::vector< proshade_double* >* CSymList, std::vector< std::vector< proshade_double* > >* ret, proshade_double axErr, proshade_double minPeakHeight );
120  void predictTetraAxes ( std::vector< proshade_double* >* CSymList, std::vector< std::vector< proshade_double* > >* ret, proshade_double axErr, proshade_double minPeakHeight );
121  std::vector< proshade_unsign > findReliableUnphasedSymmetries ( std::vector < proshade_double* >* allCs, std::vector < std::vector < proshade_double* > >* allDs, proshade_signed verbose, proshade_signed messageShift, proshade_double tolerance );
122  void allocateCentreOfMapFourierTransforms ( proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim, fftw_complex *&origMap, fftw_complex *&origCoeffs, fftw_complex *&rotMapComplex,
123  fftw_complex *&rotCoeffs, fftw_complex *&trFunc, fftw_complex *&trFuncCoeffs, fftw_plan *planForwardFourier, fftw_plan *planForwardFourierRot,
124  fftw_plan *planReverseFourierComb );
125  void releaseCentreOfMapFourierTransforms ( fftw_complex *origMap, fftw_complex *origCoeffs, fftw_complex *rotMapComplex, fftw_complex *rotCoeffs, fftw_complex *trFunc, fftw_complex *trFuncCoeffs,
126  fftw_plan planForwardFourier, fftw_plan planForwardFourierRot, fftw_plan planReverseFourierComb );
127  std::vector< proshade_double > findTranslationBetweenRotatedAndOriginalMap ( ProSHADE_internal_data::ProSHADE_data* symStr, std::vector < proshade_double > symElem, fftw_complex *origCoeffs,
128  fftw_complex* rotMapComplex, fftw_complex* rotCoeffs, fftw_plan planForwardFourierRot, fftw_complex* trFuncCoeffs,
129  fftw_complex* trFunc, fftw_plan planReverseFourierComb );
130  std::vector< proshade_double > findPointFromTranslations ( ProSHADE_internal_data::ProSHADE_data* symStr, std::vector < std::vector < proshade_double > > symElems, fftw_complex *origCoeffs,
131  fftw_complex* rotMapComplex, fftw_complex* rotCoeffs, fftw_plan planForwardFourierRot, fftw_complex* trFuncCoeffs,
132  fftw_complex* trFunc, fftw_plan planReverseFourierComb );
133 }
134 
135 #endif
ProSHADE_internal_symmetry::findTetra3C2s
void findTetra3C2s(std::vector< proshade_double * > *CSymList, std::vector< proshade_double * > *ret, proshade_double axErr, ProSHADE_internal_data::ProSHADE_data *dataObj, proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight)
This function takes the list of C symmetries and finds the 3 C2 symmetries with correct angles requir...
Definition: ProSHADE_symmetry.cpp:998
ProSHADE_internal_symmetry::detectOctahedralSymmetry
bool detectOctahedralSymmetry(std::vector< proshade_double * > *CSymList, proshade_double axErr, proshade_double minPeakHeight)
This function takes the list of C symmetries and decides whether basic requirements for octahhedral s...
Definition: ProSHADE_symmetry.cpp:1127
ProSHADE_distances.hpp
This is the header file containing declarations of functions required for computation of shape distan...
ProSHADE_internal_symmetry::findMissingAxes
bool findMissingAxes(std::vector< std::vector< proshade_unsign > > *possibilities, std::vector< proshade_double * > *CSymList, proshade_unsign requiredNoAxes, proshade_double axErr, proshade_double angle, proshade_unsign fold, ProSHADE_internal_data::ProSHADE_data *dataObj, proshade_double minPeakHeight)
This function tries to find an axis which would complete a particular group of axes for polyhedral sy...
Definition: ProSHADE_symmetry.cpp:598
ProSHADE_internal_symmetry::sortArrVecHlp
bool sortArrVecHlp(const proshade_double *a, const proshade_double *b)
This function compares two arrays of two based on the first number.
Definition: ProSHADE_symmetry.cpp:658
ProSHADE_internal_symmetry::findIcos6C5s
void findIcos6C5s(std::vector< proshade_double * > *CSymList, std::vector< proshade_double * > *ret, proshade_double axErr, ProSHADE_internal_data::ProSHADE_data *dataObj, proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight)
This function takes the list of C symmetries and finds the six C5 symmetries with given angles requir...
Definition: ProSHADE_symmetry.cpp:1912
ProSHADE_internal_symmetry::testGroupAgainstGroup
bool testGroupAgainstGroup(std::vector< proshade_double * > *CSymList, std::vector< proshade_unsign > *grp1, std::vector< proshade_double * > *RetList, std::vector< proshade_unsign > *grp2, proshade_double angle, proshade_double axErr)
This function compares two groups of axes for a single pair having the required angle.
Definition: ProSHADE_symmetry.cpp:1083
ProSHADE_internal_symmetry::addAxisUnlessSame
proshade_signed addAxisUnlessSame(proshade_unsign fold, proshade_double axX, proshade_double axY, proshade_double axZ, proshade_double axHeight, proshade_double averageFSC, std::vector< proshade_double * > *prosp, proshade_double axErr)
This function simply creates a new axis from information in aruments and tests if no such axis alread...
Definition: ProSHADE_symmetry.cpp:1511
ProSHADE_internal_symmetry::findPredictedAxesHeights
void findPredictedAxesHeights(std::vector< proshade_double * > *ret, ProSHADE_internal_data::ProSHADE_data *dataObj, ProSHADE_settings *settings)
This function finds the rotation function value for all axes supplied in the ret parameter.
Definition: ProSHADE_symmetry.cpp:3063
ProSHADE_internal_symmetry::findPointFromTranslations
std::vector< proshade_double > findPointFromTranslations(ProSHADE_internal_data::ProSHADE_data *symStr, std::vector< std::vector< proshade_double > > symElems, fftw_complex *origCoeffs, fftw_complex *rotMapComplex, fftw_complex *rotCoeffs, fftw_plan planForwardFourierRot, fftw_complex *trFuncCoeffs, fftw_complex *trFunc, fftw_plan planReverseFourierComb)
This function computes the average of optimal translations for a cyclic point group.
Definition: ProSHADE_symmetry.cpp:3985
ProSHADE_internal_symmetry::allocateCentreOfMapFourierTransforms
void allocateCentreOfMapFourierTransforms(proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim, fftw_complex *&origMap, fftw_complex *&origCoeffs, fftw_complex *&rotMapComplex, fftw_complex *&rotCoeffs, fftw_complex *&trFunc, fftw_complex *&trFuncCoeffs, fftw_plan *planForwardFourier, fftw_plan *planForwardFourierRot, fftw_plan *planReverseFourierComb)
This function allocates the required memory for the Fourier transforms required to find the centre of...
Definition: ProSHADE_symmetry.cpp:3841
ProSHADE_internal_data::ProSHADE_data
This class contains all inputed and derived data for a single structure.
Definition: ProSHADE_data.hpp:49
ProSHADE_internal_symmetry::missingAxisHeight
proshade_double missingAxisHeight(proshade_double xVal, proshade_double yVal, proshade_double zVal, ProSHADE_internal_data::ProSHADE_data *dataObj, proshade_unsign fold, proshade_double axErr)
This function searches for the highest peaks average that would produce the required axis and fold.
Definition: ProSHADE_symmetry.cpp:679
ProSHADE_internal_symmetry::searchMissingSymmetrySpace
void searchMissingSymmetrySpace(ProSHADE_internal_data::ProSHADE_data *dataObj, std::vector< proshade_double * > *CSymList, std::vector< proshade_unsign > *grp, std::vector< proshade_double * > *hlpVec, proshade_double axErr, proshade_double angle, proshade_unsign fold, proshade_double minPeakHeight)
This function tests feasible axes against the missing axis criteria, returning a set of matching axes...
Definition: ProSHADE_symmetry.cpp:890
ProSHADE_internal_symmetry::findMissingAxesDual
bool findMissingAxesDual(std::vector< proshade_unsign > *possibilities, std::vector< proshade_double * > *CSymList, std::vector< proshade_double * > *ret, std::vector< proshade_unsign > *retGroup, proshade_unsign requiredNoAxes, proshade_double axErr, proshade_unsign noMatchesG1, proshade_double angle1, proshade_unsign noMatchesG2, proshade_double angle2, proshade_unsign fold, ProSHADE_internal_data::ProSHADE_data *dataObj)
This function tries to find a particular symmetry axes which would complete a group of symmetries wit...
Definition: ProSHADE_symmetry.cpp:1416
ProSHADE_internal_symmetry::optimiseDGroupAngleFromAxesHeights
void optimiseDGroupAngleFromAxesHeights(std::vector< std::vector< proshade_double > > *ret, ProSHADE_internal_data::ProSHADE_data *dataObj, ProSHADE_settings *settings)
This function takes two axes with almost dihedral angle and optimises their relative positions as wel...
Definition: ProSHADE_symmetry.cpp:3269
ProSHADE_internal_symmetry::isSymmetrySame
bool isSymmetrySame(std::vector< proshade_double * > *ret, proshade_double *sym, proshade_double simThres, proshade_signed *matchedPos)
This function checks if a very similar symmetry is not already saved.
Definition: ProSHADE_symmetry.cpp:215
ProSHADE_internal_symmetry::saveMissingAxisNewOnly
void saveMissingAxisNewOnly(std::vector< proshade_double * > *axVec, proshade_double axX, proshade_double axY, proshade_double axZ, proshade_double height, proshade_unsign fold, proshade_double axErr)
This function saves the recovered information about missing axis into a full symmetry,...
Definition: ProSHADE_symmetry.cpp:821
ProSHADE_internal_symmetry::findIcos10C3s
void findIcos10C3s(std::vector< proshade_double * > *CSymList, std::vector< proshade_double * > *ret, proshade_double axErr, ProSHADE_internal_data::ProSHADE_data *dataObj, proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight)
This function takes the list of C symmetries and finds the ten C3 symmetries with correct angles requ...
Definition: ProSHADE_symmetry.cpp:2337
ProSHADE_internal_symmetry
This namespace contains the symmetry detection related code.
Definition: ProSHADE_data.cpp:75
ProSHADE_internal_symmetry::detectIcosahedralSymmetry
bool detectIcosahedralSymmetry(std::vector< proshade_double * > *CSymList, proshade_double axErr, proshade_double minPeakHeight)
This function takes the list of C symmetries and decides whether basic requirements for isosahedral s...
Definition: ProSHADE_symmetry.cpp:1854
ProSHADE_settings
This class stores all the settings and is passed to the executive classes instead of a multitude of p...
Definition: ProSHADE_settings.hpp:37
ProSHADE_internal_symmetry::findOcta4C3s
void findOcta4C3s(std::vector< proshade_double * > *CSymList, std::vector< proshade_double * > *ret, proshade_double axErr, ProSHADE_internal_data::ProSHADE_data *dataObj, proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight)
This function takes the list of C symmetries and finds the four C3 symmetries with correct angles req...
Definition: ProSHADE_symmetry.cpp:1252
ProSHADE_internal_symmetry::findReliableUnphasedSymmetries
std::vector< proshade_unsign > findReliableUnphasedSymmetries(std::vector< proshade_double * > *allCs, std::vector< std::vector< proshade_double * > > *allDs, proshade_signed verbose, proshade_signed messageShift, proshade_double tolerance)
This function checks the list of detected axes (presumably from phaseless symmetry detection) and ret...
Definition: ProSHADE_symmetry.cpp:3735
ProSHADE_internal_symmetry::findTranslationBetweenRotatedAndOriginalMap
std::vector< proshade_double > findTranslationBetweenRotatedAndOriginalMap(ProSHADE_internal_data::ProSHADE_data *symStr, std::vector< proshade_double > symElem, fftw_complex *origCoeffs, fftw_complex *rotMapComplex, fftw_complex *rotCoeffs, fftw_plan planForwardFourierRot, fftw_complex *trFuncCoeffs, fftw_complex *trFunc, fftw_plan planReverseFourierComb)
This function takes a single rotation matrix and procceds to compute the optimal translation between ...
Definition: ProSHADE_symmetry.cpp:3919
ProSHADE_internal_symmetry::releaseCentreOfMapFourierTransforms
void releaseCentreOfMapFourierTransforms(fftw_complex *origMap, fftw_complex *origCoeffs, fftw_complex *rotMapComplex, fftw_complex *rotCoeffs, fftw_complex *trFunc, fftw_complex *trFuncCoeffs, fftw_plan planForwardFourier, fftw_plan planForwardFourierRot, fftw_plan planReverseFourierComb)
This function releases the allocated memory for the Fourier transforms used to find the centre of the...
Definition: ProSHADE_symmetry.cpp:3881
ProSHADE_internal_symmetry::findMissingAxesTriple
bool findMissingAxesTriple(std::vector< proshade_unsign > *possibilities, std::vector< proshade_double * > *CSymList, std::vector< proshade_double * > *ret, std::vector< proshade_unsign > *retGroup, proshade_unsign requiredNoAxes, proshade_double axErr, proshade_unsign noMatchesG1, proshade_double angle1, proshade_unsign noMatchesG2, proshade_double angle2, proshade_unsign noMatchesG3, proshade_double angle3, proshade_unsign fold, ProSHADE_internal_data::ProSHADE_data *dataObj)
This function tries to find a particular symmetry axis which would complete a group of symmetries wit...
Definition: ProSHADE_symmetry.cpp:2493
ProSHADE_internal_symmetry::predictOctaAxes
void predictOctaAxes(std::vector< proshade_double * > *CSymList, std::vector< std::vector< proshade_double * > > *ret, proshade_double axErr, proshade_double minPeakHeight)
This function predicts all octahedral point group symmetry axes from the cyclic point groups list.
Definition: ProSHADE_symmetry.cpp:2214
ProSHADE_internal_symmetry::findOcta6C2s
void findOcta6C2s(std::vector< proshade_double * > *CSymList, std::vector< proshade_double * > *ret, proshade_double axErr, ProSHADE_internal_data::ProSHADE_data *dataObj, proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight)
This function takes the list of C symmetries and finds the six C2 symmetries with correct angles requ...
Definition: ProSHADE_symmetry.cpp:1335
ProSHADE_internal_symmetry::testGroupAgainstSymmetry
bool testGroupAgainstSymmetry(std::vector< proshade_double * > *CSymList, std::vector< proshade_unsign > *grp, proshade_double *sym, proshade_double axErr, proshade_double angle, bool improve, proshade_unsign pos=0)
This function tests whether a symmetry has particular angle to all members of a group.
Definition: ProSHADE_symmetry.cpp:531
ProSHADE_internal_symmetry::findOcta3C4s
void findOcta3C4s(std::vector< proshade_double * > *CSymList, std::vector< proshade_double * > *ret, proshade_double axErr, ProSHADE_internal_data::ProSHADE_data *dataObj, proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight)
This function takes the list of C symmetries and finds the 3 C4 symmetries with perpendicular angles ...
Definition: ProSHADE_symmetry.cpp:1185
ProSHADE_internal_symmetry::predictTetraAxes
void predictTetraAxes(std::vector< proshade_double * > *CSymList, std::vector< std::vector< proshade_double * > > *ret, proshade_double axErr, proshade_double minPeakHeight)
This function predicts all tetrahedral point group symmetry axes from the cyclic point groups list.
Definition: ProSHADE_symmetry.cpp:3609
ProSHADE_internal_symmetry::findTetra4C3s
void findTetra4C3s(std::vector< proshade_double * > *CSymList, std::vector< proshade_double * > *ret, proshade_double axErr, ProSHADE_internal_data::ProSHADE_data *dataObj, proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight)
This function takes the list of C symmetries and finds the 4 C3 symmetries with correct angles requir...
Definition: ProSHADE_symmetry.cpp:461
ProSHADE_internal_symmetry::findPredictedSingleAxisHeight
proshade_double findPredictedSingleAxisHeight(proshade_double *axis, proshade_double fold, ProSHADE_internal_data::ProSHADE_data *dataObj, ProSHADE_settings *settings)
This function finds the rotation function value for a single axis.
Definition: ProSHADE_symmetry.cpp:3409
ProSHADE_internal_symmetry::checkFittingAxisTripleAndSave
void checkFittingAxisTripleAndSave(std::vector< proshade_unsign > *retGroup, std::vector< proshade_double * > *ret, proshade_unsign fold, proshade_double axX, proshade_double axY, proshade_double axZ, std::vector< proshade_double * > *prosp, proshade_double axErr, proshade_unsign noMatchesG1, proshade_double angle1, proshade_unsign noMatchesG2, proshade_double angle2, proshade_unsign noMatchesG3, proshade_double angle3, ProSHADE_internal_data::ProSHADE_data *dataObj)
This function takes a newly detected "missing" axis and tests it for belonging to the group,...
Definition: ProSHADE_symmetry.cpp:2604
ProSHADE_internal_symmetry::findIcos15C2s
void findIcos15C2s(std::vector< proshade_double * > *CSymList, std::vector< proshade_double * > *ret, proshade_double axErr, ProSHADE_internal_data::ProSHADE_data *dataObj, proshade_signed verbose, proshade_signed messageShift, proshade_double minPeakHeight)
This function takes the list of C symmetries and finds the fifteen C3 symmetries with correct angles ...
Definition: ProSHADE_symmetry.cpp:2410
ProSHADE_internal_symmetry::predictIcosAxes
void predictIcosAxes(std::vector< proshade_double * > *CSymList, std::vector< std::vector< proshade_double * > > *ret, proshade_double axErr, proshade_double minPeakHeight)
This function predicts all possible icosahedral point groups symmetry axes from the cyclic point grou...
Definition: ProSHADE_symmetry.cpp:2034
ProSHADE_internal_symmetry::detectTetrahedralSymmetry
bool detectTetrahedralSymmetry(std::vector< proshade_double * > *CSymList, proshade_double axErr, proshade_double minPeakHeight)
This function takes the list of C symmetries and decides whether basic requirements for tetrahedral s...
Definition: ProSHADE_symmetry.cpp:410
ProSHADE_internal_symmetry::checkFittingAxisDualAndSave
bool checkFittingAxisDualAndSave(std::vector< proshade_unsign > *retGroup, std::vector< proshade_double * > *ret, proshade_unsign fold, proshade_double axX, proshade_double axY, proshade_double axZ, std::vector< proshade_double * > *prosp, proshade_double axErr, proshade_unsign noMatchesG1, proshade_double angle1, proshade_unsign noMatchesG2, proshade_double angle2, ProSHADE_internal_data::ProSHADE_data *dataObj)
This function takes a newly detected "missing" axis and tests it for belonging to the group,...
Definition: ProSHADE_symmetry.cpp:1610
ProSHADE_internal_symmetry::findMissingAxisPoints
std::vector< proshade_double * > findMissingAxisPoints(proshade_double xVal, proshade_double yVal, proshade_double zVal, ProSHADE_internal_data::ProSHADE_data *dataObj, proshade_double axErr)
This function searches for all the self-rotation map points conforming to the axis,...
Definition: ProSHADE_symmetry.cpp:743