ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/applications/staticProps/BondOrderParameter.cpp
Revision: 3017
Committed: Fri Sep 22 01:36:27 2006 UTC (17 years, 10 months ago) by gezelter
File size: 9933 byte(s)
Log Message:
Following Rein ten Wolde article

File Contents

# User Rev Content
1 chuckv 3005 /*
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. Acknowledgement of the program authors must be made in any
10     * publication of scientific results based in part on use of the
11     * program. An acceptable form of acknowledgement is citation of
12     * the article in which the program was described (Matthew
13     * A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14     * J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15     * Parallel Simulation Engine for Molecular Dynamics,"
16     * J. Comput. Chem. 26, pp. 252-271 (2005))
17     *
18     * 2. Redistributions of source code must retain the above copyright
19     * notice, this list of conditions and the following disclaimer.
20     *
21     * 3. Redistributions in binary form must reproduce the above copyright
22     * notice, this list of conditions and the following disclaimer in the
23     * documentation and/or other materials provided with the
24     * distribution.
25     *
26     * This software is provided "AS IS," without a warranty of any
27     * kind. All express or implied conditions, representations and
28     * warranties, including any implied warranty of merchantability,
29     * fitness for a particular purpose or non-infringement, are hereby
30     * excluded. The University of Notre Dame and its licensors shall not
31     * be liable for any damages suffered by licensee as a result of
32     * using, modifying or distributing the software or its
33     * derivatives. In no event will the University of Notre Dame or its
34     * licensors be liable for any lost revenue, profit or data, or for
35     * direct, indirect, special, consequential, incidental or punitive
36     * damages, however caused and regardless of the theory of liability,
37     * arising out of the use of or inability to use software, even if the
38     * University of Notre Dame has been advised of the possibility of
39     * such damages.
40     */
41    
42    
43     /* Creates orientational bond order parameters as outlined by
44     * Bond-orientaional order in liquids and glasses, Steinhart,Nelson,Ronchetti
45     * Phys Rev B, 28,784,1983
46     *
47     */
48 gezelter 3007
49 chuckv 3005 #include "applications/staticProps/BondOrderParameter.hpp"
50     #include "utils/simError.h"
51     #include "io/DumpReader.hpp"
52     #include "primitives/Molecule.hpp"
53     #include "utils/NumericConstant.hpp"
54 gezelter 3010 #include "math/SphericalHarmonic.hpp"
55 gezelter 3009
56 chuckv 3005 namespace oopse {
57    
58 gezelter 3007 BondOrderParameter::BondOrderParameter(SimInfo* info,
59     const std::string& filename,
60     const std::string& sele,
61 gezelter 3009 double rCut, int lNumber, int nbins) : StaticAnalyser(info, filename), selectionScript_(sele), evaluator_(info), seleMan_(info){
62 gezelter 3007
63 gezelter 3009 setOutputName(getPrefix(filename) + ".bo");
64 chuckv 3005
65 gezelter 3007 evaluator_.loadScriptString(sele);
66     if (!evaluator_.isDynamic()) {
67     seleMan_.setSelectionSet(evaluator_.evaluate());
68 chuckv 3005 }
69    
70 gezelter 3007 // Set up cutoff radius and order of the Legendre Polynomial:
71    
72 chuckv 3005 lNumber_ = lNumber;
73     rCut_ = rCut;
74 gezelter 3017 mSize_ = 2*lNumber_+1;
75    
76     // Q can take values from 0 to 1
77    
78     MinQ_ = 0.0;
79     MaxQ_ = 3.0;
80     deltaQ_ = (MaxQ_ - MinQ_) / nbins;
81     Q_histogram_.resize(nbins);
82    
83     // W_6 for icosahedral clusters is 11 / sqrt(4199) = 0.169754, so we'll
84     // use values for MinW_ and MaxW_ that are slightly larger than this:
85    
86     MinW_ = -0.18;
87     MaxW_ = 0.18;
88     deltaW_ = (MaxW_ - MinW_) / nbins;
89     W_histogram_.resize(nbins);
90    
91 gezelter 3007 }
92 chuckv 3005
93 gezelter 3009 BondOrderParameter::~BondOrderParameter() {
94 gezelter 3017 Q_histogram_.clear();
95     W_histogram_.clear();
96 gezelter 3009 }
97    
98 gezelter 3017 void BondOrderParameter::initalizeHistogram() {
99     std::fill(Q_histogram_.begin(), Q_histogram_.end(), 0);
100     std::fill(W_histogram_.begin(), W_histogram_.end(), 0);
101     }
102    
103 gezelter 3007 void BondOrderParameter::process() {
104     Molecule* mol;
105     Atom* atom;
106     RigidBody* rb;
107 gezelter 3011 int myIndex;
108 gezelter 3007 SimInfo::MoleculeIterator mi;
109     Molecule::RigidBodyIterator rbIter;
110     Molecule::AtomIterator ai;
111     StuntDouble* sd;
112 gezelter 3011 Vector3d vec;
113 gezelter 3010 RealType costheta;
114 gezelter 3007 RealType phi;
115     RealType r;
116     RealType dist;
117 gezelter 3017 std::map<int,ComplexType> q_lm;
118 gezelter 3014 std::map<int,ComplexType> QBar_lm;
119 gezelter 3007 RealType QSq_l;
120     RealType Q_l;
121 gezelter 3011 ComplexType W_l;
122     ComplexType W_l_hat;
123 gezelter 3017 int nBonds, Nbonds;
124 gezelter 3010 SphericalHarmonic sphericalHarmonic;
125 gezelter 3007 int i, j;
126 gezelter 3011 // Make arrays for Wigner3jm
127     double* THRCOF = new double[mSize_];
128     // Variables for Wigner routine
129     double l_ = (double)lNumber_;
130     double m1Pass, m2Min, m2Max;
131     int error, m1, m2, m3;
132    
133 gezelter 3008 // Set the l for the spherical harmonic, it doesn't change
134     sphericalHarmonic.setL(lNumber_);
135    
136 gezelter 3007 DumpReader reader(info_, dumpFilename_);
137     int nFrames = reader.getNFrames();
138 gezelter 3009 frameCounter_ = 0;
139 chuckv 3005
140 gezelter 3007 for (int istep = 0; istep < nFrames; istep += step_) {
141     reader.readFrame(istep);
142 gezelter 3009 frameCounter_++;
143 gezelter 3007 currentSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot();
144    
145     if (evaluator_.isDynamic()) {
146     seleMan_.setSelectionSet(evaluator_.evaluate());
147     }
148 chuckv 3005
149 gezelter 3007 // update the positions of atoms which belong to the rigidbodies
150 chuckv 3005
151 gezelter 3007 for (mol = info_->beginMolecule(mi); mol != NULL;
152     mol = info_->nextMolecule(mi)) {
153     for (rb = mol->beginRigidBody(rbIter); rb != NULL;
154     rb = mol->nextRigidBody(rbIter)) {
155     rb->updateAtoms();
156     }
157 gezelter 3017 }
158    
159 gezelter 3007 // outer loop is over the selected StuntDoubles:
160 chuckv 3005
161 gezelter 3007 for (sd = seleMan_.beginSelected(i); sd != NULL;
162     sd = seleMan_.nextSelected(i)) {
163 chuckv 3005
164 gezelter 3011 myIndex = sd->getGlobalIndex();
165 gezelter 3017 nBonds = 0;
166     for (int m = -lNumber_; m <= lNumber_; m++) {
167     q_lm[m] = 0.0;
168     }
169 gezelter 3007
170     // inner loop is over all other atoms in the system:
171    
172     for (mol = info_->beginMolecule(mi); mol != NULL;
173     mol = info_->nextMolecule(mi)) {
174     for (atom = mol->beginAtom(ai); atom != NULL;
175     atom = mol->nextAtom(ai)) {
176 chuckv 3005
177 gezelter 3011 if (atom->getGlobalIndex() != myIndex) {
178 chuckv 3005
179 gezelter 3011 vec = sd->getPos() - atom->getPos();
180     currentSnapshot_->wrapVector(vec);
181 gezelter 3010
182 gezelter 3011 // Calculate "bonds" and build Q_lm(r) where
183     // Q_lm = Y_lm(theta(r),phi(r))
184     // The spherical harmonics are wrt any arbitrary coordinate
185     // system, we choose standard spherical coordinates
186    
187     r = vec.length();
188    
189     // Check to see if neighbor is in bond cutoff
190    
191     if (r < rCut_) {
192     costheta = vec.z() / r;
193     phi = atan2(vec.y(), vec.x());
194    
195     for(int m = -lNumber_; m <= lNumber_; m++){
196     sphericalHarmonic.setM(m);
197 gezelter 3017 q_lm[m] += sphericalHarmonic.getValueAt(costheta, phi);
198 gezelter 3011 }
199     nBonds++;
200     }
201     }
202 gezelter 3007 }
203     }
204 gezelter 3017 RealType ql = 0.0;
205     for(int m=-lNumber_; m<=lNumber_; m++) {
206     ql += norm(QBar_lm[m]);
207     }
208     ql *= 4.0*NumericConstant::PI/(RealType)(2*lNumber_+1);
209     collectHistogram(sqrt(ql)/(RealType)nBonds);
210    
211     Nbonds += nBonds;
212     for (int m=-lNumber_; m<=lNumber_; m++) {
213     QBar_lm[m] += q_lm[m];
214     }
215     }
216 gezelter 3015 }
217 chuckv 3005
218 gezelter 3015 // Normalize Qbar2
219     for (int m = -lNumber_;m <= lNumber_; m++){
220 gezelter 3017 QBar_lm[m] /= Nbonds;
221 gezelter 3007 }
222    
223 gezelter 3015 // Find second order invariant Q_l
224 gezelter 3007
225 gezelter 3015 QSq_l = 0.0;
226     for (int m = -lNumber_; m <= lNumber_; m++){
227     QSq_l += norm(QBar_lm[m]);
228 gezelter 3007 }
229 gezelter 3015
230     std::cout << "qsl = " << QSq_l << "\n";
231     Q_l = sqrt(QSq_l * 4.0 * NumericConstant::PI / (RealType)(2*lNumber_ + 1));
232    
233     // Find Third Order Invariant W_l
234    
235     W_l = 0.0;
236     for (int m1 = -lNumber_; m1 <= lNumber_; m1++) {
237     // Zero work array
238     for (int ii = 0; ii < mSize_; ii++){
239     THRCOF[ii] = 0.0;
240     }
241     // Get Wigner coefficients
242     m1Pass = (double)m1;
243 gezelter 3009
244 gezelter 3015 Wigner3jm(&l_, &l_, &l_,
245     &m1Pass, &m2Min, &m2Max,
246     THRCOF, &mSize_, &error);
247 gezelter 3009
248 gezelter 3015 for (int mmm = 0; mmm < (int)(m2Max - m2Min); mmm++) {
249     m2 = (int)floor(m2Min) + mmm;
250     m3 = -m1-m2;
251     W_l += THRCOF[mmm] * QBar_lm[m1] * QBar_lm[m2] * QBar_lm[m3];
252 gezelter 3009 }
253 gezelter 3007 }
254 gezelter 3015
255     W_l_hat = W_l / pow(QSq_l, 1.5);
256    
257     writeOrderParameter(Q_l, real(W_l_hat));
258     }
259 chuckv 3005
260 gezelter 3017 void BondOrderParameter::collectHistogram(RealType q_l) {
261 chuckv 3005
262 gezelter 3017 if (q_l >= MinQ_ && q_l < MaxQ_) {
263     int qbin = (q_l - MinQ_) / deltaQ_;
264     Q_histogram_[qbin] += 1;
265     Qcount_++;
266     sumQ_ += q_l;
267     sumQ2_ += q_l * q_l;
268     } else {
269     sprintf( painCave.errMsg,
270     "q_l value outside reasonable range\n");
271     painCave.severity = OOPSE_ERROR;
272     painCave.isFatal = 1;
273     simError();
274     }
275    
276     }
277    
278    
279 gezelter 3015 void BondOrderParameter::writeOrderParameter(RealType ql, RealType Wlhat) {
280 chuckv 3005
281 gezelter 3015 std::ofstream os(getOutputFileName().c_str());
282    
283     if (os.is_open()) {
284 gezelter 3009
285 gezelter 3015 os << "# Bond Order Parameters\n";
286     os << "# selection: (" << selectionScript_ << ")\n";
287     os << "# \n";
288     os << "# <Q_" << lNumber_ << ">: " << ql << "\n";
289     os << "# <W_" << lNumber_ << ">: " << Wlhat << "\n";
290 gezelter 3017 // Normalize by number of frames and write it out:
291     for (int i = 0; i < Q_histogram_.size(); ++i) {
292     RealType Qval = MinQ_ + (i + 0.5) * deltaQ_;
293     osq << Qval << "\t" << (RealType)Q_histogram_[i] / (RealType)Qcount_ << "\n";
294     }
295    
296 gezelter 3015 os.close();
297    
298 gezelter 3009 } else {
299     sprintf(painCave.errMsg, "BondOrderParameter: unable to open %s\n",
300 gezelter 3015 getOutputFileName().c_str());
301 gezelter 3009 painCave.isFatal = 1;
302     simError();
303     }
304     }
305 gezelter 3007 }

Properties

Name Value
svn:executable *