ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/io/DumpWriter.cpp
Revision: 2008
Committed: Sun Feb 13 19:10:25 2005 UTC (19 years, 5 months ago) by tim
File size: 18925 byte(s)
Log Message:
dynamicProps get built

File Contents

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