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

Comparing trunk/src/io/RestWriter.cpp (file contents):
Revision 1969 by gezelter, Wed Feb 26 14:14:50 2014 UTC vs.
Revision 2021 by gezelter, Tue Sep 23 15:25:08 2014 UTC

# Line 56 | Line 56 | namespace OpenMD {
56    RestWriter::RestWriter(SimInfo* info, const std::string& filename,
57                           std::vector<Restraint*> restraints ) :
58      info_(info){
59 <    createRestFile_ = true;
59 >
60 >    std::vector<Restraint*>::const_iterator resti;
61 >
62 >    createRestFile_ = false;  
63 >
64 > #ifdef IS_MPI    
65 >    MPI_Status* istatus;
66 > #endif
67      
68 +    int printAny = 0;
69 +    for(resti=restraints.begin(); resti != restraints.end(); ++resti){
70 +      if ((*resti)->getPrintRestraint()) {
71 +        printAny = 1;
72 +      }
73 +    }
74 +    
75   #ifdef IS_MPI
76 +    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      if(worldRank == 0){
84   #endif
85    
86 <      output_ = new std::ofstream(filename.c_str());
87 <
88 <      if(!output_){
89 <        sprintf( painCave.errMsg,
90 <                 "Could not open %s for restraint output.\n",
91 <                 filename.c_str());
92 <        painCave.isFatal = 1;
93 <        simError();
86 >      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        }
97 <
97 >        
98   #ifdef IS_MPI
99      }
100   #endif // is_mpi
101  
102  
80 #ifdef IS_MPI
81    MPI_Status* istatus;
82 #endif
83    
103   #ifndef IS_MPI
104          
105 <    (*output_) << "#time\t";
105 >    if (createRestFile_) (*output_) << "#time\t";
106  
88    std::vector<Restraint*>::const_iterator resti;
89
107      for(resti=restraints.begin(); resti != restraints.end(); ++resti){
108        if ((*resti)->getPrintRestraint()) {
92        
109          std::string myName = (*resti)->getRestraintName();
110          int myType = (*resti)->getRestraintType();
111          
# Line 109 | Line 125 | namespace OpenMD {
125        }
126      }
127  
128 <    (*output_) << "\n";
129 <    (*output_).flush();
128 >    if (createRestFile_) (*output_) << "\n";
129 >    if (createRestFile_) (*output_).flush();
130      
131   #else
132      
133      std::string buffer;
134  
119    std::vector<Restraint*>::const_iterator resti;
120
135      for(resti=restraints.begin(); resti != restraints.end(); ++resti){
136        if ((*resti)->getPrintRestraint()) {
123        
137          std::string myName = (*resti)->getRestraintName();
138          int myType = (*resti)->getRestraintType();
139  
# Line 145 | Line 158 | namespace OpenMD {
158      const int masterNode = 0;
159      
160      if (worldRank == masterNode) {
161 <      (*output_) << "#time\t";
162 <      (*output_) << buffer;
161 >      if (createRestFile_) (*output_) << "#time\t";
162 >      if (createRestFile_) (*output_) << buffer;
163        
164        int nProc;
165        MPI_Comm_size( MPI_COMM_WORLD, &nProc);
# Line 163 | Line 176 | namespace OpenMD {
176          } else {
177            MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD,
178                     istatus);
179 <          (*output_) << recvBuffer;
179 >          if (createRestFile_) (*output_) << recvBuffer;
180            delete [] recvBuffer;
181          }
182        }
183 <      (*output_).flush();
183 >       if (createRestFile_) (*output_).flush();
184      } else {
185        int sendBufferLength = buffer.size() + 1;
186        MPI_Send(&sendBufferLength, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
# Line 186 | Line 199 | namespace OpenMD {
199   #endif
200      
201   #ifndef IS_MPI
202 <    (*output_) << info_->getSnapshotManager()->getCurrentSnapshot()->getTime();
202 >     if (createRestFile_)  (*output_) << info_->getSnapshotManager()->getCurrentSnapshot()->getTime();
203      
204      // 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      
208 <    for( i = restInfo.begin(); i != restInfo.end(); ++i){
209 <      for(j = (*i).begin(); j != (*i).end(); ++j){                
210 <        (*output_) << "\t" << (j->second).first << "\t" << (j->second).second;
211 <      }
212 <      (*output_) << std::endl;
208 >    if ( createRestFile_ ) {
209 >      
210 >      for( i = restInfo.begin(); i != restInfo.end(); ++i){        
211 >        for(j = (*i).begin(); j != (*i).end(); ++j){                
212 >          (*output_) << "\t" << (j->second).first << "\t" << (j->second).second;
213 >        }
214 >        (*output_) << std::endl;
215 >      }      
216 >      (*output_).flush();
217      }
201    (*output_).flush();
218   #else
219      std::string buffer, first, second;
220      std::stringstream ss;
# Line 206 | Line 222 | namespace OpenMD {
222      std::vector<std::map<int, Restraint::RealPair> >::const_iterator i;
223      std::map<int, Restraint::RealPair>::const_iterator j;
224      
225 <    for( i = restInfo.begin(); i != restInfo.end(); ++i){
226 <      for(j = (*i).begin(); j != (*i).end(); ++j){
227 <        ss.clear();
228 <        ss << (j->second).first;
229 <        ss >> first;
230 <        ss.clear();
231 <        ss << (j->second).second;
232 <        ss >> second;
233 <        buffer += ("\t" + first + "\t" + second);      
225 >    if ( createRestFile_ ) {
226 >      for( i = restInfo.begin(); i != restInfo.end(); ++i){
227 >        
228 >        for(j = (*i).begin(); j != (*i).end(); ++j){
229 >          ss.clear();
230 >          ss << (j->second).first;
231 >          ss >> first;
232 >          ss.clear();
233 >          ss << (j->second).second;
234 >          ss >> second;
235 >          buffer += ("\t" + first + "\t" + second);      
236 >        }
237 >        buffer += "\n";    
238        }
219      buffer += "\n";    
239      }
240      
241      const int masterNode = 0;
242      
243 <    if (worldRank == masterNode) {
244 <      (*output_) << info_->getSnapshotManager()->getCurrentSnapshot()->getTime();
226 <      (*output_) << buffer;
227 <      
228 <      int nProc;
229 <      MPI_Comm_size( MPI_COMM_WORLD, &nProc);
230 <      for (int i = 1; i < nProc; ++i) {
243 >    if (createRestFile_) {
244 >      if (worldRank == masterNode) {
245          
246 <        // receive the length of the string buffer that was
247 <        // prepared by processor i
248 <        
249 <        int recvLength;
250 <        MPI_Recv(&recvLength, 1, MPI_INT, i, 0, MPI_COMM_WORLD, istatus);
251 <        char* recvBuffer = new char[recvLength];
238 <        if (recvBuffer == NULL) {
239 <        } else {
240 <          MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD,
241 <                   istatus);
242 <          (*output_) << recvBuffer;
246 >        (*output_) << info_->getSnapshotManager()->getCurrentSnapshot()->getTime();
247 >        (*output_) << buffer;
248 >      
249 >        int nProc;
250 >        MPI_Comm_size( MPI_COMM_WORLD, &nProc);
251 >        for (int i = 1; i < nProc; ++i) {
252            
253 <          delete [] recvBuffer;
254 <        }
255 <      }
256 <      (*output_).flush();
257 <    } else {
258 <      int sendBufferLength = buffer.size() + 1;
259 <      MPI_Send(&sendBufferLength, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
260 <      MPI_Send((void *)buffer.c_str(), sendBufferLength,
261 <               MPI_CHAR, masterNode, 0, MPI_COMM_WORLD);
253 >          // receive the length of the string buffer that was
254 >          // prepared by processor i
255 >          
256 >          int recvLength;
257 >          MPI_Recv(&recvLength, 1, MPI_INT, i, 0, MPI_COMM_WORLD, istatus);
258 >          char* recvBuffer = new char[recvLength];
259 >          if (recvBuffer == NULL) {
260 >          } else {
261 >            MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD,
262 >                     istatus);
263 >            if (createRestFile_) (*output_) << recvBuffer;
264 >            
265 >            delete [] recvBuffer;
266 >          }
267 >        }      
268 >        (*output_).flush();
269 >      } else {
270 >        int sendBufferLength = buffer.size() + 1;
271 >        MPI_Send(&sendBufferLength, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
272 >        MPI_Send((void *)buffer.c_str(), sendBufferLength,
273 >                 MPI_CHAR, masterNode, 0, MPI_COMM_WORLD);
274 >      }
275      }
276   #endif // is_mpi
277    }

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines