1 |
< |
#define _LARGEFILE_SOURCE64 |
2 |
< |
#define _FILE_OFFSET_BITS 64 |
1 |
> |
/* |
2 |
> |
* Copyright (c) 2005 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 |
19 |
> |
* notice, this list of conditions and the following disclaimer. |
20 |
> |
* |
21 |
> |
* 3. Redistributions in binary form must reproduce the above copyright |
22 |
> |
* notice, this list of conditions and the following disclaimer in the |
23 |
> |
* documentation and/or other materials provided with the |
24 |
> |
* distribution. |
25 |
> |
* |
26 |
> |
* This software is provided "AS IS," without a warranty of any |
27 |
> |
* kind. All express or implied conditions, representations and |
28 |
> |
* warranties, including any implied warranty of merchantability, |
29 |
> |
* fitness for a particular purpose or non-infringement, are hereby |
30 |
> |
* excluded. The University of Notre Dame and its licensors shall not |
31 |
> |
* be liable for any damages suffered by licensee as a result of |
32 |
> |
* using, modifying or distributing the software or its |
33 |
> |
* derivatives. In no event will the University of Notre Dame or its |
34 |
> |
* licensors be liable for any lost revenue, profit or data, or for |
35 |
> |
* direct, indirect, special, consequential, incidental or punitive |
36 |
> |
* damages, however caused and regardless of the theory of liability, |
37 |
> |
* arising out of the use of or inability to use software, even if the |
38 |
> |
* University of Notre Dame has been advised of the possibility of |
39 |
> |
* such damages. |
40 |
> |
*/ |
41 |
> |
|
42 |
> |
#include "io/DumpWriter.hpp" |
43 |
> |
#include "primitives/Molecule.hpp" |
44 |
> |
#include "utils/simError.h" |
45 |
> |
#include "io/basic_teebuf.hpp" |
46 |
> |
#include "io/gzstream.hpp" |
47 |
> |
#include "io/Globals.hpp" |
48 |
|
|
4 |
– |
#include <string.h> |
5 |
– |
#include <iostream> |
6 |
– |
#include <fstream> |
7 |
– |
#include <algorithm> |
8 |
– |
#include <utility> |
9 |
– |
|
49 |
|
#ifdef IS_MPI |
50 |
|
#include <mpi.h> |
12 |
– |
#include "mpiSimulation.hpp" |
13 |
– |
|
14 |
– |
namespace dWrite{ |
15 |
– |
void DieDieDie( void ); |
16 |
– |
} |
17 |
– |
|
18 |
– |
using namespace dWrite; |
51 |
|
#endif //is_mpi |
52 |
|
|
53 |
< |
#include "ReadWrite.hpp" |
22 |
< |
#include "simError.h" |
53 |
> |
namespace oopse { |
54 |
|
|
55 |
< |
DumpWriter::DumpWriter( SimInfo* the_entry_plug ){ |
55 |
> |
DumpWriter::DumpWriter(SimInfo* info) |
56 |
> |
: info_(info), filename_(info->getDumpFileName()), eorFilename_(info->getFinalConfigFileName()){ |
57 |
|
|
58 |
< |
entry_plug = the_entry_plug; |
58 |
> |
Globals* simParams = info->getSimParams(); |
59 |
> |
needCompression_ = simParams->getCompressDumpFile(); |
60 |
|
|
61 |
+ |
#ifdef HAVE_LIBZ |
62 |
+ |
if (needCompression_) { |
63 |
+ |
filename_ += ".gz"; |
64 |
+ |
eorFilename_ += ".gz"; |
65 |
+ |
} |
66 |
+ |
#endif |
67 |
+ |
|
68 |
|
#ifdef IS_MPI |
69 |
< |
if(worldRank == 0 ){ |
69 |
> |
|
70 |
> |
if (worldRank == 0) { |
71 |
|
#endif // is_mpi |
72 |
|
|
32 |
– |
dumpFile.open(entry_plug->sampleName.c_str(), ios::out | ios::trunc ); |
73 |
|
|
74 |
< |
if( !dumpFile ){ |
74 |
> |
dumpFile_ = createOStream(filename_); |
75 |
|
|
76 |
< |
sprintf( painCave.errMsg, |
77 |
< |
"Could not open \"%s\" for dump output.\n", |
78 |
< |
entry_plug->sampleName.c_str()); |
79 |
< |
painCave.isFatal = 1; |
80 |
< |
simError(); |
81 |
< |
} |
76 |
> |
if (!dumpFile_) { |
77 |
> |
sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n", |
78 |
> |
filename_.c_str()); |
79 |
> |
painCave.isFatal = 1; |
80 |
> |
simError(); |
81 |
> |
} |
82 |
|
|
83 |
|
#ifdef IS_MPI |
44 |
– |
} |
84 |
|
|
85 |
< |
//sort the local atoms by global index |
47 |
< |
sortByGlobalIndex(); |
48 |
< |
|
49 |
< |
sprintf( checkPointMsg, |
50 |
< |
"Sucessfully opened output file for dumping.\n"); |
51 |
< |
MPIcheckPoint(); |
52 |
< |
#endif // is_mpi |
53 |
< |
} |
85 |
> |
} |
86 |
|
|
87 |
< |
DumpWriter::~DumpWriter( ){ |
87 |
> |
sprintf(checkPointMsg, "Sucessfully opened output file for dumping.\n"); |
88 |
> |
MPIcheckPoint(); |
89 |
|
|
57 |
– |
#ifdef IS_MPI |
58 |
– |
if(worldRank == 0 ){ |
90 |
|
#endif // is_mpi |
91 |
|
|
92 |
< |
dumpFile.close(); |
92 |
> |
} |
93 |
|
|
63 |
– |
#ifdef IS_MPI |
64 |
– |
} |
65 |
– |
#endif // is_mpi |
66 |
– |
} |
94 |
|
|
95 |
< |
#ifdef IS_MPI |
95 |
> |
DumpWriter::DumpWriter(SimInfo* info, const std::string& filename) |
96 |
> |
: info_(info), filename_(filename){ |
97 |
|
|
98 |
< |
/** |
99 |
< |
* A hook function to load balancing |
72 |
< |
*/ |
98 |
> |
Globals* simParams = info->getSimParams(); |
99 |
> |
eorFilename_ = filename_.substr(0, filename_.rfind(".")) + ".eor"; |
100 |
|
|
101 |
< |
void DumpWriter::update(){ |
75 |
< |
sortByGlobalIndex(); |
76 |
< |
} |
77 |
< |
|
78 |
< |
/** |
79 |
< |
* Auxiliary sorting function |
80 |
< |
*/ |
81 |
< |
|
82 |
< |
bool indexSortingCriterion(const pair<int, int>& p1, const pair<int, int>& p2){ |
83 |
< |
return p1.second < p2.second; |
84 |
< |
} |
101 |
> |
needCompression_ = simParams->getCompressDumpFile(); |
102 |
|
|
103 |
< |
/** |
104 |
< |
* Sorting the local index by global index |
105 |
< |
*/ |
106 |
< |
|
107 |
< |
void DumpWriter::sortByGlobalIndex(){ |
91 |
< |
Molecule* mols = entry_plug->molecules; |
92 |
< |
indexArray.clear(); |
93 |
< |
|
94 |
< |
for(int i = 0; i < entry_plug->n_mol;i++) |
95 |
< |
indexArray.push_back(make_pair(i, mols[i].getGlobalIndex())); |
96 |
< |
|
97 |
< |
sort(indexArray.begin(), indexArray.end(), indexSortingCriterion); |
98 |
< |
} |
99 |
< |
|
103 |
> |
#ifdef HAVE_LIBZ |
104 |
> |
if (needCompression_) { |
105 |
> |
filename_ += ".gz"; |
106 |
> |
eorFilename_ += ".gz"; |
107 |
> |
} |
108 |
|
#endif |
109 |
< |
|
102 |
< |
void DumpWriter::writeDump(double currentTime){ |
103 |
< |
|
104 |
< |
ofstream finalOut; |
105 |
< |
vector<ofstream*> fileStreams; |
106 |
< |
|
109 |
> |
|
110 |
|
#ifdef IS_MPI |
111 |
< |
if(worldRank == 0 ){ |
112 |
< |
#endif |
110 |
< |
finalOut.open( entry_plug->finalName.c_str(), ios::out | ios::trunc ); |
111 |
< |
if( !finalOut ){ |
112 |
< |
sprintf( painCave.errMsg, |
113 |
< |
"Could not open \"%s\" for final dump output.\n", |
114 |
< |
entry_plug->finalName.c_str() ); |
115 |
< |
painCave.isFatal = 1; |
116 |
< |
simError(); |
117 |
< |
} |
118 |
< |
#ifdef IS_MPI |
119 |
< |
} |
111 |
> |
|
112 |
> |
if (worldRank == 0) { |
113 |
|
#endif // is_mpi |
114 |
|
|
122 |
– |
fileStreams.push_back(&finalOut); |
123 |
– |
fileStreams.push_back(&dumpFile); |
115 |
|
|
116 |
< |
writeFrame(fileStreams, currentTime); |
116 |
> |
dumpFile_ = createOStream(filename_); |
117 |
|
|
118 |
+ |
if (!dumpFile_) { |
119 |
+ |
sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n", |
120 |
+ |
filename_.c_str()); |
121 |
+ |
painCave.isFatal = 1; |
122 |
+ |
simError(); |
123 |
+ |
} |
124 |
+ |
|
125 |
|
#ifdef IS_MPI |
128 |
– |
finalOut.close(); |
129 |
– |
#endif |
130 |
– |
|
131 |
– |
} |
126 |
|
|
127 |
< |
void DumpWriter::writeFinal(double currentTime){ |
127 |
> |
} |
128 |
|
|
129 |
< |
ofstream finalOut; |
130 |
< |
vector<ofstream*> fileStreams; |
129 |
> |
sprintf(checkPointMsg, "Sucessfully opened output file for dumping.\n"); |
130 |
> |
MPIcheckPoint(); |
131 |
|
|
138 |
– |
#ifdef IS_MPI |
139 |
– |
if(worldRank == 0 ){ |
132 |
|
#endif // is_mpi |
133 |
|
|
142 |
– |
finalOut.open( entry_plug->finalName.c_str(), ios::out | ios::trunc ); |
143 |
– |
|
144 |
– |
if( !finalOut ){ |
145 |
– |
sprintf( painCave.errMsg, |
146 |
– |
"Could not open \"%s\" for final dump output.\n", |
147 |
– |
entry_plug->finalName.c_str() ); |
148 |
– |
painCave.isFatal = 1; |
149 |
– |
simError(); |
134 |
|
} |
135 |
|
|
136 |
< |
#ifdef IS_MPI |
153 |
< |
} |
154 |
< |
#endif // is_mpi |
155 |
< |
|
156 |
< |
fileStreams.push_back(&finalOut); |
157 |
< |
writeFrame(fileStreams, currentTime); |
136 |
> |
DumpWriter::~DumpWriter() { |
137 |
|
|
138 |
|
#ifdef IS_MPI |
160 |
– |
finalOut.close(); |
161 |
– |
#endif |
162 |
– |
|
163 |
– |
} |
139 |
|
|
140 |
< |
void DumpWriter::writeFrame( vector<ofstream*>& outFile, double currentTime ){ |
140 |
> |
if (worldRank == 0) { |
141 |
> |
#endif // is_mpi |
142 |
|
|
143 |
< |
const int BUFFERSIZE = 2000; |
168 |
< |
const int MINIBUFFERSIZE = 100; |
143 |
> |
delete dumpFile_; |
144 |
|
|
170 |
– |
char tempBuffer[BUFFERSIZE]; |
171 |
– |
char writeLine[BUFFERSIZE]; |
172 |
– |
|
173 |
– |
int i; |
174 |
– |
unsigned int k; |
175 |
– |
|
145 |
|
#ifdef IS_MPI |
177 |
– |
|
178 |
– |
/********************************************************************* |
179 |
– |
* Documentation? You want DOCUMENTATION? |
180 |
– |
* |
181 |
– |
* Why all the potatoes below? |
182 |
– |
* |
183 |
– |
* To make a long story short, the original version of DumpWriter |
184 |
– |
* worked in the most inefficient way possible. Node 0 would |
185 |
– |
* poke each of the node for an individual atom's formatted data |
186 |
– |
* as node 0 worked its way down the global index. This was particularly |
187 |
– |
* inefficient since the method blocked all processors at every atom |
188 |
– |
* (and did it twice!). |
189 |
– |
* |
190 |
– |
* An intermediate version of DumpWriter could be described from Node |
191 |
– |
* zero's perspective as follows: |
192 |
– |
* |
193 |
– |
* 1) Have 100 of your friends stand in a circle. |
194 |
– |
* 2) When you say go, have all of them start tossing potatoes at |
195 |
– |
* you (one at a time). |
196 |
– |
* 3) Catch the potatoes. |
197 |
– |
* |
198 |
– |
* It was an improvement, but MPI has buffers and caches that could |
199 |
– |
* best be described in this analogy as "potato nets", so there's no |
200 |
– |
* need to block the processors atom-by-atom. |
201 |
– |
* |
202 |
– |
* This new and improved DumpWriter works in an even more efficient |
203 |
– |
* way: |
204 |
– |
* |
205 |
– |
* 1) Have 100 of your friend stand in a circle. |
206 |
– |
* 2) When you say go, have them start tossing 5-pound bags of |
207 |
– |
* potatoes at you. |
208 |
– |
* 3) Once you've caught a friend's bag of potatoes, |
209 |
– |
* toss them a spud to let them know they can toss another bag. |
210 |
– |
* |
211 |
– |
* How's THAT for documentation? |
212 |
– |
* |
213 |
– |
*********************************************************************/ |
146 |
|
|
147 |
< |
int *potatoes; |
216 |
< |
int myPotato; |
147 |
> |
} |
148 |
|
|
149 |
< |
int nProc; |
219 |
< |
int j, which_node, done, which_atom, local_index, currentIndex; |
220 |
< |
double atomData[13]; |
221 |
< |
int isDirectional; |
222 |
< |
char* atomTypeString; |
223 |
< |
char MPIatomTypeString[MINIBUFFERSIZE]; |
224 |
< |
int nObjects; |
225 |
< |
int msgLen; // the length of message actually recieved at master nodes |
226 |
< |
#endif //is_mpi |
149 |
> |
#endif // is_mpi |
150 |
|
|
151 |
< |
double q[4], ji[3]; |
229 |
< |
DirectionalAtom* dAtom; |
230 |
< |
double pos[3], vel[3]; |
231 |
< |
int nTotObjects; |
232 |
< |
StuntDouble* sd; |
233 |
< |
char* molName; |
234 |
< |
vector<StuntDouble*> integrableObjects; |
235 |
< |
vector<StuntDouble*>::iterator iter; |
236 |
< |
nTotObjects = entry_plug->getTotIntegrableObjects(); |
237 |
< |
#ifndef IS_MPI |
238 |
< |
|
239 |
< |
for(k = 0; k < outFile.size(); k++){ |
240 |
< |
*outFile[k] << nTotObjects << "\n"; |
151 |
> |
} |
152 |
|
|
153 |
< |
*outFile[k] << currentTime << ";\t" |
243 |
< |
<< entry_plug->Hmat[0][0] << "\t" |
244 |
< |
<< entry_plug->Hmat[1][0] << "\t" |
245 |
< |
<< entry_plug->Hmat[2][0] << ";\t" |
246 |
< |
|
247 |
< |
<< entry_plug->Hmat[0][1] << "\t" |
248 |
< |
<< entry_plug->Hmat[1][1] << "\t" |
249 |
< |
<< entry_plug->Hmat[2][1] << ";\t" |
153 |
> |
void DumpWriter::writeCommentLine(std::ostream& os, Snapshot* s) { |
154 |
|
|
155 |
< |
<< entry_plug->Hmat[0][2] << "\t" |
156 |
< |
<< entry_plug->Hmat[1][2] << "\t" |
157 |
< |
<< entry_plug->Hmat[2][2] << ";"; |
155 |
> |
double currentTime; |
156 |
> |
Mat3x3d hmat; |
157 |
> |
double chi; |
158 |
> |
double integralOfChiDt; |
159 |
> |
Mat3x3d eta; |
160 |
> |
|
161 |
> |
currentTime = s->getTime(); |
162 |
> |
hmat = s->getHmat(); |
163 |
> |
chi = s->getChi(); |
164 |
> |
integralOfChiDt = s->getIntegralOfChiDt(); |
165 |
> |
eta = s->getEta(); |
166 |
> |
|
167 |
> |
os << currentTime << ";\t" |
168 |
> |
<< hmat(0, 0) << "\t" << hmat(1, 0) << "\t" << hmat(2, 0) << ";\t" |
169 |
> |
<< hmat(0, 1) << "\t" << hmat(1, 1) << "\t" << hmat(2, 1) << ";\t" |
170 |
> |
<< hmat(0, 2) << "\t" << hmat(1, 2) << "\t" << hmat(2, 2) << ";\t"; |
171 |
|
|
172 |
|
//write out additional parameters, such as chi and eta |
173 |
< |
*outFile[k] << entry_plug->the_integrator->getAdditionalParameters() << endl; |
173 |
> |
|
174 |
> |
os << chi << "\t" << integralOfChiDt << "\t;"; |
175 |
> |
|
176 |
> |
os << eta(0, 0) << "\t" << eta(1, 0) << "\t" << eta(2, 0) << ";\t" |
177 |
> |
<< eta(0, 1) << "\t" << eta(1, 1) << "\t" << eta(2, 1) << ";\t" |
178 |
> |
<< eta(0, 2) << "\t" << eta(1, 2) << "\t" << eta(2, 2) << ";"; |
179 |
> |
|
180 |
> |
os << "\n"; |
181 |
|
} |
258 |
– |
|
259 |
– |
for( i=0; i< entry_plug->n_mol; i++ ){ |
182 |
|
|
183 |
< |
integrableObjects = entry_plug->molecules[i].getIntegrableObjects(); |
184 |
< |
molName = (entry_plug->compStamps[entry_plug->molecules[i].getStampID()])->getID(); |
185 |
< |
|
264 |
< |
for( iter = integrableObjects.begin();iter != integrableObjects.end(); ++iter){ |
265 |
< |
sd = *iter; |
266 |
< |
sd->getPos(pos); |
267 |
< |
sd->getVel(vel); |
183 |
> |
void DumpWriter::writeFrame(std::ostream& os) { |
184 |
> |
const int BUFFERSIZE = 2000; |
185 |
> |
const int MINIBUFFERSIZE = 100; |
186 |
|
|
187 |
< |
sprintf( tempBuffer, |
188 |
< |
"%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t", |
271 |
< |
sd->getType(), |
272 |
< |
pos[0], |
273 |
< |
pos[1], |
274 |
< |
pos[2], |
275 |
< |
vel[0], |
276 |
< |
vel[1], |
277 |
< |
vel[2]); |
278 |
< |
strcpy( writeLine, tempBuffer ); |
187 |
> |
char tempBuffer[BUFFERSIZE]; |
188 |
> |
char writeLine[BUFFERSIZE]; |
189 |
|
|
190 |
< |
if( sd->isDirectional() ){ |
190 |
> |
Quat4d q; |
191 |
> |
Vector3d ji; |
192 |
> |
Vector3d pos; |
193 |
> |
Vector3d vel; |
194 |
|
|
195 |
< |
sd->getQ( q ); |
196 |
< |
sd->getJ( ji ); |
195 |
> |
Molecule* mol; |
196 |
> |
StuntDouble* integrableObject; |
197 |
> |
SimInfo::MoleculeIterator mi; |
198 |
> |
Molecule::IntegrableObjectIterator ii; |
199 |
> |
|
200 |
> |
int nTotObjects; |
201 |
> |
nTotObjects = info_->getNGlobalIntegrableObjects(); |
202 |
|
|
203 |
< |
sprintf( tempBuffer, |
286 |
< |
"%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n", |
287 |
< |
q[0], |
288 |
< |
q[1], |
289 |
< |
q[2], |
290 |
< |
q[3], |
291 |
< |
ji[0], |
292 |
< |
ji[1], |
293 |
< |
ji[2]); |
294 |
< |
strcat( writeLine, tempBuffer ); |
295 |
< |
} |
296 |
< |
else |
297 |
< |
strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" ); |
298 |
< |
|
299 |
< |
for(k = 0; k < outFile.size(); k++) |
300 |
< |
*outFile[k] << writeLine; |
301 |
< |
} |
203 |
> |
#ifndef IS_MPI |
204 |
|
|
303 |
– |
} |
205 |
|
|
206 |
< |
#else // is_mpi |
206 |
> |
os << nTotObjects << "\n"; |
207 |
> |
|
208 |
> |
writeCommentLine(os, info_->getSnapshotManager()->getCurrentSnapshot()); |
209 |
|
|
210 |
< |
/* code to find maximum tag value */ |
308 |
< |
|
309 |
< |
int *tagub, flag, MAXTAG; |
310 |
< |
MPI_Attr_get(MPI_COMM_WORLD, MPI_TAG_UB, &tagub, &flag); |
311 |
< |
if (flag) { |
312 |
< |
MAXTAG = *tagub; |
313 |
< |
} else { |
314 |
< |
MAXTAG = 32767; |
315 |
< |
} |
210 |
> |
for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) { |
211 |
|
|
212 |
< |
int haveError; |
212 |
> |
for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
213 |
> |
integrableObject = mol->nextIntegrableObject(ii)) { |
214 |
> |
|
215 |
|
|
216 |
< |
MPI_Status istatus; |
217 |
< |
int nCurObj; |
321 |
< |
int *MolToProcMap = mpiSim->getMolToProcMap(); |
216 |
> |
pos = integrableObject->getPos(); |
217 |
> |
vel = integrableObject->getVel(); |
218 |
|
|
219 |
< |
// write out header and node 0's coordinates |
219 |
> |
sprintf(tempBuffer, "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t", |
220 |
> |
integrableObject->getType().c_str(), |
221 |
> |
pos[0], pos[1], pos[2], |
222 |
> |
vel[0], vel[1], vel[2]); |
223 |
|
|
224 |
< |
if( worldRank == 0 ){ |
224 |
> |
strcpy(writeLine, tempBuffer); |
225 |
|
|
226 |
< |
// Node 0 needs a list of the magic potatoes for each processor; |
226 |
> |
if (integrableObject->isDirectional()) { |
227 |
> |
q = integrableObject->getQ(); |
228 |
> |
ji = integrableObject->getJ(); |
229 |
|
|
230 |
< |
nProc = mpiSim->getNProcessors(); |
231 |
< |
potatoes = new int[nProc]; |
230 |
> |
sprintf(tempBuffer, "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n", |
231 |
> |
q[0], q[1], q[2], q[3], |
232 |
> |
ji[0], ji[1], ji[2]); |
233 |
> |
strcat(writeLine, tempBuffer); |
234 |
> |
} else { |
235 |
> |
strcat(writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n"); |
236 |
> |
} |
237 |
|
|
238 |
< |
//write out the comment lines |
333 |
< |
for (i = 0; i < nProc; i++) |
334 |
< |
potatoes[i] = 0; |
335 |
< |
|
336 |
< |
for(k = 0; k < outFile.size(); k++){ |
337 |
< |
*outFile[k] << nTotObjects << "\n"; |
238 |
> |
os << writeLine; |
239 |
|
|
240 |
< |
*outFile[k] << currentTime << ";\t" |
241 |
< |
<< entry_plug->Hmat[0][0] << "\t" |
341 |
< |
<< entry_plug->Hmat[1][0] << "\t" |
342 |
< |
<< entry_plug->Hmat[2][0] << ";\t" |
240 |
> |
} |
241 |
> |
} |
242 |
|
|
243 |
< |
<< entry_plug->Hmat[0][1] << "\t" |
244 |
< |
<< entry_plug->Hmat[1][1] << "\t" |
245 |
< |
<< entry_plug->Hmat[2][1] << ";\t" |
243 |
> |
os.flush(); |
244 |
> |
#else // is_mpi |
245 |
> |
/********************************************************************* |
246 |
> |
* Documentation? You want DOCUMENTATION? |
247 |
> |
* |
248 |
> |
* Why all the potatoes below? |
249 |
> |
* |
250 |
> |
* To make a long story short, the original version of DumpWriter |
251 |
> |
* worked in the most inefficient way possible. Node 0 would |
252 |
> |
* poke each of the node for an individual atom's formatted data |
253 |
> |
* as node 0 worked its way down the global index. This was particularly |
254 |
> |
* inefficient since the method blocked all processors at every atom |
255 |
> |
* (and did it twice!). |
256 |
> |
* |
257 |
> |
* An intermediate version of DumpWriter could be described from Node |
258 |
> |
* zero's perspective as follows: |
259 |
> |
* |
260 |
> |
* 1) Have 100 of your friends stand in a circle. |
261 |
> |
* 2) When you say go, have all of them start tossing potatoes at |
262 |
> |
* you (one at a time). |
263 |
> |
* 3) Catch the potatoes. |
264 |
> |
* |
265 |
> |
* It was an improvement, but MPI has buffers and caches that could |
266 |
> |
* best be described in this analogy as "potato nets", so there's no |
267 |
> |
* need to block the processors atom-by-atom. |
268 |
> |
* |
269 |
> |
* This new and improved DumpWriter works in an even more efficient |
270 |
> |
* way: |
271 |
> |
* |
272 |
> |
* 1) Have 100 of your friend stand in a circle. |
273 |
> |
* 2) When you say go, have them start tossing 5-pound bags of |
274 |
> |
* potatoes at you. |
275 |
> |
* 3) Once you've caught a friend's bag of potatoes, |
276 |
> |
* toss them a spud to let them know they can toss another bag. |
277 |
> |
* |
278 |
> |
* How's THAT for documentation? |
279 |
> |
* |
280 |
> |
*********************************************************************/ |
281 |
> |
const int masterNode = 0; |
282 |
|
|
283 |
< |
<< entry_plug->Hmat[0][2] << "\t" |
284 |
< |
<< entry_plug->Hmat[1][2] << "\t" |
285 |
< |
<< entry_plug->Hmat[2][2] << ";"; |
286 |
< |
|
287 |
< |
*outFile[k] << entry_plug->the_integrator->getAdditionalParameters() << endl; |
283 |
> |
int * potatoes; |
284 |
> |
int myPotato; |
285 |
> |
int nProc; |
286 |
> |
int which_node; |
287 |
> |
double atomData[13]; |
288 |
> |
int isDirectional; |
289 |
> |
char MPIatomTypeString[MINIBUFFERSIZE]; |
290 |
> |
int msgLen; // the length of message actually recieved at master nodes |
291 |
> |
int haveError; |
292 |
> |
MPI_Status istatus; |
293 |
> |
int nCurObj; |
294 |
> |
|
295 |
> |
// code to find maximum tag value |
296 |
> |
int * tagub; |
297 |
> |
int flag; |
298 |
> |
int MAXTAG; |
299 |
> |
MPI_Attr_get(MPI_COMM_WORLD, MPI_TAG_UB, &tagub, &flag); |
300 |
> |
|
301 |
> |
if (flag) { |
302 |
> |
MAXTAG = *tagub; |
303 |
> |
} else { |
304 |
> |
MAXTAG = 32767; |
305 |
|
} |
306 |
|
|
307 |
< |
currentIndex = 0; |
307 |
> |
if (worldRank == masterNode) { //master node (node 0) is responsible for writing the dump file |
308 |
|
|
309 |
< |
for (i = 0 ; i < mpiSim->getNMolGlobal(); i++ ) { |
358 |
< |
|
359 |
< |
// Get the Node number which has this atom; |
360 |
< |
|
361 |
< |
which_node = MolToProcMap[i]; |
362 |
< |
|
363 |
< |
if (which_node != 0) { |
364 |
< |
|
365 |
< |
if (potatoes[which_node] + 1 >= MAXTAG) { |
366 |
< |
// The potato was going to exceed the maximum value, |
367 |
< |
// so wrap this processor potato back to 0: |
309 |
> |
// Node 0 needs a list of the magic potatoes for each processor; |
310 |
|
|
311 |
< |
potatoes[which_node] = 0; |
312 |
< |
MPI_Send(&potatoes[which_node], 1, MPI_INT, which_node, 0, MPI_COMM_WORLD); |
371 |
< |
|
372 |
< |
} |
311 |
> |
MPI_Comm_size(MPI_COMM_WORLD, &nProc); |
312 |
> |
potatoes = new int[nProc]; |
313 |
|
|
314 |
< |
myPotato = potatoes[which_node]; |
314 |
> |
//write out the comment lines |
315 |
> |
for(int i = 0; i < nProc; i++) { |
316 |
> |
potatoes[i] = 0; |
317 |
> |
} |
318 |
|
|
376 |
– |
//recieve the number of integrableObject in current molecule |
377 |
– |
MPI_Recv(&nCurObj, 1, MPI_INT, which_node, |
378 |
– |
myPotato, MPI_COMM_WORLD, &istatus); |
379 |
– |
myPotato++; |
380 |
– |
|
381 |
– |
for(int l = 0; l < nCurObj; l++){ |
319 |
|
|
320 |
< |
if (potatoes[which_node] + 2 >= MAXTAG) { |
321 |
< |
// The potato was going to exceed the maximum value, |
385 |
< |
// so wrap this processor potato back to 0: |
320 |
> |
os << nTotObjects << "\n"; |
321 |
> |
writeCommentLine(os, info_->getSnapshotManager()->getCurrentSnapshot()); |
322 |
|
|
323 |
< |
potatoes[which_node] = 0; |
388 |
< |
MPI_Send(&potatoes[which_node], 1, MPI_INT, which_node, 0, MPI_COMM_WORLD); |
389 |
< |
|
390 |
< |
} |
323 |
> |
for(int i = 0; i < info_->getNGlobalMolecules(); i++) { |
324 |
|
|
325 |
< |
MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, which_node, |
393 |
< |
myPotato, MPI_COMM_WORLD, &istatus); |
325 |
> |
// Get the Node number which has this atom; |
326 |
|
|
327 |
< |
atomTypeString = MPIatomTypeString; |
327 |
> |
which_node = info_->getMolToProc(i); |
328 |
|
|
329 |
< |
myPotato++; |
329 |
> |
if (which_node != masterNode) { //current molecule is in slave node |
330 |
> |
if (potatoes[which_node] + 1 >= MAXTAG) { |
331 |
> |
// The potato was going to exceed the maximum value, |
332 |
> |
// so wrap this processor potato back to 0: |
333 |
|
|
334 |
< |
MPI_Recv(atomData, 13, MPI_DOUBLE, which_node, myPotato, MPI_COMM_WORLD, &istatus); |
335 |
< |
myPotato++; |
334 |
> |
potatoes[which_node] = 0; |
335 |
> |
MPI_Send(&potatoes[which_node], 1, MPI_INT, which_node, 0, |
336 |
> |
MPI_COMM_WORLD); |
337 |
> |
} |
338 |
|
|
339 |
< |
MPI_Get_count(&istatus, MPI_DOUBLE, &msgLen); |
339 |
> |
myPotato = potatoes[which_node]; |
340 |
|
|
341 |
< |
if(msgLen == 13) |
342 |
< |
isDirectional = 1; |
343 |
< |
else |
344 |
< |
isDirectional = 0; |
345 |
< |
|
346 |
< |
// If we've survived to here, format the line: |
347 |
< |
|
348 |
< |
if (!isDirectional) { |
349 |
< |
|
350 |
< |
sprintf( writeLine, |
351 |
< |
"%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t", |
352 |
< |
atomTypeString, |
353 |
< |
atomData[0], |
354 |
< |
atomData[1], |
355 |
< |
atomData[2], |
356 |
< |
atomData[3], |
357 |
< |
atomData[4], |
358 |
< |
atomData[5]); |
422 |
< |
|
423 |
< |
strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" ); |
424 |
< |
|
425 |
< |
} |
426 |
< |
else { |
427 |
< |
|
428 |
< |
sprintf( writeLine, |
429 |
< |
"%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n", |
430 |
< |
atomTypeString, |
431 |
< |
atomData[0], |
432 |
< |
atomData[1], |
433 |
< |
atomData[2], |
434 |
< |
atomData[3], |
435 |
< |
atomData[4], |
436 |
< |
atomData[5], |
437 |
< |
atomData[6], |
438 |
< |
atomData[7], |
439 |
< |
atomData[8], |
440 |
< |
atomData[9], |
441 |
< |
atomData[10], |
442 |
< |
atomData[11], |
443 |
< |
atomData[12]); |
444 |
< |
|
445 |
< |
} |
446 |
< |
|
447 |
< |
for(k = 0; k < outFile.size(); k++) |
448 |
< |
*outFile[k] << writeLine; |
341 |
> |
//recieve the number of integrableObject in current molecule |
342 |
> |
MPI_Recv(&nCurObj, 1, MPI_INT, which_node, myPotato, |
343 |
> |
MPI_COMM_WORLD, &istatus); |
344 |
> |
myPotato++; |
345 |
> |
|
346 |
> |
for(int l = 0; l < nCurObj; l++) { |
347 |
> |
if (potatoes[which_node] + 2 >= MAXTAG) { |
348 |
> |
// The potato was going to exceed the maximum value, |
349 |
> |
// so wrap this processor potato back to 0: |
350 |
> |
|
351 |
> |
potatoes[which_node] = 0; |
352 |
> |
MPI_Send(&potatoes[which_node], 1, MPI_INT, which_node, |
353 |
> |
0, MPI_COMM_WORLD); |
354 |
> |
} |
355 |
> |
|
356 |
> |
MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, |
357 |
> |
which_node, myPotato, MPI_COMM_WORLD, |
358 |
> |
&istatus); |
359 |
|
|
360 |
< |
}// end for(int l =0) |
451 |
< |
potatoes[which_node] = myPotato; |
360 |
> |
myPotato++; |
361 |
|
|
362 |
< |
} |
363 |
< |
else { |
364 |
< |
|
456 |
< |
haveError = 0; |
457 |
< |
|
458 |
< |
local_index = indexArray[currentIndex].first; |
362 |
> |
MPI_Recv(atomData, 13, MPI_DOUBLE, which_node, myPotato, |
363 |
> |
MPI_COMM_WORLD, &istatus); |
364 |
> |
myPotato++; |
365 |
|
|
366 |
< |
integrableObjects = (entry_plug->molecules[local_index]).getIntegrableObjects(); |
366 |
> |
MPI_Get_count(&istatus, MPI_DOUBLE, &msgLen); |
367 |
|
|
368 |
< |
for(iter= integrableObjects.begin(); iter != integrableObjects.end(); ++iter){ |
369 |
< |
sd = *iter; |
370 |
< |
atomTypeString = sd->getType(); |
371 |
< |
|
466 |
< |
sd->getPos(pos); |
467 |
< |
sd->getVel(vel); |
468 |
< |
|
469 |
< |
atomData[0] = pos[0]; |
470 |
< |
atomData[1] = pos[1]; |
471 |
< |
atomData[2] = pos[2]; |
368 |
> |
if (msgLen == 13) |
369 |
> |
isDirectional = 1; |
370 |
> |
else |
371 |
> |
isDirectional = 0; |
372 |
|
|
373 |
< |
atomData[3] = vel[0]; |
474 |
< |
atomData[4] = vel[1]; |
475 |
< |
atomData[5] = vel[2]; |
476 |
< |
|
477 |
< |
isDirectional = 0; |
373 |
> |
// If we've survived to here, format the line: |
374 |
|
|
375 |
< |
if( sd->isDirectional() ){ |
375 |
> |
if (!isDirectional) { |
376 |
> |
sprintf(writeLine, "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t", |
377 |
> |
MPIatomTypeString, atomData[0], |
378 |
> |
atomData[1], atomData[2], |
379 |
> |
atomData[3], atomData[4], |
380 |
> |
atomData[5]); |
381 |
|
|
382 |
< |
isDirectional = 1; |
383 |
< |
|
384 |
< |
sd->getQ( q ); |
385 |
< |
sd->getJ( ji ); |
382 |
> |
strcat(writeLine, |
383 |
> |
"0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n"); |
384 |
> |
} else { |
385 |
> |
sprintf(writeLine, |
386 |
> |
"%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n", |
387 |
> |
MPIatomTypeString, |
388 |
> |
atomData[0], |
389 |
> |
atomData[1], |
390 |
> |
atomData[2], |
391 |
> |
atomData[3], |
392 |
> |
atomData[4], |
393 |
> |
atomData[5], |
394 |
> |
atomData[6], |
395 |
> |
atomData[7], |
396 |
> |
atomData[8], |
397 |
> |
atomData[9], |
398 |
> |
atomData[10], |
399 |
> |
atomData[11], |
400 |
> |
atomData[12]); |
401 |
> |
} |
402 |
|
|
403 |
< |
for (int j = 0; j < 6 ; j++) |
487 |
< |
atomData[j] = atomData[j]; |
488 |
< |
|
489 |
< |
atomData[6] = q[0]; |
490 |
< |
atomData[7] = q[1]; |
491 |
< |
atomData[8] = q[2]; |
492 |
< |
atomData[9] = q[3]; |
493 |
< |
|
494 |
< |
atomData[10] = ji[0]; |
495 |
< |
atomData[11] = ji[1]; |
496 |
< |
atomData[12] = ji[2]; |
497 |
< |
} |
498 |
< |
|
499 |
< |
// If we've survived to here, format the line: |
500 |
< |
|
501 |
< |
if (!isDirectional) { |
502 |
< |
|
503 |
< |
sprintf( writeLine, |
504 |
< |
"%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t", |
505 |
< |
atomTypeString, |
506 |
< |
atomData[0], |
507 |
< |
atomData[1], |
508 |
< |
atomData[2], |
509 |
< |
atomData[3], |
510 |
< |
atomData[4], |
511 |
< |
atomData[5]); |
512 |
< |
|
513 |
< |
strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" ); |
514 |
< |
|
515 |
< |
} |
516 |
< |
else { |
517 |
< |
|
518 |
< |
sprintf( writeLine, |
519 |
< |
"%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n", |
520 |
< |
atomTypeString, |
521 |
< |
atomData[0], |
522 |
< |
atomData[1], |
523 |
< |
atomData[2], |
524 |
< |
atomData[3], |
525 |
< |
atomData[4], |
526 |
< |
atomData[5], |
527 |
< |
atomData[6], |
528 |
< |
atomData[7], |
529 |
< |
atomData[8], |
530 |
< |
atomData[9], |
531 |
< |
atomData[10], |
532 |
< |
atomData[11], |
533 |
< |
atomData[12]); |
534 |
< |
|
535 |
< |
} |
536 |
< |
|
537 |
< |
for(k = 0; k < outFile.size(); k++) |
538 |
< |
*outFile[k] << writeLine; |
539 |
< |
|
540 |
< |
|
541 |
< |
}//end for(iter = integrableObject.begin()) |
542 |
< |
|
543 |
< |
currentIndex++; |
544 |
< |
} |
403 |
> |
os << writeLine; |
404 |
|
|
405 |
< |
}//end for(i = 0; i < mpiSim->getNmol()) |
547 |
< |
|
548 |
< |
for(k = 0; k < outFile.size(); k++) |
549 |
< |
outFile[k]->flush(); |
550 |
< |
|
551 |
< |
sprintf( checkPointMsg, |
552 |
< |
"Sucessfully took a dump.\n"); |
553 |
< |
|
554 |
< |
MPIcheckPoint(); |
555 |
< |
|
556 |
< |
delete[] potatoes; |
557 |
< |
|
558 |
< |
} else { |
405 |
> |
} // end for(int l =0) |
406 |
|
|
407 |
< |
// worldRank != 0, so I'm a remote node. |
407 |
> |
potatoes[which_node] = myPotato; |
408 |
> |
} else { //master node has current molecule |
409 |
|
|
410 |
< |
// Set my magic potato to 0: |
410 |
> |
mol = info_->getMoleculeByGlobalIndex(i); |
411 |
|
|
412 |
< |
myPotato = 0; |
413 |
< |
currentIndex = 0; |
414 |
< |
|
415 |
< |
for (i = 0 ; i < mpiSim->getNMolGlobal(); i++ ) { |
416 |
< |
|
417 |
< |
// Am I the node which has this integrableObject? |
418 |
< |
|
419 |
< |
if (MolToProcMap[i] == worldRank) { |
412 |
> |
if (mol == NULL) { |
413 |
> |
sprintf(painCave.errMsg, "Molecule not found on node %d!", worldRank); |
414 |
> |
painCave.isFatal = 1; |
415 |
> |
simError(); |
416 |
> |
} |
417 |
> |
|
418 |
> |
for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
419 |
> |
integrableObject = mol->nextIntegrableObject(ii)) { |
420 |
|
|
421 |
+ |
pos = integrableObject->getPos(); |
422 |
+ |
vel = integrableObject->getVel(); |
423 |
|
|
424 |
< |
if (myPotato + 1 >= MAXTAG) { |
425 |
< |
|
426 |
< |
// The potato was going to exceed the maximum value, |
577 |
< |
// so wrap this processor potato back to 0 (and block until |
578 |
< |
// node 0 says we can go: |
579 |
< |
|
580 |
< |
MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD, &istatus); |
581 |
< |
|
582 |
< |
} |
424 |
> |
atomData[0] = pos[0]; |
425 |
> |
atomData[1] = pos[1]; |
426 |
> |
atomData[2] = pos[2]; |
427 |
|
|
428 |
< |
local_index = indexArray[currentIndex].first; |
429 |
< |
integrableObjects = entry_plug->molecules[local_index].getIntegrableObjects(); |
430 |
< |
|
587 |
< |
nCurObj = integrableObjects.size(); |
588 |
< |
|
589 |
< |
MPI_Send(&nCurObj, 1, MPI_INT, 0, |
590 |
< |
myPotato, MPI_COMM_WORLD); |
591 |
< |
myPotato++; |
428 |
> |
atomData[3] = vel[0]; |
429 |
> |
atomData[4] = vel[1]; |
430 |
> |
atomData[5] = vel[2]; |
431 |
|
|
432 |
< |
for( iter = integrableObjects.begin(); iter != integrableObjects.end(); iter++){ |
432 |
> |
isDirectional = 0; |
433 |
|
|
434 |
< |
if (myPotato + 2 >= MAXTAG) { |
435 |
< |
|
597 |
< |
// The potato was going to exceed the maximum value, |
598 |
< |
// so wrap this processor potato back to 0 (and block until |
599 |
< |
// node 0 says we can go: |
600 |
< |
|
601 |
< |
MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD, &istatus); |
602 |
< |
|
603 |
< |
} |
604 |
< |
|
605 |
< |
sd = *iter; |
606 |
< |
|
607 |
< |
atomTypeString = sd->getType(); |
434 |
> |
if (integrableObject->isDirectional()) { |
435 |
> |
isDirectional = 1; |
436 |
|
|
437 |
< |
sd->getPos(pos); |
438 |
< |
sd->getVel(vel); |
437 |
> |
q = integrableObject->getQ(); |
438 |
> |
ji = integrableObject->getJ(); |
439 |
|
|
440 |
< |
atomData[0] = pos[0]; |
441 |
< |
atomData[1] = pos[1]; |
442 |
< |
atomData[2] = pos[2]; |
440 |
> |
for(int j = 0; j < 6; j++) { |
441 |
> |
atomData[j] = atomData[j]; |
442 |
> |
} |
443 |
|
|
444 |
< |
atomData[3] = vel[0]; |
445 |
< |
atomData[4] = vel[1]; |
446 |
< |
atomData[5] = vel[2]; |
447 |
< |
|
620 |
< |
isDirectional = 0; |
444 |
> |
atomData[6] = q[0]; |
445 |
> |
atomData[7] = q[1]; |
446 |
> |
atomData[8] = q[2]; |
447 |
> |
atomData[9] = q[3]; |
448 |
|
|
449 |
< |
if( sd->isDirectional() ){ |
449 |
> |
atomData[10] = ji[0]; |
450 |
> |
atomData[11] = ji[1]; |
451 |
> |
atomData[12] = ji[2]; |
452 |
> |
} |
453 |
|
|
454 |
< |
isDirectional = 1; |
455 |
< |
|
456 |
< |
sd->getQ( q ); |
457 |
< |
sd->getJ( ji ); |
458 |
< |
|
454 |
> |
// If we've survived to here, format the line: |
455 |
> |
|
456 |
> |
if (!isDirectional) { |
457 |
> |
sprintf(writeLine, "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t", |
458 |
> |
integrableObject->getType().c_str(), atomData[0], |
459 |
> |
atomData[1], atomData[2], |
460 |
> |
atomData[3], atomData[4], |
461 |
> |
atomData[5]); |
462 |
> |
|
463 |
> |
strcat(writeLine, |
464 |
> |
"0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n"); |
465 |
> |
} else { |
466 |
> |
sprintf(writeLine, |
467 |
> |
"%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n", |
468 |
> |
integrableObject->getType().c_str(), |
469 |
> |
atomData[0], |
470 |
> |
atomData[1], |
471 |
> |
atomData[2], |
472 |
> |
atomData[3], |
473 |
> |
atomData[4], |
474 |
> |
atomData[5], |
475 |
> |
atomData[6], |
476 |
> |
atomData[7], |
477 |
> |
atomData[8], |
478 |
> |
atomData[9], |
479 |
> |
atomData[10], |
480 |
> |
atomData[11], |
481 |
> |
atomData[12]); |
482 |
> |
} |
483 |
> |
|
484 |
> |
|
485 |
> |
os << writeLine; |
486 |
> |
|
487 |
> |
} //end for(iter = integrableObject.begin()) |
488 |
> |
} |
489 |
> |
} //end for(i = 0; i < mpiSim->getNmol()) |
490 |
> |
|
491 |
> |
os.flush(); |
492 |
> |
|
493 |
> |
sprintf(checkPointMsg, "Sucessfully took a dump.\n"); |
494 |
> |
MPIcheckPoint(); |
495 |
> |
|
496 |
> |
delete [] potatoes; |
497 |
> |
} else { |
498 |
> |
|
499 |
> |
// worldRank != 0, so I'm a remote node. |
500 |
> |
|
501 |
> |
// Set my magic potato to 0: |
502 |
> |
|
503 |
> |
myPotato = 0; |
504 |
> |
|
505 |
> |
for(int i = 0; i < info_->getNGlobalMolecules(); i++) { |
506 |
> |
|
507 |
> |
// Am I the node which has this integrableObject? |
508 |
> |
int whichNode = info_->getMolToProc(i); |
509 |
> |
if (whichNode == worldRank) { |
510 |
> |
if (myPotato + 1 >= MAXTAG) { |
511 |
> |
|
512 |
> |
// The potato was going to exceed the maximum value, |
513 |
> |
// so wrap this processor potato back to 0 (and block until |
514 |
> |
// node 0 says we can go: |
515 |
> |
|
516 |
> |
MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD, |
517 |
> |
&istatus); |
518 |
> |
} |
519 |
> |
|
520 |
> |
mol = info_->getMoleculeByGlobalIndex(i); |
521 |
> |
|
522 |
|
|
523 |
< |
atomData[6] = q[0]; |
631 |
< |
atomData[7] = q[1]; |
632 |
< |
atomData[8] = q[2]; |
633 |
< |
atomData[9] = q[3]; |
634 |
< |
|
635 |
< |
atomData[10] = ji[0]; |
636 |
< |
atomData[11] = ji[1]; |
637 |
< |
atomData[12] = ji[2]; |
638 |
< |
} |
523 |
> |
nCurObj = mol->getNIntegrableObjects(); |
524 |
|
|
525 |
< |
|
526 |
< |
strncpy(MPIatomTypeString, atomTypeString, MINIBUFFERSIZE); |
525 |
> |
MPI_Send(&nCurObj, 1, MPI_INT, 0, myPotato, MPI_COMM_WORLD); |
526 |
> |
myPotato++; |
527 |
|
|
528 |
< |
// null terminate the string before sending (just in case): |
529 |
< |
MPIatomTypeString[MINIBUFFERSIZE-1] = '\0'; |
528 |
> |
for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
529 |
> |
integrableObject = mol->nextIntegrableObject(ii)) { |
530 |
|
|
531 |
< |
MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0, |
647 |
< |
myPotato, MPI_COMM_WORLD); |
648 |
< |
|
649 |
< |
myPotato++; |
650 |
< |
|
651 |
< |
if (isDirectional) { |
531 |
> |
if (myPotato + 2 >= MAXTAG) { |
532 |
|
|
533 |
< |
MPI_Send(atomData, 13, MPI_DOUBLE, 0, |
534 |
< |
myPotato, MPI_COMM_WORLD); |
535 |
< |
|
656 |
< |
} else { |
533 |
> |
// The potato was going to exceed the maximum value, |
534 |
> |
// so wrap this processor potato back to 0 (and block until |
535 |
> |
// node 0 says we can go: |
536 |
|
|
537 |
< |
MPI_Send(atomData, 6, MPI_DOUBLE, 0, |
538 |
< |
myPotato, MPI_COMM_WORLD); |
539 |
< |
} |
537 |
> |
MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD, |
538 |
> |
&istatus); |
539 |
> |
} |
540 |
|
|
541 |
< |
myPotato++; |
541 |
> |
pos = integrableObject->getPos(); |
542 |
> |
vel = integrableObject->getVel(); |
543 |
|
|
544 |
< |
} |
544 |
> |
atomData[0] = pos[0]; |
545 |
> |
atomData[1] = pos[1]; |
546 |
> |
atomData[2] = pos[2]; |
547 |
|
|
548 |
< |
currentIndex++; |
549 |
< |
|
550 |
< |
} |
551 |
< |
|
548 |
> |
atomData[3] = vel[0]; |
549 |
> |
atomData[4] = vel[1]; |
550 |
> |
atomData[5] = vel[2]; |
551 |
> |
|
552 |
> |
isDirectional = 0; |
553 |
> |
|
554 |
> |
if (integrableObject->isDirectional()) { |
555 |
> |
isDirectional = 1; |
556 |
> |
|
557 |
> |
q = integrableObject->getQ(); |
558 |
> |
ji = integrableObject->getJ(); |
559 |
> |
|
560 |
> |
atomData[6] = q[0]; |
561 |
> |
atomData[7] = q[1]; |
562 |
> |
atomData[8] = q[2]; |
563 |
> |
atomData[9] = q[3]; |
564 |
> |
|
565 |
> |
atomData[10] = ji[0]; |
566 |
> |
atomData[11] = ji[1]; |
567 |
> |
atomData[12] = ji[2]; |
568 |
> |
} |
569 |
> |
|
570 |
> |
strncpy(MPIatomTypeString, integrableObject->getType().c_str(), MINIBUFFERSIZE); |
571 |
> |
|
572 |
> |
// null terminate the std::string before sending (just in case): |
573 |
> |
MPIatomTypeString[MINIBUFFERSIZE - 1] = '\0'; |
574 |
> |
|
575 |
> |
MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0, |
576 |
> |
myPotato, MPI_COMM_WORLD); |
577 |
> |
|
578 |
> |
myPotato++; |
579 |
> |
|
580 |
> |
if (isDirectional) { |
581 |
> |
MPI_Send(atomData, 13, MPI_DOUBLE, 0, myPotato, |
582 |
> |
MPI_COMM_WORLD); |
583 |
> |
} else { |
584 |
> |
MPI_Send(atomData, 6, MPI_DOUBLE, 0, myPotato, |
585 |
> |
MPI_COMM_WORLD); |
586 |
> |
} |
587 |
> |
|
588 |
> |
myPotato++; |
589 |
> |
} |
590 |
> |
|
591 |
> |
} |
592 |
> |
|
593 |
|
} |
594 |
+ |
sprintf(checkPointMsg, "Sucessfully took a dump.\n"); |
595 |
+ |
MPIcheckPoint(); |
596 |
+ |
} |
597 |
|
|
598 |
< |
sprintf( checkPointMsg, |
599 |
< |
"Sucessfully took a dump.\n"); |
600 |
< |
MPIcheckPoint(); |
598 |
> |
#endif // is_mpi |
599 |
> |
|
600 |
> |
} |
601 |
> |
|
602 |
> |
void DumpWriter::writeDump() { |
603 |
> |
writeFrame(*dumpFile_); |
604 |
> |
} |
605 |
> |
|
606 |
> |
void DumpWriter::writeEor() { |
607 |
> |
std::ostream* eorStream; |
608 |
|
|
609 |
+ |
#ifdef IS_MPI |
610 |
+ |
if (worldRank == 0) { |
611 |
+ |
#endif // is_mpi |
612 |
+ |
|
613 |
+ |
eorStream = createOStream(eorFilename_); |
614 |
+ |
|
615 |
+ |
#ifdef IS_MPI |
616 |
|
} |
617 |
+ |
#endif // is_mpi |
618 |
|
|
619 |
+ |
writeFrame(*eorStream); |
620 |
|
|
621 |
< |
|
621 |
> |
#ifdef IS_MPI |
622 |
> |
if (worldRank == 0) { |
623 |
|
#endif // is_mpi |
624 |
< |
} |
624 |
> |
delete eorStream; |
625 |
|
|
626 |
|
#ifdef IS_MPI |
627 |
+ |
} |
628 |
+ |
#endif // is_mpi |
629 |
|
|
630 |
< |
// a couple of functions to let us escape the write loop |
630 |
> |
} |
631 |
|
|
687 |
– |
void dWrite::DieDieDie( void ){ |
632 |
|
|
633 |
< |
MPI_Finalize(); |
634 |
< |
exit (0); |
633 |
> |
void DumpWriter::writeDumpAndEor() { |
634 |
> |
std::vector<std::streambuf*> buffers; |
635 |
> |
std::ostream* eorStream; |
636 |
> |
#ifdef IS_MPI |
637 |
> |
if (worldRank == 0) { |
638 |
> |
#endif // is_mpi |
639 |
> |
|
640 |
> |
buffers.push_back(dumpFile_->rdbuf()); |
641 |
> |
|
642 |
> |
eorStream = createOStream(eorFilename_); |
643 |
> |
|
644 |
> |
buffers.push_back(eorStream->rdbuf()); |
645 |
> |
|
646 |
> |
#ifdef IS_MPI |
647 |
> |
} |
648 |
> |
#endif // is_mpi |
649 |
> |
|
650 |
> |
TeeBuf tbuf(buffers.begin(), buffers.end()); |
651 |
> |
std::ostream os(&tbuf); |
652 |
> |
|
653 |
> |
writeFrame(os); |
654 |
> |
|
655 |
> |
#ifdef IS_MPI |
656 |
> |
if (worldRank == 0) { |
657 |
> |
#endif // is_mpi |
658 |
> |
delete eorStream; |
659 |
> |
|
660 |
> |
#ifdef IS_MPI |
661 |
> |
} |
662 |
> |
#endif // is_mpi |
663 |
> |
|
664 |
> |
} |
665 |
> |
|
666 |
> |
std::ostream* DumpWriter::createOStream(const std::string& filename) { |
667 |
> |
|
668 |
> |
std::ostream* newOStream; |
669 |
> |
#ifdef HAVE_LIBZ |
670 |
> |
if (needCompression_) { |
671 |
> |
newOStream = new ogzstream(filename.c_str()); |
672 |
> |
} else { |
673 |
> |
newOStream = new std::ofstream(filename.c_str()); |
674 |
> |
} |
675 |
> |
#else |
676 |
> |
newOStream = new std::ofstream(filename.c_str()); |
677 |
> |
#endif |
678 |
> |
return newOStream; |
679 |
|
} |
680 |
|
|
681 |
< |
#endif //is_mpi |
681 |
> |
}//end namespace oopse |