OpenMD 3.0
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
ZConsWriter.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#include "io/ZConsWriter.hpp"
46
47#include <algorithm>
48#include <iostream>
49#include <vector>
50
51#ifdef IS_MPI
52#include <mpi.h>
53#endif
54
55#include "utils/simError.h"
56
57namespace OpenMD {
58 ZConsWriter::ZConsWriter(SimInfo* info, const std::string& filename) :
59 info_(info) {
60 // use a primary - secondary model, only the primary node writes
61 // to disk
62#ifdef IS_MPI
63 if (worldRank == 0) {
64#endif
65 output_.open(filename.c_str());
66
67 if (!output_) {
68 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
69 "Could not open %s for z constrain output_ \n",
70 filename.c_str());
71 painCave.isFatal = 1;
72 simError();
73 }
74
75 output_ << "//time(fs)" << std::endl;
76 output_ << "//number of fixed z-constrain molecules" << std::endl;
77 output_ << "//global Index of molecule\tzconstrain force\tcurrentZPos"
78 << std::endl;
79#ifdef IS_MPI
80 }
81#endif
82 }
83
84 ZConsWriter::~ZConsWriter() {
85#ifdef IS_MPI
86 if (worldRank == 0) {
87#endif
88 output_.close();
89#ifdef IS_MPI
90 }
91#endif
92 }
93
94 void ZConsWriter::writeFZ(const std::list<ZconstraintMol>& fixedZmols) {
95#ifndef IS_MPI
96 output_ << info_->getSnapshotManager()->getCurrentSnapshot()->getTime()
97 << std::endl;
98 output_ << fixedZmols.size() << std::endl;
99
100 std::list<ZconstraintMol>::const_iterator i;
101 for (i = fixedZmols.begin(); i != fixedZmols.end(); ++i) {
102 output_ << i->mol->getGlobalIndex() << "\t" << i->fz << "\t" << i->zpos
103 << "\t" << i->param.zTargetPos << std::endl;
104 }
105#else
106
107 const int primaryNode = 0;
108 int nproc;
109 int myNode;
110 MPI_Comm_size(MPI_COMM_WORLD, &nproc);
111 MPI_Comm_rank(MPI_COMM_WORLD, &myNode);
112
113 std::vector<int> tmpNFixedZmols(nproc, 0);
114 std::vector<int> nFixedZmolsInProc(nproc, 0);
115 tmpNFixedZmols[myNode] = fixedZmols.size();
116
117 // do MPI_ALLREDUCE to exchange the total number of atoms,
118 // rigidbodies and cutoff groups
119 MPI_Allreduce(&tmpNFixedZmols[0], &nFixedZmolsInProc[0], nproc, MPI_INT,
120 MPI_SUM, MPI_COMM_WORLD);
121
122 MPI_Status* ierr = NULL;
123 int zmolIndex;
124 RealType data[3];
125
126 if (myNode == primaryNode) {
127 std::vector<ZconsData> zconsData;
128 ZconsData tmpData;
129 for (int i = 0; i < nproc; ++i) {
130 if (i == primaryNode) {
131 std::list<ZconstraintMol>::const_iterator j;
132 for (j = fixedZmols.begin(); j != fixedZmols.end(); ++j) {
133 tmpData.zmolIndex = j->mol->getGlobalIndex();
134 tmpData.zforce = j->fz;
135 tmpData.zpos = j->zpos;
136 tmpData.zconsPos = j->param.zTargetPos;
137 zconsData.push_back(tmpData);
138 }
139
140 } else {
141 for (int k = 0; k < nFixedZmolsInProc[i]; ++k) {
142 MPI_Recv(&zmolIndex, 1, MPI_INT, i, 0, MPI_COMM_WORLD, ierr);
143 MPI_Recv(data, 3, MPI_REALTYPE, i, 0, MPI_COMM_WORLD, ierr);
144 tmpData.zmolIndex = zmolIndex;
145 tmpData.zforce = data[0];
146 tmpData.zpos = data[1];
147 tmpData.zconsPos = data[2];
148 zconsData.push_back(tmpData);
149 }
150 }
151 }
152
153 output_ << info_->getSnapshotManager()->getCurrentSnapshot()->getTime()
154 << std::endl;
155 output_ << zconsData.size() << std::endl;
156
157 std::vector<ZconsData>::iterator l;
158 for (l = zconsData.begin(); l != zconsData.end(); ++l) {
159 output_ << l->zmolIndex << "\t" << l->zforce << "\t" << l->zpos << "\t"
160 << l->zconsPos << std::endl;
161 }
162
163 } else {
164 std::list<ZconstraintMol>::const_iterator j;
165 for (j = fixedZmols.begin(); j != fixedZmols.end(); ++j) {
166 zmolIndex = j->mol->getGlobalIndex();
167 data[0] = j->fz;
168 data[1] = j->zpos;
169 data[2] = j->param.zTargetPos;
170 MPI_Send(&zmolIndex, 1, MPI_INT, primaryNode, 0, MPI_COMM_WORLD);
171 MPI_Send(data, 3, MPI_REALTYPE, primaryNode, 0, MPI_COMM_WORLD);
172 }
173 }
174#endif
175 }
176} // namespace OpenMD
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.