ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/io/RestWriter.cpp
(Generate patch)

Comparing:
trunk/src/io/RestWriter.cpp (file contents), Revision 417 by chrisfen, Thu Mar 10 15:10:24 2005 UTC vs.
branches/development/src/io/RestWriter.cpp (file contents), Revision 1665 by gezelter, Tue Nov 22 20:38:56 2011 UTC

# Line 1 | Line 1
1 < /*
2 < * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
1 > /*
2 > * Copyright (c) 2009 The University of Notre Dame. All Rights Reserved.
3   *
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 < * 1. Acknowledgement of the program authors must be made in any
10 < *    publication of scientific results based in part on use of the
11 < *    program.  An acceptable form of acknowledgement is citation of
12 < *    the article in which the program was described (Matthew
13 < *    A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 < *    J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 < *    Parallel Simulation Engine for Molecular Dynamics,"
16 < *    J. Comput. Chem. 26, pp. 252-271 (2005))
17 < *
18 < * 2. Redistributions of source code must retain the above copyright
9 > * 1. Redistributions of source code must retain the above copyright
10   *    notice, this list of conditions and the following disclaimer.
11   *
12 < * 3. Redistributions in binary form must reproduce the above copyright
12 > * 2. Redistributions in binary form must reproduce the above copyright
13   *    notice, this list of conditions and the following disclaimer in the
14   *    documentation and/or other materials provided with the
15   *    distribution.
# Line 37 | Line 28
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 + *
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 + * [4]  Kuang & Gezelter,  J. Chem. Phys. 133, 164101 (2010).
40 + * [5]  Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41   */
42 <
43 < #include <algorithm>
42 >
43 >
44 > #include <string>
45 > #include <sstream>
46   #include <iostream>
44 #include <map>
47  
46 #include "primitives/Molecule.hpp"
48   #include "io/RestWriter.hpp"
49   #include "utils/simError.h"
50 + #include "brains/SnapshotManager.hpp"
51 + #ifdef IS_MPI
52 + #include <mpi.h>
53 + #endif
54  
55 <
56 < namespace oopse {
57 <  RestWriter::RestWriter(SimInfo* info) :
58 <  info_(info) {
55 > namespace OpenMD {
56 >  RestWriter::RestWriter(SimInfo* info, const std::string& filename,
57 >                         std::vector<Restraint*> restraints ) :
58 >    info_(info){
59 >    createRestFile_ = true;
60      
55    //we use master - slave mode, only master node writes to disk
56    outName = info_->getRestFileName();
57  }
58  
59  RestWriter::~RestWriter() {}
60  
61  void RestWriter::writeZangle(){
62    const int BUFFERSIZE = 2000;
63    char tempBuffer[BUFFERSIZE];
64    char writeLine[BUFFERSIZE];
65    
66    std::ofstream finalOut;
67    
68    Molecule* mol;
69    StuntDouble* integrableObject;
70    SimInfo::MoleculeIterator mi;
71    Molecule::IntegrableObjectIterator ii;
72    
61   #ifdef IS_MPI
62 <    if(worldRank == 0 ){
63 < #endif    
64 <      finalOut.open( outName.c_str(), std::ios::out | std::ios::trunc );
65 <      if( !finalOut ){
62 >    if(worldRank == 0){
63 > #endif
64 >  
65 >      output_ = new std::ofstream(filename.c_str());
66 >
67 >      if(!output_){
68          sprintf( painCave.errMsg,
69 <                 "Could not open \"%s\" for zAngle output.\n",
70 <                 outName );
69 >                 "Could not open %s for restraint output.\n",
70 >                 filename.c_str());
71          painCave.isFatal = 1;
72          simError();
73        }
74 +
75   #ifdef IS_MPI
76      }
77   #endif // is_mpi
78 +
79 +
80 + #ifdef IS_MPI
81 +    MPI_Status istatus;
82 + #endif
83      
84   #ifndef IS_MPI
85 <    // first we do output for the single processor version
86 <    finalOut
87 <      << info_->getSnapshotManager()->getCurrentSnapshot()->getTime()
88 <      << " : omega values at this time\n";
89 <    
90 <    for (mol = info_->beginMolecule(mi); mol != NULL;
91 <         mol = info_->nextMolecule(mi)) {
96 <      
97 <      for (integrableObject = mol->beginIntegrableObject(ii);
98 <           integrableObject != NULL;
99 <           integrableObject = mol->nextIntegrableObject(ii)) {    
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 <        sprintf( tempBuffer,
94 <                 "%14.10lf\n",
103 <                 integrableObject->getZangle());
104 <        strcpy( writeLine, tempBuffer );
93 >        std::string myName = (*resti)->getRestraintName();
94 >        int myType = (*resti)->getRestraintType();
95          
96 <        finalOut << writeLine;      
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 >          
107 >        if (myType & Restraint::rtSwingY)
108 >          (*output_) << "\tSwingYAngle(radians)\tEnergy(kcal/mol)";
109        }
108      
110      }
111 +
112 +    (*output_) << "\n";
113 +    (*output_).flush();
114      
115   #else
112    int nproc;
113    MPI_Comm_size(MPI_COMM_WORLD, &nproc);
114    const int masterNode = 0;
115    int myNode = worldRank;
116    std::vector<int> tmpNIntObjects(nproc, 0);
117    std::vector<int> nIntObjectsInProc(nproc, 0);
118    tmpNIntObjects[myNode] = info_->getNGlobalIntegrableObjects();
116      
117 <    //do MPI_ALLREDUCE to exchange the total number of atoms, rigidbodies and cutoff groups
118 <    MPI_Allreduce(&tmpNIntObjects[0], &nIntObjectsInProc[0], nproc, MPI_INT,
119 <                  MPI_SUM, MPI_COMM_WORLD);
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 <    MPI_Status ierr;
125 <    int intObIndex;
126 <    double zAngle;
145 >    const int masterNode = 0;
146      
147 <    if (masterNode == 0) {
148 <      std::map<int, double> zAngData;
149 <      for(int i = 0 ; i < nproc; ++i) {
150 <        if (i == masterNode) {
151 <          for (mol = info_->beginMolecule(mi); mol != NULL;
152 <               mol = info_->nextMolecule(mi)) {
153 <            
154 <            for (integrableObject = mol->beginIntegrableObject(ii);
155 <                 integrableObject != NULL;
156 <                 integrableObject = mol->nextIntegrableObject(ii)) {
157 <              
158 <              intObIndex = integrableObject->getGlobalIndex() ;
159 <              zAngle = integrableObject->getZangle();
160 <              zAngData.insert(pair<int, double>(intObIndex, zAngle));
161 <            }      
143 <          }
144 <          
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 <          for(int k = 0; k < nIntObjectsInProc[i]; ++k) {
164 <            MPI_Recv(&intObIndex, 1, MPI_INT, i, 0, MPI_COMM_WORLD,&ierr);
165 <            MPI_Recv(&zAngle, 1, MPI_DOUBLE, i, 0, MPI_COMM_WORLD,&ierr);
166 <            zAngData.insert(pair<int, double>(intObIndex, zAngle));
150 <          }
163 >          MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD,
164 >                   &istatus);
165 >          (*output_) << recvBuffer;
166 >          delete [] recvBuffer;
167          }
168 <        
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 >  
181 >  void RestWriter::writeRest(std::vector<std::map<int, Restraint::RealPair> > restInfo) {
182 >    
183 > #ifdef IS_MPI
184 >    MPI_Status istatus;
185 > #endif
186 >    
187 > #ifndef IS_MPI
188 >    (*output_) << info_->getSnapshotManager()->getCurrentSnapshot()->getTime();
189 >    
190 >    // 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 >    
194 >    for( i = restInfo.begin(); i != restInfo.end(); ++i){
195 >      for(j = (*i).begin(); j != (*i).end(); ++j){                
196 >        (*output_) << "\t" << (j->second).first << "\t" << (j->second).second;
197        }
198 <      
199 <      finalOut
200 <        << info_->getSnapshotManager()->getCurrentSnapshot()->getTime()
201 <        << " : omega values at this time\n";
202 <      
203 <      std::map<int, double>::iterator l;
204 <      for (l = zAngData.begin(); l != zAngData.end(); ++l) {
205 <        finalOut << l->second << "\n";
198 >      (*output_) << std::endl;
199 >    }
200 >    (*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 <    } else {
228 <      
229 <      for (mol = info_->beginMolecule(mi); mol != NULL;
167 <           mol = info_->nextMolecule(mi)) {
227 >      int nProc;
228 >      MPI_Comm_size(MPI_COMM_WORLD, &nProc);
229 >      for (int i = 1; i < nProc; ++i) {
230          
231 <        for (integrableObject = mol->beginIntegrableObject(ii);
232 <             integrableObject != NULL;
233 <             integrableObject = mol->nextIntegrableObject(ii)) {
234 <          intObIndex = integrableObject->getGlobalIndex();            
235 <          zAngle = integrableObject->getZangle();
236 <          MPI_Send(&intObIndex, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
237 <          MPI_Send(&zAngle, 1, MPI_DOUBLE, masterNode, 0, MPI_COMM_WORLD);
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 <      }
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
253 > #endif // is_mpi
254 >  }
255 >  
256 >  
257 >  RestWriter::~RestWriter() {
258      
259   #ifdef IS_MPI
182    finalOut.close();
183 #endif
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 < }
272 >  void RestWriter::writeClosing(std::ostream& os) {
273 >    os.flush();
274 >  }
275 >  
276 > }// end namespace OpenMD
277 >

Comparing:
trunk/src/io/RestWriter.cpp (property svn:keywords), Revision 417 by chrisfen, Thu Mar 10 15:10:24 2005 UTC vs.
branches/development/src/io/RestWriter.cpp (property svn:keywords), Revision 1665 by gezelter, Tue Nov 22 20:38:56 2011 UTC

# Line 0 | Line 1
1 + Author Id Revision Date

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines