ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/io/RestWriter.cpp
Revision: 1782
Committed: Wed Aug 22 02:28:28 2012 UTC (12 years, 8 months ago) by gezelter
File size: 8589 byte(s)
Log Message:
MERGE OpenMD development branch 1465:1781 into trunk

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     * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (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 cli2 1360
44 cli2 1407 #include <string>
45     #include <sstream>
46 cli2 1360 #include <iostream>
47    
48 chrisfen 993 #include "io/RestWriter.hpp"
49 chrisfen 417 #include "utils/simError.h"
50 cli2 1360 #include "brains/SnapshotManager.hpp"
51 chrisfen 993 #ifdef IS_MPI
52     #include <mpi.h>
53 cli2 1360 #endif
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 cli2 1407 createRestFile_ = true;
60    
61 chrisfen 990 #ifdef IS_MPI
62 cli2 1360 if(worldRank == 0){
63     #endif
64 cli2 1407
65     output_ = new std::ofstream(filename.c_str());
66    
67 cli2 1360 if(!output_){
68     sprintf( painCave.errMsg,
69     "Could not open %s for restraint output.\n",
70     filename.c_str());
71     painCave.isFatal = 1;
72     simError();
73     }
74    
75 cli2 1407 #ifdef IS_MPI
76     }
77     #endif // is_mpi
78 cli2 1364
79 cli2 1407
80     #ifdef IS_MPI
81     MPI_Status istatus;
82     #endif
83    
84     #ifndef IS_MPI
85    
86     (*output_) << "#time\t";
87    
88     std::vector<Restraint*>::const_iterator resti;
89    
90     for(resti=restraints.begin(); resti != restraints.end(); ++resti){
91     if ((*resti)->getPrintRestraint()) {
92    
93     std::string myName = (*resti)->getRestraintName();
94     int myType = (*resti)->getRestraintType();
95    
96     (*output_) << myName << ":";
97    
98     if (myType & Restraint::rtDisplacement)
99     (*output_) << "\tPosition(angstroms)\tEnergy(kcal/mol)";
100    
101     if (myType & Restraint::rtTwist)
102     (*output_) << "\tTwistAngle(radians)\tEnergy(kcal/mol)";
103    
104     if (myType & Restraint::rtSwingX)
105     (*output_) << "\tSwingXAngle(radians)\tEnergy(kcal/mol)";
106 cli2 1364
107 cli2 1407 if (myType & Restraint::rtSwingY)
108     (*output_) << "\tSwingYAngle(radians)\tEnergy(kcal/mol)";
109 cli2 1360 }
110 gezelter 507 }
111 cli2 1407
112     (*output_) << "\n";
113     (*output_).flush();
114    
115     #else
116    
117     std::string buffer;
118    
119     std::vector<Restraint*>::const_iterator resti;
120    
121     for(resti=restraints.begin(); resti != restraints.end(); ++resti){
122     if ((*resti)->getPrintRestraint()) {
123    
124     std::string myName = (*resti)->getRestraintName();
125     int myType = (*resti)->getRestraintType();
126    
127     buffer += (myName + ":");
128    
129     if (myType & Restraint::rtDisplacement)
130     buffer += "\tPosition(angstroms)\tEnergy(kcal/mol)";
131    
132     if (myType & Restraint::rtTwist)
133     buffer += "\tTwistAngle(radians)\tEnergy(kcal/mol)";
134    
135     if (myType & Restraint::rtSwingX)
136     buffer += "\tSwingXAngle(radians)\tEnergy(kcal/mol)";
137    
138     if (myType & Restraint::rtSwingY)
139     buffer += "\tSwingYAngle(radians)\tEnergy(kcal/mol)";
140    
141     buffer += "\n";
142     }
143     }
144    
145     const int masterNode = 0;
146    
147     if (worldRank == masterNode) {
148     (*output_) << "#time\t";
149     (*output_) << buffer;
150    
151     int nProc;
152     MPI_Comm_size(MPI_COMM_WORLD, &nProc);
153     for (int i = 1; i < nProc; ++i) {
154    
155     // receive the length of the string buffer that was
156     // prepared by processor i
157    
158     int recvLength;
159     MPI_Recv(&recvLength, 1, MPI_INT, i, 0, MPI_COMM_WORLD, &istatus);
160     char* recvBuffer = new char[recvLength];
161     if (recvBuffer == NULL) {
162     } else {
163     MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD,
164     &istatus);
165     (*output_) << recvBuffer;
166     delete [] recvBuffer;
167     }
168     }
169     (*output_).flush();
170     } else {
171     int sendBufferLength = buffer.size() + 1;
172     MPI_Send(&sendBufferLength, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
173     MPI_Send((void *)buffer.c_str(), sendBufferLength, MPI_CHAR, masterNode,
174     0, MPI_COMM_WORLD);
175     }
176    
177     #endif // is_mpi
178    
179     }
180 cli2 1360
181 cli2 1407 void RestWriter::writeRest(std::vector<std::map<int, Restraint::RealPair> > restInfo) {
182    
183 chrisfen 993 #ifdef IS_MPI
184 cli2 1407 MPI_Status istatus;
185 cli2 1360 #endif
186 chrisfen 417
187 cli2 1407 #ifndef IS_MPI
188     (*output_) << info_->getSnapshotManager()->getCurrentSnapshot()->getTime();
189    
190 cli2 1360 // output some information about the molecules
191     std::vector<std::map<int, Restraint::RealPair> >::const_iterator i;
192     std::map<int, Restraint::RealPair>::const_iterator j;
193 cli2 1407
194 cli2 1360 for( i = restInfo.begin(); i != restInfo.end(); ++i){
195     for(j = (*i).begin(); j != (*i).end(); ++j){
196 cli2 1407 (*output_) << "\t" << (j->second).first << "\t" << (j->second).second;
197 chrisfen 417 }
198 cli2 1407 (*output_) << std::endl;
199 chrisfen 417 }
200 cli2 1407 (*output_).flush();
201     #else
202     std::string buffer, first, second;
203     std::stringstream ss;
204    
205     std::vector<std::map<int, Restraint::RealPair> >::const_iterator i;
206     std::map<int, Restraint::RealPair>::const_iterator j;
207    
208     for( i = restInfo.begin(); i != restInfo.end(); ++i){
209     for(j = (*i).begin(); j != (*i).end(); ++j){
210     ss.clear();
211     ss << (j->second).first;
212     ss >> first;
213     ss.clear();
214     ss << (j->second).second;
215     ss >> second;
216     buffer += ("\t" + first + "\t" + second);
217     }
218     buffer += "\n";
219     }
220    
221     const int masterNode = 0;
222    
223     if (worldRank == masterNode) {
224     (*output_) << info_->getSnapshotManager()->getCurrentSnapshot()->getTime();
225     (*output_) << buffer;
226    
227     int nProc;
228     MPI_Comm_size(MPI_COMM_WORLD, &nProc);
229     for (int i = 1; i < nProc; ++i) {
230    
231     // receive the length of the string buffer that was
232     // prepared by processor i
233    
234     int recvLength;
235     MPI_Recv(&recvLength, 1, MPI_INT, i, 0, MPI_COMM_WORLD, &istatus);
236     char* recvBuffer = new char[recvLength];
237     if (recvBuffer == NULL) {
238     } else {
239     MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD,
240     &istatus);
241     (*output_) << recvBuffer;
242    
243     delete [] recvBuffer;
244     }
245     }
246     (*output_).flush();
247     } else {
248     int sendBufferLength = buffer.size() + 1;
249     MPI_Send(&sendBufferLength, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
250     MPI_Send((void *)buffer.c_str(), sendBufferLength, MPI_CHAR, masterNode,
251     0, MPI_COMM_WORLD);
252     }
253     #endif // is_mpi
254 chrisfen 417 }
255    
256 cli2 1407
257     RestWriter::~RestWriter() {
258    
259     #ifdef IS_MPI
260    
261     if (worldRank == 0) {
262     #endif // is_mpi
263     if (createRestFile_){
264     writeClosing(*output_);
265     delete output_;
266     }
267     #ifdef IS_MPI
268     }
269     #endif // is_mpi
270     }
271    
272     void RestWriter::writeClosing(std::ostream& os) {
273     os.flush();
274     }
275    
276 gezelter 1390 }// end namespace OpenMD
277 cli2 1407

Properties

Name Value
svn:keywords Author Id Revision Date