ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/applications/staticProps/Hxy.cpp
Revision: 1465
Committed: Fri Jul 9 23:08:25 2010 UTC (14 years, 9 months ago) by chuckv
File size: 12968 byte(s)
Log Message:
Creating busticated version of OpenMD

File Contents

# User Rev Content
1 xsun 955 /*
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 gezelter 1390 * 1. Redistributions of source code must retain the above copyright
10 xsun 955 * notice, this list of conditions and the following disclaimer.
11     *
12 gezelter 1390 * 2. Redistributions in binary form must reproduce the above copyright
13 xsun 955 * 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 gezelter 1390 * 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, 24107 (2008).
39     * [4] Vardeman & Gezelter, in progress (2009).
40 xsun 955 *
41     *
42     * Created by Xiuquan Sun on 05/09/06.
43     * @author Xiuquan Sun
44 gezelter 1442 * @version $Id$
45 xsun 955 *
46     */
47    
48     /* Calculates the undulation spectrum of the lipid membrance. */
49    
50     #include <algorithm>
51     #include <fstream>
52     #include "applications/staticProps/Hxy.hpp"
53     #include "utils/simError.h"
54     #include "io/DumpReader.hpp"
55     #include "primitives/Molecule.hpp"
56     #include<stdio.h>
57     #include<string.h>
58     #include<stdlib.h>
59     #include<math.h>
60    
61 gezelter 1390 namespace OpenMD {
62 xsun 955
63     Hxy::Hxy(SimInfo* info, const std::string& filename, const std::string& sele, int nbins_x, int nbins_y, int nrbins)
64     : StaticAnalyser(info, filename), selectionScript_(sele), evaluator_(info), seleMan_(info), nBinsX_(nbins_x), nBinsY_(nbins_y), nbins_(nrbins){
65    
66     evaluator_.loadScriptString(sele);
67     if (!evaluator_.isDynamic()) {
68     seleMan_.setSelectionSet(evaluator_.evaluate());
69     }
70    
71     gridsample_.resize(nBinsX_*nBinsY_);
72     gridZ_.resize(nBinsX_*nBinsY_);
73 xsun 967 mag.resize(nBinsX_*nBinsY_);
74     newmag.resize(nBinsX_*nBinsY_);
75 xsun 955
76 gezelter 956 sum_bin.resize(nbins_);
77     avg_bin.resize(nbins_);
78     errbin_sum.resize(nbins_);
79     errbin.resize(nbins_);
80     sum_bin_sq.resize(nbins_);
81     avg_bin_sq.resize(nbins_);
82     errbin_sum_sq.resize(nbins_);
83     errbin_sq.resize(nbins_);
84 xsun 955
85 xsun 967 bin.resize(nbins_);
86     samples.resize(nbins_);
87    
88 xsun 955 setOutputName(getPrefix(filename) + ".Hxy");
89     }
90    
91 xsun 967 Hxy::~Hxy(){
92     gridsample_.clear();
93     gridZ_.clear();
94     sum_bin.clear();
95     avg_bin.clear();
96     errbin_sum.clear();
97     errbin.clear();
98     sum_bin_sq.clear();
99     avg_bin_sq.clear();
100     errbin_sum_sq.clear();
101     errbin_sq.clear();
102    
103     for(int i=0; i < bin.size(); i++)
104     bin[i].clear();
105     for(int i=0; i < samples.size(); i++)
106     samples[i].clear();
107    
108     mag.clear();
109     newmag.clear();
110     }
111    
112 xsun 955 void Hxy::process() {
113 gezelter 957 #if defined(HAVE_FFTW_H) || defined(HAVE_DFFTW_H) || defined(HAVE_FFTW3_H)
114 xsun 955 DumpReader reader(info_, dumpFilename_);
115     int nFrames = reader.getNFrames();
116     nProcessed_ = nFrames/step_;
117 gezelter 956
118 xsun 967 for(int k=0; k < bin.size(); k++)
119     bin[k].resize(nFrames);
120     for(int k=0; k < samples.size(); k++)
121     samples[k].resize(nFrames);
122    
123 tim 963 RealType lenX_, lenY_;
124     RealType gridX_, gridY_;
125     RealType halfBoxX_, halfBoxY_;
126 xsun 967
127 xsun 955 int binNoX, binNoY;
128 tim 963 RealType interpsum, value;
129 xsun 955 int ninterp, px, py, newp;
130     int newx, newy, newindex, index;
131     int new_i, new_j, new_index;
132 xsun 967
133 tim 963 RealType freq_x, freq_y, zero_freq_x, zero_freq_y, freq;
134 xsun 967 RealType maxfreqx, maxfreqy, maxfreq;
135    
136 xsun 955 int whichbin;
137 gezelter 956 int nMolecules;
138 xsun 967
139 xsun 969 std::fill(sum_bin.begin(), sum_bin.end(), 0.0);
140     std::fill(avg_bin.begin(), avg_bin.end(), 0.0);
141     std::fill(errbin_sum.begin(), errbin_sum.end(), 0.0);
142     std::fill(errbin.begin(), errbin.end(), 0.0);
143     std::fill(sum_bin_sq.begin(), sum_bin_sq.end(), 0.0);
144     std::fill(avg_bin_sq.begin(), avg_bin_sq.end(), 0.0);
145     std::fill(errbin_sum_sq.begin(), errbin_sum_sq.end(), 0.0);
146     std::fill(errbin_sq.begin(), errbin_sq.end(), 0.0);
147    
148     for(int i=0; i < bin.size(); i++)
149     std::fill(bin[i].begin(), bin[i].end(), 0.0);
150    
151     for(int i=0; i < samples.size(); i++)
152     std::fill(samples[i].begin(), samples[i].end(), 0);
153    
154 xsun 955 for (int istep = 0; istep < nFrames; istep += step_) {
155    
156 gezelter 956 reader.readFrame(istep);
157     currentSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot();
158     nMolecules = info_->getNGlobalMolecules();
159 gezelter 957
160 gezelter 956 Mat3x3d hmat = currentSnapshot_->getHmat();
161    
162 gezelter 957 #ifdef HAVE_FFTW3_H
163     fftw_plan p;
164     #else
165     fftwnd_plan p;
166     #endif
167 xsun 955 fftw_complex *in, *out;
168 gezelter 956
169     in = (fftw_complex*) fftw_malloc(sizeof(fftw_complex) * (nBinsX_*nBinsY_));
170     out = (fftw_complex*) fftw_malloc(sizeof(fftw_complex) *(nBinsX_*nBinsY_));
171 gezelter 957
172     #ifdef HAVE_FFTW3_H
173     p = fftw_plan_dft_2d(nBinsX_, nBinsY_, in, out,
174     FFTW_FORWARD, FFTW_ESTIMATE);
175     #else
176     p = fftw2d_create_plan(nBinsX_, nBinsY_, FFTW_FORWARD, FFTW_ESTIMATE);
177     #endif
178 xsun 969
179 xsun 967 std::fill(gridsample_.begin(), gridsample_.end(), 0);
180     std::fill(gridZ_.begin(), gridZ_.end(), 0.0);
181     std::fill(mag.begin(), mag.end(), 0.0);
182     std::fill(newmag.begin(), newmag.end(), 0.0);
183    
184 xsun 969 int i, j;
185 gezelter 956
186 xsun 955 StuntDouble* sd;
187 gezelter 956
188 xsun 955 lenX_ = hmat(0,0);
189     lenY_ = hmat(1,1);
190 gezelter 956
191 xsun 955 gridX_ = lenX_ /(nBinsX_);
192     gridY_ = lenY_ /(nBinsY_);
193 gezelter 956
194 xsun 967 halfBoxX_ = lenX_ / 2.0;
195     halfBoxY_ = lenY_ / 2.0;
196 gezelter 956
197 xsun 955 if (evaluator_.isDynamic()) {
198     seleMan_.setSelectionSet(evaluator_.evaluate());
199     }
200    
201     //wrap the stuntdoubles into a cell
202     for (sd = seleMan_.beginSelected(i); sd != NULL; sd = seleMan_.nextSelected(i)) {
203     Vector3d pos = sd->getPos();
204 gezelter 1078 if (usePeriodicBoundaryConditions_)
205     currentSnapshot_->wrapVector(pos);
206 xsun 955 sd->setPos(pos);
207 xsun 967 }
208 xsun 955
209     //determine which atom belongs to which grid
210     for (sd = seleMan_.beginSelected(i); sd != NULL; sd = seleMan_.nextSelected(i)) {
211     Vector3d pos = sd->getPos();
212     //int binNo = (pos.z() /deltaR_) - 1;
213 xsun 967 int binNoX = (int) ((pos.x() + halfBoxX_) / gridX_);
214     int binNoY = (int) ((pos.y() + halfBoxY_) / gridY_);
215 xsun 955 //std::cout << "pos.z = " << pos.z() << " halfBoxZ_ = " << halfBoxZ_ << " deltaR_ = " << deltaR_ << " binNo = " << binNo << "\n";
216     gridZ_[binNoX*nBinsY_+binNoY] += pos.z();
217     gridsample_[binNoX*nBinsY_+binNoY]++;
218     }
219    
220     // FFT stuff depends on nx and ny, so delay allocation until we have
221     // that information
222    
223     for(i = 0; i < nBinsX_; i++){
224     for(j = 0; j < nBinsY_; j++){
225     newindex = i * nBinsY_ + j;
226     if(gridsample_[newindex] > 0){
227 tim 963 gridZ_[newindex] = gridZ_[newindex] / (RealType)gridsample_[newindex];
228 xsun 955 }
229     }
230     }
231 gezelter 956
232 xsun 955 for (i=0; i< nBinsX_; i++) {
233     for(j=0; j< nBinsY_; j++) {
234     newindex = i*nBinsY_ + j;
235     if (gridsample_[newindex] == 0) {
236     // interpolate from surrounding points:
237    
238     interpsum = 0.0;
239     ninterp = 0;
240    
241     //point1 = bottom;
242    
243     px = i;
244     py = j - 1;
245     newp = px*nBinsY_ + py;
246     if ((py >= 0) && (gridsample_[newp] > 0)) {
247     interpsum += gridZ_[newp];
248     ninterp++;
249     }
250    
251     //point2 = top;
252    
253     px = i;
254     py = j + 1;
255     newp = px*nBinsY_ + py;
256     if ((py < nBinsY_) && (gridsample_[newp] > 0)) {
257     interpsum += gridZ_[newp];
258     ninterp++;
259     }
260    
261     //point3 = left;
262    
263     px = i - 1;
264     py = j;
265     newp = px*nBinsY_ + py;
266     if ((px >= 0) && (gridsample_[newp] > 0)) {
267     interpsum += gridZ_[newp];
268     ninterp++;
269     }
270    
271     //point4 = right;
272    
273     px = i + 1;
274     py = j;
275     newp = px*nBinsY_ + py;
276 gezelter 956 if ( (px < nBinsX_ ) && ( gridsample_[newp] > 0 )) {
277     interpsum += gridZ_[newp];
278     ninterp++;
279     }
280    
281 tim 963 value = interpsum / (RealType)ninterp;
282 xsun 955
283     gridZ_[newindex] = value;
284     }
285     }
286     }
287    
288     for (i=0; i < nBinsX_; i++) {
289     for (j=0; j < nBinsY_; j++) {
290     newindex = i*nBinsY_ + j;
291    
292     c_re(in[newindex]) = gridZ_[newindex];
293     c_im(in[newindex]) = 0.0;
294     }
295     }
296 xsun 967
297 gezelter 957 #ifdef HAVE_FFTW3_H
298     fftw_execute(p);
299     #else
300     fftwnd_one(p, in, out);
301     #endif
302 xsun 955
303     for (i=0; i< nBinsX_; i++) {
304     for(j=0; j< nBinsY_; j++) {
305     newindex = i*nBinsY_ + j;
306     mag[newindex] = pow(c_re(out[newindex]),2) + pow(c_im(out[newindex]),2);
307     }
308     }
309 xsun 967
310 gezelter 957 #ifdef HAVE_FFTW3_H
311 xsun 955 fftw_destroy_plan(p);
312 gezelter 957 #else
313     fftwnd_destroy_plan(p);
314     #endif
315 xsun 955 fftw_free(out);
316     fftw_free(in);
317    
318     for (i=0; i< (nBinsX_/2); i++) {
319     for(j=0; j< (nBinsY_/2); j++) {
320     index = i*nBinsY_ + j;
321     new_i = i + (nBinsX_/2);
322     new_j = j + (nBinsY_/2);
323     new_index = new_i*nBinsY_ + new_j;
324     newmag[new_index] = mag[index];
325     }
326     }
327    
328     for (i=(nBinsX_/2); i< nBinsX_; i++) {
329     for(j=0; j< (nBinsY_/2); j++) {
330     index = i*nBinsY_ + j;
331     new_i = i - (nBinsX_/2);
332     new_j = j + (nBinsY_/2);
333     new_index = new_i*nBinsY_ + new_j;
334     newmag[new_index] = mag[index];
335     }
336     }
337    
338     for (i=0; i< (nBinsX_/2); i++) {
339     for(j=(nBinsY_/2); j< nBinsY_; j++) {
340     index = i*nBinsY_ + j;
341     new_i = i + (nBinsX_/2);
342     new_j = j - (nBinsY_/2);
343     new_index = new_i*nBinsY_ + new_j;
344     newmag[new_index] = mag[index];
345     }
346     }
347    
348     for (i=(nBinsX_/2); i< nBinsX_; i++) {
349     for(j=(nBinsY_/2); j< nBinsY_; j++) {
350     index = i*nBinsY_ + j;
351     new_i = i - (nBinsX_/2);
352     new_j = j - (nBinsY_/2);
353     new_index = new_i*nBinsY_ + new_j;
354     newmag[new_index] = mag[index];
355     }
356     }
357    
358     maxfreqx = 1.0 / gridX_;
359     maxfreqy = 1.0 / gridY_;
360    
361     // printf("%lf\t%lf\t%lf\t%lf\n", dx, dy, maxfreqx, maxfreqy);
362    
363     maxfreq = sqrt(maxfreqx*maxfreqx + maxfreqy*maxfreqy);
364 tim 963 dfreq = maxfreq/(RealType)(nbins_-1);
365 xsun 955
366     //printf("%lf\n", dfreq);
367    
368     zero_freq_x = nBinsX_/2;
369     zero_freq_y = nBinsY_/2;
370    
371     for (i=0; i< nBinsX_; i++) {
372     for(j=0; j< nBinsY_; j++) {
373    
374 tim 963 freq_x = (RealType)(i - zero_freq_x)*maxfreqx*2 / nBinsX_;
375     freq_y = (RealType)(j - zero_freq_y)*maxfreqy*2 / nBinsY_;
376 xsun 955
377     freq = sqrt(freq_x*freq_x + freq_y*freq_y);
378    
379     whichbin = (int) (freq / dfreq);
380     newindex = i*nBinsY_ + j;
381     // printf("%d %d %lf %lf\n", whichbin, newindex, freq, dfreq);
382     bin[whichbin][istep] += newmag[newindex];
383     samples[whichbin][istep]++;
384     }
385     }
386    
387 gezelter 956 for ( i = 0; i < nbins_; i++) {
388 xsun 955 if ( samples[i][istep] > 0) {
389 xsun 969 bin[i][istep] = 4.0 * sqrt(bin[i][istep] / (RealType)samples[i][istep]) / (RealType)nBinsX_ / (RealType)nBinsY_;
390 xsun 955 }
391     }
392     }
393 xsun 967
394 gezelter 956 for (int i = 0; i < nbins_; i++) {
395     for (int j = 0; j < nFrames; j++) {
396 xsun 955 sum_bin[i] += bin[i][j];
397     sum_bin_sq[i] += bin[i][j] * bin[i][j];
398     }
399 tim 963 avg_bin[i] = sum_bin[i] / (RealType)nFrames;
400     avg_bin_sq[i] = sum_bin_sq[i] / (RealType)nFrames;
401 gezelter 956 for (int j = 0; j < nFrames; j++) {
402 xsun 955 errbin_sum[i] += pow((bin[i][j] - avg_bin[i]), 2);
403     errbin_sum_sq[i] += pow((bin[i][j] * bin[i][j] - avg_bin_sq[i]), 2);
404     }
405 tim 963 errbin[i] = sqrt( errbin_sum[i] / (RealType)nFrames );
406     errbin_sq[i] = sqrt( errbin_sum_sq[i] / (RealType)nFrames );
407 xsun 955 }
408 xsun 967
409 gezelter 956 printSpectrum();
410 xsun 967
411 gezelter 956 #else
412     sprintf(painCave.errMsg, "Hxy: FFTW support was not compiled in!\n");
413     painCave.isFatal = 1;
414     simError();
415 xsun 955
416 gezelter 956 #endif
417 xsun 967 }
418 xsun 955
419     void Hxy::printSpectrum() {
420     std::ofstream rdfStream(outputFilename_.c_str());
421     if (rdfStream.is_open()) {
422 xsun 967
423     for (int i = 0; i < nbins_; ++i) {
424 xsun 955 if ( avg_bin[i] > 0 ){
425 xsun 967 rdfStream << (RealType)i * dfreq << "\t"
426 xsun 955 <<pow(avg_bin[i], 2)<<"\t"
427     <<errbin_sq[i]<<"\t"
428     <<avg_bin[i]<<"\t"
429     <<errbin[i]<<"\n";
430     }
431     }
432     } else {
433    
434     sprintf(painCave.errMsg, "Hxy: unable to open %s\n", outputFilename_.c_str());
435     painCave.isFatal = 1;
436     simError();
437     }
438 xsun 967
439 xsun 955 rdfStream.close();
440 xsun 967
441 xsun 955 }
442    
443     }

Properties

Name Value
svn:keywords Author Id Revision Date