ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/applications/staticProps/Hxy.cpp
Revision: 2753
Committed: Tue May 16 20:38:23 2006 UTC (18 years, 2 months ago) by gezelter
File size: 11922 byte(s)
Log Message:
Autoconf fixes for FFTW.  Multiple FFTW version support.

File Contents

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