ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/branches/new_design/OOPSE-2.0/src/io/DumpWriter.cpp
Revision: 1927
Committed: Wed Jan 12 17:14:35 2005 UTC (20 years, 6 months ago) by tim
File size: 18903 byte(s)
Log Message:
change license

File Contents

# User Rev Content
1 tim 1927 /*
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 tim 1819 #include "io/DumpWriter.hpp"
43 tim 1808 #include "primitives/Molecule.hpp"
44 tim 1739 #include "utils/simError.h"
45 tim 1727
46     #ifdef IS_MPI
47     #include <mpi.h>
48     #endif //is_mpi
49 tim 1808
50 tim 1739 namespace oopse {
51 tim 1808
52 tim 1739 DumpWriter::DumpWriter(SimInfo* info, const std::string& filename)
53     : info_(info), filename_(filename){
54 tim 1727 #ifdef IS_MPI
55    
56     if (worldRank == 0) {
57     #endif // is_mpi
58    
59 tim 1832 dumpFile_.open(filename_.c_str(), std::ios::out | std::ios::trunc);
60 tim 1727
61 tim 1808 if (!dumpFile_) {
62 tim 1727 sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n",
63 tim 1808 filename_.c_str());
64 tim 1727 painCave.isFatal = 1;
65     simError();
66     }
67    
68     #ifdef IS_MPI
69    
70     }
71    
72     sprintf(checkPointMsg, "Sucessfully opened output file for dumping.\n");
73     MPIcheckPoint();
74    
75     #endif // is_mpi
76    
77     }
78    
79     DumpWriter::~DumpWriter() {
80    
81     #ifdef IS_MPI
82    
83     if (worldRank == 0) {
84     #endif // is_mpi
85    
86 tim 1808 dumpFile_.close();
87 tim 1727
88     #ifdef IS_MPI
89    
90     }
91    
92     #endif // is_mpi
93    
94     }
95    
96 tim 1808 void DumpWriter::writeCommentLine(std::ostream& os, Snapshot* s) {
97 tim 1727
98 tim 1739 double currentTime;
99     Mat3x3d hmat;
100     double chi;
101     double integralOfChiDt;
102     Mat3x3d eta;
103    
104     currentTime = s->getTime();
105     hmat = s->getHmat();
106     chi = s->getChi();
107     integralOfChiDt = s->getIntegralOfChiDt();
108     eta = s->getEta();
109    
110     os << currentTime << ";\t"
111 tim 1808 << hmat(0, 0) << "\t" << hmat(1, 0) << "\t" << hmat(2, 0) << ";\t"
112     << hmat(0, 1) << "\t" << hmat(1, 1) << "\t" << hmat(2, 1) << ";\t"
113 tim 1867 << hmat(0, 2) << "\t" << hmat(1, 2) << "\t" << hmat(2, 2) << ";\t";
114 tim 1727
115 tim 1739 //write out additional parameters, such as chi and eta
116 tim 1727
117 tim 1808 os << chi << "\t" << integralOfChiDt << "\t;";
118 tim 1727
119 tim 1808 os << eta(0, 0) << "\t" << eta(1, 0) << "\t" << eta(2, 0) << ";\t"
120     << eta(0, 1) << "\t" << eta(1, 1) << "\t" << eta(2, 1) << ";\t"
121     << eta(0, 2) << "\t" << eta(1, 2) << "\t" << eta(2, 2) << ";";
122 tim 1739
123 tim 1830 os << std::endl;
124 tim 1727 }
125    
126 tim 1739 void DumpWriter::writeFrame(std::ostream& os) {
127     const int BUFFERSIZE = 2000;
128     const int MINIBUFFERSIZE = 100;
129 tim 1727
130 tim 1739 char tempBuffer[BUFFERSIZE];
131     char writeLine[BUFFERSIZE];
132 tim 1727
133 tim 1739 Quat4d q;
134     Vector3d ji;
135     Vector3d pos;
136     Vector3d vel;
137 tim 1727
138 tim 1739 Molecule* mol;
139     StuntDouble* integrableObject;
140 tim 1808 SimInfo::MoleculeIterator mi;
141     Molecule::IntegrableObjectIterator ii;
142 tim 1739
143     int nTotObjects;
144     nTotObjects = info_->getNGlobalIntegrableObjects();
145 tim 1727
146 tim 1739 #ifndef IS_MPI
147 tim 1727
148    
149 tim 1739 os << nTotObjects << "\n";
150    
151 tim 1808 writeCommentLine(os, info_->getSnapshotManager()->getCurrentSnapshot());
152 tim 1727
153 tim 1739 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
154 tim 1727
155 tim 1739 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
156     integrableObject = mol->nextIntegrableObject(ii)) {
157    
158 tim 1727
159 tim 1739 pos = integrableObject->getPos();
160     vel = integrableObject->getVel();
161 tim 1727
162 tim 1739 sprintf(tempBuffer, "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
163 tim 1808 integrableObject->getType().c_str(),
164 tim 1739 pos[0], pos[1], pos[2],
165     vel[0], vel[1], vel[2]);
166 tim 1727
167 tim 1739 strcpy(writeLine, tempBuffer);
168 tim 1727
169 tim 1739 if (integrableObject->isDirectional()) {
170     q = integrableObject->getQ();
171     ji = integrableObject->getJ();
172 tim 1727
173 tim 1739 sprintf(tempBuffer, "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
174     q[0], q[1], q[2], q[3],
175     ji[0], ji[1], ji[2]);
176     strcat(writeLine, tempBuffer);
177     } else {
178     strcat(writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n");
179     }
180 tim 1727
181 tim 1739 os << writeLine;
182 tim 1727
183     }
184     }
185    
186 tim 1739 #else // is_mpi
187 tim 1727 /*********************************************************************
188     * Documentation? You want DOCUMENTATION?
189     *
190     * Why all the potatoes below?
191     *
192     * To make a long story short, the original version of DumpWriter
193     * worked in the most inefficient way possible. Node 0 would
194     * poke each of the node for an individual atom's formatted data
195     * as node 0 worked its way down the global index. This was particularly
196     * inefficient since the method blocked all processors at every atom
197     * (and did it twice!).
198     *
199     * An intermediate version of DumpWriter could be described from Node
200     * zero's perspective as follows:
201     *
202     * 1) Have 100 of your friends stand in a circle.
203     * 2) When you say go, have all of them start tossing potatoes at
204     * you (one at a time).
205     * 3) Catch the potatoes.
206     *
207     * It was an improvement, but MPI has buffers and caches that could
208     * best be described in this analogy as "potato nets", so there's no
209     * need to block the processors atom-by-atom.
210     *
211     * This new and improved DumpWriter works in an even more efficient
212     * way:
213     *
214     * 1) Have 100 of your friend stand in a circle.
215     * 2) When you say go, have them start tossing 5-pound bags of
216     * potatoes at you.
217     * 3) Once you've caught a friend's bag of potatoes,
218     * toss them a spud to let them know they can toss another bag.
219     *
220     * How's THAT for documentation?
221     *
222     *********************************************************************/
223 tim 1739 const int masterNode = 0;
224 tim 1727
225     int * potatoes;
226     int myPotato;
227     int nProc;
228     int which_node;
229     double atomData[13];
230     int isDirectional;
231 tim 1739 const char * atomTypeString;
232 tim 1727 char MPIatomTypeString[MINIBUFFERSIZE];
233     int msgLen; // the length of message actually recieved at master nodes
234 tim 1739 int haveError;
235     MPI_Status istatus;
236     int nCurObj;
237    
238     // code to find maximum tag value
239     int * tagub;
240     int flag;
241     int MAXTAG;
242 tim 1727 MPI_Attr_get(MPI_COMM_WORLD, MPI_TAG_UB, &tagub, &flag);
243    
244     if (flag) {
245     MAXTAG = *tagub;
246     } else {
247     MAXTAG = 32767;
248     }
249    
250 tim 1739 if (worldRank == masterNode) { //master node (node 0) is responsible for writing the dump file
251 tim 1727
252     // Node 0 needs a list of the magic potatoes for each processor;
253    
254 tim 1739 MPI_Comm_size(MPI_COMM_WORLD, &nProc);
255 tim 1727 potatoes = new int[nProc];
256    
257     //write out the comment lines
258 tim 1883 for(int i = 0; i < nProc; i++) {
259 tim 1727 potatoes[i] = 0;
260     }
261    
262    
263 tim 1739 os << nTotObjects << "\n";
264 tim 1883 writeCommentLine(os, info_->getSnapshotManager()->getCurrentSnapshot());
265 tim 1727
266 tim 1883 for(int i = 0; i < info_->getNGlobalMolecules(); i++) {
267 tim 1727
268     // Get the Node number which has this atom;
269    
270 tim 1739 which_node = info_->getMolToProc(i);
271 tim 1727
272 tim 1739 if (which_node != masterNode) { //current molecule is in slave node
273 tim 1727 if (potatoes[which_node] + 1 >= MAXTAG) {
274     // The potato was going to exceed the maximum value,
275     // so wrap this processor potato back to 0:
276    
277     potatoes[which_node] = 0;
278     MPI_Send(&potatoes[which_node], 1, MPI_INT, which_node, 0,
279     MPI_COMM_WORLD);
280     }
281    
282     myPotato = potatoes[which_node];
283    
284     //recieve the number of integrableObject in current molecule
285     MPI_Recv(&nCurObj, 1, MPI_INT, which_node, myPotato,
286     MPI_COMM_WORLD, &istatus);
287     myPotato++;
288    
289     for(int l = 0; l < nCurObj; l++) {
290     if (potatoes[which_node] + 2 >= MAXTAG) {
291     // The potato was going to exceed the maximum value,
292     // so wrap this processor potato back to 0:
293    
294     potatoes[which_node] = 0;
295     MPI_Send(&potatoes[which_node], 1, MPI_INT, which_node,
296     0, MPI_COMM_WORLD);
297     }
298    
299     MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR,
300     which_node, myPotato, MPI_COMM_WORLD,
301     &istatus);
302    
303     atomTypeString = MPIatomTypeString;
304    
305     myPotato++;
306    
307     MPI_Recv(atomData, 13, MPI_DOUBLE, which_node, myPotato,
308     MPI_COMM_WORLD, &istatus);
309     myPotato++;
310    
311     MPI_Get_count(&istatus, MPI_DOUBLE, &msgLen);
312    
313     if (msgLen == 13)
314     isDirectional = 1;
315     else
316     isDirectional = 0;
317    
318     // If we've survived to here, format the line:
319    
320     if (!isDirectional) {
321     sprintf(writeLine, "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
322     atomTypeString, atomData[0],
323     atomData[1], atomData[2],
324     atomData[3], atomData[4],
325     atomData[5]);
326    
327     strcat(writeLine,
328     "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n");
329     } else {
330     sprintf(writeLine,
331     "%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",
332     atomTypeString,
333     atomData[0],
334     atomData[1],
335     atomData[2],
336     atomData[3],
337     atomData[4],
338     atomData[5],
339     atomData[6],
340     atomData[7],
341     atomData[8],
342     atomData[9],
343     atomData[10],
344     atomData[11],
345     atomData[12]);
346     }
347    
348 tim 1739 os << writeLine;
349    
350 tim 1727 } // end for(int l =0)
351    
352     potatoes[which_node] = myPotato;
353 tim 1739 } else { //master node has current molecule
354 tim 1727
355 tim 1739 mol = info_->getMoleculeByGlobalIndex(i);
356 tim 1727
357 tim 1739 if (mol == NULL) {
358 tim 1884 sprintf(painCave.errMsg, "Molecule not found on node %d!", worldRank);
359 tim 1739 painCave.isFatal = 1;
360     simError();
361     }
362    
363     for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
364     integrableObject = mol->nextIntegrableObject(ii)) {
365    
366     atomTypeString = integrableObject->getType().c_str();
367 tim 1727
368 tim 1739 pos = integrableObject->getPos();
369     vel = integrableObject->getVel();
370 tim 1727
371     atomData[0] = pos[0];
372     atomData[1] = pos[1];
373     atomData[2] = pos[2];
374    
375     atomData[3] = vel[0];
376     atomData[4] = vel[1];
377     atomData[5] = vel[2];
378    
379     isDirectional = 0;
380    
381 tim 1739 if (integrableObject->isDirectional()) {
382 tim 1727 isDirectional = 1;
383    
384 tim 1739 q = integrableObject->getQ();
385     ji = integrableObject->getJ();
386 tim 1727
387     for(int j = 0; j < 6; j++) {
388     atomData[j] = atomData[j];
389     }
390    
391     atomData[6] = q[0];
392     atomData[7] = q[1];
393     atomData[8] = q[2];
394     atomData[9] = q[3];
395    
396     atomData[10] = ji[0];
397     atomData[11] = ji[1];
398     atomData[12] = ji[2];
399     }
400    
401     // If we've survived to here, format the line:
402    
403     if (!isDirectional) {
404     sprintf(writeLine, "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
405     atomTypeString, atomData[0],
406     atomData[1], atomData[2],
407     atomData[3], atomData[4],
408     atomData[5]);
409    
410     strcat(writeLine,
411     "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n");
412     } else {
413     sprintf(writeLine,
414     "%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",
415     atomTypeString,
416     atomData[0],
417     atomData[1],
418     atomData[2],
419     atomData[3],
420     atomData[4],
421     atomData[5],
422     atomData[6],
423     atomData[7],
424     atomData[8],
425     atomData[9],
426     atomData[10],
427     atomData[11],
428     atomData[12]);
429     }
430    
431 tim 1739
432     os << writeLine;
433    
434 tim 1727 } //end for(iter = integrableObject.begin())
435     }
436     } //end for(i = 0; i < mpiSim->getNmol())
437    
438 tim 1739 os.flush();
439 tim 1727
440     sprintf(checkPointMsg, "Sucessfully took a dump.\n");
441    
442     MPIcheckPoint();
443    
444     delete [] potatoes;
445     } else {
446    
447     // worldRank != 0, so I'm a remote node.
448    
449     // Set my magic potato to 0:
450    
451     myPotato = 0;
452    
453 tim 1883 for(int i = 0; i < info_->getNGlobalMolecules(); i++) {
454 tim 1727
455     // Am I the node which has this integrableObject?
456 tim 1883 int whichNode = info_->getMolToProc(i);
457 tim 1739 if (whichNode == worldRank) {
458 tim 1727 if (myPotato + 1 >= MAXTAG) {
459    
460     // The potato was going to exceed the maximum value,
461     // so wrap this processor potato back to 0 (and block until
462     // node 0 says we can go:
463    
464     MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD,
465     &istatus);
466     }
467    
468 tim 1739 mol = info_->getMoleculeByGlobalIndex(i);
469 tim 1727
470 tim 1739
471     nCurObj = mol->getNIntegrableObjects();
472 tim 1727
473     MPI_Send(&nCurObj, 1, MPI_INT, 0, myPotato, MPI_COMM_WORLD);
474     myPotato++;
475    
476 tim 1739 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
477     integrableObject = mol->nextIntegrableObject(ii)) {
478    
479 tim 1727 if (myPotato + 2 >= MAXTAG) {
480    
481     // The potato was going to exceed the maximum value,
482     // so wrap this processor potato back to 0 (and block until
483     // node 0 says we can go:
484    
485     MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD,
486     &istatus);
487     }
488    
489 tim 1739 atomTypeString = integrableObject->getType().c_str();
490 tim 1727
491 tim 1739 pos = integrableObject->getPos();
492     vel = integrableObject->getVel();
493 tim 1727
494     atomData[0] = pos[0];
495     atomData[1] = pos[1];
496     atomData[2] = pos[2];
497    
498     atomData[3] = vel[0];
499     atomData[4] = vel[1];
500     atomData[5] = vel[2];
501    
502     isDirectional = 0;
503    
504 tim 1739 if (integrableObject->isDirectional()) {
505 tim 1727 isDirectional = 1;
506    
507 tim 1739 q = integrableObject->getQ();
508     ji = integrableObject->getJ();
509 tim 1727
510     atomData[6] = q[0];
511     atomData[7] = q[1];
512     atomData[8] = q[2];
513     atomData[9] = q[3];
514    
515     atomData[10] = ji[0];
516     atomData[11] = ji[1];
517     atomData[12] = ji[2];
518     }
519    
520     strncpy(MPIatomTypeString, atomTypeString, MINIBUFFERSIZE);
521    
522 tim 1829 // null terminate the std::string before sending (just in case):
523 tim 1727 MPIatomTypeString[MINIBUFFERSIZE - 1] = '\0';
524    
525     MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
526     myPotato, MPI_COMM_WORLD);
527    
528     myPotato++;
529    
530     if (isDirectional) {
531     MPI_Send(atomData, 13, MPI_DOUBLE, 0, myPotato,
532     MPI_COMM_WORLD);
533     } else {
534     MPI_Send(atomData, 6, MPI_DOUBLE, 0, myPotato,
535     MPI_COMM_WORLD);
536     }
537    
538     myPotato++;
539     }
540 tim 1739
541 tim 1727 }
542 tim 1739
543 tim 1727 }
544    
545     sprintf(checkPointMsg, "Sucessfully took a dump.\n");
546     MPIcheckPoint();
547     }
548    
549     #endif // is_mpi
550    
551     }
552    
553 tim 1739 }//end namespace oopse