ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/applications/staticProps/Hxy.cpp
Revision: 2752
Committed: Tue May 16 02:06:37 2006 UTC (18 years, 3 months ago) by gezelter
File size: 11632 byte(s)
Log Message:
*** empty log message ***

File Contents

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