ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/io/DumpWriter.cpp
Revision: 2553
Committed: Thu Jan 12 19:56:30 2006 UTC (18 years, 6 months ago) by chuckv
File size: 21400 byte(s)
Log Message:
Changed formating.

File Contents

# User Rev Content
1 gezelter 2204 /*
2 gezelter 1930 * 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 tim 2060 #include "io/basic_teebuf.hpp"
46 tim 2314 #include "io/gzstream.hpp"
47     #include "io/Globals.hpp"
48    
49 gezelter 1490 #ifdef IS_MPI
50     #include <mpi.h>
51 gezelter 1930 #endif //is_mpi
52 gezelter 1490
53 gezelter 1930 namespace oopse {
54 gezelter 1490
55 gezelter 2204 DumpWriter::DumpWriter(SimInfo* info)
56     : info_(info), filename_(info->getDumpFileName()), eorFilename_(info->getFinalConfigFileName()){
57 tim 2314
58     Globals* simParams = info->getSimParams();
59     needCompression_ = simParams->getCompressDumpFile();
60 chrisfen 2425 needForceVector_ = simParams->getOutputForceVector();
61 chuckv 2490 createDumpFile_ = true;
62 tim 2318 #ifdef HAVE_LIBZ
63 tim 2314 if (needCompression_) {
64     filename_ += ".gz";
65     eorFilename_ += ".gz";
66     }
67 tim 2318 #endif
68 tim 2314
69 tim 2060 #ifdef IS_MPI
70    
71 gezelter 2204 if (worldRank == 0) {
72 tim 2060 #endif // is_mpi
73    
74 chuckv 2490
75 tim 2314 dumpFile_ = createOStream(filename_);
76    
77 tim 2060 if (!dumpFile_) {
78 gezelter 2204 sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n",
79     filename_.c_str());
80     painCave.isFatal = 1;
81     simError();
82 tim 2060 }
83    
84     #ifdef IS_MPI
85    
86 gezelter 2204 }
87 tim 2060
88 gezelter 2204 sprintf(checkPointMsg, "Sucessfully opened output file for dumping.\n");
89     MPIcheckPoint();
90 tim 2060
91     #endif // is_mpi
92    
93 gezelter 2204 }
94 tim 2060
95    
96 gezelter 2204 DumpWriter::DumpWriter(SimInfo* info, const std::string& filename)
97     : info_(info), filename_(filename){
98 tim 2314
99     Globals* simParams = info->getSimParams();
100     eorFilename_ = filename_.substr(0, filename_.rfind(".")) + ".eor";
101    
102     needCompression_ = simParams->getCompressDumpFile();
103 chrisfen 2425 needForceVector_ = simParams->getOutputForceVector();
104 chuckv 2490 createDumpFile_ = true;
105 tim 2318 #ifdef HAVE_LIBZ
106 tim 2314 if (needCompression_) {
107     filename_ += ".gz";
108     eorFilename_ += ".gz";
109     }
110 tim 2318 #endif
111 tim 2314
112 gezelter 1930 #ifdef IS_MPI
113 gezelter 1490
114 gezelter 2204 if (worldRank == 0) {
115 gezelter 1930 #endif // is_mpi
116 gezelter 1490
117 chuckv 2490
118 tim 2314 dumpFile_ = createOStream(filename_);
119    
120 gezelter 1930 if (!dumpFile_) {
121 gezelter 2204 sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n",
122     filename_.c_str());
123     painCave.isFatal = 1;
124     simError();
125 gezelter 1930 }
126 gezelter 1490
127     #ifdef IS_MPI
128    
129 gezelter 2204 }
130 gezelter 1490
131 gezelter 2204 sprintf(checkPointMsg, "Sucessfully opened output file for dumping.\n");
132     MPIcheckPoint();
133 gezelter 1490
134     #endif // is_mpi
135 gezelter 1930
136 gezelter 2204 }
137 chuckv 2490
138     DumpWriter::DumpWriter(SimInfo* info, const std::string& filename, bool writeDumpFile)
139     : info_(info), filename_(filename){
140    
141     Globals* simParams = info->getSimParams();
142     eorFilename_ = filename_.substr(0, filename_.rfind(".")) + ".eor";
143    
144     needCompression_ = simParams->getCompressDumpFile();
145     needForceVector_ = simParams->getOutputForceVector();
146    
147     #ifdef HAVE_LIBZ
148     if (needCompression_) {
149     filename_ += ".gz";
150     eorFilename_ += ".gz";
151     }
152     #endif
153    
154     #ifdef IS_MPI
155    
156     if (worldRank == 0) {
157     #endif // is_mpi
158    
159     createDumpFile_ = writeDumpFile;
160     if (createDumpFile_) {
161     dumpFile_ = createOStream(filename_);
162    
163     if (!dumpFile_) {
164     sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n",
165     filename_.c_str());
166     painCave.isFatal = 1;
167     simError();
168     }
169     }
170     #ifdef IS_MPI
171    
172     }
173    
174     sprintf(checkPointMsg, "Sucessfully opened output file for dumping.\n");
175     MPIcheckPoint();
176    
177     #endif // is_mpi
178    
179     }
180    
181    
182    
183    
184    
185 gezelter 1490
186 gezelter 2204 DumpWriter::~DumpWriter() {
187 gezelter 1490
188     #ifdef IS_MPI
189 gezelter 1930
190     if (worldRank == 0) {
191 gezelter 1490 #endif // is_mpi
192 chuckv 2490 if (createDumpFile_){
193     delete dumpFile_;
194     }
195 gezelter 1490 #ifdef IS_MPI
196 gezelter 1930
197     }
198    
199 gezelter 1490 #endif // is_mpi
200 gezelter 1930
201 gezelter 2204 }
202 gezelter 1490
203 gezelter 2204 void DumpWriter::writeCommentLine(std::ostream& os, Snapshot* s) {
204 gezelter 1490
205 gezelter 1930 double currentTime;
206     Mat3x3d hmat;
207     double chi;
208     double integralOfChiDt;
209     Mat3x3d eta;
210    
211     currentTime = s->getTime();
212     hmat = s->getHmat();
213     chi = s->getChi();
214     integralOfChiDt = s->getIntegralOfChiDt();
215     eta = s->getEta();
216    
217     os << currentTime << ";\t"
218 gezelter 2204 << hmat(0, 0) << "\t" << hmat(1, 0) << "\t" << hmat(2, 0) << ";\t"
219     << hmat(0, 1) << "\t" << hmat(1, 1) << "\t" << hmat(2, 1) << ";\t"
220     << hmat(0, 2) << "\t" << hmat(1, 2) << "\t" << hmat(2, 2) << ";\t";
221 gezelter 1490
222 gezelter 1930 //write out additional parameters, such as chi and eta
223    
224 chuckv 2553 os << chi << "\t" << integralOfChiDt << ";\t";
225 gezelter 1930
226     os << eta(0, 0) << "\t" << eta(1, 0) << "\t" << eta(2, 0) << ";\t"
227 gezelter 2204 << eta(0, 1) << "\t" << eta(1, 1) << "\t" << eta(2, 1) << ";\t"
228     << eta(0, 2) << "\t" << eta(1, 2) << "\t" << eta(2, 2) << ";";
229 gezelter 1930
230 tim 2060 os << "\n";
231 gezelter 2204 }
232 gezelter 1490
233 gezelter 2204 void DumpWriter::writeFrame(std::ostream& os) {
234 gezelter 1930 const int BUFFERSIZE = 2000;
235     const int MINIBUFFERSIZE = 100;
236    
237     char tempBuffer[BUFFERSIZE];
238     char writeLine[BUFFERSIZE];
239    
240     Quat4d q;
241     Vector3d ji;
242     Vector3d pos;
243     Vector3d vel;
244 chrisfen 2420 Vector3d frc;
245     Vector3d trq;
246 gezelter 1930
247     Molecule* mol;
248     StuntDouble* integrableObject;
249     SimInfo::MoleculeIterator mi;
250     Molecule::IntegrableObjectIterator ii;
251 gezelter 1490
252 gezelter 1930 int nTotObjects;
253     nTotObjects = info_->getNGlobalIntegrableObjects();
254 gezelter 1490
255 gezelter 1930 #ifndef IS_MPI
256 gezelter 1490
257    
258 gezelter 1930 os << nTotObjects << "\n";
259    
260     writeCommentLine(os, info_->getSnapshotManager()->getCurrentSnapshot());
261 gezelter 1490
262 gezelter 1930 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
263 gezelter 1490
264 gezelter 2204 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
265     integrableObject = mol->nextIntegrableObject(ii)) {
266 gezelter 1930
267 gezelter 1490
268 gezelter 2204 pos = integrableObject->getPos();
269     vel = integrableObject->getVel();
270 gezelter 1490
271 gezelter 2204 sprintf(tempBuffer, "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
272     integrableObject->getType().c_str(),
273     pos[0], pos[1], pos[2],
274     vel[0], vel[1], vel[2]);
275 gezelter 1490
276 gezelter 2204 strcpy(writeLine, tempBuffer);
277 gezelter 1490
278 gezelter 2204 if (integrableObject->isDirectional()) {
279     q = integrableObject->getQ();
280     ji = integrableObject->getJ();
281 gezelter 1490
282 chrisfen 2420 sprintf(tempBuffer, "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf",
283 gezelter 2204 q[0], q[1], q[2], q[3],
284     ji[0], ji[1], ji[2]);
285     strcat(writeLine, tempBuffer);
286     } else {
287 chrisfen 2420 strcat(writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0");
288 gezelter 2204 }
289 gezelter 1490
290 chrisfen 2420 if (needForceVector_) {
291     frc = integrableObject->getFrc();
292     trq = integrableObject->getTrq();
293    
294     sprintf(tempBuffer, "\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf",
295     frc[0], frc[1], frc[2],
296     trq[0], trq[1], trq[2]);
297     strcat(writeLine, tempBuffer);
298     }
299    
300     strcat(writeLine, "\n");
301 gezelter 2204 os << writeLine;
302 gezelter 1490
303 gezelter 2204 }
304 gezelter 1490 }
305    
306 tim 2008 os.flush();
307 gezelter 1930 #else // is_mpi
308     /*********************************************************************
309     * Documentation? You want DOCUMENTATION?
310     *
311     * Why all the potatoes below?
312     *
313     * To make a long story short, the original version of DumpWriter
314     * worked in the most inefficient way possible. Node 0 would
315     * poke each of the node for an individual atom's formatted data
316     * as node 0 worked its way down the global index. This was particularly
317     * inefficient since the method blocked all processors at every atom
318     * (and did it twice!).
319     *
320     * An intermediate version of DumpWriter could be described from Node
321     * zero's perspective as follows:
322     *
323     * 1) Have 100 of your friends stand in a circle.
324     * 2) When you say go, have all of them start tossing potatoes at
325     * you (one at a time).
326     * 3) Catch the potatoes.
327     *
328     * It was an improvement, but MPI has buffers and caches that could
329     * best be described in this analogy as "potato nets", so there's no
330     * need to block the processors atom-by-atom.
331     *
332     * This new and improved DumpWriter works in an even more efficient
333     * way:
334     *
335     * 1) Have 100 of your friend stand in a circle.
336     * 2) When you say go, have them start tossing 5-pound bags of
337     * potatoes at you.
338     * 3) Once you've caught a friend's bag of potatoes,
339     * toss them a spud to let them know they can toss another bag.
340     *
341     * How's THAT for documentation?
342     *
343     *********************************************************************/
344     const int masterNode = 0;
345 gezelter 1490
346 gezelter 1930 int * potatoes;
347     int myPotato;
348     int nProc;
349     int which_node;
350 chrisfen 2420 double atomData[19];
351 gezelter 1930 int isDirectional;
352     char MPIatomTypeString[MINIBUFFERSIZE];
353     int msgLen; // the length of message actually recieved at master nodes
354     int haveError;
355     MPI_Status istatus;
356     int nCurObj;
357    
358     // code to find maximum tag value
359     int * tagub;
360     int flag;
361     int MAXTAG;
362     MPI_Attr_get(MPI_COMM_WORLD, MPI_TAG_UB, &tagub, &flag);
363 gezelter 1490
364 gezelter 1930 if (flag) {
365 gezelter 2204 MAXTAG = *tagub;
366 gezelter 1930 } else {
367 gezelter 2204 MAXTAG = 32767;
368 gezelter 1930 }
369 gezelter 1490
370 gezelter 1930 if (worldRank == masterNode) { //master node (node 0) is responsible for writing the dump file
371 gezelter 1490
372 gezelter 2204 // Node 0 needs a list of the magic potatoes for each processor;
373 gezelter 1490
374 gezelter 2204 MPI_Comm_size(MPI_COMM_WORLD, &nProc);
375     potatoes = new int[nProc];
376 gezelter 1490
377 gezelter 2204 //write out the comment lines
378     for(int i = 0; i < nProc; i++) {
379     potatoes[i] = 0;
380     }
381 gezelter 1490
382    
383 gezelter 2204 os << nTotObjects << "\n";
384     writeCommentLine(os, info_->getSnapshotManager()->getCurrentSnapshot());
385 gezelter 1490
386 gezelter 2204 for(int i = 0; i < info_->getNGlobalMolecules(); i++) {
387 gezelter 1490
388 gezelter 2204 // Get the Node number which has this atom;
389 gezelter 1490
390 gezelter 2204 which_node = info_->getMolToProc(i);
391 gezelter 1490
392 gezelter 2204 if (which_node != masterNode) { //current molecule is in slave node
393     if (potatoes[which_node] + 1 >= MAXTAG) {
394     // The potato was going to exceed the maximum value,
395     // so wrap this processor potato back to 0:
396 gezelter 1490
397 gezelter 2204 potatoes[which_node] = 0;
398     MPI_Send(&potatoes[which_node], 1, MPI_INT, which_node, 0,
399     MPI_COMM_WORLD);
400     }
401 gezelter 1490
402 gezelter 2204 myPotato = potatoes[which_node];
403 gezelter 1490
404 gezelter 2204 //recieve the number of integrableObject in current molecule
405     MPI_Recv(&nCurObj, 1, MPI_INT, which_node, myPotato,
406     MPI_COMM_WORLD, &istatus);
407     myPotato++;
408 gezelter 1490
409 gezelter 2204 for(int l = 0; l < nCurObj; l++) {
410     if (potatoes[which_node] + 2 >= MAXTAG) {
411     // The potato was going to exceed the maximum value,
412     // so wrap this processor potato back to 0:
413 gezelter 1490
414 gezelter 2204 potatoes[which_node] = 0;
415     MPI_Send(&potatoes[which_node], 1, MPI_INT, which_node,
416     0, MPI_COMM_WORLD);
417     }
418 gezelter 1490
419 gezelter 2204 MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR,
420     which_node, myPotato, MPI_COMM_WORLD,
421     &istatus);
422 gezelter 1490
423 gezelter 2204 myPotato++;
424 gezelter 1490
425 chrisfen 2420 MPI_Recv(atomData, 19, MPI_DOUBLE, which_node, myPotato,
426 gezelter 2204 MPI_COMM_WORLD, &istatus);
427     myPotato++;
428 gezelter 1490
429 gezelter 2204 MPI_Get_count(&istatus, MPI_DOUBLE, &msgLen);
430 gezelter 1490
431 chrisfen 2420 if (msgLen == 13 || msgLen == 19)
432 gezelter 2204 isDirectional = 1;
433     else
434     isDirectional = 0;
435 gezelter 1490
436 gezelter 2204 // If we've survived to here, format the line:
437 gezelter 1490
438 gezelter 2204 if (!isDirectional) {
439     sprintf(writeLine, "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
440     MPIatomTypeString, atomData[0],
441     atomData[1], atomData[2],
442     atomData[3], atomData[4],
443     atomData[5]);
444 gezelter 1490
445 gezelter 2204 strcat(writeLine,
446 chrisfen 2420 "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0");
447 gezelter 2204 } else {
448     sprintf(writeLine,
449 chrisfen 2420 "%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",
450 gezelter 2204 MPIatomTypeString,
451     atomData[0],
452     atomData[1],
453     atomData[2],
454     atomData[3],
455     atomData[4],
456     atomData[5],
457     atomData[6],
458     atomData[7],
459     atomData[8],
460     atomData[9],
461     atomData[10],
462     atomData[11],
463     atomData[12]);
464     }
465 chrisfen 2420
466     if (needForceVector_) {
467     if (!isDirectional) {
468     sprintf(writeLine, "\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf",
469     atomData[6],
470     atomData[7],
471     atomData[8],
472     atomData[9],
473     atomData[10],
474     atomData[11]);
475     } else {
476     sprintf(writeLine, "\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf",
477     atomData[13],
478     atomData[14],
479     atomData[15],
480     atomData[16],
481     atomData[17],
482     atomData[18]);
483     }
484     }
485 gezelter 1490
486 chrisfen 2420 sprintf(writeLine, "\n");
487 gezelter 2204 os << writeLine;
488 gezelter 1490
489 gezelter 2204 } // end for(int l =0)
490 gezelter 1490
491 gezelter 2204 potatoes[which_node] = myPotato;
492     } else { //master node has current molecule
493 gezelter 1490
494 gezelter 2204 mol = info_->getMoleculeByGlobalIndex(i);
495 gezelter 1490
496 gezelter 2204 if (mol == NULL) {
497     sprintf(painCave.errMsg, "Molecule not found on node %d!", worldRank);
498     painCave.isFatal = 1;
499     simError();
500     }
501 gezelter 1930
502 gezelter 2204 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
503     integrableObject = mol->nextIntegrableObject(ii)) {
504 gezelter 1490
505 gezelter 2204 pos = integrableObject->getPos();
506     vel = integrableObject->getVel();
507 gezelter 1490
508 gezelter 2204 atomData[0] = pos[0];
509     atomData[1] = pos[1];
510     atomData[2] = pos[2];
511 gezelter 1490
512 gezelter 2204 atomData[3] = vel[0];
513     atomData[4] = vel[1];
514     atomData[5] = vel[2];
515 gezelter 1490
516 gezelter 2204 isDirectional = 0;
517 gezelter 1490
518 gezelter 2204 if (integrableObject->isDirectional()) {
519     isDirectional = 1;
520 gezelter 1490
521 gezelter 2204 q = integrableObject->getQ();
522     ji = integrableObject->getJ();
523 gezelter 1490
524 gezelter 2204 for(int j = 0; j < 6; j++) {
525     atomData[j] = atomData[j];
526     }
527 gezelter 1490
528 gezelter 2204 atomData[6] = q[0];
529     atomData[7] = q[1];
530     atomData[8] = q[2];
531     atomData[9] = q[3];
532 gezelter 1490
533 gezelter 2204 atomData[10] = ji[0];
534     atomData[11] = ji[1];
535     atomData[12] = ji[2];
536     }
537 gezelter 1490
538 chrisfen 2420 if (needForceVector_) {
539     frc = integrableObject->getFrc();
540     trq = integrableObject->getTrq();
541    
542     if (!isDirectional) {
543     atomData[6] = frc[0];
544     atomData[7] = frc[1];
545     atomData[8] = frc[2];
546     atomData[9] = trq[0];
547     atomData[10] = trq[1];
548     atomData[11] = trq[2];
549     } else {
550     atomData[13] = frc[0];
551     atomData[14] = frc[1];
552     atomData[15] = frc[2];
553     atomData[16] = trq[0];
554     atomData[17] = trq[1];
555     atomData[18] = trq[2];
556     }
557     }
558    
559 gezelter 2204 // If we've survived to here, format the line:
560 gezelter 1490
561 gezelter 2204 if (!isDirectional) {
562     sprintf(writeLine, "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
563     integrableObject->getType().c_str(), atomData[0],
564     atomData[1], atomData[2],
565     atomData[3], atomData[4],
566     atomData[5]);
567 gezelter 1490
568 gezelter 2204 strcat(writeLine,
569 chrisfen 2420 "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0");
570 gezelter 2204 } else {
571     sprintf(writeLine,
572 chrisfen 2420 "%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",
573 gezelter 2204 integrableObject->getType().c_str(),
574     atomData[0],
575     atomData[1],
576     atomData[2],
577     atomData[3],
578     atomData[4],
579     atomData[5],
580     atomData[6],
581     atomData[7],
582     atomData[8],
583     atomData[9],
584     atomData[10],
585     atomData[11],
586     atomData[12]);
587     }
588 gezelter 1490
589 chrisfen 2420 if (needForceVector_) {
590     if (!isDirectional) {
591     sprintf(writeLine, "\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf",
592     atomData[6],
593     atomData[7],
594     atomData[8],
595     atomData[9],
596     atomData[10],
597     atomData[11]);
598     } else {
599     sprintf(writeLine, "\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf",
600     atomData[13],
601     atomData[14],
602     atomData[15],
603     atomData[16],
604     atomData[17],
605     atomData[18]);
606     }
607     }
608 gezelter 1490
609 chrisfen 2420 sprintf(writeLine, "\n");
610 gezelter 2204 os << writeLine;
611 gezelter 1490
612 gezelter 2204 } //end for(iter = integrableObject.begin())
613     }
614     } //end for(i = 0; i < mpiSim->getNmol())
615 gezelter 1490
616 gezelter 2204 os.flush();
617 tim 1935
618 gezelter 2204 sprintf(checkPointMsg, "Sucessfully took a dump.\n");
619     MPIcheckPoint();
620 gezelter 1490
621 gezelter 2204 delete [] potatoes;
622 gezelter 1930 } else {
623 gezelter 1490
624 gezelter 2204 // worldRank != 0, so I'm a remote node.
625 gezelter 1490
626 gezelter 2204 // Set my magic potato to 0:
627 gezelter 1490
628 gezelter 2204 myPotato = 0;
629 gezelter 1490
630 gezelter 2204 for(int i = 0; i < info_->getNGlobalMolecules(); i++) {
631 gezelter 1490
632 gezelter 2204 // Am I the node which has this integrableObject?
633     int whichNode = info_->getMolToProc(i);
634     if (whichNode == worldRank) {
635     if (myPotato + 1 >= MAXTAG) {
636 gezelter 1490
637 gezelter 2204 // The potato was going to exceed the maximum value,
638     // so wrap this processor potato back to 0 (and block until
639     // node 0 says we can go:
640 gezelter 1490
641 gezelter 2204 MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD,
642     &istatus);
643     }
644 gezelter 1490
645 gezelter 2204 mol = info_->getMoleculeByGlobalIndex(i);
646 gezelter 1930
647 gezelter 1490
648 gezelter 2204 nCurObj = mol->getNIntegrableObjects();
649 gezelter 1490
650 gezelter 2204 MPI_Send(&nCurObj, 1, MPI_INT, 0, myPotato, MPI_COMM_WORLD);
651     myPotato++;
652 gezelter 1490
653 gezelter 2204 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
654     integrableObject = mol->nextIntegrableObject(ii)) {
655 gezelter 1490
656 gezelter 2204 if (myPotato + 2 >= MAXTAG) {
657 gezelter 1490
658 gezelter 2204 // The potato was going to exceed the maximum value,
659     // so wrap this processor potato back to 0 (and block until
660     // node 0 says we can go:
661 gezelter 1490
662 gezelter 2204 MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD,
663     &istatus);
664     }
665 gezelter 1490
666 gezelter 2204 pos = integrableObject->getPos();
667     vel = integrableObject->getVel();
668 gezelter 1490
669 gezelter 2204 atomData[0] = pos[0];
670     atomData[1] = pos[1];
671     atomData[2] = pos[2];
672 gezelter 1490
673 gezelter 2204 atomData[3] = vel[0];
674     atomData[4] = vel[1];
675     atomData[5] = vel[2];
676 gezelter 1490
677 gezelter 2204 isDirectional = 0;
678 gezelter 1490
679 gezelter 2204 if (integrableObject->isDirectional()) {
680     isDirectional = 1;
681 gezelter 1490
682 gezelter 2204 q = integrableObject->getQ();
683     ji = integrableObject->getJ();
684 gezelter 1490
685 gezelter 2204 atomData[6] = q[0];
686     atomData[7] = q[1];
687     atomData[8] = q[2];
688     atomData[9] = q[3];
689 gezelter 1930
690 gezelter 2204 atomData[10] = ji[0];
691     atomData[11] = ji[1];
692     atomData[12] = ji[2];
693     }
694 gezelter 1930
695 chrisfen 2420 if (needForceVector_) {
696     frc = integrableObject->getFrc();
697     trq = integrableObject->getTrq();
698    
699     if (!isDirectional) {
700     atomData[6] = frc[0];
701     atomData[7] = frc[1];
702     atomData[8] = frc[2];
703    
704     atomData[9] = trq[0];
705     atomData[10] = trq[1];
706     atomData[11] = trq[2];
707     } else {
708     atomData[13] = frc[0];
709     atomData[14] = frc[1];
710     atomData[15] = frc[2];
711    
712     atomData[16] = trq[0];
713     atomData[17] = trq[1];
714     atomData[18] = trq[2];
715     }
716     }
717    
718 gezelter 2204 strncpy(MPIatomTypeString, integrableObject->getType().c_str(), MINIBUFFERSIZE);
719 gezelter 1930
720 gezelter 2204 // null terminate the std::string before sending (just in case):
721     MPIatomTypeString[MINIBUFFERSIZE - 1] = '\0';
722 gezelter 1930
723 gezelter 2204 MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
724     myPotato, MPI_COMM_WORLD);
725 gezelter 1930
726 gezelter 2204 myPotato++;
727 gezelter 1930
728 chrisfen 2420 if (isDirectional && needForceVector_) {
729     MPI_Send(atomData, 19, MPI_DOUBLE, 0, myPotato,
730     MPI_COMM_WORLD);
731     } else if (isDirectional) {
732 gezelter 2204 MPI_Send(atomData, 13, MPI_DOUBLE, 0, myPotato,
733     MPI_COMM_WORLD);
734 chrisfen 2420 } else if (needForceVector_) {
735     MPI_Send(atomData, 12, MPI_DOUBLE, 0, myPotato,
736     MPI_COMM_WORLD);
737 gezelter 2204 } else {
738     MPI_Send(atomData, 6, MPI_DOUBLE, 0, myPotato,
739     MPI_COMM_WORLD);
740     }
741 gezelter 1930
742 gezelter 2204 myPotato++;
743     }
744 gezelter 1930
745 gezelter 2204 }
746 gezelter 1930
747 gezelter 2204 }
748     sprintf(checkPointMsg, "Sucessfully took a dump.\n");
749     MPIcheckPoint();
750 gezelter 1930 }
751    
752     #endif // is_mpi
753    
754 gezelter 2204 }
755 gezelter 1490
756 gezelter 2204 void DumpWriter::writeDump() {
757 tim 2314 writeFrame(*dumpFile_);
758 gezelter 2204 }
759 tim 2060
760 gezelter 2204 void DumpWriter::writeEor() {
761 tim 2314 std::ostream* eorStream;
762 tim 2060
763     #ifdef IS_MPI
764     if (worldRank == 0) {
765     #endif // is_mpi
766    
767 tim 2314 eorStream = createOStream(eorFilename_);
768 tim 2060
769     #ifdef IS_MPI
770     }
771     #endif // is_mpi
772    
773 tim 2314 writeFrame(*eorStream);
774    
775     #ifdef IS_MPI
776     if (worldRank == 0) {
777     #endif // is_mpi
778     delete eorStream;
779    
780     #ifdef IS_MPI
781     }
782     #endif // is_mpi
783    
784 gezelter 2204 }
785 tim 2060
786    
787 gezelter 2204 void DumpWriter::writeDumpAndEor() {
788 tim 2060 std::vector<std::streambuf*> buffers;
789 tim 2314 std::ostream* eorStream;
790 tim 2060 #ifdef IS_MPI
791     if (worldRank == 0) {
792     #endif // is_mpi
793    
794 tim 2314 buffers.push_back(dumpFile_->rdbuf());
795 tim 2060
796 tim 2314 eorStream = createOStream(eorFilename_);
797 tim 2060
798 tim 2314 buffers.push_back(eorStream->rdbuf());
799 tim 2060
800     #ifdef IS_MPI
801     }
802     #endif // is_mpi
803    
804     TeeBuf tbuf(buffers.begin(), buffers.end());
805     std::ostream os(&tbuf);
806    
807     writeFrame(os);
808 tim 2314
809     #ifdef IS_MPI
810     if (worldRank == 0) {
811     #endif // is_mpi
812     delete eorStream;
813    
814     #ifdef IS_MPI
815     }
816     #endif // is_mpi
817 tim 2060
818 gezelter 2204 }
819 tim 2060
820 tim 2314 std::ostream* DumpWriter::createOStream(const std::string& filename) {
821 tim 2318
822 tim 2314 std::ostream* newOStream;
823 tim 2318 #ifdef HAVE_LIBZ
824 tim 2314 if (needCompression_) {
825     newOStream = new ogzstream(filename.c_str());
826     } else {
827     newOStream = new std::ofstream(filename.c_str());
828     }
829 tim 2318 #else
830     newOStream = new std::ofstream(filename.c_str());
831     #endif
832 tim 2314 return newOStream;
833     }
834 tim 2060
835 gezelter 1930 }//end namespace oopse