53 std::string path, std::vector<std::vector<int>> rings,
54 std::vector<std::vector<int>> nList,
56 bool printClusters,
bool onlyTetrahedral) {
60 std::vector<ring::strucType>
64 std::vector<cage::Cage> cageList;
65 std::vector<std::vector<int>>
69 std::vector<cage::iceType>
72 int numHC, numDDC, mixedRings, prismaticRings, basalRings;
75 std::vector<double> quat;
76 std::vector<std::vector<double>> quatList;
78 Eigen::MatrixXd refPntsDDC(14,
80 Eigen::MatrixXd refPntsHC(12,
84 std::vector<double> rmsdPerAtom;
89 if (onlyTetrahedral) {
96 atomTypes.resize(yCloud->nop);
98 rmsdPerAtom.resize(yCloud->nop, -1);
99 noOfCommonElements.resize(yCloud->nop);
104 for (
int ringSize = initRingSize; ringSize <= maxRingSize; ringSize++) {
110 if (ringsOneType.size() == 0) {
116 ringsOneType.size());
123 firstFrame, &numHC, &numDDC, &ringType);
145 std::string filePathHC =
"templates/hc.xyz" ;
146 std::string filePathDDC =
"templates/ddc.xyz" ;
153 for (
int i = 0; i < yCloud->nop; i++) {
154 if (rmsdPerAtom[i] != -1) {
155 noOfCommonElements[i] = 1;
164 for (
int icage = 0; icage < numHC; icage++) {
169 quatList.push_back(quat);
172 &noOfCommonElements, atomTypes);
176 for (
int icage = numHC; icage < cageList.size(); icage++) {
181 quatList.push_back(quat);
184 &noOfCommonElements, atomTypes);
201 std::vector<int> dumpAtomTypes;
202 dumpAtomTypes.resize(atomTypes.size());
204 for (
int i = 0; i < atomTypes.size(); i++) {
206 dumpAtomTypes[i] = 1;
209 dumpAtomTypes[i] = 2;
212 dumpAtomTypes[i] = 3;
215 dumpAtomTypes[i] = 4;
218 dumpAtomTypes[i] = 5;
221 dumpAtomTypes[i] = 0;
253 const Eigen::MatrixXd &refPoints,
cage::Cage cageUnit,
254 std::vector<std::vector<int>> rings, std::vector<std::vector<int>> nList,
255 std::vector<double> *quat,
double *rmsd) {
259 std::vector<int> basal1,
262 Eigen::MatrixXd targetPointSet(12, 3);
264 std::vector<double> rmsdList;
267 std::vector<double> currentQuat;
274 iring = cageUnit.
rings [0];
275 jring = cageUnit.
rings [1];
276 rmsdList.resize(yCloud->nop);
285 for (
int i = 0; i < ringSize; i++) {
291 ¤tRmsd, &rmsdList, ¤tScale);
295 scale = currentScale;
298 if (currentRmsd < *rmsd) {
301 scale = currentScale;
325 const Eigen::MatrixXd &refPoints, std::vector<cage::Cage> cageList,
326 int cageIndex, std::vector<std::vector<int>> rings,
327 std::vector<double> *quat,
double *rmsd) {
329 std::vector<int> ddcOrder;
331 Eigen::MatrixXd targetPointSet(14, 3);
333 std::vector<double> rmsdList;
336 std::vector<double> currentQuat;
348 for (
int i = 0; i < ringSize; i++) {
354 ¤tRmsd, &rmsdList, ¤tScale);
358 scale = currentScale;
361 if (currentRmsd < *rmsd) {
364 scale = currentScale;
379 Eigen::MatrixXd refPnts(12, 3);
384 std::vector<std::vector<int>> nList;
385 std::vector<std::vector<int>> rings;
386 std::vector<ring::strucType>
388 std::vector<int> listHC;
390 std::vector<cage::Cage> cageList;
397 for (
int i = 0; i < 3; i++) {
398 setCloud.box[i] = 50;
408 ringType.resize(rings.size());
410 listHC =
ring::findHC (rings, &ringType, nList, &cageList);
412 iring = cageList[0].rings[0];
413 jring = cageList[0].rings[1];
415 std::vector<int> matchedBasal1,
423 &matchedBasal1, &matchedBasal2);
437 Eigen::MatrixXd refPnts(14, 3);
442 std::vector<std::vector<int>> nList;
443 std::vector<std::vector<int>> rings;
444 std::vector<ring::strucType>
446 std::vector<int> listDDC,
448 std::vector<int> ddcOrder;
450 std::vector<cage::Cage> cageList;
458 for (
int i = 0; i < 3; i++) {
459 setCloud.
box [i] = 50;
469 ringType.resize(rings.size());
471 listDDC =
ring::findDDC (rings, &ringType, listHC, &cageList);
488 std::vector<double> *rmsdPerAtom,
489 std::vector<int> *noOfCommonAtoms,
490 std::vector<cage::iceType> atomTypes) {
492 int nRings = cageUnit.
rings .size();
496 int ringSize = rings[0].size();
499 for (
int i = 0; i < nRings; i++) {
500 iring = cageUnit.
rings [i];
503 for (
int j = 0; j < ringSize; j++) {
504 iatom = rings[iring][j];
512 if ((*rmsdPerAtom)[iatom] == -1) {
513 (*rmsdPerAtom)[iatom] = rmsd;
514 (*noOfCommonAtoms)[iatom] = 1;
516 (*rmsdPerAtom)[iatom] += rmsd;
517 (*noOfCommonAtoms)[iatom] += 1;
530 std::vector<int> *noOfCommonAtoms) {
532 int nop = (*rmsdPerAtom).size();
534 for (
int iatom = 0; iatom < nop; iatom++) {
536 if ((*noOfCommonAtoms)[iatom] == 0) {
537 (*rmsdPerAtom)[iatom] = -1;
539 (*rmsdPerAtom)[iatom] /= (*noOfCommonAtoms)[iatom];
582 std::string path, std::vector<std::vector<int>> rings,
583 std::vector<std::vector<int>> nList,
585 int *numHC,
int *numDDC, std::vector<ring::strucType> *ringType) {
588 std::vector<int> listDDC;
589 std::vector<int> listHC;
594 std::vector<cage::Cage> cageList;
596 int mixedRings, prismaticRings, basalRings;
604 if (rings.size() == 0) {
609 (*ringType).resize(rings.size());
614 listHC =
ring::findHC (rings, ringType, nList, &cageList);
626 &prismaticRings, &basalRings);
630 mixedRings, basalRings, prismaticRings, firstFrame);
646 std::vector<std::vector<int>> rings, std::vector<cage::Cage> cageList,
647 int numHC,
int numDDC) {
649 std::vector<int> linkedList;
650 int nCages = numHC + numDDC;
655 int singleDDCs, singleHCs;
657 double volDDC = 72.7;
658 double volHC = 40.63;
663 std::vector<int> nClusters;
665 std::vector<int> atoms;
670 linkedList.resize(nCages);
672 for (
int icage = 0; icage < nCages; icage++) {
674 linkedList[icage] = icage;
680 for (
int i = 0; i < nCages - 1; i++) {
682 if (linkedList[i] != i) {
691 for (
int k = i + 1; k < nCages; k++) {
693 if (linkedList[k] != k) {
698 if (cageList[k].type != cageList[j].type) {
706 temp = linkedList[j];
707 linkedList[j] = linkedList[k];
708 linkedList[k] = temp;
719 visited.resize(nCages);
723 for (
int i = 0; i < nCages; i++) {
731 if (linkedList[i] == i) {
745 type = cageList[i].type;
747 nextElement = linkedList[currentIndex];
752 while (nextElement != index) {
754 currentIndex = nextElement;
755 visited[currentIndex] =
true ;
757 nextElement = linkedList[currentIndex];
760 nClusters.push_back(iClusterNumber);
764 int clusterID = nClusters.size();
773 std::ofstream outputFile;
776 std::string outputDirName = path +
"bulkTopo" ;
778 outputDirName = path +
"bulkTopo/clusterXYZ/" ;
780 outputDirName = path +
"bulkTopo/clusterXYZ/frame-" +
781 std::to_string(yCloud->currentFrame);
785 outputFile.open(path +
"bulkTopo/clusterXYZ/frame-" +
786 std::to_string(yCloud->currentFrame) +
"/info.dat" );
787 outputFile <<
"# volDDC volHC in Angstrom^3\n" ;
788 outputFile << singleDDCs * volDDC <<
" " << singleHCs * volHC <<
"\n" ;
789 outputFile <<
"There are " << singleDDCs <<
" single DDCs\n" ;
790 outputFile <<
"There are " << singleHCs <<
" single HCs\n" ;
802 std::vector<cage::Cage> cageList,
803 std::vector<int> clusterCages) {
805 std::vector<int> atoms;
806 int ringSize = rings[0].size();
815 for (
int j = 0; j < cageList[icage].rings.size(); j++) {
816 iring = cageList[icage].rings[j];
818 for (
int k = 0; k < ringSize; k++) {
819 iatom = rings[iring][k];
820 atoms.push_back(iatom);
828 sort(atoms.begin(), atoms.end());
829 auto ip = std::unique(atoms.begin(), atoms.end());
831 atoms.resize(std::distance(atoms.begin(), ip));
File containing functions used specific to bulk topological unit matching (TUM) criterion.
std::vector< int > rings
type of the cage : can be DDC or HC
std::vector< T > box
Number of atoms.
std::vector< std::vector< int > > neighbourListByIndex(molSys::PointCloud< molSys::Point< double >, double > *yCloud, std::vector< std::vector< int > > nList)
std::vector< std::vector< int > > neighListO(double rcutoff, molSys::PointCloud< molSys::Point< double >, double > *yCloud, int typeI)
std::vector< std::vector< int > > ringNetwork(std::vector< std::vector< int > > nList, int maxDepth)
int getAtomTypesTopoBulk(std::vector< std::vector< int > > rings, std::vector< ring::strucType > ringType, std::vector< cage::iceType > *atomTypes)
int findBulkPrisms(std::vector< std::vector< int > > rings, std::vector< ring::strucType > *ringType, std::vector< std::vector< int > > nList, molSys::PointCloud< molSys::Point< double >, double > *yCloud, std::vector< double > *rmsdPerAtom, double heightCutoff=8)
Find out which rings are prisms.
int getStrucNumbers(std::vector< ring::strucType > ringType, std::vector< cage::Cage > cageList, int *numHC, int *numDDC, int *mixedRings, int *prismaticRings, int *basalRings)
Determines the number of HCs, DDCs, Mixed rings, prismatic and basal rings.
std::vector< int > findDDC(std::vector< std::vector< int > > rings, std::vector< strucType > *ringType, std::vector< int > listHC, std::vector< cage::Cage > *cageList)
std::vector< int > findHC(std::vector< std::vector< int > > rings, std::vector< strucType > *ringType, std::vector< std::vector< int > > nList, std::vector< cage::Cage > *cageList)
std::vector< int > findMixedRings(std::vector< std::vector< int > > rings, std::vector< strucType > *ringType, std::vector< int > *listDDC, std::vector< int > *listHC)
Eigen::MatrixXd buildRefHC(std::string fileName)
Build a reference Hexagonal cage, reading in from a template XYZ file.
std::vector< cage::Cage > topoBulkCriteria(std::string path, std::vector< std::vector< int > > rings, std::vector< std::vector< int > > nList, molSys::PointCloud< molSys::Point< double >, double > *yCloud, int firstFrame, int *numHC, int *numDDC, std::vector< ring::strucType > *ringType)
int topoUnitMatchingBulk(std::string path, std::vector< std::vector< int > > rings, std::vector< std::vector< int > > nList, molSys::PointCloud< molSys::Point< double >, double > *yCloud, int firstFrame, bool printClusters, bool onlyTetrahedral)
int averageRMSDatom(std::vector< double > *rmsdPerAtom, std::vector< int > *noOfCommonAtoms)
Average the RMSD per atom.
std::vector< int > atomsFromCages(std::vector< std::vector< int > > rings, std::vector< cage::Cage > cageList, std::vector< int > clusterCages)
Gets the atoms in the cages of a given cluster.
int shapeMatchDDC(molSys::PointCloud< molSys::Point< double >, double > *yCloud, const Eigen::MatrixXd &refPoints, std::vector< cage::Cage > cageList, int cageIndex, std::vector< std::vector< int > > rings, std::vector< double > *quat, double *rmsd)
Shape-matching for a target DDC.
int updateRMSDatom(std::vector< std::vector< int > > rings, cage::Cage cageUnit, double rmsd, std::vector< double > *rmsdPerAtom, std::vector< int > *noOfCommonAtoms, std::vector< cage::iceType > atomTypes)
int shapeMatchHC(molSys::PointCloud< molSys::Point< double >, double > *yCloud, const Eigen::MatrixXd &refPoints, cage::Cage cageUnit, std::vector< std::vector< int > > rings, std::vector< std::vector< int > > nList, std::vector< double > *quat, double *rmsd)
Shape-matching for a target HC.
int clearRingList(std::vector< std::vector< int > > &rings)
Erases memory for a vector of vectors for a list of rings.
Eigen::MatrixXd buildRefDDC(std::string fileName)
Build a reference Double-Diamond cage, reading in from a template XYZ file.
bool hasCommonElements(std::vector< int > ring1, std::vector< int > ring2)
std::vector< std::vector< int > > getSingleRingSize(std::vector< std::vector< int > > rings, int ringSize)
Returns a vector of vectors of rings of a single size.
int clusterCages(molSys::PointCloud< molSys::Point< double >, double > *yCloud, std::string path, std::vector< std::vector< int > > rings, std::vector< cage::Cage > cageList, int numHC, int numDDC)
molSys::PointCloud< molSys::Point< double >, double > readXYZ(std::string filename)
Function for reading in atom coordinates from an XYZ file.
int hornAbsOrientation(const Eigen::MatrixXd &refPoints, const Eigen::MatrixXd &targetPoints, std::vector< double > *quat, double *rmsd, std::vector< double > *rmsdList, double *scale)
Get the absolute orientation using Horn's algorithm (with quaternions)
Eigen::MatrixXd changeDiaCageOrder(molSys::PointCloud< molSys::Point< double >, double > *yCloud, std::vector< int > ddcOrder, int startingIndex=0)
std::vector< int > relOrderDDC(int index, std::vector< std::vector< int > > rings, std::vector< cage::Cage > cageList)
Matches the order of the basal rings of an DDC or a potential HC.
Eigen::MatrixXd changeHexCageOrder(molSys::PointCloud< molSys::Point< double >, double > *yCloud, std::vector< int > basal1, std::vector< int > basal2, int startingIndex=0)
int relOrderHC(molSys::PointCloud< molSys::Point< double >, double > *yCloud, std::vector< int > basal1, std::vector< int > basal2, std::vector< std::vector< int > > nList, std::vector< int > *matchedBasal1, std::vector< int > *matchedBasal2)
Matches the order of the basal rings of an HC or a potential HC.
int makePath(const std::string &path)
int writeTopoBulkData(std::string path, int currentFrame, int numHC, int numDDC, int mixedRings, int basalRings, int prismaticRings, int firstFrame)
int writeLAMMPSdataTopoBulk(molSys::PointCloud< molSys::Point< double >, double > *yCloud, std::vector< std::vector< int > > nList, std::vector< cage::iceType > atomTypes, std::string path, bool bondsBetweenDummy=false)
int writeLAMMPSdumpCages(molSys::PointCloud< molSys::Point< double >, double > *yCloud, std::vector< double > rmsdPerAtom, std::vector< int > atomTypes, std::string path, int firstFrame)
Write out a LAMMPS dump file containing the RMSD per atom for bulk ice.
int writeXYZcluster(std::string path, molSys::PointCloud< molSys::Point< double >, double > *yCloud, std::vector< int > atoms, int clusterID, cage::cageType type)
Function for writing out the XYZ files for each cluster.
This contains a cage, with the constituent rings.
This contains a collection of points; contains information for a particular frame.
This contains per-particle information.
Copy