ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/applications/staticProps/Hxy.cpp
Revision: 3054
Committed: Wed Oct 18 21:58:48 2006 UTC (17 years, 8 months ago) by gezelter
File size: 13004 byte(s)
Log Message:
fixing a wrapVector problem in staticProps, also making Shifted force
and electrostatic damping the default behavior

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