1 |
/* |
2 |
* Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
3 |
* |
4 |
* The University of Notre Dame grants you ("Licensee") a |
5 |
* non-exclusive, royalty free, license to use, modify and |
6 |
* redistribute this software in source and binary code form, provided |
7 |
* that the following conditions are met: |
8 |
* |
9 |
* 1. Redistributions of source code must retain the above copyright |
10 |
* notice, this list of conditions and the following disclaimer. |
11 |
* |
12 |
* 2. Redistributions in binary form must reproduce the above copyright |
13 |
* notice, this list of conditions and the following disclaimer in the |
14 |
* documentation and/or other materials provided with the |
15 |
* distribution. |
16 |
* |
17 |
* This software is provided "AS IS," without a warranty of any |
18 |
* kind. All express or implied conditions, representations and |
19 |
* warranties, including any implied warranty of merchantability, |
20 |
* fitness for a particular purpose or non-infringement, are hereby |
21 |
* excluded. The University of Notre Dame and its licensors shall not |
22 |
* be liable for any damages suffered by licensee as a result of |
23 |
* using, modifying or distributing the software or its |
24 |
* derivatives. In no event will the University of Notre Dame or its |
25 |
* licensors be liable for any lost revenue, profit or data, or for |
26 |
* direct, indirect, special, consequential, incidental or punitive |
27 |
* damages, however caused and regardless of the theory of liability, |
28 |
* arising out of the use of or inability to use software, even if the |
29 |
* University of Notre Dame has been advised of the possibility of |
30 |
* such damages. |
31 |
* |
32 |
* SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your |
33 |
* research, please cite the appropriate papers when you publish your |
34 |
* work. Good starting points are: |
35 |
* |
36 |
* [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005). |
37 |
* [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006). |
38 |
* [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008). |
39 |
* [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010). |
40 |
* [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011). |
41 |
*/ |
42 |
|
43 |
#include <algorithm> |
44 |
#include <functional> |
45 |
#include "applications/sequentialProps/ContactAngle2.hpp" |
46 |
#include "utils/simError.h" |
47 |
#include "io/DumpReader.hpp" |
48 |
#include "primitives/Molecule.hpp" |
49 |
#include "utils/NumericConstant.hpp" |
50 |
#include "utils/PhysicalConstants.hpp" |
51 |
#include "math/Eigenvalue.hpp" |
52 |
|
53 |
namespace OpenMD { |
54 |
|
55 |
ContactAngle2::ContactAngle2(SimInfo* info, const std::string& filename, |
56 |
const std::string& sele, RealType solidZ, |
57 |
RealType threshDens, RealType bufferLength, |
58 |
int nrbins, int nzbins) |
59 |
: SequentialAnalyzer(info, filename), solidZ_(solidZ), |
60 |
threshDens_(threshDens), bufferLength_(bufferLength), nRBins_(nrbins), |
61 |
nZBins_(nzbins), selectionScript_(sele), seleMan_(info), |
62 |
evaluator_(info) { |
63 |
|
64 |
setOutputName(getPrefix(filename) + ".ca2"); |
65 |
|
66 |
evaluator_.loadScriptString(sele); |
67 |
|
68 |
if (!evaluator_.isDynamic()) { |
69 |
seleMan_.setSelectionSet(evaluator_.evaluate()); |
70 |
} |
71 |
} |
72 |
|
73 |
void ContactAngle2::doFrame() { |
74 |
StuntDouble* sd; |
75 |
int i; |
76 |
|
77 |
// set up the bins for density analysis |
78 |
|
79 |
Mat3x3d hmat = info_->getSnapshotManager()->getCurrentSnapshot()->getHmat(); |
80 |
RealType len = std::min(hmat(0, 0), hmat(1, 1)); |
81 |
RealType zLen = hmat(2,2); |
82 |
|
83 |
RealType dr = len / (RealType) nRBins_; |
84 |
RealType dz = zLen / (RealType) nZBins_; |
85 |
|
86 |
std::vector<std::vector<RealType> > histo; |
87 |
histo.resize(nRBins_); |
88 |
for (unsigned int i = 0; i < histo.size(); ++i){ |
89 |
histo[i].resize(nZBins_); |
90 |
std::fill(histo[i].begin(), histo[i].end(), 0.0); |
91 |
} |
92 |
|
93 |
if (evaluator_.isDynamic()) { |
94 |
seleMan_.setSelectionSet(evaluator_.evaluate()); |
95 |
} |
96 |
|
97 |
|
98 |
RealType mtot = 0.0; |
99 |
Vector3d com(V3Zero); |
100 |
RealType mass; |
101 |
|
102 |
for (sd = seleMan_.beginSelected(i); sd != NULL; |
103 |
sd = seleMan_.nextSelected(i)) { |
104 |
mass = sd->getMass(); |
105 |
mtot += mass; |
106 |
com += sd->getPos() * mass; |
107 |
} |
108 |
|
109 |
com /= mtot; |
110 |
|
111 |
// now that we have the centroid, we can make cylindrical density maps |
112 |
Vector3d pos; |
113 |
RealType r; |
114 |
RealType z; |
115 |
|
116 |
for (sd = seleMan_.beginSelected(i); sd != NULL; |
117 |
sd = seleMan_.nextSelected(i)) { |
118 |
pos = sd->getPos() - com; |
119 |
|
120 |
// r goes from zero upwards |
121 |
r = sqrt(pow(pos.x(), 2) + pow(pos.y(), 2)); |
122 |
// z is possibly symmetric around 0 |
123 |
z = pos.z(); |
124 |
|
125 |
int whichRBin = int(r / dr); |
126 |
int whichZBin = int( (zLen/2.0 + z) / dz); |
127 |
|
128 |
if ((whichRBin < nRBins_) && (whichZBin >= 0) && (whichZBin < nZBins_)) { |
129 |
std::size_t i = static_cast<std::size_t>(whichRBin); |
130 |
std::size_t j = static_cast<std::size_t>(whichZBin); |
131 |
histo[i][j] += sd->getMass(); |
132 |
} |
133 |
|
134 |
} |
135 |
|
136 |
for(unsigned int i = 0 ; i < histo.size(); ++i){ |
137 |
|
138 |
RealType rL = i * dr; |
139 |
RealType rU = rL + dr; |
140 |
RealType volSlice = NumericConstant::PI * dz * (( rU*rU ) - ( rL*rL )); |
141 |
|
142 |
for (unsigned int j = 0; j < histo[i].size(); ++j) { |
143 |
histo[i][j] *= PhysicalConstants::densityConvert / volSlice; |
144 |
} |
145 |
} |
146 |
|
147 |
std::vector<Vector<RealType, 2> > points; |
148 |
points.clear(); |
149 |
|
150 |
for (unsigned int j = 0; j < nZBins_; ++j) { |
151 |
|
152 |
// The z coordinates were measured relative to the selection |
153 |
// center of mass. However, we're interested in the elevation |
154 |
// above the solid surface. Also, the binning was done around |
155 |
// zero with enough bins to cover the zLength of the box: |
156 |
|
157 |
RealType thez = com.z() - solidZ_ - zLen/2.0 + dz * (j + 0.5); |
158 |
bool aboveThresh = false; |
159 |
bool foundThresh = false; |
160 |
int rloc = 0; |
161 |
|
162 |
for (std::size_t i = 0; i < nRBins_; ++i) { |
163 |
|
164 |
if (histo[i][j] >= threshDens_) aboveThresh = true; |
165 |
|
166 |
if (aboveThresh && (histo[i][j] <= threshDens_)) { |
167 |
rloc = i; |
168 |
foundThresh = true; |
169 |
aboveThresh = false; |
170 |
} |
171 |
|
172 |
} |
173 |
if (foundThresh) { |
174 |
Vector<RealType,2> point; |
175 |
point[0] = dr*(rloc+0.5); |
176 |
point[1] = thez; |
177 |
|
178 |
if (thez > bufferLength_) { |
179 |
points.push_back( point ); |
180 |
} |
181 |
} |
182 |
} |
183 |
|
184 |
int numPoints = points.size(); |
185 |
|
186 |
// Compute the average of the data points. |
187 |
Vector<RealType, 2> average = points[0]; |
188 |
int i0; |
189 |
for (i0 = 1; i0 < numPoints; ++i0) { |
190 |
average += points[i0]; |
191 |
} |
192 |
RealType invNumPoints = ((RealType)1)/(RealType)numPoints; |
193 |
average *= invNumPoints; |
194 |
|
195 |
DynamicRectMatrix<RealType> mat(4, 4); |
196 |
int row, col; |
197 |
for (row = 0; row < 4; ++row) { |
198 |
for (col = 0; col < 4; ++col){ |
199 |
mat(row,col) = 0.0; |
200 |
} |
201 |
} |
202 |
for (int i = 0; i < numPoints; ++i) { |
203 |
RealType x = points[i][0]; |
204 |
RealType y = points[i][1]; |
205 |
RealType x2 = x*x; |
206 |
RealType y2 = y*y; |
207 |
RealType xy = x*y; |
208 |
RealType r2 = x2+y2; |
209 |
RealType xr2 = x*r2; |
210 |
RealType yr2 = y*r2; |
211 |
RealType r4 = r2*r2; |
212 |
|
213 |
mat(0,1) += x; |
214 |
mat(0,2) += y; |
215 |
mat(0,3) += r2; |
216 |
mat(1,1) += x2; |
217 |
mat(1,2) += xy; |
218 |
mat(1,3) += xr2; |
219 |
mat(2,2) += y2; |
220 |
mat(2,3) += yr2; |
221 |
mat(3,3) += r4; |
222 |
} |
223 |
mat(0,0) = (RealType)numPoints; |
224 |
|
225 |
for (row = 0; row < 4; ++row) { |
226 |
for (col = 0; col < row; ++col) { |
227 |
mat(row,col) = mat(col,row); |
228 |
} |
229 |
} |
230 |
|
231 |
for (row = 0; row < 4; ++row) { |
232 |
for (col = 0; col < 4; ++col) { |
233 |
mat(row,col) *= invNumPoints; |
234 |
} |
235 |
} |
236 |
|
237 |
JAMA::Eigenvalue<RealType> eigensystem(mat); |
238 |
DynamicRectMatrix<RealType> evects(4, 4); |
239 |
DynamicVector<RealType> evals(4); |
240 |
|
241 |
eigensystem.getRealEigenvalues(evals); |
242 |
eigensystem.getV(evects); |
243 |
|
244 |
DynamicVector<RealType> evector = evects.getColumn(0); |
245 |
RealType inv = ((RealType)1)/evector[3]; // beware zero divide |
246 |
RealType coeff[3]; |
247 |
for (row = 0; row < 3; ++row) { |
248 |
coeff[row] = inv*evector[row]; |
249 |
} |
250 |
|
251 |
Vector<RealType, 2> center; |
252 |
|
253 |
center[0] = -((RealType)0.5)*coeff[1]; |
254 |
center[1] = -((RealType)0.5)*coeff[2]; |
255 |
RealType radius = sqrt(fabs(center[0]*center[0] + center[1]*center[1] |
256 |
- coeff[0])); |
257 |
|
258 |
int i1; |
259 |
for (i1 = 0; i1 < 100; ++i1) { |
260 |
// Update the iterates. |
261 |
Vector<RealType, 2> current = center; |
262 |
|
263 |
// Compute average L, dL/da, dL/db. |
264 |
RealType lenAverage = (RealType)0; |
265 |
Vector<RealType, 2> derLenAverage = Vector<RealType, 2>(0.0); |
266 |
for (i0 = 0; i0 < numPoints; ++i0) { |
267 |
Vector<RealType, 2> diff = points[i0] - center; |
268 |
RealType length = diff.length(); |
269 |
if (length > 1e-6) { |
270 |
lenAverage += length; |
271 |
RealType invLength = ((RealType)1)/length; |
272 |
derLenAverage -= invLength*diff; |
273 |
} |
274 |
} |
275 |
lenAverage *= invNumPoints; |
276 |
derLenAverage *= invNumPoints; |
277 |
|
278 |
center = average + lenAverage*derLenAverage; |
279 |
radius = lenAverage; |
280 |
|
281 |
Vector<RealType, 2> diff = center - current; |
282 |
if (fabs(diff[0]) <= 1e-6 && fabs(diff[1]) <= 1e-6) { |
283 |
break; |
284 |
} |
285 |
} |
286 |
|
287 |
RealType zCen = center[1]; |
288 |
RealType rDrop = radius; |
289 |
RealType ca; |
290 |
|
291 |
if (fabs(zCen) > rDrop) { |
292 |
ca = 180.0; |
293 |
} else { |
294 |
ca = 90.0 + asin(zCen/rDrop)*(180.0/M_PI); |
295 |
} |
296 |
|
297 |
values_.push_back( ca ); |
298 |
|
299 |
} |
300 |
} |
301 |
|
302 |
|