ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/io/RestWriter.cpp
Revision: 2020
Committed: Mon Sep 22 19:18:35 2014 UTC (10 years, 7 months ago) by gezelter
File size: 9417 byte(s)
Log Message:
Fixes for restraints, renaming of UniformField

File Contents

# User Rev Content
1 gezelter 507 /*
2 cli2 1360 * Copyright (c) 2009 The University of Notre Dame. All Rights Reserved.
3 chrisfen 417 *
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 gezelter 1390 * 1. Redistributions of source code must retain the above copyright
10 chrisfen 417 * notice, this list of conditions and the following disclaimer.
11     *
12 gezelter 1390 * 2. Redistributions in binary form must reproduce the above copyright
13 chrisfen 417 * notice, this list of conditions and the following disclaimer in the
14     * documentation and/or other materials provided with the
15     * distribution.
16     *
17     * This software is provided "AS IS," without a warranty of any
18     * kind. All express or implied conditions, representations and
19     * warranties, including any implied warranty of merchantability,
20     * fitness for a particular purpose or non-infringement, are hereby
21     * excluded. The University of Notre Dame and its licensors shall not
22     * be liable for any damages suffered by licensee as a result of
23     * using, modifying or distributing the software or its
24     * derivatives. In no event will the University of Notre Dame or its
25     * licensors be liable for any lost revenue, profit or data, or for
26     * direct, indirect, special, consequential, incidental or punitive
27     * damages, however caused and regardless of the theory of liability,
28     * arising out of the use of or inability to use software, even if the
29     * University of Notre Dame has been advised of the possibility of
30     * such damages.
31 gezelter 1390 *
32     * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33     * research, please cite the appropriate papers when you publish your
34     * work. Good starting points are:
35     *
36     * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37     * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38 gezelter 1879 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
39 gezelter 1782 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40     * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 chrisfen 417 */
42    
43 gezelter 1938 #ifdef IS_MPI
44     #include <mpi.h>
45     #endif
46 cli2 1360
47 cli2 1407 #include <string>
48     #include <sstream>
49 cli2 1360 #include <iostream>
50    
51 chrisfen 993 #include "io/RestWriter.hpp"
52 chrisfen 417 #include "utils/simError.h"
53 cli2 1360 #include "brains/SnapshotManager.hpp"
54 chrisfen 417
55 gezelter 1390 namespace OpenMD {
56 cli2 1360 RestWriter::RestWriter(SimInfo* info, const std::string& filename,
57     std::vector<Restraint*> restraints ) :
58     info_(info){
59 gezelter 2020
60     std::vector<Restraint*>::const_iterator resti;
61    
62     createRestFile_ = false;
63    
64     #ifdef IS_MPI
65     MPI_Status* istatus;
66     #endif
67 cli2 1407
68 gezelter 2020 int printAny = 0;
69     for(resti=restraints.begin(); resti != restraints.end(); ++resti){
70     if ((*resti)->getPrintRestraint()) {
71     printAny = 1;
72     }
73     }
74    
75 chrisfen 990 #ifdef IS_MPI
76 gezelter 2020 MPI_Allreduce(MPI_IN_PLACE, &printAny, 1, MPI_INT, MPI_SUM,
77     MPI_COMM_WORLD);
78     #endif
79    
80     if (printAny) createRestFile_ = true;
81    
82     #ifdef IS_MPI
83 cli2 1360 if(worldRank == 0){
84     #endif
85 cli2 1407
86 gezelter 2020 if (createRestFile_) {
87     output_ = new std::ofstream(filename.c_str());
88    
89     if(!output_){
90     sprintf( painCave.errMsg,
91     "Could not open %s for restraint output.\n",
92     filename.c_str());
93     painCave.isFatal = 1;
94     simError();
95     }
96 cli2 1360 }
97 gezelter 2020
98 cli2 1407 #ifdef IS_MPI
99     }
100     #endif // is_mpi
101 cli2 1364
102 cli2 1407
103     #ifndef IS_MPI
104    
105 gezelter 2020 if (createRestFile_) (*output_) << "#time\t";
106 cli2 1407
107     for(resti=restraints.begin(); resti != restraints.end(); ++resti){
108     if ((*resti)->getPrintRestraint()) {
109     std::string myName = (*resti)->getRestraintName();
110     int myType = (*resti)->getRestraintType();
111    
112     (*output_) << myName << ":";
113    
114     if (myType & Restraint::rtDisplacement)
115     (*output_) << "\tPosition(angstroms)\tEnergy(kcal/mol)";
116    
117     if (myType & Restraint::rtTwist)
118     (*output_) << "\tTwistAngle(radians)\tEnergy(kcal/mol)";
119    
120     if (myType & Restraint::rtSwingX)
121     (*output_) << "\tSwingXAngle(radians)\tEnergy(kcal/mol)";
122 cli2 1364
123 cli2 1407 if (myType & Restraint::rtSwingY)
124     (*output_) << "\tSwingYAngle(radians)\tEnergy(kcal/mol)";
125 cli2 1360 }
126 gezelter 507 }
127 cli2 1407
128 gezelter 2020 if (createRestFile_) (*output_) << "\n";
129     if (createRestFile_) (*output_).flush();
130 cli2 1407
131     #else
132    
133     std::string buffer;
134    
135     for(resti=restraints.begin(); resti != restraints.end(); ++resti){
136     if ((*resti)->getPrintRestraint()) {
137     std::string myName = (*resti)->getRestraintName();
138     int myType = (*resti)->getRestraintType();
139    
140     buffer += (myName + ":");
141    
142     if (myType & Restraint::rtDisplacement)
143     buffer += "\tPosition(angstroms)\tEnergy(kcal/mol)";
144    
145     if (myType & Restraint::rtTwist)
146     buffer += "\tTwistAngle(radians)\tEnergy(kcal/mol)";
147    
148     if (myType & Restraint::rtSwingX)
149     buffer += "\tSwingXAngle(radians)\tEnergy(kcal/mol)";
150    
151     if (myType & Restraint::rtSwingY)
152     buffer += "\tSwingYAngle(radians)\tEnergy(kcal/mol)";
153    
154     buffer += "\n";
155     }
156     }
157    
158     const int masterNode = 0;
159    
160     if (worldRank == masterNode) {
161 gezelter 2020 if (createRestFile_) (*output_) << "#time\t";
162     if (createRestFile_) (*output_) << buffer;
163 cli2 1407
164 gezelter 1969 int nProc;
165     MPI_Comm_size( MPI_COMM_WORLD, &nProc);
166    
167 cli2 1407 for (int i = 1; i < nProc; ++i) {
168    
169     // receive the length of the string buffer that was
170     // prepared by processor i
171    
172     int recvLength;
173 gezelter 1969 MPI_Recv(&recvLength, 1, MPI_INT, i, 0, MPI_COMM_WORLD, istatus);
174 cli2 1407 char* recvBuffer = new char[recvLength];
175     if (recvBuffer == NULL) {
176     } else {
177 gezelter 1969 MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD,
178     istatus);
179 gezelter 2020 if (createRestFile_) (*output_) << recvBuffer;
180 cli2 1407 delete [] recvBuffer;
181     }
182     }
183 gezelter 2020 if (createRestFile_) (*output_).flush();
184 cli2 1407 } else {
185     int sendBufferLength = buffer.size() + 1;
186 gezelter 1969 MPI_Send(&sendBufferLength, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
187     MPI_Send((void *)buffer.c_str(), sendBufferLength, MPI_CHAR,
188     masterNode, 0, MPI_COMM_WORLD);
189 cli2 1407 }
190    
191     #endif // is_mpi
192    
193     }
194 cli2 1360
195 cli2 1407 void RestWriter::writeRest(std::vector<std::map<int, Restraint::RealPair> > restInfo) {
196    
197 chrisfen 993 #ifdef IS_MPI
198 gezelter 1969 MPI_Status* istatus;
199 cli2 1360 #endif
200 chrisfen 417
201 cli2 1407 #ifndef IS_MPI
202 gezelter 2020 if (createRestFile_) (*output_) << info_->getSnapshotManager()->getCurrentSnapshot()->getTime();
203 cli2 1407
204 cli2 1360 // output some information about the molecules
205     std::vector<std::map<int, Restraint::RealPair> >::const_iterator i;
206     std::map<int, Restraint::RealPair>::const_iterator j;
207 cli2 1407
208 gezelter 2020 cerr << "risize = " << restInfo.size() << "\n";
209    
210     if ( createRestFile_ ) {
211    
212     for( i = restInfo.begin(); i != restInfo.end(); ++i){
213     for(j = (*i).begin(); j != (*i).end(); ++j){
214     (*output_) << "\t" << (j->second).first << "\t" << (j->second).second;
215     }
216     (*output_) << std::endl;
217     }
218     (*output_).flush();
219 chrisfen 417 }
220 cli2 1407 #else
221     std::string buffer, first, second;
222     std::stringstream ss;
223    
224     std::vector<std::map<int, Restraint::RealPair> >::const_iterator i;
225     std::map<int, Restraint::RealPair>::const_iterator j;
226    
227 gezelter 2020 if ( createRestFile_ ) {
228     for( i = restInfo.begin(); i != restInfo.end(); ++i){
229    
230     for(j = (*i).begin(); j != (*i).end(); ++j){
231     ss.clear();
232     ss << (j->second).first;
233     ss >> first;
234     ss.clear();
235     ss << (j->second).second;
236     ss >> second;
237     buffer += ("\t" + first + "\t" + second);
238     }
239     buffer += "\n";
240 cli2 1407 }
241     }
242    
243     const int masterNode = 0;
244    
245 gezelter 2020 if (createRestFile_) {
246     if (worldRank == masterNode) {
247    
248     (*output_) << info_->getSnapshotManager()->getCurrentSnapshot()->getTime();
249     (*output_) << buffer;
250 cli2 1407
251 gezelter 2020 int nProc;
252     MPI_Comm_size( MPI_COMM_WORLD, &nProc);
253     for (int i = 1; i < nProc; ++i) {
254 cli2 1407
255 gezelter 2020 // receive the length of the string buffer that was
256     // prepared by processor i
257    
258     int recvLength;
259     MPI_Recv(&recvLength, 1, MPI_INT, i, 0, MPI_COMM_WORLD, istatus);
260     char* recvBuffer = new char[recvLength];
261     if (recvBuffer == NULL) {
262     } else {
263     MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD,
264     istatus);
265     if (createRestFile_) (*output_) << recvBuffer;
266    
267     delete [] recvBuffer;
268     }
269     }
270     (*output_).flush();
271     } else {
272     int sendBufferLength = buffer.size() + 1;
273     MPI_Send(&sendBufferLength, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
274     MPI_Send((void *)buffer.c_str(), sendBufferLength,
275     MPI_CHAR, masterNode, 0, MPI_COMM_WORLD);
276     }
277 cli2 1407 }
278     #endif // is_mpi
279 chrisfen 417 }
280    
281 cli2 1407
282     RestWriter::~RestWriter() {
283    
284     #ifdef IS_MPI
285    
286     if (worldRank == 0) {
287     #endif // is_mpi
288     if (createRestFile_){
289     writeClosing(*output_);
290     delete output_;
291     }
292     #ifdef IS_MPI
293     }
294     #endif // is_mpi
295     }
296    
297     void RestWriter::writeClosing(std::ostream& os) {
298     os.flush();
299     }
300    
301 gezelter 1390 }// end namespace OpenMD
302 cli2 1407

Properties

Name Value
svn:keywords Author Id Revision Date