Gen
+ Collaboration diagram for Gen:

Modules

 Ring
 

Functions

double gen::radDeg (double angle)
 
double gen::eigenVecAngle (std::vector< double > OO, std::vector< double > OH)
 Eigen function for getting the angle (in radians) between the O–O and O-H vectors. More...
 
double gen::getAverageWithoutOutliers (std::vector< double > inpVec)
 Get the average, after excluding the outliers, using quartiles. More...
 
double gen::calcMedian (std::vector< double > *input)
 Inline generic function for calculating the median given a vector of the values. More...
 
double gen::periodicDist (molSys::PointCloud< molSys::Point< double >, double > *yCloud, int iatom, int jatom)
 Inline generic function for obtaining the unwrapped periodic distance between two particles, whose indices (not IDs) have been given. More...
 
double gen::unWrappedDistFromPoint (molSys::PointCloud< molSys::Point< double >, double > *yCloud, int iatom, std::vector< double > singlePoint)
 
double gen::distance (molSys::PointCloud< molSys::Point< double >, double > *yCloud, int iatom, int jatom)
 Inline generic function for obtaining the wrapped distance between two particles WITHOUT applying PBCs, whose indices (not IDs) have been given. More...
 
std::array< double, 3 > gen::relDist (molSys::PointCloud< molSys::Point< double >, double > *yCloud, int iatom, int jatom)
 
bool gen::compareByAtomID (const molSys::Point< double > &a, const molSys::Point< double > &b)
 
int gen::prettyPrintYoda (molSys::PointCloud< molSys::Point< double >, double > *yCloud, std::string outFile)
 Generic function for printing all the struct information. More...
 
int gen::unwrappedCoordShift (molSys::PointCloud< molSys::Point< double >, double > *yCloud, int iatomIndex, int jatomIndex, double *x_i, double *y_i, double *z_i, double *x_j, double *y_j, double *z_j)
 Shift particles (unwrapped coordinates) More...
 
double gen::angDistDegQuaternions (std::vector< double > quat1, std::vector< double > quat2)
 
std::vector< std::stringgen::tokenizer (std::string line)
 Function for tokenizing line strings into words (strings) delimited by whitespace. This returns a vector with the words in it. More...
 
std::vector< double > gen::tokenizerDouble (std::string line)
 Function for tokenizing line strings into a vector of doubles. More...
 
std::vector< int > gen::tokenizerInt (std::string line)
 Function for tokenizing line strings into a vector of ints. More...
 
bool gen::file_exists (const std::string &name)
 Function for checking if a file exists or not. More...
 
std::vector< std::complex< double > > gen::avgVector (std::vector< std::complex< double >> v, int l, int neigh)
 

Variables

const double gen::pi = boost::math::constants::pi<double>()
 

Detailed Description

Function Documentation

◆ angDistDegQuaternions()

double gen::angDistDegQuaternions ( std::vector< double >  quat1,
std::vector< double >  quat2 
)

Function for getting the angular distance between two quaternions. Returns the result in degrees

Function for getting the angular distance between two quaternions. Returns the result in degrees.

Parameters
[in]quat1The first quaternion
[in]quat2The second quaternion
Returns
The output angle between the input quaternions, in degrees

Definition at line 267 of file generic.cpp.

268  {
269  //
270  double prod; // Product of quat1 and conjugate of quat2
271  // The angular distance is
272  // angularDistance = 2*cosInverse(quat1*conj(quat2))
273  prod = quat1[0] * quat2[0] - quat1[1] * quat2[1] - quat1[2] * quat2[2] -
274  quat1[3] * quat2[3];
275  // The angular distance is:
276  double angDist = 2 * acos(prod) * 180.0 / (gen::pi);
277  // Return the angular distance
278  return angDist;
279 }
T acos(T... args)
const double pi
Definition: generic.hpp:54

◆ avgVector()

std::vector<std::complex<double> > gen::avgVector ( std::vector< std::complex< double >>  v,
int  l,
int  neigh 
)
inline

Calculates the complex vector, normalized by the number of nearest neighbours, of length \(2l+1\).

Parameters
[in]vThe complex vector to be normalized, of length \(2l+1\)
[in]lA free integer parameter
[in]neighThe number of nearest neighbours
Returns
length \(2l+1\), normalized by the number of nearest neighbours

Definition at line 313 of file generic.hpp.

313  {
314  if (neigh == 0) {
315  return v;
316  }
317  for (int m = 0; m < 2 * l + 1; m++) {
318  v[m] = (1.0 / (double)neigh) * v[m];
319  }
320 
321  return v;
322 }

◆ calcMedian()

double gen::calcMedian ( std::vector< double > *  input)
inline

Inline generic function for calculating the median given a vector of the values.

Parameters
[in]yCloudThe input PointCloud, which contains the particle coordinates, simulation box lengths etc.
[in]inputThe input vector with the values
Returns
The median value

Definition at line 77 of file generic.hpp.

77  {
78  int n = (*input).size(); // Number of elements
79  double median; // Output median value
80 
81  // Sort the vector
82  std::sort((*input).begin(), (*input).end());
83 
84  // Calculate the median
85  // For even values, the median is the average of the two middle values
86  if (n % 2 == 0) {
87  median = 0.5 * ((*input)[n / 2] + (*input)[n / 2 - 1]); // n/2+n/2-1
88  } // median is average of middle values
89  else {
90  median = (*input)[(n + 1) / 2 -
91  1]; // middle value of 7 elements is the 4th element
92  } // if odd, it is the middle value
93 
94  return median;
95 }
T sort(T... args)

◆ compareByAtomID()

bool gen::compareByAtomID ( const molSys::Point< double > &  a,
const molSys::Point< double > &  b 
)
inline

Inline generic function for sorting or comparing two particles, according to the atom ID when the entire Point objects have been passed.

Parameters
[in]aThe input Point for A.
[in]bThe input Point for B.
Returns
True if the atom ID of A is less than the atom ID of B

Definition at line 234 of file generic.hpp.

235  {
236  return a.atomID < b.atomID;
237 }

◆ distance()

double gen::distance ( molSys::PointCloud< molSys::Point< double >, double > *  yCloud,
int  iatom,
int  jatom 
)
inline

Inline generic function for obtaining the wrapped distance between two particles WITHOUT applying PBCs, whose indices (not IDs) have been given.

Parameters
[in]yCloudThe input PointCloud, which contains the particle coordinates, simulation box lengths etc.
[in]iatomThe index of the \( i^{th} \) atom.
[in]jatomThe index of the \( j^{th} \) atom.
Returns
The wrapped distance.

Definition at line 170 of file generic.hpp.

171  {
173  double r2 = 0.0; // Squared absolute distance
174 
175  // Get x1-x2 etc
176  dr[0] = fabs(yCloud->pts[iatom].x - yCloud->pts[jatom].x);
177  dr[1] = fabs(yCloud->pts[iatom].y - yCloud->pts[jatom].y);
178  dr[2] = fabs(yCloud->pts[iatom].z - yCloud->pts[jatom].z);
179 
180  // Get the squared absolute distance
181  for (int k = 0; k < 3; k++) {
182  r2 += pow(dr[k], 2.0);
183  }
184 
185  return sqrt(r2);
186 }
T fabs(T... args)
std::vector< S > pts
Definition: mol_sys.hpp:171
T pow(T... args)
T sqrt(T... args)

◆ eigenVecAngle()

double gen::eigenVecAngle ( std::vector< double >  OO,
std::vector< double >  OH 
)

Eigen function for getting the angle (in radians) between the O–O and O-H vectors.

Function for obtaining the angle between two input vectors (std::vector). Internally, the vectors are converted to Eigen vectors. The dot product between the input vectors is used to calculate the angle between them.

Parameters
[in]OOThe O–O vector (but can be any vector, in general).
[in]OHThe O-H vector (but can be any vector, in general).
Returns
The output angle between the input vectors, in radians

Definition at line 155 of file generic.cpp.

155  {
156  Eigen::Vector3d eigOO = Eigen::Map<Eigen::Vector3d>(OO.data(), OO.size());
157  Eigen::Vector3d eigOH = Eigen::Map<Eigen::Vector3d>(OH.data(), OH.size());
158  double angle;
159  angle = acos(eigOO.dot(eigOH) / (eigOO.norm() * eigOH.norm()));
160  return angle;
161 }
T data(T... args)
T size(T... args)

◆ file_exists()

bool gen::file_exists ( const std::string name)
inline

Function for checking if a file exists or not.

Parameters
[in]nameThe name of the file

Definition at line 298 of file generic.hpp.

298  {
299  // Replace by boost function later
300  struct stat buffer;
301  return (stat(name.c_str(), &buffer) == 0);
302 }
T c_str(T... args)

◆ getAverageWithoutOutliers()

double gen::getAverageWithoutOutliers ( std::vector< double >  inpVec)

Get the average, after excluding the outliers, using quartiles.

Get the average, after excluding the outliers, using quartiles

Parameters
[in]inpVecThe vector containing values whose average is desired
Returns
The desired average value

Definition at line 169 of file generic.cpp.

169  {
170  //
171  double avgVal = 0.0; // Average value, excluding the outliers
172  double median; // Median value
173  int n = inpVec.size(); // Number of values
174  std::vector<double> lowerRange,
175  upperRange; // n/2 smallest and n/2 largest numbers
176  double firstQuartile, thirdQuartile; // First and third quartiles
177  double iqr; // Interquartile range
178  double outlierLimLow, outlierLimHigh; // Outliers limit
179  int numOfObservations = 0; // Number of observations used for the average
180  // ----------------------
181  // Calculate the median (the vector is sorted inside the function)
182  median = calcMedian(&inpVec);
183  // ----------------------
184  // Get the n/2 smallest and largest numbers
185  //
186  if (n % 2 == 0) {
187  for (int i = 0; i < n / 2; i++) {
188  // n/2 smallest numbers
189  lowerRange.push_back(inpVec[i]);
190  // n/2 largest numbers
191  upperRange.push_back(inpVec[n / 2 + i]);
192  } // end of loop to fill up the n/2 smallest and n/2 largest
193  } // even
194  else {
195  //
196  int halfN = (n + 1) / 2;
197  // Exclude the median
198  for (int i = 0; i < halfN; i++) {
199  // (n+1)/2 smallest numbers
200  lowerRange.push_back(inpVec[i]);
201  // (n+1)/2 largest numbers
202  upperRange.push_back(inpVec[halfN + i]);
203  } // end of filling up the smallest and largest half-ranges
204  } // for odd numbers
205  // ----------------------
206  // Calculate the first and third quartiles, and interquartile range
207  //
208  // First quartile
209  firstQuartile = calcMedian(&lowerRange);
210  // Third quartile
211  thirdQuartile = calcMedian(&upperRange);
212  // Interquartile range
213  iqr = thirdQuartile - firstQuartile;
214  // ----------------------
215  // Calculate the average without outliers
216  // Outliers are defined as values which
217  // are less than Q1 - 1.5IQR
218  // or greater than Q3 + 1.5IQR
219  //
220  // Get the limits for the outliers
221  outlierLimLow = firstQuartile - 1.5 * iqr;
222  outlierLimHigh = thirdQuartile + 1.5 * iqr;
223  //
224  // Loop through the values in inpVec to get the average, excluding outliers
225  for (int i = 0; i < n; i++) {
226  //
227  if (inpVec[i] < outlierLimLow) {
228  continue;
229  } // lower limit outlier
230  else if (inpVec[i] > outlierLimHigh) {
231  continue;
232  } // higher limit outlier
233  else {
234  // Number of observations added
235  numOfObservations++;
236  // Add to the average
237  avgVal += inpVec[i];
238  } // take the average
239  } // end of loop for getting the average
240  //
241  // Divide by the number of observations used
242  avgVal /= numOfObservations;
243  // ----------------------
244  // This fails if there are not enough observations (ring size = 3)
245  if (numOfObservations == 0) {
246  double sumVal = 0.0;
247  // Loop through all the values and sum
248  for (int i = 0; i <= n; i++) {
249  sumVal += inpVec[i];
250  } // end of sum
251  // Normalize
252  avgVal = sumVal / n;
253  return avgVal;
254  } // for triangular prism blocks
255  // ----------------------
256 
257  return avgVal;
258 }
double calcMedian(std::vector< double > *input)
Inline generic function for calculating the median given a vector of the values.
Definition: generic.hpp:77
T push_back(T... args)

◆ periodicDist()

double gen::periodicDist ( molSys::PointCloud< molSys::Point< double >, double > *  yCloud,
int  iatom,
int  jatom 
)
inline

Inline generic function for obtaining the unwrapped periodic distance between two particles, whose indices (not IDs) have been given.

Parameters
[in]yCloudThe input PointCloud, which contains the particle coordinates, simulation box lengths etc.
[in]iatomThe index of the \( i^{th} \) atom.
[in]jatomThe index of the \( j^{th} \) atom.
Returns
The unwrapped periodic distance.

Definition at line 108 of file generic.hpp.

109  {
111  double r2 = 0.0; // Squared absolute distance
112 
113  // Get x1-x2 etc
114  dr[0] = fabs(yCloud->pts[iatom].x - yCloud->pts[jatom].x);
115  dr[1] = fabs(yCloud->pts[iatom].y - yCloud->pts[jatom].y);
116  dr[2] = fabs(yCloud->pts[iatom].z - yCloud->pts[jatom].z);
117 
118  // Get the squared absolute distance
119  for (int k = 0; k < 3; k++) {
120  // Correct for periodicity
121  dr[k] -= yCloud->box[k] * round(dr[k] / yCloud->box[k]);
122  r2 += pow(dr[k], 2.0);
123  }
124 
125  return sqrt(r2);
126 }
std::vector< T > box
Number of atoms.
Definition: mol_sys.hpp:174
T round(T... args)

◆ prettyPrintYoda()

int gen::prettyPrintYoda ( molSys::PointCloud< molSys::Point< double >, double > *  yCloud,
std::string  outFile 
)

Generic function for printing all the struct information.

Function for printing out info in a PointCloud object.

Parameters
[in]yCloudThe input PointCloud to be printed.
[in]outFileThe name of the output file to which the information will be printed.

Definition at line 25 of file generic.cpp.

27  {
28  std::ofstream outputFile;
29  // Create a new file in the output directory
30  outputFile.open(outFile);
31 
32  if (outputFile.is_open()) {
33  // First line
34  outputFile << "# Frame\tAtomID\tx\ty\tz\tcij\ticeType\n";
35  // Write out all the information out line by line
36  for (int i = 0; i < yCloud->nop; i++) {
37  outputFile << yCloud->currentFrame << "\t" << yCloud->pts[i].atomID
38  << "\t" << yCloud->pts[i].x << "\t" << yCloud->pts[i].y << "\t"
39  << yCloud->pts[i].z << "\t";
40  // Print out cij
41  // for(int c=0; c<yCloud->pts[i].c_ij.size(); c++){outputFile <<
42  // yCloud->pts[i].c_ij[c]<<"\t";} Print out the classifier
43  // TODO: Should print string representation
44  outputFile << static_cast<int>(yCloud->pts[i].iceType) << "\n";
45  }
46  }
47  // Close the file
48  outputFile.close();
49  return 0;
50 }
T close(T... args)
int nop
Current frame number.
Definition: mol_sys.hpp:173
int currentFrame
Collection of points.
Definition: mol_sys.hpp:172
T is_open(T... args)
T open(T... args)

◆ radDeg()

double gen::radDeg ( double  angle)
inline

Inline function for converting radians->degrees.

Parameters
[in]angleThe input angle, in radians
Returns
The input angle, in degrees

Definition at line 61 of file generic.hpp.

61 { return (angle * 180) / gen::pi; }

◆ relDist()

std::array<double, 3> gen::relDist ( molSys::PointCloud< molSys::Point< double >, double > *  yCloud,
int  iatom,
int  jatom 
)
inline

Inline generic function for getting the relative unwrapped distance between two particles for each dimension. The indices (not IDs) of the particles have been given.

Parameters
[in]yCloudThe input PointCloud, which contains the particle coordinates, simulation box lengths etc.
[in]iatomThe index of the \( i^{th} \) atom.
[in]jatomThe index of the \( j^{th} \) atom.
Returns
The unwrapped relative distances for each dimension.

Definition at line 200 of file generic.hpp.

201  {
203  std::array<double, 3> box = {yCloud->box[0], yCloud->box[1], yCloud->box[2]};
204  double r2 = 0.0; // Squared absolute distance
205 
206  // Get x1-x2 etc
207  dr[0] = yCloud->pts[iatom].x - yCloud->pts[jatom].x;
208  dr[1] = yCloud->pts[iatom].y - yCloud->pts[jatom].y;
209  dr[2] = yCloud->pts[iatom].z - yCloud->pts[jatom].z;
210 
211  // Get the relative distance
212  for (int k = 0; k < 3; k++) {
213  //
214  if (dr[k] < -box[k] * 0.5) {
215  dr[k] = dr[k] + box[k];
216  }
217  if (dr[k] >= box[k] * 0.5) {
218  dr[k] = dr[k] - box[k];
219  }
220  }
221 
222  return dr;
223 }

◆ tokenizer()

std::vector<std::string> gen::tokenizer ( std::string  line)
inline

Function for tokenizing line strings into words (strings) delimited by whitespace. This returns a vector with the words in it.

Parameters
[in]lineThe string containing the line to be tokenized

Definition at line 259 of file generic.hpp.

◆ tokenizerDouble()

std::vector<double> gen::tokenizerDouble ( std::string  line)
inline

Function for tokenizing line strings into a vector of doubles.

Parameters
[in]lineThe string containing the line to be tokenized

Definition at line 270 of file generic.hpp.

270  {
271  std::istringstream iss(line);
272  std::vector<double> tokens;
273  double number; // Each number being read in from the line
274  while (iss >> number) {
275  tokens.push_back(number);
276  }
277  return tokens;
278 }

◆ tokenizerInt()

std::vector<int> gen::tokenizerInt ( std::string  line)
inline

Function for tokenizing line strings into a vector of ints.

Parameters
[in]lineThe string containing the line to be tokenized

Definition at line 284 of file generic.hpp.

284  {
285  std::istringstream iss(line);
286  std::vector<int> tokens;
287  int number; // Each number being read in from the line
288  while (iss >> number) {
289  tokens.push_back(number);
290  }
291  return tokens;
292 }

◆ unwrappedCoordShift()

int gen::unwrappedCoordShift ( molSys::PointCloud< molSys::Point< double >, double > *  yCloud,
int  iatomIndex,
int  jatomIndex,
double *  x_i,
double *  y_i,
double *  z_i,
double *  x_j,
double *  y_j,
double *  z_j 
)

Shift particles (unwrapped coordinates)

Function for getting the unwrapped coordinates of a pair of atoms.

Parameters
[in]yCloudThe input PointCloud to be printed.
[in]iatomIndexIndex of the \( i^{th} \) atom.
[in]jatomIndexIndex of the \( j^{th} \) atom.
[in,out]x_iX Coordinate of the \( i^{th} \) atom corresponding to the unwrapped distance.
[in,out]y_iY Coordinate of the \( i^{th} \) atom corresponding to the unwrapped distance.
[in,out]z_iZ Coordinate of the \( i^{th} \) atom corresponding to the unwrapped distance.
[in,out]x_jX Coordinate of the \( j^{th} \) atom corresponding to the unwrapped distance.
[in,out]y_jY Coordinate of the \( j^{th} \) atom corresponding to the unwrapped distance.
[in,out]z_jZ Coordinate of the \( j^{th} \) atom corresponding to the unwrapped distance.

Definition at line 71 of file generic.cpp.

74  {
75  //
76  double x_iatom, y_iatom, z_iatom;
77  double x_jatom, y_jatom, z_jatom;
78  double x_ij, y_ij, z_ij; // Relative distance
79  std::vector<double> box = yCloud->box;
80  double xPBC, yPBC, zPBC; // Actual unwrapped distance
81 
82  // ----------------------------------------------------------------------
83  // INIT
84  // iatom
85  x_iatom = yCloud->pts[iatomIndex].x;
86  y_iatom = yCloud->pts[iatomIndex].y;
87  z_iatom = yCloud->pts[iatomIndex].z;
88  // jatom
89  x_jatom = yCloud->pts[jatomIndex].x;
90  y_jatom = yCloud->pts[jatomIndex].y;
91  z_jatom = yCloud->pts[jatomIndex].z;
92  // ----------------------------------------------------------------------
93  // GET RELATIVE DISTANCE
94  x_ij = x_iatom - x_jatom;
95  y_ij = y_iatom - y_jatom;
96  z_ij = z_iatom - z_jatom;
97  // ----------------------------------------------------------------------
98  // SHIFT COORDINATES IF REQUIRED
99  // Shift x
100  if (fabs(x_ij) > 0.5 * box[0]) {
101  // Get the actual distance
102  xPBC = box[0] - fabs(x_ij);
103  if (x_ij < 0) {
104  x_jatom = x_iatom - xPBC;
105  } // To the -x side of currentIndex
106  else {
107  x_jatom = x_iatom + xPBC;
108  } // Add to the + side
109  } // Shift nextElement
110  //
111  // Shift y
112  if (fabs(y_ij) > 0.5 * box[1]) {
113  // Get the actual distance
114  yPBC = box[1] - fabs(y_ij);
115  if (y_ij < 0) {
116  y_jatom = y_iatom - yPBC;
117  } // To the -y side of currentIndex
118  else {
119  y_jatom = y_iatom + yPBC;
120  } // Add to the + side
121  } // Shift nextElement
122  //
123  // Shift z
124  if (fabs(z_ij) > 0.5 * box[2]) {
125  // Get the actual distance
126  zPBC = box[2] - fabs(z_ij);
127  if (z_ij < 0) {
128  z_jatom = z_iatom - zPBC;
129  } // To the -z side of currentIndex
130  else {
131  z_jatom = z_iatom + zPBC;
132  } // Add to the + side
133  } // Shift nextElement
134  // ----------------------------------------------------------------------
135  // Assign values
136  *x_i = x_iatom;
137  *y_i = y_iatom;
138  *z_i = z_iatom;
139  *x_j = x_jatom;
140  *y_j = y_jatom;
141  *z_j = z_jatom;
142 
143  return 0;
144 }

◆ unWrappedDistFromPoint()

double gen::unWrappedDistFromPoint ( molSys::PointCloud< molSys::Point< double >, double > *  yCloud,
int  iatom,
std::vector< double >  singlePoint 
)
inline

Inline generic function for obtaining the unwrapped periodic distance between one particle and another point, whose index has been given.

Parameters
[in]yCloudThe input PointCloud, which contains the particle coordinates, simulation box lengths etc.
[in]iatomThe index of the \( i^{th} \) atom.
[in]singlePointVector containing coordinate values
Returns
The unwrapped periodic distance.

Definition at line 138 of file generic.hpp.

140  {
142  double r2 = 0.0; // Squared absolute distance
143 
144  // Get x1-x2 etc
145  dr[0] = fabs(yCloud->pts[iatom].x - singlePoint[0]);
146  dr[1] = fabs(yCloud->pts[iatom].y - singlePoint[1]);
147  dr[2] = fabs(yCloud->pts[iatom].z - singlePoint[2]);
148 
149  // Get the squared absolute distance
150  for (int k = 0; k < 3; k++) {
151  // Correct for periodicity
152  dr[k] -= yCloud->box[k] * round(dr[k] / yCloud->box[k]);
153  r2 += pow(dr[k], 2.0);
154  }
155 
156  return sqrt(r2);
157 }

Variable Documentation

◆ pi

const double gen::pi = boost::math::constants::pi<double>()

Uses Boost to get the value of pi.

Definition at line 54 of file generic.hpp.