OpenMD 3.0
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
Hxy.cpp
1/*
2 * Copyright (c) 2004-present, The University of Notre Dame. All rights
3 * reserved.
4 *
5 * Redistribution and use in source and binary forms, with or without
6 * modification, are permitted provided that the following conditions are met:
7 *
8 * 1. Redistributions of source code must retain the above copyright notice,
9 * this list of conditions and the following disclaimer.
10 *
11 * 2. Redistributions in binary form must reproduce the above copyright notice,
12 * this list of conditions and the following disclaimer in the documentation
13 * and/or other materials provided with the distribution.
14 *
15 * 3. Neither the name of the copyright holder nor the names of its
16 * contributors may be used to endorse or promote products derived from
17 * this software without specific prior written permission.
18 *
19 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
20 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
21 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
22 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
23 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
24 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
25 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
26 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
27 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
28 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
29 * POSSIBILITY OF SUCH DAMAGE.
30 *
31 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
32 * research, please cite the appropriate papers when you publish your
33 * work. Good starting points are:
34 *
35 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
36 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
37 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
38 * [4] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
39 * [5] Kuang & Gezelter, Mol. Phys., 110, 691-701 (2012).
40 * [6] Lamichhane, Gezelter & Newman, J. Chem. Phys. 141, 134109 (2014).
41 * [7] Lamichhane, Newman & Gezelter, J. Chem. Phys. 141, 134110 (2014).
42 * [8] Bhattarai, Newman & Gezelter, Phys. Rev. B 99, 094106 (2019).
43 */
44
45/* Calculates the undulation spectrum of the lipid membrance. */
46
47#include "applications/staticProps/Hxy.hpp"
48
49#include <algorithm>
50#include <cmath>
51#include <cstdio>
52#include <cstdlib>
53#include <cstring>
54#include <fstream>
55
56#include "io/DumpReader.hpp"
58#include "types/LennardJonesAdapter.hpp"
59#include "utils/Accumulator.hpp"
60#include "utils/AccumulatorView.hpp"
61#include "utils/BaseAccumulator.hpp"
62#include "utils/Utility.hpp"
63#include "utils/simError.h"
64
65using namespace OpenMD::Utils;
66
67namespace OpenMD {
68
69 Hxy::Hxy(SimInfo* info, const std::string& filename, const std::string& sele,
70 int nbins_x, int nbins_y, int nbins_z, int nrbins) :
71 StaticAnalyser(info, filename, nrbins),
72 selectionScript_(sele), evaluator_(info), seleMan_(info),
73 nBinsX_(nbins_x), nBinsY_(nbins_y), nBinsZ_(nbins_z) {
74 evaluator_.loadScriptString(sele);
75 if (!evaluator_.isDynamic()) {
76 seleMan_.setSelectionSet(evaluator_.evaluate());
77 }
78
79 // dens_ stores the local density, rho(x,y,z) on a 3-D grid
80 dens_.resize(nBinsX_);
81 // bin stores the upper and lower surface cutoff locations (z) for
82 // a column through grid location x,y
83 minHeight_.resize(nBinsX_);
84 maxHeight_.resize(nBinsX_);
85
86 for (unsigned int i = 0; i < nBinsX_; i++) {
87 dens_[i].resize(nBinsY_);
88 minHeight_[i].resize(nBinsY_);
89 maxHeight_[i].resize(nBinsY_);
90 for (unsigned int j = 0; j < nBinsY_; j++) {
91 dens_[i][j].resize(nBinsZ_);
92 }
93 }
94
95 mag1.resize(nBinsX_ * nBinsY_);
96 newmag1.resize(nBinsX_ * nBinsY_);
97 mag2.resize(nBinsX_ * nBinsY_);
98 newmag2.resize(nBinsX_ * nBinsY_);
99
100 // Pre-load the OutputData
101 data_.resize(Hxy::ENDINDEX);
102
103 OutputData freq;
104 freq.units = "Angstroms^-1";
105 freq.title = "Spatial Frequency";
106 freq.dataHandling = DataHandling::Average;
107 for (unsigned int i = 0; i < nBins_; i++)
108 freq.accumulator.push_back(
109 std::make_unique<AccumulatorView<RealAccumulator>>());
110 data_[FREQUECY] = std::move(freq);
111
112 OutputData top;
113 top.units = "Angstroms";
114 top.title = "Hxy (Upper surface)";
115 freq.dataHandling = DataHandling::Max;
116 for (unsigned int i = 0; i < nBins_; i++)
117 top.accumulator.push_back(
118 std::make_unique<AccumulatorView<RealAccumulator>>());
119 data_[TOP] = std::move(top);
120
121 OutputData bottom;
122 bottom.units = "Angstroms";
123 bottom.title = "Hxy (Lower surface)";
124 freq.dataHandling = DataHandling::Max;
125 for (unsigned int i = 0; i < nBins_; i++)
126 bottom.accumulator.push_back(
127 std::make_unique<AccumulatorView<RealAccumulator>>());
128 data_[BOTTOM] = std::move(bottom);
129
130 setOutputName(getPrefix(filename) + ".Hxy");
131 }
132
133 void Hxy::process() {
134#if defined(HAVE_FFTW_H) || defined(HAVE_DFFTW_H) || defined(HAVE_FFTW3_H)
135 StuntDouble* sd;
136 int ii;
137 bool usePeriodicBoundaryConditions_ =
138 info_->getSimParams()->getUsePeriodicBoundaryConditions();
139 std::cerr << "usePeriodicBoundaryConditions_ = "
140 << usePeriodicBoundaryConditions_ << "\n";
141
142 DumpReader reader(info_, dumpFilename_);
143 int nFrames = reader.getNFrames();
144 nProcessed_ = nFrames / step_;
145
146 for (int istep = 0; istep < nFrames; istep += step_) {
147 for (unsigned int i = 0; i < nBinsX_; i++) {
148 std::fill(minHeight_[i].begin(), minHeight_[i].end(), 0.0);
149 std::fill(maxHeight_[i].begin(), maxHeight_[i].end(), 0.0);
150 for (unsigned int j = 0; j < nBinsY_; j++) {
151 std::fill(dens_[i][j].begin(), dens_[i][j].end(), 0.0);
152 }
153 }
154
155 reader.readFrame(istep);
156 currentSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot();
157 if (evaluator_.isDynamic()) {
158 seleMan_.setSelectionSet(evaluator_.evaluate());
159 }
160
161#ifdef HAVE_FFTW3_H
162 fftw_plan p1, p2;
163#else
164 fftwnd_plan p1, p2;
165#endif
166 fftw_complex *in1, *in2, *out1, *out2;
167
168 in1 = (fftw_complex*)fftw_malloc(sizeof(fftw_complex) *
169 (nBinsX_ * nBinsY_));
170 out1 = (fftw_complex*)fftw_malloc(sizeof(fftw_complex) *
171 (nBinsX_ * nBinsY_));
172 in2 = (fftw_complex*)fftw_malloc(sizeof(fftw_complex) *
173 (nBinsX_ * nBinsY_));
174 out2 = (fftw_complex*)fftw_malloc(sizeof(fftw_complex) *
175 (nBinsX_ * nBinsY_));
176
177#ifdef HAVE_FFTW3_H
178 p1 = fftw_plan_dft_2d(nBinsX_, nBinsY_, in1, out1, FFTW_FORWARD,
179 FFTW_ESTIMATE);
180 p2 = fftw_plan_dft_2d(nBinsX_, nBinsY_, in2, out2, FFTW_FORWARD,
181 FFTW_ESTIMATE);
182#else
183 p1 = fftw2d_create_plan(nBinsX_, nBinsY_, FFTW_FORWARD, FFTW_ESTIMATE);
184 p2 = fftw2d_create_plan(nBinsX_, nBinsY_, FFTW_FORWARD, FFTW_ESTIMATE);
185#endif
186
187 Mat3x3d hmat = currentSnapshot_->getHmat();
188 Mat3x3d invBox = currentSnapshot_->getInvHmat();
189 RealType lenX_ = hmat(0, 0);
190 RealType lenY_ = hmat(1, 1);
191 RealType lenZ_ = hmat(2, 2);
192
193 RealType x, y, z, dx, dy, dz;
194 RealType sigma, rcut;
195 int di, dj, dk, ibin, jbin, kbin;
196 int igrid, jgrid, kgrid;
197 Vector3d scaled;
198
199 dx = lenX_ / nBinsX_;
200 dy = lenY_ / nBinsY_;
201 dz = lenZ_ / nBinsZ_;
202
203 for (sd = seleMan_.beginSelected(ii); sd != NULL;
204 sd = seleMan_.nextSelected(ii)) {
205 if (sd->isAtom()) {
206 Atom* atom = static_cast<Atom*>(sd);
207 Vector3d pos = sd->getPos();
208 LennardJonesAdapter lja = LennardJonesAdapter(atom->getAtomType());
209 // For SPC/E water, this yields the Willard-Chandler
210 // distance of 2.4 Angstroms:
211 sigma = lja.getSigma() * 0.758176459;
212 rcut = 3.0 * sigma;
213
214 // scaled positions relative to the box vectors
215 // -> the atom's position in numbers of box lengths (more accurately
216 // box vectors)
217 scaled = invBox * pos;
218
219 // wrap the vector back into the unit box by subtracting
220 // integer box numbers
221 for (int j = 0; j < 3; j++) {
222 scaled[j] -= roundMe(scaled[j]);
223 scaled[j] += 0.5;
224 // Handle the special case when an object is exactly on
225 // the boundary (a scaled coordinate of 1.0 is the same as
226 // scaled coordinate of 0.0)
227 if (scaled[j] >= 1.0) scaled[j] -= 1.0;
228 }
229 // find ijk-indices of voxel that atom is in.
230 ibin = nBinsX_ * scaled.x();
231 jbin = nBinsY_ * scaled.y();
232 kbin = nBinsZ_ * scaled.z();
233
234 di = (int)(rcut / dx);
235 dj = (int)(rcut / dy);
236 dk = (int)(rcut / dz);
237
238 for (int i = -di; i <= di; i++) {
239 igrid = ibin + i;
240 while (igrid >= int(nBinsX_)) {
241 igrid -= int(nBinsX_);
242 }
243 while (igrid < 0) {
244 igrid += int(nBinsX_);
245 }
246
247 x = lenX_ * (RealType(i) / RealType(nBinsX_));
248
249 for (int j = -dj; j <= dj; j++) {
250 jgrid = jbin + j;
251 while (jgrid >= int(nBinsY_)) {
252 jgrid -= int(nBinsY_);
253 }
254 while (jgrid < 0) {
255 jgrid += int(nBinsY_);
256 }
257
258 y = lenY_ * (RealType(j) / RealType(nBinsY_));
259
260 for (int k = -dk; k <= dk; k++) {
261 kgrid = kbin + k;
262 while (kgrid >= int(nBinsZ_)) {
263 kgrid -= int(nBinsZ_);
264 }
265 while (kgrid < 0) {
266 kgrid += int(nBinsZ_);
267 }
268
269 z = lenZ_ * (RealType(k) / RealType(nBinsZ_));
270
271 RealType dist = sqrt(x * x + y * y + z * z);
272
273 dens_[igrid][jgrid][kgrid] += getDensity(dist, sigma, rcut);
274 }
275 }
276 }
277 }
278 }
279
280 RealType maxDens(0.0);
281 for (unsigned int i = 0; i < nBinsX_; i++) {
282 for (unsigned int j = 0; j < nBinsY_; j++) {
283 for (unsigned int k = 0; k < nBinsZ_; k++) {
284 if (dens_[i][j][k] > maxDens) maxDens = dens_[i][j][k];
285 }
286 }
287 }
288
289 RealType threshold = maxDens / 2.0;
290 RealType z0, z1, h0, h1;
291 std::cerr << "maxDens = "
292 << "\t" << maxDens << "\n";
293 std::cerr << "threshold value = "
294 << "\t" << threshold << "\n";
295
296 for (unsigned int i = 0; i < nBinsX_; i++) {
297 for (unsigned int j = 0; j < nBinsY_; j++) {
298 // There are two cases if we are periodic in z and the
299 // density is localized in z. Either we're starting below
300 // the isodensity, or above it:
301 // ______ _______ ______
302 // ____/ \_____ or: \______/
303 //
304 // In either case, there are two crossings of the
305 // isodensity.
306
307 bool minFound = false;
308 bool maxFound = false;
309
310 if (dens_[i][j][0] < threshold) {
311 for (unsigned int k = 0; k < nBinsZ_ - 1; k++) {
312 z0 = lenZ_ * (RealType(k) / RealType(nBinsZ_));
313 z1 = lenZ_ * (RealType(k + 1) / RealType(nBinsZ_));
314 h0 = dens_[i][j][k];
315 h1 = dens_[i][j][k + 1];
316
317 if (h0 < threshold && h1 > threshold && !minFound) {
318 // simple linear interpolation to find the height:
319 minHeight_[i][j] =
320 z0 + (z1 - z0) * (threshold - h0) / (h1 - h0);
321 minFound = true;
322 }
323 if (h0 > threshold && h1 < threshold && minFound) {
324 // simple linear interpolation to find the height:
325 maxHeight_[i][j] =
326 z0 + (z1 - z0) * (threshold - h0) / (h1 - h0);
327 maxFound = true;
328 }
329 }
330
331 } else {
332 for (unsigned int k = 0; k < nBinsZ_ - 1; k++) {
333 z0 = lenZ_ * (RealType(k) / RealType(nBinsZ_));
334 z1 = lenZ_ * (RealType(k + 1) / RealType(nBinsZ_));
335 h0 = dens_[i][j][k];
336 h1 = dens_[i][j][k + 1];
337
338 if (h0 > threshold && h1 < threshold && !maxFound) {
339 // simple linear interpolation to find the height:
340 maxHeight_[i][j] =
341 z0 + (z1 - z0) * (threshold - h0) / (h1 - h0);
342 maxFound = true;
343 }
344 if (h0 < threshold && h1 > threshold && maxFound) {
345 // simple linear interpolation to find the height:
346 minHeight_[i][j] =
347 z0 + (z1 - z0) * (threshold - h0) / (h1 - h0);
348 }
349 }
350 }
351 }
352 }
353
354 RealType minBar = 0.0;
355 RealType maxBar = 0.0;
356 int count = 0;
357 for (unsigned int i = 0; i < nBinsX_; i++) {
358 for (unsigned int j = 0; j < nBinsY_; j++) {
359 // if (minHeight_[i][j] < 0.0) std::cerr << "minHeight[i][j] = " <<
360 // "\t"
361 // << minHeight_[i][j] << "\n";
362 minBar += minHeight_[i][j];
363 maxBar += maxHeight_[i][j];
364 count++;
365 }
366 }
367 minBar /= count;
368 maxBar /= count;
369
370 std::cerr << "bottomSurf = " << minBar << "\ttopSurf = " << maxBar
371 << "\n";
372 int newindex;
373 // RealType Lx = 10.0;
374 // RealType Ly = 10.0;
375 for (unsigned int i = 0; i < nBinsX_; i++) {
376 for (unsigned int j = 0; j < nBinsY_; j++) {
377 newindex = i * nBinsY_ + j;
378 // if(minHeight_[i][j] < 0.0) std::cerr << minHeight_[i][j] << "\t";
379 c_re(in1[newindex]) = maxHeight_[i][j] - maxBar;
380 // c_re(in1[newindex]) = 2.0*cos(2.0*Constants::PI*i/Lx/Ly);
381 c_im(in1[newindex]) = 0.0;
382 c_re(in2[newindex]) =
383 minHeight_[i][j] - minBar; // this was the original line.
384 // c_re(in2[newindex]) = 1.0*cos(2.0*Constants::PI*i/Lx/Ly);
385 c_im(in2[newindex]) = 0.0;
386 }
387 }
388
389#ifdef HAVE_FFTW3_H
390 fftw_execute(p1);
391 fftw_execute(p2);
392#else
393 fftwnd_one(p1, in1, out1);
394 fftwnd_one(p2, in2, out2);
395#endif
396
397 for (unsigned int i = 0; i < nBinsX_; i++) {
398 for (unsigned int j = 0; j < nBinsY_; j++) {
399 newindex = i * nBinsY_ + j;
400 mag1[newindex] = sqrt(pow(c_re(out1[newindex]), 2) +
401 pow(c_im(out1[newindex]), 2)) /
402 (RealType(nBinsX_ * nBinsY_));
403 mag2[newindex] = sqrt(pow(c_re(out2[newindex]), 2) +
404 pow(c_im(out2[newindex]), 2)) /
405 (RealType(nBinsX_ * nBinsY_));
406 }
407 }
408
409#ifdef HAVE_FFTW3_H
410 fftw_destroy_plan(p1);
411 fftw_destroy_plan(p2);
412#else
413 fftwnd_destroy_plan(p1);
414 fftwnd_destroy_plan(p2);
415#endif
416 fftw_free(out1);
417 fftw_free(in1);
418 fftw_free(out2);
419 fftw_free(in2);
420
421 int index, new_i, new_j, new_index;
422 for (unsigned int i = 0; i < (nBinsX_ / 2); i++) {
423 for (unsigned int j = 0; j < (nBinsY_ / 2); j++) {
424 index = i * nBinsY_ + j;
425 new_i = i + (nBinsX_ / 2);
426 new_j = j + (nBinsY_ / 2);
427 new_index = new_i * nBinsY_ + new_j;
428 newmag1[new_index] = mag1[index];
429 newmag2[new_index] = mag2[index];
430 }
431 }
432
433 for (unsigned int i = (nBinsX_ / 2); i < nBinsX_; i++) {
434 for (unsigned int j = 0; j < (nBinsY_ / 2); j++) {
435 index = i * nBinsY_ + j;
436 new_i = i - (nBinsX_ / 2);
437 new_j = j + (nBinsY_ / 2);
438 new_index = new_i * nBinsY_ + new_j;
439 newmag1[new_index] = mag1[index];
440 newmag2[new_index] = mag2[index];
441 }
442 }
443
444 for (unsigned int i = 0; i < (nBinsX_ / 2); i++) {
445 for (unsigned int j = (nBinsY_ / 2); j < nBinsY_; j++) {
446 index = i * nBinsY_ + j;
447 new_i = i + (nBinsX_ / 2);
448 new_j = j - (nBinsY_ / 2);
449 new_index = new_i * nBinsY_ + new_j;
450 newmag1[new_index] = mag1[index];
451 newmag2[new_index] = mag2[index];
452 }
453 }
454
455 for (unsigned int i = (nBinsX_ / 2); i < nBinsX_; i++) {
456 for (unsigned int j = (nBinsY_ / 2); j < nBinsY_; j++) {
457 index = i * nBinsY_ + j;
458 new_i = i - (nBinsX_ / 2);
459 new_j = j - (nBinsY_ / 2);
460 new_index = new_i * nBinsY_ + new_j;
461 newmag1[new_index] = mag1[index];
462 newmag2[new_index] = mag2[index];
463 }
464 }
465
466 /*
467 for (unsigned int i=0; i< nBinsX_; i++) {
468 for(unsigned int j=0; j< nBinsY_; j++) {
469 newindex = i*nBinsY_ + j;
470 std::cout << newmag1[newindex] << "\t";
471 }
472 std::cout << "\n";
473 }
474 */
475
476 RealType maxfreqx = RealType(nBinsX_) / lenX_;
477 RealType maxfreqy = RealType(nBinsY_) / lenY_;
478
479 RealType maxfreq = sqrt(maxfreqx * maxfreqx + maxfreqy * maxfreqy);
480 dfreq_ = maxfreq / (RealType)(nBins_ - 1);
481
482 int zero_freq_x = nBinsX_ / 2;
483 int zero_freq_y = nBinsY_ / 2;
484
485 for (int i = 0; i < (int)nBinsX_; i++) {
486 for (int j = 0; j < (int)nBinsY_; j++) {
487 RealType freq_x =
488 (RealType)(i - zero_freq_x) * maxfreqx / (RealType)nBinsX_;
489 RealType freq_y =
490 (RealType)(j - zero_freq_y) * maxfreqy / (RealType)nBinsY_;
491
492 RealType freq = sqrt(freq_x * freq_x + freq_y * freq_y);
493
494 unsigned int whichbin = (unsigned int)(freq / dfreq_);
495
496 newindex = i * nBinsY_ + j;
497
498 data_[FREQUECY].accumulator[whichbin]->add(freq);
499 data_[TOP].accumulator[whichbin]->add(newmag1[newindex]);
500 data_[BOTTOM].accumulator[whichbin]->add(newmag2[newindex]);
501 }
502 }
503 }
504
505 writeOutput();
506
507#else
508 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
509 "Hxy: FFTW support was not compiled in!\n");
510 painCave.isFatal = 1;
511 simError();
512
513#endif
514 }
515
516 RealType Hxy::getDensity(RealType r, RealType sigma, RealType rcut) {
517 RealType sigma2 = sigma * sigma;
518 RealType dens =
519 exp(-r * r / (sigma2 * 2.0)) / (pow(2.0 * Constants::PI * sigma2, 3));
520 RealType dcut = exp(-rcut * rcut / (sigma2 * 2.0)) /
521 (pow(2.0 * Constants::PI * sigma2, 3));
522 if (r < rcut)
523 return dens - dcut;
524 else
525 return 0.0;
526 }
527} // namespace OpenMD
One of the heavy-weight classes of OpenMD, SimInfo maintains objects and variables relating to the cu...
Definition SimInfo.hpp:93
"Don't move, or you're dead! Stand up! Captain, we've got them!"
Vector3d getPos()
Returns the current position of this stuntDouble.
bool isAtom()
Tests if this stuntDouble is an atom.
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
std::string getPrefix(const std::string &str)