topo_one_dim.cpp
Go to the documentation of this file.
1 //-----------------------------------------------------------------------------------
2 // d-SEAMS - Deferred Structural Elucidation Analysis for Molecular Simulations
3 //
4 // Copyright (c) 2018--present d-SEAMS core team
5 //
6 // This program is free software: you can redistribute it and/or modify
7 // it under the terms of the MIT License as published by
8 // the Open Source Initiative.
9 //
10 // A copy of the MIT License is included in the LICENSE file of this repository.
11 // You should have received a copy of the MIT License along with this program.
12 // If not, see <https://opensource.org/licenses/MIT>.
13 //-----------------------------------------------------------------------------------
14 
15 #include <topo_one_dim.hpp>
16 
17 // -----------------------------------------------------------------------------------------------------
18 // PRISM ALGORITHMS
19 // -----------------------------------------------------------------------------------------------------
20 
52  molSys::PointCloud<molSys::Point<double>, double> *yCloud, int maxDepth,
53  int *atomID, int firstFrame, int currentFrame, bool doShapeMatching) {
54  //
56  ringsOneType; // Vector of vectors of rings of a single size
57  std::vector<int> listPrism; // Vector for ring indices of n-sided prism
59  ringType; // This vector will have a value for each ring inside
61  atomState; // This vector will have a value for each atom, depending on
62  // the ring type
63  int nPerfectPrisms; // Number of perfect prisms of each type
64  int nImperfectPrisms; // Number of deformed prisms of each type
65  std::vector<int> nPrismList; // Vector of the values of the number of perfect
66  // prisms for a particular frame
67  std::vector<int> nDefPrismList; // Vector of the values of the number of
68  // deformed prisms for a particular frame
70  heightPercent; // Height percent for a particular n and frame
72  atomTypes; // contains int values for each prism type considered
73  double avgPrismHeight = 2.845; // A value of 2.7-2.85 Angstrom is reasonable
74  // Qualifier for the RMSD per atom:
75  std::vector<double> rmsdPerAtom;
76  // -------------------------------------------------------------------------------
77  // Init
78  nPrismList.resize(
79  maxDepth -
80  2); // Has a value for every value of ringSize from 3, upto maxDepth
81  nDefPrismList.resize(maxDepth - 2);
82  heightPercent.resize(maxDepth - 2);
83  // The atomTypes vector is the same size as the pointCloud atoms
84  atomTypes.resize(yCloud->nop, 1); // The dummy or unclassified value is 1
85  // The rmsdPerAtom vector is the same size as the pointCloud atoms and has an
86  // RMSD value for every atom
87  rmsdPerAtom.resize(yCloud->nop, -1);
88  // Resize the atom state vector
89  atomState.resize(yCloud->nop); // Dummy or unclassified
90  // -------------------------------------------------------------------------------
91  // Run this loop for rings of sizes upto maxDepth
92  // The smallest possible ring is of size 3
93  for (int ringSize = 3; ringSize <= maxDepth; ringSize++) {
94  // Clear ringsOneType
95  ring::clearRingList(ringsOneType);
96  // Get rings of the current ring size
97  ringsOneType = ring::getSingleRingSize(rings, ringSize);
98  //
99  // Continue if there are zero rings of ringSize
100  if (ringsOneType.size() == 0) {
101  nPrismList[ringSize - 3] = 0; // Update the number of prisms
102  nDefPrismList[ringSize - 3] = 0; // Update the number of deformed prisms
103  heightPercent[ringSize - 3] = 0.0; // Update the height%
104  continue;
105  } // skip if there are no rings
106  //
107  // -------------
108  // Init of variables specific to ringSize prisms
109  listPrism.resize(0);
110  ringType.resize(0);
111  nPerfectPrisms = 0;
112  nImperfectPrisms = 0;
113  ringType.resize(
114  ringsOneType.size()); // Has a value for each ring. init to zero.
115  // -------------
116  // Now that you have rings of a certain size:
117  // Find prisms, saving the ring indices to listPrism
118  listPrism = ring::findPrisms(ringsOneType, &ringType, &nPerfectPrisms,
119  &nImperfectPrisms, nList, yCloud, &rmsdPerAtom,
120  doShapeMatching);
121  // -------------
122  nPrismList[ringSize - 3] =
123  nPerfectPrisms; // Update the number of perfect prisms
124  nDefPrismList[ringSize - 3] =
125  nImperfectPrisms; // Update the number of defromed prisms
126  int totalPrisms = nPerfectPrisms + nImperfectPrisms;
127  // Update the height% for the phase
128  heightPercent[ringSize - 3] =
129  topoparam::normHeightPercent(yCloud, totalPrisms, avgPrismHeight);
130  // Continue if there are no prism units
131  if (nPerfectPrisms + nImperfectPrisms == 0) {
132  continue;
133  } // skip for no prisms
134  // Do a bunch of write-outs and calculations
135  // TODO: Write out each individual prism as data files (maybe with an
136  // option)
137  // Get the atom types for a particular prism type
138  ring::assignPrismType(ringsOneType, listPrism, ringSize, ringType,
139  &atomTypes, &atomState);
140  // -------------
141  } // end of loop through every possible ringSize
142 
143  // Calculate the height%
144 
145  // Write out the prism information
146  sout::writePrismNum(path, nPrismList, nDefPrismList, heightPercent, maxDepth,
147  yCloud->currentFrame, firstFrame);
148  // Reassign the prism block types for the deformed prisms
149  if (doShapeMatching) {
150  ring::deformedPrismTypes(atomState, &atomTypes, maxDepth);
151  } // reassign prism block types for deformed prisms
152 
153  // Get rid of translations along the axial dimension
154  ring::rmAxialTranslations(yCloud, atomID, firstFrame, currentFrame);
155 
156  // Write out dump files with the RMSD per atom for the shape-matching
157  // criterion
158  if (doShapeMatching) {
159  sout::writeLAMMPSdumpINT(yCloud, rmsdPerAtom, atomTypes, maxDepth, path);
160  } // reassign prism block types for deformed prisms
161 
162  // Write out the lammps data file for the particular frame
163  sout::writeLAMMPSdataAllPrisms(yCloud, nList, atomTypes, maxDepth, path,
164  doShapeMatching);
165 
166  return 0;
167 }
168 
194  std::vector<ring::strucType> *ringType, int *nPerfectPrisms,
195  int *nImperfectPrisms, std::vector<std::vector<int>> nList,
196  molSys::PointCloud<molSys::Point<double>, double> *yCloud,
197  std::vector<double> *rmsdPerAtom, bool doShapeMatching) {
198  std::vector<int> listPrism;
199  int totalRingNum = rings.size(); // Total number of rings
200  std::vector<int> basal1; // First basal ring
201  std::vector<int> basal2; // Second basal ring
202  bool cond1, cond2; // Conditions for rings to be basal (true) or not (false)
203  bool relaxedCond; // Condition so that at least one bond exists between the
204  // two basal rings
205  bool isAxialPair; // Basal rings should be parallel in one dimension to
206  // prevent overcounting
207  int ringSize = rings[0].size(); // Number of nodes in each ring
208  *nImperfectPrisms = 0; // Number of undeformed prisms
209  *nPerfectPrisms = 0; // Number of undeformed prisms
210  // Matrix for the reference ring for a given ringSize.
211  Eigen::MatrixXd refPointSet(ringSize, 3);
212 
213  // Get the reference ring point set for a given ring size.
214  // Get the axial dimension
215  int axialDim = std::max_element(yCloud->box.begin(), yCloud->box.end()) -
216  yCloud->box.begin();
217  refPointSet = pntToPnt::getPointSetRefRing(ringSize, axialDim);
218  //
219 
220  // Two loops through all the rings are required to find pairs of basal rings
221  for (int iring = 0; iring < totalRingNum - 1; iring++) {
222  cond1 = false;
223  cond2 = false;
224  basal1 = rings[iring]; // Assign iring to basal1
225  // Loop through the other rings to get a pair
226  for (int jring = iring + 1; jring < totalRingNum; jring++) {
227  basal2 = rings[jring]; // Assign jring to basal2
228  // ------------
229  // Put extra check for axial basal rings if shapeMatching is being done
230  if (doShapeMatching == true || ringSize == 4) {
231  isAxialPair = false; // init
232  isAxialPair =
233  ring::discardExtraTetragonBlocks(&basal1, &basal2, yCloud);
234  if (isAxialPair == false) {
235  continue;
236  }
237  } // end of check for tetragonal prism blocks
238  // ------------
239  // Step one: Check to see if basal1 and basal2 have common
240  // elements or not. If they don't, then they cannot be basal rings
241  cond1 = ring::hasCommonElements(basal1, basal2);
242  if (cond1 == true) {
243  continue;
244  }
245  // -----------
246  // Step two and three: One of the elements of basal2 must be the nearest
247  // neighbour of the first (index0; l1) If m_k is the nearest neighbour of
248  // l1, m_{k+1} ... m_{k+(n-1)} must be neighbours of l_i+1 etc or l_i-1
249  cond2 = ring::basalPrismConditions(nList, &basal1, &basal2);
250  // If cond2 is false, the strict criteria for prisms has not been met
251  if (cond2 == false) {
252  // Skip if shape-matching is not desired
253  if (!doShapeMatching) {
254  continue;
255  } // shape-matching not desired
256  //
257  // If shape-matching is to be done:
258  // Check for the reduced criteria fulfilment
259  relaxedCond = ring::relaxedPrismConditions(nList, &basal1, &basal2);
260  // Skip if relaxed criteria are not met
261  if (relaxedCond == false) {
262  continue;
263  } // end of skipping if the prisms do not fulfil relaxed criteria
264 
265  // Do shape matching here
266  bool isDeformedPrism = match::matchPrism(
267  yCloud, nList, refPointSet, &basal1, &basal2, rmsdPerAtom, false);
268 
269  // Success! The rings are basal rings of a deformed prism!
270  if (isDeformedPrism) {
271  //
272  // // Update the number of prism blocks
273  *nImperfectPrisms += 1;
274  // Update iring
275  if ((*ringType)[iring] == ring::strucType::unclassified) {
276  (*ringType)[iring] = ring::strucType::deformedPrism;
277  listPrism.push_back(iring);
278  } else if ((*ringType)[iring] == ring::strucType::Prism) {
279  (*ringType)[iring] = ring::strucType::mixedPrismRing;
280  } // if it is deformed
281  // Update jring
282  if ((*ringType)[jring] == ring::strucType::unclassified) {
283  (*ringType)[jring] = ring::strucType::deformedPrism;
284  listPrism.push_back(jring);
285  } else if ((*ringType)[jring] == ring::strucType::Prism) {
286  (*ringType)[jring] = ring::strucType::mixedPrismRing;
287  } // if it is deformed
288  } // end of update of ring types
289 
290  // // Write outs
291  // // Now write out axial basal rings
292  // sout::writeBasalRingsPrism(&basal1, &basal2, nDeformedPrisms, nList,
293  // yCloud, true);
294  } // end of reduced criteria
295  // Strict criteria
296  else {
297  // Update the number of prism blocks
298  *nPerfectPrisms += 1;
299  // Update iring
300  if ((*ringType)[iring] == ring::strucType::unclassified) {
301  (*ringType)[iring] = ring::strucType::Prism;
302  listPrism.push_back(iring);
303  } else if ((*ringType)[iring] == ring::strucType::deformedPrism) {
304  (*ringType)[iring] = ring::strucType::mixedPrismRing;
305  } // if it is deformed
306  // Update jring
307  if ((*ringType)[jring] == ring::strucType::unclassified) {
308  (*ringType)[jring] = ring::strucType::Prism;
309  listPrism.push_back(jring);
310  } else if ((*ringType)[jring] == ring::strucType::deformedPrism) {
311  (*ringType)[jring] = ring::strucType::mixedPrismRing;
312  } // if it is deformed
313  //
314  // Shape-matching to get the RMSD (if shape-matching is desired)
315  if (doShapeMatching) {
316  bool isKnownPrism = match::matchPrism(
317  yCloud, nList, refPointSet, &basal1, &basal2, rmsdPerAtom, true);
318  } // end of shape-matching to get rmsd
319  //
320  // // Now write out axial basal rings for convex hull calculations
321  // sout::writePrisms(&basal1, &basal2, *nPrisms, yCloud);
322  // // Write out prisms for shape-matching
323  // sout::writeBasalRingsPrism(&basal1, &basal2, *nPrisms, nList, yCloud,
324  // false);
325  // -----------
326  } // end of strict criteria
327 
328  } // end of loop through rest of the rings to get the second basal ring
329  } // end of loop through all rings for first basal ring
330 
331  sort(listPrism.begin(), listPrism.end());
332  auto ip = std::unique(listPrism.begin(), listPrism.end());
333  // Resize peripheral rings to remove undefined terms
334  listPrism.resize(std::distance(listPrism.begin(), ip));
335 
336  return listPrism;
337 }
338 
352  std::vector<int> *basal1,
353  std::vector<int> *basal2) {
354  int l1 = (*basal1)[0]; // first element of basal1 ring
355  int ringSize =
356  (*basal1).size(); // Size of the ring; each ring contains n elements
357  int m_k; // Atom ID of element in basal2
358  bool l1_neighbour; // m_k is a neighbour of l1(true) or not (false)
359 
360  // isNeighbour is initialized to false for all basal2 elements; indication if
361  // basal2 elements are neighbours of basal1
362  std::vector<bool> isNeighbour(ringSize, false);
363  int kIndex; // m_k index
364  int lAtomID; // atomID of the current element of basal1
365  int kAtomID; // atomID of the current element of basal2
366 
367  // ---------------------------------------------
368  // COMPARISON OF basal2 ELEMENTS WITH l1
369  for (int k = 0; k < ringSize; k++) {
370  l1_neighbour = false;
371  m_k = (*basal2)[k];
372  // =================================
373  // Checking to seee if m_k is be a neighbour of l1
374  // Find m_k inside l1 neighbour list
375  auto it = std::find(nList[l1].begin() + 1, nList[l1].end(), m_k);
376 
377  // If the element has been found, for l1
378  if (it != nList[l1].end()) {
379  l1_neighbour = true;
380  kIndex = k;
381  break;
382  }
383  } // l1 is a neighbour of m_k
384 
385  // If there is no nearest neighbour, then the two rings are not part of the
386  // prism
387  if (!l1_neighbour) {
388  return false;
389  }
390 
391  // ---------------------------------------------
392  // NEIGHBOURS of basal1 in basal2
393  isNeighbour[kIndex] = true;
394 
395  // All elements of basal1 must be neighbours of basal2
396  for (int i = 1; i < ringSize; i++) {
397  lAtomID = (*basal1)[i]; // element of basal1 ring
398  for (int k = 0; k < ringSize; k++) {
399  // Skip if already a neighbour
400  if (isNeighbour[k]) {
401  continue;
402  }
403  // Get the comparison basal2 element
404  kAtomID = (*basal2)[k]; // element of basal2 ring;
405 
406  // Checking to see if kAtomID is a neighbour of lAtomID
407  // Find kAtomID inside lAtomID neighbour list
408  auto it1 =
409  std::find(nList[lAtomID].begin() + 1, nList[lAtomID].end(), kAtomID);
410 
411  // If the element has been found, for l1
412  if (it1 != nList[lAtomID].end()) {
413  isNeighbour[k] = true;
414  }
415  } // Loop through basal2
416  } // end of check for neighbours of basal1
417 
418  // ---------------------------------------------
419 
420  // They should all be neighbours
421  for (int k = 0; k < ringSize; k++) {
422  // Check to see if any element is false
423  if (!isNeighbour[k]) {
424  return false;
425  }
426  }
427 
428  // Everything works out!
429  return true;
430 }
431 
437  std::vector<int> *basal1,
438  std::vector<int> *basal2) {
439  int ringSize =
440  (*basal1).size(); // Size of the ring; each ring contains n elements
441  int m_k; // Atom ID of element in basal2
442  bool isNeighbour = false; // This is true if there is at least one bond
443  // between the basal rings
444  int l_k; // Atom ID of element in basal1
445 
446  // ---------------------------------------------
447  // COMPARISON OF basal2 ELEMENTS (m_k) WITH basal1 ELEMENTS (l_k)
448  // Loop through all the elements of basal1
449  for (int l = 0; l < ringSize; l++) {
450  l_k = (*basal1)[l];
451  // Search for the nearest neighbour of l_k in basal2
452  // Loop through basal2 elements
453  for (int m = 0; m < ringSize; m++) {
454  m_k = (*basal2)[m];
455  // Find m_k inside l_k neighbour list
456  auto it = std::find(nList[l_k].begin() + 1, nList[l_k].end(), m_k);
457 
458  // If the element has been found, for l1
459  if (it != nList[l_k].end()) {
460  isNeighbour = true;
461  break;
462  } // found element
463  } // end of loop through all the elements of basal2
464 
465  // If a neighbour has been found then
466  if (isNeighbour) {
467  return true;
468  }
469  } // end of loop through all the elements of basal1
470 
471  // If a neighbour has not been found, return false
472  return false;
473 }
474 
488  std::vector<int> *basal1, std::vector<int> *basal2,
489  molSys::PointCloud<molSys::Point<double>, double> *yCloud) {
490  int ringSize =
491  (*basal1).size(); // Size of the ring; each ring contains n elements
492  int iatomIndex,
493  jatomIndex; // Indices of the elements in basal1 and basal2 respectively
494  double r_i, r_j; // Coordinates in the axial dimension of iatom and jatom of
495  // basal1 and basal2 respectively
496  int axialDim; // 0 for x, 1 for y and 2 for z dimensions respectively
497  // Variables for getting the projected area
498  bool axialBasal1, axialBasal2; // bools for checking if basal1 and basal2 are
499  // axial (true) respectively
500  double areaXY, areaXZ,
501  areaYZ; // Projected area on the XY, XZ and YZ planes respectively
502  // ----------------------------------------
503  // Find the axial dimension for a quasi-one-dimensional ice nanotube
504  // The axial dimension will have the largest box length
505  // Index -> axial dimension
506  // 0 -> x dim
507  // 1 -> y dim
508  // 2 -> z dim
509  axialDim = std::max_element(yCloud->box.begin(), yCloud->box.end()) -
510  yCloud->box.begin();
511  // ----------------------------------------
512  // Calculate projected area onto the XY, YZ and XZ planes for basal1
513  axialBasal1 = false; // Init to false
514  axialBasal2 = false; // Init
515 
516  // Init the projected area
517  areaXY = 0.0;
518  areaXZ = 0.0;
519  areaYZ = 0.0;
520 
521  jatomIndex = (*basal1)[0];
522 
523  // All points except the first pair
524  for (int k = 1; k < ringSize; k++) {
525  iatomIndex = (*basal1)[k]; // Current vertex
526 
527  // Add to the polygon area
528  // ------
529  // XY plane
530  areaXY += (yCloud->pts[jatomIndex].x + yCloud->pts[iatomIndex].x) *
531  (yCloud->pts[jatomIndex].y - yCloud->pts[iatomIndex].y);
532  // ------
533  // XZ plane
534  areaXZ += (yCloud->pts[jatomIndex].x + yCloud->pts[iatomIndex].x) *
535  (yCloud->pts[jatomIndex].z - yCloud->pts[iatomIndex].z);
536  // ------
537  // YZ plane
538  areaYZ += (yCloud->pts[jatomIndex].y + yCloud->pts[iatomIndex].y) *
539  (yCloud->pts[jatomIndex].z - yCloud->pts[iatomIndex].z);
540  // ------
541  jatomIndex = iatomIndex;
542  }
543 
544  // Closure point
545  iatomIndex = (*basal1)[0];
546  // ------
547  // XY plane
548  areaXY += (yCloud->pts[jatomIndex].x + yCloud->pts[iatomIndex].x) *
549  (yCloud->pts[jatomIndex].y - yCloud->pts[iatomIndex].y);
550  // ------
551  // XZ plane
552  areaXZ += (yCloud->pts[jatomIndex].x + yCloud->pts[iatomIndex].x) *
553  (yCloud->pts[jatomIndex].z - yCloud->pts[iatomIndex].z);
554  // ------
555  // YZ plane
556  areaYZ += (yCloud->pts[jatomIndex].y + yCloud->pts[iatomIndex].y) *
557  (yCloud->pts[jatomIndex].z - yCloud->pts[iatomIndex].z);
558  // ------
559  // The actual projected area is half of this
560  areaXY *= 0.5;
561  areaXZ *= 0.5;
562  areaYZ *= 0.5;
563  // Get the absolute value
564  areaXY = fabs(areaXY);
565  areaXZ = fabs(areaXZ);
566  areaYZ = fabs(areaYZ);
567 
568  // If the axial dimension is x, y, or z:
569  // then the maximum basal area should be in the YZ, XZ and XY dimensions
570  // respectively
571  // x dim
572  if (axialDim == 0) {
573  if (areaYZ > areaXY && areaYZ > areaXZ) {
574  axialBasal1 = true;
575  } // end of check for axial ring for basal1
576  } // x dim
577  // y dim
578  else if (axialDim == 1) {
579  if (areaXZ > areaXY && areaXZ > areaYZ) {
580  axialBasal1 = true;
581  } // end of check for axial ring for basal1
582  } // x dim
583  // z dim
584  else if (axialDim == 2) {
585  if (areaXY > areaXZ && areaXY > areaYZ) {
586  axialBasal1 = true;
587  } // end of check for axial ring for basal1
588  } // x dim
589  else {
590  std::cerr << "Could not find the axial dimension.\n";
591  return false;
592  }
593  // ----------------------------------------
594  // Calculate projected area onto the XY, YZ and XZ planes for basal2
595 
596  // Init the projected area
597  areaXY = 0.0;
598  areaXZ = 0.0;
599  areaYZ = 0.0;
600 
601  jatomIndex = (*basal2)[0];
602 
603  // All points except the first pair
604  for (int k = 1; k < ringSize; k++) {
605  iatomIndex = (*basal2)[k]; // Current vertex
606 
607  // Add to the polygon area
608  // ------
609  // XY plane
610  areaXY += (yCloud->pts[jatomIndex].x + yCloud->pts[iatomIndex].x) *
611  (yCloud->pts[jatomIndex].y - yCloud->pts[iatomIndex].y);
612  // ------
613  // XZ plane
614  areaXZ += (yCloud->pts[jatomIndex].x + yCloud->pts[iatomIndex].x) *
615  (yCloud->pts[jatomIndex].z - yCloud->pts[iatomIndex].z);
616  // ------
617  // YZ plane
618  areaYZ += (yCloud->pts[jatomIndex].y + yCloud->pts[iatomIndex].y) *
619  (yCloud->pts[jatomIndex].z - yCloud->pts[iatomIndex].z);
620  // ------
621  jatomIndex = iatomIndex;
622  }
623 
624  // Closure point
625  iatomIndex = (*basal2)[0];
626  // ------
627  // XY plane
628  areaXY += (yCloud->pts[jatomIndex].x + yCloud->pts[iatomIndex].x) *
629  (yCloud->pts[jatomIndex].y - yCloud->pts[iatomIndex].y);
630  // ------
631  // XZ plane
632  areaXZ += (yCloud->pts[jatomIndex].x + yCloud->pts[iatomIndex].x) *
633  (yCloud->pts[jatomIndex].z - yCloud->pts[iatomIndex].z);
634  // ------
635  // YZ plane
636  areaYZ += (yCloud->pts[jatomIndex].y + yCloud->pts[iatomIndex].y) *
637  (yCloud->pts[jatomIndex].z - yCloud->pts[iatomIndex].z);
638  // ------
639  // The actual projected area is half of this
640  areaXY *= 0.5;
641  areaXZ *= 0.5;
642  areaYZ *= 0.5;
643  // Get the absolute value
644  areaXY = fabs(areaXY);
645  areaXZ = fabs(areaXZ);
646  areaYZ = fabs(areaYZ);
647 
648  // Check if xy projected area is the greatest
649  // If the axial dimension is x, y, or z:
650  // then the maximum basal area should be in the YZ, XZ and XY dimensions
651  // respectively
652  // x dim
653  if (axialDim == 0) {
654  if (areaYZ > areaXY && areaYZ > areaXZ) {
655  axialBasal2 = true;
656  } // end of check for axial ring for basal1
657  } // x dim
658  // y dim
659  else if (axialDim == 1) {
660  if (areaXZ > areaXY && areaXZ > areaYZ) {
661  axialBasal2 = true;
662  } // end of check for axial ring for basal1
663  } // x dim
664  // z dim
665  else if (axialDim == 2) {
666  if (areaXY > areaXZ && areaXY > areaYZ) {
667  axialBasal2 = true;
668  } // end of check for axial ring for basal1
669  } // x dim
670  else {
671  std::cerr << "Could not find the axial dimension.\n";
672  return false;
673  }
674  // ----------------------------------------
675 
676  // Now check if basal1 and basal2 are axial or not
677  if (axialBasal1 == true && axialBasal2 == true) {
678  return true;
679  } else {
680  return false;
681  } // Check for basal1 and basal2
682 }
683 
697  std::vector<int> listPrism, int ringSize,
699  std::vector<int> *atomTypes,
700  std::vector<ring::strucType> *atomState) {
701  // Every value in listPrism corresponds to an index in rings.
702  // Every ring contains atom indices, corresponding to the indices (not atom
703  // IDs) in rings
704  int iring; // Index of current ring
705  int iatom; // Index of current atom
706  ring::strucType currentState; // Current state of the ring being filled
707 
708  // Dummy value corresponds to a value of 1.
709  // Each value is initialized to the value of 1.
710 
711  // Loop through every ring in rings
712  for (int i = 0; i < listPrism.size(); i++) {
713  iring = listPrism[i];
714  // Get the state of all atoms in the ring
715  currentState = ringType[iring];
716  // Loop through every element in iring
717  for (int j = 0; j < ringSize; j++) {
718  iatom = rings[iring][j]; // Atom index
719  // Update the atom type
720  (*atomTypes)[iatom] = ringSize;
721  //
722  // Update the state of the atom
723  if ((*atomState)[iatom] == ring::strucType::unclassified) {
724  (*atomState)[iatom] = currentState;
725  } // Update the unclassified atom
726  else {
727  if ((*atomState)[iatom] != currentState) {
728  // For mixed, there is a preference
729  if (currentState == ring::strucType::mixedPrismRing) {
730  (*atomState)[iatom] = currentState;
731  } // fill
732  else if ((*atomState)[iatom] == ring::strucType::deformedPrism &&
733  currentState == ring::strucType::Prism) {
734  (*atomState)[iatom] = ring::strucType::mixedPrismRing;
735  } else if ((*atomState)[iatom] == ring::strucType::Prism &&
736  currentState == ring::strucType::deformedPrism) {
737  (*atomState)[iatom] = ring::strucType::mixedPrismRing;
738  }
739  } //
740  } // already filled?
741 
742  } // end of loop through every atom in iring
743  } // end of loop through every ring
744 
745  return 0;
746 } // end of function
747 
748 // Get the atom type values for deformed prisms
760  std::vector<int> *atomTypes, int maxDepth) {
761  //
762  int nop = atomState.size(); // Number of particles
763 
764  // Loop through all particles
765  for (int iatom = 0; iatom < nop; iatom++) {
766  // Check the atom state
767  // Deformed
768  if (atomState[iatom] == ring::strucType::deformedPrism) {
769  (*atomTypes)[iatom] += maxDepth - 2;
770  } // type for a deformed prism atom
771  else if (atomState[iatom] == ring::strucType::mixedPrismRing) {
772  (*atomTypes)[iatom] = 2;
773  } // type for a mixed prism ring
774  } // end of reassignation
775 
776  //
777  return 0;
778 }
779 
780 // Shift the entire ice nanotube and remove axial translations
792  molSys::PointCloud<molSys::Point<double>, double> *yCloud, int *atomID,
793  int firstFrame, int currentFrame) {
794  //
795  int atomIndex; // Index of the atom to be shifted first
796  double shiftDistance; // Value by which all atoms will be shifted
797  double distFrmBoundary; // Distance from either the upper or lower boundary
798  // along the axial dimension
799  // Get the axial dimension
800  int axialDim = std::max_element(yCloud->box.begin(), yCloud->box.end()) -
801  yCloud->box.begin();
802  // Check: (default z)
803  if (axialDim < 0 || axialDim > 2) {
804  axialDim = 2;
805  } // end of check to set the axial dimension
806  // Lower and higher limits of the box in the axial dimension
807  double boxLowAxial = yCloud->boxLow[axialDim];
808  double boxHighAxial = boxLowAxial + yCloud->box[axialDim];
809  //
810  // -----------------------------------
811  // Update the atomID of the 'first' atom in yCloud if the current frame is the
812  // first frame.
813  // Get the atom index of the first atom to be shifted down or up
814  if (currentFrame == firstFrame) {
815  *atomID = yCloud->pts[0].atomID;
816  atomIndex = 0;
817  } // end of update of the atom ID to be shifted for the first frame
818  else {
819  //
820  int iatomID = *atomID; // Atom ID of the first particle to be moved down
821  // Find the index given the atom ID
822  auto it = yCloud->idIndexMap.find(iatomID);
823 
824  if (it != yCloud->idIndexMap.end()) {
825  atomIndex = it->second;
826  } // found jatom
827  else {
828  std::cerr << "Lost atoms.\n";
829  return 1;
830  } // error handling
831  //
832  } // Update of atomIndex for all other frames
833  // -----------------------------------
834  // Calculate the distance by which the atoms must be shifted (negative value)
835  if (axialDim == 0) {
836  shiftDistance = boxLowAxial - yCloud->pts[atomIndex].x;
837  } // x coordinate
838  else if (axialDim == 1) {
839  shiftDistance = boxLowAxial - yCloud->pts[atomIndex].y;
840  } // y coordinate
841  else {
842  shiftDistance = boxLowAxial - yCloud->pts[atomIndex].z;
843  } // z coordinate
844  // -----------------------------------
845  // Shift all the particles
846  for (int iatom = 0; iatom < yCloud->nop; iatom++) {
847  // Shift the particles along the axial dimension only
848  if (axialDim == 0) {
849  yCloud->pts[iatom].x += shiftDistance;
850  // Wrap-around
851  if (yCloud->pts[iatom].x < boxLowAxial) {
852  distFrmBoundary = boxLowAxial - yCloud->pts[iatom].x; // positive value
853  yCloud->pts[iatom].x = boxHighAxial - distFrmBoundary;
854  } // end of wrap-around
855  } // x coordinate
856  else if (axialDim == 1) {
857  yCloud->pts[iatom].y += shiftDistance;
858  // Wrap-around
859  if (yCloud->pts[iatom].y < boxLowAxial) {
860  distFrmBoundary = boxLowAxial - yCloud->pts[iatom].y; // positive value
861  yCloud->pts[iatom].y = boxHighAxial - distFrmBoundary;
862  } // end of wrap-around
863  } // y coordinate
864  else {
865  yCloud->pts[iatom].z += shiftDistance;
866  // Wrap-around
867  if (yCloud->pts[iatom].z < boxLowAxial) {
868  distFrmBoundary = boxLowAxial - yCloud->pts[iatom].z; // positive value
869  yCloud->pts[iatom].z = boxHighAxial - distFrmBoundary;
870  } // end of wrap-around
871  } // z coordinate
872  } // end of shifting all paritcles downward
873  // -----------------------------------
874  return 0;
875 } // end of function
T begin(T... args)
T distance(T... args)
T end(T... args)
T find(T... args)
bool relaxedPrismConditions(std::vector< std::vector< int >> nList, std::vector< int > *basal1, std::vector< int > *basal2)
int prismAnalysis(std::string path, std::vector< std::vector< int >> rings, std::vector< std::vector< int >> nList, molSys::PointCloud< molSys::Point< double >, double > *yCloud, int maxDepth, int *atomID, int firstFrame, int currentFrame, bool doShapeMatching=false)
std::vector< int > findPrisms(std::vector< std::vector< int >> rings, std::vector< strucType > *ringType, int *nPerfectPrisms, int *nImperfectPrisms, std::vector< std::vector< int >> nList, molSys::PointCloud< molSys::Point< double >, double > *yCloud, std::vector< double > *rmsdPerAtom, bool doShapeMatching=false)
strucType
Definition: ring.hpp:116
int deformedPrismTypes(std::vector< ring::strucType > atomState, std::vector< int > *atomTypes, int maxDepth)
Get the atom type values for deformed prisms.
int clearRingList(std::vector< std::vector< int >> &rings)
Erases memory for a vector of vectors for a list of rings.
Definition: ring.cpp:22
bool discardExtraTetragonBlocks(std::vector< int > *basal1, std::vector< int > *basal2, molSys::PointCloud< molSys::Point< double >, double > *yCloud)
int rmAxialTranslations(molSys::PointCloud< molSys::Point< double >, double > *yCloud, int *atomID, int firstFrame, int currentFrame)
Shift the entire ice nanotube and remove axial translations.
int assignPrismType(std::vector< std::vector< int >> rings, std::vector< int > listPrism, int ringSize, std::vector< ring::strucType > ringType, std::vector< int > *atomTypes, std::vector< ring::strucType > *atomState)
bool hasCommonElements(std::vector< int > ring1, std::vector< int > ring2)
Definition: ring.cpp:226
bool basalPrismConditions(std::vector< std::vector< int >> nList, std::vector< int > *basal1, std::vector< int > *basal2)
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.
Definition: ring.cpp:200
@ Prism
The ring belongs to a prism block, classified according to the prism identification scheme.
@ unclassified
The ring is unclassified, which may be either water or a deformed type which cannot be classified by ...
T max_element(T... args)
bool matchPrism(molSys::PointCloud< molSys::Point< double >, double > *yCloud, std::vector< std::vector< int >> nList, const Eigen::MatrixXd &refPoints, std::vector< int > *basal1, std::vector< int > *basal2, std::vector< double > *rmsdPerAtom, bool isPerfect=true)
Definition: shapeMatch.cpp:21
Eigen::MatrixXd getPointSetRefRing(int n, int axialDim)
int writePrismNum(std::string path, std::vector< int > nPrisms, std::vector< int > nDefPrisms, std::vector< double > heightPercent, int maxDepth, int currentFrame, int firstFrame)
int writeLAMMPSdumpINT(molSys::PointCloud< molSys::Point< double >, double > *yCloud, std::vector< double > rmsdPerAtom, std::vector< int > atomTypes, int maxDepth, std::string path)
Write out a LAMMPS dump file containing the RMSD per atom.
int writeLAMMPSdataAllPrisms(molSys::PointCloud< molSys::Point< double >, double > *yCloud, std::vector< std::vector< int >> nList, std::vector< int > atomTypes, int maxDepth, std::string path, bool doShapeMatching=false)
Write a data file for prisms of every type.
double normHeightPercent(molSys::PointCloud< molSys::Point< double >, double > *yCloud, int nPrisms, double avgPrismHeight)
T push_back(T... args)
T resize(T... args)
T size(T... args)
This contains a collection of points; contains information for a particular frame.
Definition: mol_sys.hpp:170
This contains per-particle information.
Definition: mol_sys.hpp:149
File containing functions used specific to quasi-one-dimensional topological network critera (the pri...
T unique(T... args)