ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/io/RestWriter.cpp
Revision: 1969
Committed: Wed Feb 26 14:14:50 2014 UTC (11 years, 2 months ago) by gezelter
File size: 8590 byte(s)
Log Message:
Fixes to deal with deprecation of MPI C++ bindings.  We've reverted back to the
C calls.

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 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 gezelter 1969 MPI_Status* istatus;
82 cli2 1407 #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 gezelter 1969 int nProc;
152     MPI_Comm_size( MPI_COMM_WORLD, &nProc);
153    
154 cli2 1407 for (int i = 1; i < nProc; ++i) {
155    
156     // receive the length of the string buffer that was
157     // prepared by processor i
158    
159     int recvLength;
160 gezelter 1969 MPI_Recv(&recvLength, 1, MPI_INT, i, 0, MPI_COMM_WORLD, istatus);
161 cli2 1407 char* recvBuffer = new char[recvLength];
162     if (recvBuffer == NULL) {
163     } else {
164 gezelter 1969 MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD,
165     istatus);
166 cli2 1407 (*output_) << recvBuffer;
167     delete [] recvBuffer;
168     }
169     }
170     (*output_).flush();
171     } else {
172     int sendBufferLength = buffer.size() + 1;
173 gezelter 1969 MPI_Send(&sendBufferLength, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
174     MPI_Send((void *)buffer.c_str(), sendBufferLength, MPI_CHAR,
175     masterNode, 0, MPI_COMM_WORLD);
176 cli2 1407 }
177    
178     #endif // is_mpi
179    
180     }
181 cli2 1360
182 cli2 1407 void RestWriter::writeRest(std::vector<std::map<int, Restraint::RealPair> > restInfo) {
183    
184 chrisfen 993 #ifdef IS_MPI
185 gezelter 1969 MPI_Status* istatus;
186 cli2 1360 #endif
187 chrisfen 417
188 cli2 1407 #ifndef IS_MPI
189     (*output_) << info_->getSnapshotManager()->getCurrentSnapshot()->getTime();
190    
191 cli2 1360 // output some information about the molecules
192     std::vector<std::map<int, Restraint::RealPair> >::const_iterator i;
193     std::map<int, Restraint::RealPair>::const_iterator j;
194 cli2 1407
195 cli2 1360 for( i = restInfo.begin(); i != restInfo.end(); ++i){
196     for(j = (*i).begin(); j != (*i).end(); ++j){
197 cli2 1407 (*output_) << "\t" << (j->second).first << "\t" << (j->second).second;
198 chrisfen 417 }
199 cli2 1407 (*output_) << std::endl;
200 chrisfen 417 }
201 cli2 1407 (*output_).flush();
202     #else
203     std::string buffer, first, second;
204     std::stringstream ss;
205    
206     std::vector<std::map<int, Restraint::RealPair> >::const_iterator i;
207     std::map<int, Restraint::RealPair>::const_iterator j;
208    
209     for( i = restInfo.begin(); i != restInfo.end(); ++i){
210     for(j = (*i).begin(); j != (*i).end(); ++j){
211     ss.clear();
212     ss << (j->second).first;
213     ss >> first;
214     ss.clear();
215     ss << (j->second).second;
216     ss >> second;
217     buffer += ("\t" + first + "\t" + second);
218     }
219     buffer += "\n";
220     }
221    
222     const int masterNode = 0;
223    
224     if (worldRank == masterNode) {
225     (*output_) << info_->getSnapshotManager()->getCurrentSnapshot()->getTime();
226     (*output_) << buffer;
227    
228     int nProc;
229 gezelter 1969 MPI_Comm_size( MPI_COMM_WORLD, &nProc);
230 cli2 1407 for (int i = 1; i < nProc; ++i) {
231    
232     // receive the length of the string buffer that was
233     // prepared by processor i
234    
235     int recvLength;
236 gezelter 1969 MPI_Recv(&recvLength, 1, MPI_INT, i, 0, MPI_COMM_WORLD, istatus);
237 cli2 1407 char* recvBuffer = new char[recvLength];
238     if (recvBuffer == NULL) {
239     } else {
240 gezelter 1969 MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD,
241     istatus);
242 cli2 1407 (*output_) << recvBuffer;
243    
244     delete [] recvBuffer;
245     }
246     }
247     (*output_).flush();
248     } else {
249     int sendBufferLength = buffer.size() + 1;
250 gezelter 1969 MPI_Send(&sendBufferLength, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
251     MPI_Send((void *)buffer.c_str(), sendBufferLength,
252     MPI_CHAR, masterNode, 0, MPI_COMM_WORLD);
253 cli2 1407 }
254     #endif // is_mpi
255 chrisfen 417 }
256    
257 cli2 1407
258     RestWriter::~RestWriter() {
259    
260     #ifdef IS_MPI
261    
262     if (worldRank == 0) {
263     #endif // is_mpi
264     if (createRestFile_){
265     writeClosing(*output_);
266     delete output_;
267     }
268     #ifdef IS_MPI
269     }
270     #endif // is_mpi
271     }
272    
273     void RestWriter::writeClosing(std::ostream& os) {
274     os.flush();
275     }
276    
277 gezelter 1390 }// end namespace OpenMD
278 cli2 1407

Properties

Name Value
svn:keywords Author Id Revision Date