| 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.6 2006-05-23 21:12:45 xsun 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 | currentSnapshot_->wrapVector(pos); | 
| 208 | sd->setPos(pos); | 
| 209 | } | 
| 210 |  | 
| 211 | //determine which atom belongs to which grid | 
| 212 | for (sd = seleMan_.beginSelected(i); sd != NULL; sd = seleMan_.nextSelected(i)) { | 
| 213 | Vector3d pos = sd->getPos(); | 
| 214 | //int binNo = (pos.z() /deltaR_) - 1; | 
| 215 | int binNoX = (int) ((pos.x() + halfBoxX_) / gridX_); | 
| 216 | int binNoY = (int) ((pos.y() + halfBoxY_) / gridY_); | 
| 217 | //std::cout << "pos.z = " << pos.z() << " halfBoxZ_ = " << halfBoxZ_ << " deltaR_ = "  << deltaR_ << " binNo = " << binNo << "\n"; | 
| 218 | gridZ_[binNoX*nBinsY_+binNoY] += pos.z(); | 
| 219 | gridsample_[binNoX*nBinsY_+binNoY]++; | 
| 220 | } | 
| 221 |  | 
| 222 | // FFT stuff depends on nx and ny, so delay allocation until we have | 
| 223 | // that information | 
| 224 |  | 
| 225 | for(i = 0; i < nBinsX_; i++){ | 
| 226 | for(j = 0; j < nBinsY_; j++){ | 
| 227 | newindex = i * nBinsY_ + j; | 
| 228 | if(gridsample_[newindex] > 0){ | 
| 229 | gridZ_[newindex] = gridZ_[newindex] / (RealType)gridsample_[newindex]; | 
| 230 | } | 
| 231 | } | 
| 232 | } | 
| 233 |  | 
| 234 | for (i=0; i< nBinsX_; i++) { | 
| 235 | for(j=0; j< nBinsY_; j++) { | 
| 236 | newindex = i*nBinsY_ + j; | 
| 237 | if (gridsample_[newindex] == 0) { | 
| 238 | // interpolate from surrounding points: | 
| 239 |  | 
| 240 | interpsum = 0.0; | 
| 241 | ninterp = 0; | 
| 242 |  | 
| 243 | //point1 = bottom; | 
| 244 |  | 
| 245 | px = i; | 
| 246 | py = j - 1; | 
| 247 | newp = px*nBinsY_ + py; | 
| 248 | if ((py >= 0) && (gridsample_[newp] > 0)) { | 
| 249 | interpsum += gridZ_[newp]; | 
| 250 | ninterp++; | 
| 251 | } | 
| 252 |  | 
| 253 | //point2 = top; | 
| 254 |  | 
| 255 | px = i; | 
| 256 | py = j + 1; | 
| 257 | newp = px*nBinsY_ + py; | 
| 258 | if ((py < nBinsY_) && (gridsample_[newp] > 0)) { | 
| 259 | interpsum += gridZ_[newp]; | 
| 260 | ninterp++; | 
| 261 | } | 
| 262 |  | 
| 263 | //point3 = left; | 
| 264 |  | 
| 265 | px = i - 1; | 
| 266 | py = j; | 
| 267 | newp = px*nBinsY_ + py; | 
| 268 | if ((px >= 0) && (gridsample_[newp] > 0)) { | 
| 269 | interpsum += gridZ_[newp]; | 
| 270 | ninterp++; | 
| 271 | } | 
| 272 |  | 
| 273 | //point4 = right; | 
| 274 |  | 
| 275 | px = i + 1; | 
| 276 | py = j; | 
| 277 | newp = px*nBinsY_ + py; | 
| 278 | if ( (px < nBinsX_ ) && ( gridsample_[newp] > 0 )) { | 
| 279 | interpsum += gridZ_[newp]; | 
| 280 | ninterp++; | 
| 281 | } | 
| 282 |  | 
| 283 | value = interpsum / (RealType)ninterp; | 
| 284 |  | 
| 285 | gridZ_[newindex] = value; | 
| 286 | } | 
| 287 | } | 
| 288 | } | 
| 289 |  | 
| 290 | for (i=0; i < nBinsX_; i++) { | 
| 291 | for (j=0; j < nBinsY_; j++) { | 
| 292 | newindex = i*nBinsY_ + j; | 
| 293 |  | 
| 294 | c_re(in[newindex]) = gridZ_[newindex]; | 
| 295 | c_im(in[newindex]) = 0.0; | 
| 296 | } | 
| 297 | } | 
| 298 |  | 
| 299 | #ifdef HAVE_FFTW3_H | 
| 300 | fftw_execute(p); | 
| 301 | #else | 
| 302 | fftwnd_one(p, in, out); | 
| 303 | #endif | 
| 304 |  | 
| 305 | for (i=0; i< nBinsX_; i++) { | 
| 306 | for(j=0; j< nBinsY_; j++) { | 
| 307 | newindex = i*nBinsY_ + j; | 
| 308 | mag[newindex] = pow(c_re(out[newindex]),2) + pow(c_im(out[newindex]),2); | 
| 309 | } | 
| 310 | } | 
| 311 |  | 
| 312 | #ifdef HAVE_FFTW3_H | 
| 313 | fftw_destroy_plan(p); | 
| 314 | #else | 
| 315 | fftwnd_destroy_plan(p); | 
| 316 | #endif | 
| 317 | fftw_free(out); | 
| 318 | fftw_free(in); | 
| 319 |  | 
| 320 | for (i=0; i< (nBinsX_/2); i++) { | 
| 321 | for(j=0; j< (nBinsY_/2); 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 | for (i=(nBinsX_/2); i< nBinsX_; i++) { | 
| 331 | for(j=0; j< (nBinsY_/2); j++) { | 
| 332 | index = i*nBinsY_ + j; | 
| 333 | new_i = i - (nBinsX_/2); | 
| 334 | new_j = j + (nBinsY_/2); | 
| 335 | new_index = new_i*nBinsY_ + new_j; | 
| 336 | newmag[new_index] = mag[index]; | 
| 337 | } | 
| 338 | } | 
| 339 |  | 
| 340 | for (i=0; i< (nBinsX_/2); i++) { | 
| 341 | for(j=(nBinsY_/2); j< nBinsY_; j++) { | 
| 342 | index = i*nBinsY_ + j; | 
| 343 | new_i = i + (nBinsX_/2); | 
| 344 | new_j = j - (nBinsY_/2); | 
| 345 | new_index = new_i*nBinsY_ + new_j; | 
| 346 | newmag[new_index] = mag[index]; | 
| 347 | } | 
| 348 | } | 
| 349 |  | 
| 350 | for (i=(nBinsX_/2); i< nBinsX_; i++) { | 
| 351 | for(j=(nBinsY_/2); j< nBinsY_; j++) { | 
| 352 | index = i*nBinsY_ + j; | 
| 353 | new_i = i - (nBinsX_/2); | 
| 354 | new_j = j - (nBinsY_/2); | 
| 355 | new_index = new_i*nBinsY_ + new_j; | 
| 356 | newmag[new_index] = mag[index]; | 
| 357 | } | 
| 358 | } | 
| 359 |  | 
| 360 | maxfreqx = 1.0 / gridX_; | 
| 361 | maxfreqy = 1.0 / gridY_; | 
| 362 |  | 
| 363 | //  printf("%lf\t%lf\t%lf\t%lf\n", dx, dy, maxfreqx, maxfreqy); | 
| 364 |  | 
| 365 | maxfreq = sqrt(maxfreqx*maxfreqx + maxfreqy*maxfreqy); | 
| 366 | dfreq = maxfreq/(RealType)(nbins_-1); | 
| 367 |  | 
| 368 | //printf("%lf\n", dfreq); | 
| 369 |  | 
| 370 | zero_freq_x = nBinsX_/2; | 
| 371 | zero_freq_y = nBinsY_/2; | 
| 372 |  | 
| 373 | for (i=0; i< nBinsX_; i++) { | 
| 374 | for(j=0; j< nBinsY_; j++) { | 
| 375 |  | 
| 376 | freq_x = (RealType)(i - zero_freq_x)*maxfreqx*2 / nBinsX_; | 
| 377 | freq_y = (RealType)(j - zero_freq_y)*maxfreqy*2 / nBinsY_; | 
| 378 |  | 
| 379 | freq = sqrt(freq_x*freq_x + freq_y*freq_y); | 
| 380 |  | 
| 381 | whichbin = (int) (freq / dfreq); | 
| 382 | newindex = i*nBinsY_ + j; | 
| 383 | //    printf("%d %d %lf %lf\n", whichbin, newindex, freq, dfreq); | 
| 384 | bin[whichbin][istep] += newmag[newindex]; | 
| 385 | samples[whichbin][istep]++; | 
| 386 | } | 
| 387 | } | 
| 388 |  | 
| 389 | for ( i = 0; i < nbins_; i++) { | 
| 390 | if ( samples[i][istep] > 0) { | 
| 391 | bin[i][istep] = 4.0 * sqrt(bin[i][istep] / (RealType)samples[i][istep]) / (RealType)nBinsX_ / (RealType)nBinsY_; | 
| 392 | } | 
| 393 | } | 
| 394 | } | 
| 395 |  | 
| 396 | for (int i = 0; i < nbins_; i++) { | 
| 397 | for (int j = 0; j < nFrames; j++) { | 
| 398 | sum_bin[i] += bin[i][j]; | 
| 399 | sum_bin_sq[i] += bin[i][j] * bin[i][j]; | 
| 400 | } | 
| 401 | avg_bin[i] = sum_bin[i] / (RealType)nFrames; | 
| 402 | avg_bin_sq[i] = sum_bin_sq[i] / (RealType)nFrames; | 
| 403 | for (int j = 0; j < nFrames; j++) { | 
| 404 | errbin_sum[i] += pow((bin[i][j] - avg_bin[i]), 2); | 
| 405 | errbin_sum_sq[i] += pow((bin[i][j] * bin[i][j] - avg_bin_sq[i]), 2); | 
| 406 | } | 
| 407 | errbin[i] = sqrt( errbin_sum[i] / (RealType)nFrames ); | 
| 408 | errbin_sq[i] = sqrt( errbin_sum_sq[i] / (RealType)nFrames ); | 
| 409 | } | 
| 410 |  | 
| 411 | printSpectrum(); | 
| 412 |  | 
| 413 | #else | 
| 414 | sprintf(painCave.errMsg, "Hxy: FFTW support was not compiled in!\n"); | 
| 415 | painCave.isFatal = 1; | 
| 416 | simError(); | 
| 417 |  | 
| 418 | #endif | 
| 419 | } | 
| 420 |  | 
| 421 | void Hxy::printSpectrum() { | 
| 422 | std::ofstream rdfStream(outputFilename_.c_str()); | 
| 423 | if (rdfStream.is_open()) { | 
| 424 |  | 
| 425 | for (int i = 0; i < nbins_; ++i) { | 
| 426 | if ( avg_bin[i] > 0 ){ | 
| 427 | rdfStream << (RealType)i * dfreq << "\t" | 
| 428 | <<pow(avg_bin[i], 2)<<"\t" | 
| 429 | <<errbin_sq[i]<<"\t" | 
| 430 | <<avg_bin[i]<<"\t" | 
| 431 | <<errbin[i]<<"\n"; | 
| 432 | } | 
| 433 | } | 
| 434 | } else { | 
| 435 |  | 
| 436 | sprintf(painCave.errMsg, "Hxy: unable to open %s\n", outputFilename_.c_str()); | 
| 437 | painCave.isFatal = 1; | 
| 438 | simError(); | 
| 439 | } | 
| 440 |  | 
| 441 | rdfStream.close(); | 
| 442 |  | 
| 443 | } | 
| 444 |  | 
| 445 | } |