ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/DumpWriter.cpp
Revision: 910
Committed: Thu Jan 8 18:05:37 2004 UTC (20 years, 5 months ago) by gezelter
File size: 19466 byte(s)
Log Message:
added strncpy to DumpWriter

File Contents

# User Rev Content
1 mmeineke 723 #define _FILE_OFFSET_BITS 64
2    
3 gezelter 829 #include <string.h>
4 mmeineke 377 #include <iostream>
5     #include <fstream>
6    
7     #ifdef IS_MPI
8     #include <mpi.h>
9     #include "mpiSimulation.hpp"
10 mmeineke 440
11     namespace dWrite{
12 gezelter 907 void DieDieDie( void );
13 mmeineke 440 }
14    
15     using namespace dWrite;
16 mmeineke 377 #endif //is_mpi
17    
18     #include "ReadWrite.hpp"
19     #include "simError.h"
20    
21     DumpWriter::DumpWriter( SimInfo* the_entry_plug ){
22    
23     entry_plug = the_entry_plug;
24    
25     #ifdef IS_MPI
26     if(worldRank == 0 ){
27     #endif // is_mpi
28 tim 837
29 mmeineke 377 strcpy( outName, entry_plug->sampleName );
30 tim 837
31 mmeineke 377 outFile.open(outName, ios::out | ios::trunc );
32 tim 837
33 mmeineke 377 if( !outFile ){
34 tim 837
35 mmeineke 377 sprintf( painCave.errMsg,
36     "Could not open \"%s\" for dump output.\n",
37     outName);
38     painCave.isFatal = 1;
39     simError();
40     }
41 mmeineke 469
42 mmeineke 377 //outFile.setf( ios::scientific );
43    
44     #ifdef IS_MPI
45     }
46    
47     sprintf( checkPointMsg,
48     "Sucessfully opened output file for dumping.\n");
49     MPIcheckPoint();
50     #endif // is_mpi
51     }
52    
53     DumpWriter::~DumpWriter( ){
54    
55     #ifdef IS_MPI
56     if(worldRank == 0 ){
57     #endif // is_mpi
58    
59     outFile.close();
60    
61     #ifdef IS_MPI
62     }
63     #endif // is_mpi
64     }
65    
66     void DumpWriter::writeDump( double currentTime ){
67 tim 837
68 mmeineke 377 const int BUFFERSIZE = 2000;
69 gezelter 907 const int MINIBUFFERSIZE = 10;
70    
71 mmeineke 377 char tempBuffer[BUFFERSIZE];
72     char writeLine[BUFFERSIZE];
73    
74 mmeineke 787 int i;
75     #ifdef IS_MPI
76     int j, which_node, done, which_atom, local_index;
77 gezelter 907 double atomTransData[6];
78     double atomOrientData[7];
79     int isDirectional;
80     char* atomTypeString;
81 gezelter 910 char MPIatomTypeString[MINIBUFFERSIZE];
82 gezelter 907 int me;
83     int atomTypeTag;
84     int atomIsDirectionalTag;
85     int atomTransDataTag;
86     int atomOrientDataTag;
87 mmeineke 787 #else //is_mpi
88     int nAtoms = entry_plug->n_atoms;
89     #endif //is_mpi
90    
91 mmeineke 377 double q[4];
92     DirectionalAtom* dAtom;
93     Atom** atoms = entry_plug->atoms;
94 mmeineke 670 double pos[3], vel[3];
95 tim 837
96 mmeineke 804 // write current frame to the eor file
97 mmeineke 377
98 mmeineke 804 this->writeFinal( currentTime );
99    
100 mmeineke 377 #ifndef IS_MPI
101 tim 837
102 mmeineke 377 outFile << nAtoms << "\n";
103 tim 837
104     outFile << currentTime << ";\t"
105 mmeineke 590 << entry_plug->Hmat[0][0] << "\t"
106     << entry_plug->Hmat[1][0] << "\t"
107     << entry_plug->Hmat[2][0] << ";\t"
108 mmeineke 572
109 mmeineke 590 << entry_plug->Hmat[0][1] << "\t"
110     << entry_plug->Hmat[1][1] << "\t"
111     << entry_plug->Hmat[2][1] << ";\t"
112 mmeineke 572
113 mmeineke 590 << entry_plug->Hmat[0][2] << "\t"
114     << entry_plug->Hmat[1][2] << "\t"
115 tim 837 << entry_plug->Hmat[2][2] << ";";
116     //write out additional parameters, such as chi and eta
117     outFile << entry_plug->the_integrator->getAdditionalParameters();
118     outFile << endl;
119    
120 mmeineke 377 for( i=0; i<nAtoms; i++ ){
121 tim 837
122 mmeineke 670 atoms[i]->getPos(pos);
123     atoms[i]->getVel(vel);
124 mmeineke 377
125     sprintf( tempBuffer,
126     "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
127     atoms[i]->getType(),
128 mmeineke 670 pos[0],
129     pos[1],
130     pos[2],
131     vel[0],
132     vel[1],
133     vel[2]);
134 mmeineke 377 strcpy( writeLine, tempBuffer );
135    
136     if( atoms[i]->isDirectional() ){
137 tim 837
138 mmeineke 377 dAtom = (DirectionalAtom *)atoms[i];
139     dAtom->getQ( q );
140 tim 837
141 mmeineke 377 sprintf( tempBuffer,
142     "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
143     q[0],
144     q[1],
145     q[2],
146     q[3],
147     dAtom->getJx(),
148     dAtom->getJy(),
149     dAtom->getJz());
150     strcat( writeLine, tempBuffer );
151     }
152     else
153     strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
154 tim 837
155 mmeineke 377 outFile << writeLine;
156     }
157     outFile.flush();
158    
159     #else // is_mpi
160 gezelter 416
161 mmeineke 440 // first thing first, suspend fatalities.
162     painCave.isEventLoop = 1;
163    
164     int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone
165     int haveError;
166    
167 mmeineke 447 MPI_Status istatus;
168 gezelter 416 int *AtomToProcMap = mpiSim->getAtomToProcMap();
169 tim 837
170 mmeineke 377 // write out header and node 0's coordinates
171 tim 837
172 mmeineke 377 if( worldRank == 0 ){
173     outFile << mpiSim->getTotAtoms() << "\n";
174 tim 837
175     outFile << currentTime << ";\t"
176 gezelter 591 << entry_plug->Hmat[0][0] << "\t"
177     << entry_plug->Hmat[1][0] << "\t"
178     << entry_plug->Hmat[2][0] << ";\t"
179 tim 837
180 gezelter 591 << entry_plug->Hmat[0][1] << "\t"
181     << entry_plug->Hmat[1][1] << "\t"
182     << entry_plug->Hmat[2][1] << ";\t"
183 tim 837
184 gezelter 591 << entry_plug->Hmat[0][2] << "\t"
185     << entry_plug->Hmat[1][2] << "\t"
186 tim 837 << entry_plug->Hmat[2][2] << ";";
187    
188     outFile << entry_plug->the_integrator->getAdditionalParameters();
189     outFile << endl;
190 chuckv 434 outFile.flush();
191 gezelter 416 for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
192 gezelter 417 // Get the Node number which has this atom;
193 tim 837
194     which_node = AtomToProcMap[i];
195    
196 gezelter 907 if (which_node != 0) {
197    
198     atomTypeTag = 4*i;
199     atomIsDirectionalTag = 4*i + 1;
200     atomTransDataTag = 4*i + 2;
201     atomOrientDataTag = 4*i + 3;
202 tim 837
203 gezelter 910 MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, which_node,
204 gezelter 907 atomTypeTag, MPI_COMM_WORLD, &istatus);
205    
206 gezelter 910 strncpy(atomTypeString, MPIatomTypeString, MINIBUFFERSIZE);
207    
208 gezelter 907 MPI_Recv(&isDirectional, 1, MPI_INT, which_node,
209     atomIsDirectionalTag, MPI_COMM_WORLD, &istatus);
210    
211     MPI_Recv(atomTransData, 6, MPI_DOUBLE, which_node,
212     atomTransDataTag, MPI_COMM_WORLD, &istatus);
213    
214     if (isDirectional) {
215    
216     MPI_Recv(atomOrientData, 7, MPI_DOUBLE, which_node,
217     atomOrientDataTag, MPI_COMM_WORLD, &istatus);
218    
219     }
220    
221     } else {
222    
223 mmeineke 440 haveError = 0;
224 chuckv 436 which_atom = i;
225 tim 837 local_index=-1;
226 gezelter 907
227 chuckv 436 for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
228     if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
229     }
230 gezelter 907
231 chuckv 436 if (local_index != -1) {
232 tim 837
233 gezelter 907 atomTypeString = atoms[local_index]->getType();
234    
235 mmeineke 670 atoms[local_index]->getPos(pos);
236     atoms[local_index]->getVel(vel);
237    
238 gezelter 907 atomTransData[0] = pos[0];
239     atomTransData[1] = pos[1];
240     atomTransData[2] = pos[2];
241 tim 837
242 gezelter 907 atomTransData[3] = vel[0];
243     atomTransData[4] = vel[1];
244     atomTransData[5] = vel[2];
245    
246     isDirectional = 0;
247    
248 chuckv 436 if( atoms[local_index]->isDirectional() ){
249 tim 837
250 gezelter 907 isDirectional = 1;
251    
252 chuckv 436 dAtom = (DirectionalAtom *)atoms[local_index];
253     dAtom->getQ( q );
254 gezelter 907
255     atomOrientData[0] = q[0];
256     atomOrientData[1] = q[1];
257     atomOrientData[2] = q[2];
258     atomOrientData[3] = q[3];
259 tim 837
260 gezelter 907 atomOrientData[4] = dAtom->getJx();
261     atomOrientData[5] = dAtom->getJy();
262     atomOrientData[6] = dAtom->getJz();
263     }
264 tim 837
265 gezelter 907 } else {
266 mmeineke 440 sprintf(painCave.errMsg,
267     "Atom %d not found on processor %d\n",
268     i, worldRank );
269     haveError= 1;
270     simError();
271 tim 837 }
272    
273 gezelter 907 if(haveError) DieDieDie();
274    
275     // If we've survived to here, format the line:
276    
277     sprintf( tempBuffer,
278     "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
279     atomTypeString,
280     atomTransData[0],
281     atomTransData[1],
282     atomTransData[2],
283     atomTransData[3],
284     atomTransData[4],
285     atomTransData[5]);
286 mmeineke 440
287 gezelter 907 strcpy( writeLine, tempBuffer );
288 tim 837
289 gezelter 907 if (isDirectional) {
290 mmeineke 440
291 gezelter 907 sprintf( tempBuffer,
292     "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
293     atomOrientData[0],
294     atomOrientData[1],
295     atomOrientData[2],
296     atomOrientData[3],
297     atomOrientData[4],
298     atomOrientData[5],
299     atomOrientData[6]);
300     strcat( writeLine, tempBuffer );
301 tim 837
302 gezelter 907 } else {
303     strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
304     }
305 tim 837
306 gezelter 907 outFile << writeLine;
307     outFile.flush();
308     }
309 mmeineke 377 }
310    
311 gezelter 907 outFile.flush();
312     sprintf( checkPointMsg,
313     "Sucessfully took a dump.\n");
314     MPIcheckPoint();
315    
316 gezelter 415 } else {
317 tim 837
318 gezelter 907 // worldRank != 0, so I'm a remote node.
319    
320     for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
321    
322     // Am I the node which has this atom?
323    
324     if (AtomToProcMap[i] == worldRank) {
325 tim 837
326 gezelter 907 local_index=-1;
327     for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
328     if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
329     }
330     if (local_index != -1) {
331    
332     atomTypeString = atoms[local_index]->getType();
333 mmeineke 440
334 gezelter 907 atoms[local_index]->getPos(pos);
335     atoms[local_index]->getVel(vel);
336 tim 837
337 gezelter 907 atomTransData[0] = pos[0];
338     atomTransData[1] = pos[1];
339     atomTransData[2] = pos[2];
340 mmeineke 440
341 gezelter 907 atomTransData[3] = vel[0];
342     atomTransData[4] = vel[1];
343     atomTransData[5] = vel[2];
344    
345     isDirectional = 0;
346 tim 837
347 gezelter 907 if( atoms[local_index]->isDirectional() ){
348 mmeineke 670
349 gezelter 907 isDirectional = 1;
350    
351     dAtom = (DirectionalAtom *)atoms[local_index];
352     dAtom->getQ( q );
353    
354     atomOrientData[0] = q[0];
355     atomOrientData[1] = q[1];
356     atomOrientData[2] = q[2];
357     atomOrientData[3] = q[3];
358 mmeineke 670
359 gezelter 907 atomOrientData[4] = dAtom->getJx();
360     atomOrientData[5] = dAtom->getJy();
361     atomOrientData[6] = dAtom->getJz();
362     }
363 tim 837
364 gezelter 907 } else {
365     sprintf(painCave.errMsg,
366     "Atom %d not found on processor %d\n",
367     i, worldRank );
368     haveError= 1;
369     simError();
370     }
371 tim 837
372 gezelter 907 // I've survived this far, so send off the data!
373 tim 837
374 gezelter 907 atomTypeTag = 4*i;
375     atomIsDirectionalTag = 4*i + 1;
376     atomTransDataTag = 4*i + 2;
377     atomOrientDataTag = 4*i + 3;
378 mmeineke 440
379 gezelter 910
380     strncpy(MPIatomTypeString, atomTypeString, MINIBUFFERSIZE);
381    
382     MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
383 gezelter 907 atomTypeTag, MPI_COMM_WORLD);
384    
385     MPI_Send(&isDirectional, 1, MPI_INT, 0,
386     atomIsDirectionalTag, MPI_COMM_WORLD);
387    
388     MPI_Send(atomTransData, 6, MPI_DOUBLE, 0,
389     atomTransDataTag, MPI_COMM_WORLD);
390    
391     if (isDirectional) {
392    
393     MPI_Send(atomOrientData, 7, MPI_DOUBLE, 0,
394     atomOrientDataTag, MPI_COMM_WORLD);
395    
396     }
397    
398 mmeineke 377 }
399 gezelter 907 }
400 mmeineke 440
401 gezelter 907 sprintf( checkPointMsg,
402     "Sucessfully took a dump.\n");
403     MPIcheckPoint();
404    
405 tim 837 }
406 gezelter 907
407 mmeineke 440 painCave.isEventLoop = 0;
408    
409 mmeineke 377 #endif // is_mpi
410     }
411    
412 mmeineke 572 void DumpWriter::writeFinal(double finalTime){
413 gezelter 416
414 mmeineke 377 char finalName[500];
415     ofstream finalOut;
416 gezelter 416
417     const int BUFFERSIZE = 2000;
418 gezelter 907 const int MINIBUFFERSIZE = 10;
419 gezelter 416 char tempBuffer[BUFFERSIZE];
420 tim 837 char writeLine[BUFFERSIZE];
421 gezelter 416
422     double q[4];
423     DirectionalAtom* dAtom;
424 mmeineke 787 Atom** atoms = entry_plug->atoms;
425     int i;
426     #ifdef IS_MPI
427     int j, which_node, done, which_atom, local_index;
428 gezelter 907 double atomTransData[6];
429     double atomOrientData[7];
430     int isDirectional;
431     char* atomTypeString;
432 gezelter 910 char MPIatomTypeString[MINIBUFFERSIZE];
433 gezelter 907 int atomTypeTag;
434     int atomIsDirectionalTag;
435     int atomTransDataTag;
436     int atomOrientDataTag;
437 mmeineke 787 #else //is_mpi
438 gezelter 416 int nAtoms = entry_plug->n_atoms;
439 mmeineke 787 #endif //is_mpi
440 tim 837
441 mmeineke 670 double pos[3], vel[3];
442 tim 837
443 mmeineke 377 #ifdef IS_MPI
444     if(worldRank == 0 ){
445     #endif // is_mpi
446 tim 837
447 mmeineke 377 strcpy( finalName, entry_plug->finalName );
448 tim 837
449 mmeineke 377 finalOut.open( finalName, ios::out | ios::trunc );
450     if( !finalOut ){
451     sprintf( painCave.errMsg,
452     "Could not open \"%s\" for final dump output.\n",
453     finalName );
454     painCave.isFatal = 1;
455     simError();
456     }
457 tim 837
458 mmeineke 377 // finalOut.setf( ios::scientific );
459 tim 837
460 mmeineke 377 #ifdef IS_MPI
461     }
462 tim 837
463 mmeineke 377 sprintf(checkPointMsg,"Opened file for final configuration\n");
464 tim 837 MPIcheckPoint();
465    
466 mmeineke 377 #endif //is_mpi
467    
468 tim 837
469 mmeineke 377 #ifndef IS_MPI
470 tim 837
471 mmeineke 377 finalOut << nAtoms << "\n";
472 tim 837
473 gezelter 591 finalOut << finalTime << ";\t"
474 mmeineke 590 << entry_plug->Hmat[0][0] << "\t"
475     << entry_plug->Hmat[1][0] << "\t"
476 gezelter 591 << entry_plug->Hmat[2][0] << ";\t"
477 tim 837
478 mmeineke 590 << entry_plug->Hmat[0][1] << "\t"
479     << entry_plug->Hmat[1][1] << "\t"
480 gezelter 591 << entry_plug->Hmat[2][1] << ";\t"
481 tim 837
482 mmeineke 590 << entry_plug->Hmat[0][2] << "\t"
483     << entry_plug->Hmat[1][2] << "\t"
484 tim 837 << entry_plug->Hmat[2][2] << ";";
485    
486     //write out additional parameters, such as chi and eta
487     finalOut << entry_plug->the_integrator->getAdditionalParameters();
488     finalOut << endl;
489    
490 mmeineke 377 for( i=0; i<nAtoms; i++ ){
491 tim 837
492 mmeineke 670 atoms[i]->getPos(pos);
493     atoms[i]->getVel(vel);
494 tim 837
495 mmeineke 377 sprintf( tempBuffer,
496     "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
497     atoms[i]->getType(),
498 mmeineke 670 pos[0],
499     pos[1],
500     pos[2],
501     vel[0],
502     vel[1],
503     vel[2]);
504 mmeineke 377 strcpy( writeLine, tempBuffer );
505    
506     if( atoms[i]->isDirectional() ){
507 tim 837
508 mmeineke 377 dAtom = (DirectionalAtom *)atoms[i];
509     dAtom->getQ( q );
510 tim 837
511 mmeineke 377 sprintf( tempBuffer,
512     "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
513     q[0],
514     q[1],
515     q[2],
516     q[3],
517     dAtom->getJx(),
518     dAtom->getJy(),
519     dAtom->getJz());
520     strcat( writeLine, tempBuffer );
521     }
522     else
523     strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
524 tim 837
525 mmeineke 377 finalOut << writeLine;
526     }
527     finalOut.flush();
528 gezelter 415 finalOut.close();
529 mmeineke 377
530     #else // is_mpi
531 tim 837
532 mmeineke 440 // first thing first, suspend fatalities.
533     painCave.isEventLoop = 1;
534    
535     int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone
536     int haveError;
537    
538 mmeineke 447 MPI_Status istatus;
539 gezelter 416 int *AtomToProcMap = mpiSim->getAtomToProcMap();
540    
541 mmeineke 377 // write out header and node 0's coordinates
542 tim 837
543 mmeineke 377 if( worldRank == 0 ){
544     finalOut << mpiSim->getTotAtoms() << "\n";
545 tim 837
546 gezelter 591 finalOut << finalTime << ";\t"
547 gezelter 907 << entry_plug->Hmat[0][0] << "\t"
548     << entry_plug->Hmat[1][0] << "\t"
549     << entry_plug->Hmat[2][0] << ";\t"
550 tim 837
551 gezelter 907 << entry_plug->Hmat[0][1] << "\t"
552     << entry_plug->Hmat[1][1] << "\t"
553     << entry_plug->Hmat[2][1] << ";\t"
554 tim 837
555 gezelter 907 << entry_plug->Hmat[0][2] << "\t"
556     << entry_plug->Hmat[1][2] << "\t"
557     << entry_plug->Hmat[2][2] << ";";
558 tim 837
559     finalOut << entry_plug->the_integrator->getAdditionalParameters();
560     finalOut << endl;
561 gezelter 907 finalOut.flush();
562 gezelter 416 for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
563 gezelter 907 // Get the Node number which has this atom;
564 tim 837
565     which_node = AtomToProcMap[i];
566    
567 gezelter 907 if (which_node != 0) {
568    
569     atomTypeTag = 4*i;
570     atomIsDirectionalTag = 4*i + 1;
571     atomTransDataTag = 4*i + 2;
572     atomOrientDataTag = 4*i + 3;
573 chuckv 437
574 gezelter 910 MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, which_node,
575 gezelter 907 atomTypeTag, MPI_COMM_WORLD, &istatus);
576    
577 gezelter 910 strncpy(atomTypeString, MPIatomTypeString, MINIBUFFERSIZE);
578    
579 gezelter 907 MPI_Recv(&isDirectional, 1, MPI_INT, which_node,
580     atomIsDirectionalTag, MPI_COMM_WORLD, &istatus);
581    
582     MPI_Recv(atomTransData, 6, MPI_DOUBLE, which_node,
583     atomTransDataTag, MPI_COMM_WORLD, &istatus);
584    
585     if (isDirectional) {
586    
587     MPI_Recv(atomOrientData, 7, MPI_DOUBLE, which_node,
588     atomOrientDataTag, MPI_COMM_WORLD, &istatus);
589    
590     }
591    
592     } else {
593    
594     haveError = 0;
595 chuckv 437 which_atom = i;
596 tim 837 local_index=-1;
597 gezelter 907
598 chuckv 437 for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
599     if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
600     }
601 gezelter 907
602 tim 837 if (local_index != -1) {
603 mmeineke 670
604 gezelter 907 atomTypeString = atoms[local_index]->getType();
605    
606 mmeineke 670 atoms[local_index]->getPos(pos);
607     atoms[local_index]->getVel(vel);
608 tim 837
609 gezelter 907 atomTransData[0] = pos[0];
610     atomTransData[1] = pos[1];
611     atomTransData[2] = pos[2];
612 tim 837
613 gezelter 907 atomTransData[3] = vel[0];
614     atomTransData[4] = vel[1];
615     atomTransData[5] = vel[2];
616    
617     isDirectional = 0;
618 tim 837
619 gezelter 907 if( atoms[local_index]->isDirectional() ){
620 tim 837
621 gezelter 907 isDirectional = 1;
622    
623     dAtom = (DirectionalAtom *)atoms[local_index];
624     dAtom->getQ( q );
625    
626     atomOrientData[0] = q[0];
627     atomOrientData[1] = q[1];
628     atomOrientData[2] = q[2];
629     atomOrientData[3] = q[3];
630    
631     atomOrientData[4] = dAtom->getJx();
632     atomOrientData[5] = dAtom->getJy();
633     atomOrientData[6] = dAtom->getJz();
634     }
635    
636     } else {
637 mmeineke 440 sprintf(painCave.errMsg,
638     "Atom %d not found on processor %d\n",
639     i, worldRank );
640     haveError= 1;
641     simError();
642 tim 837 }
643 mmeineke 440
644 gezelter 907 if(haveError) DieDieDie();
645    
646     // If we've survived to here, format the line:
647    
648     sprintf( tempBuffer,
649     "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
650     atomTypeString,
651     atomTransData[0],
652     atomTransData[1],
653     atomTransData[2],
654     atomTransData[3],
655     atomTransData[4],
656     atomTransData[5]);
657 tim 837
658 gezelter 907 strcpy( writeLine, tempBuffer );
659 tim 837
660 gezelter 907 if (isDirectional) {
661 tim 837
662 gezelter 907 sprintf( tempBuffer,
663     "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
664     atomOrientData[0],
665     atomOrientData[1],
666     atomOrientData[2],
667     atomOrientData[3],
668     atomOrientData[4],
669     atomOrientData[5],
670     atomOrientData[6]);
671     strcat( writeLine, tempBuffer );
672 tim 837
673 gezelter 907 } else {
674     strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
675     }
676 tim 837
677 gezelter 907 finalOut << writeLine;
678     finalOut.flush();
679     }
680 mmeineke 377 }
681    
682 gezelter 907 finalOut.flush();
683     sprintf( checkPointMsg,
684     "Sucessfully took a dump.\n");
685     MPIcheckPoint();
686    
687 gezelter 415 } else {
688 tim 837
689 gezelter 907 // worldRank != 0, so I'm a remote node.
690    
691     for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
692    
693     // Am I the node which has this atom?
694    
695     if (AtomToProcMap[i] == worldRank) {
696 mmeineke 440
697 gezelter 907 local_index=-1;
698     for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
699     if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
700     }
701     if (local_index != -1) {
702    
703     atomTypeString = atoms[local_index]->getType();
704 tim 837
705 gezelter 907 atoms[local_index]->getPos(pos);
706     atoms[local_index]->getVel(vel);
707 tim 837
708 gezelter 907 atomTransData[0] = pos[0];
709     atomTransData[1] = pos[1];
710     atomTransData[2] = pos[2];
711 tim 837
712 gezelter 907 atomTransData[3] = vel[0];
713     atomTransData[4] = vel[1];
714     atomTransData[5] = vel[2];
715    
716     isDirectional = 0;
717 tim 837
718 gezelter 907 if( atoms[local_index]->isDirectional() ){
719 mmeineke 377
720 gezelter 907 isDirectional = 1;
721    
722     dAtom = (DirectionalAtom *)atoms[local_index];
723     dAtom->getQ( q );
724    
725     atomOrientData[0] = q[0];
726     atomOrientData[1] = q[1];
727     atomOrientData[2] = q[2];
728     atomOrientData[3] = q[3];
729 mmeineke 670
730 gezelter 907 atomOrientData[4] = dAtom->getJx();
731     atomOrientData[5] = dAtom->getJy();
732     atomOrientData[6] = dAtom->getJz();
733     }
734 tim 837
735 gezelter 907 } else {
736     sprintf(painCave.errMsg,
737     "Atom %d not found on processor %d\n",
738     i, worldRank );
739     haveError= 1;
740     simError();
741     }
742 tim 837
743 gezelter 907 // I've survived this far, so send off the data!
744 tim 837
745 gezelter 907 atomTypeTag = 4*i;
746     atomIsDirectionalTag = 4*i + 1;
747     atomTransDataTag = 4*i + 2;
748     atomOrientDataTag = 4*i + 3;
749 tim 837
750 gezelter 910 strncpy(MPIatomTypeString, atomTypeString, MINIBUFFERSIZE);
751    
752     MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
753 gezelter 907 atomTypeTag, MPI_COMM_WORLD);
754    
755     MPI_Send(&isDirectional, 1, MPI_INT, 0,
756     atomIsDirectionalTag, MPI_COMM_WORLD);
757    
758     MPI_Send(atomTransData, 6, MPI_DOUBLE, 0,
759     atomTransDataTag, MPI_COMM_WORLD);
760    
761     if (isDirectional) {
762    
763     MPI_Send(atomOrientData, 7, MPI_DOUBLE, 0,
764     atomOrientDataTag, MPI_COMM_WORLD);
765    
766     }
767    
768 mmeineke 377 }
769 gezelter 907 }
770 mmeineke 440
771 gezelter 907 sprintf( checkPointMsg,
772     "Sucessfully wrote final file.\n");
773     MPIcheckPoint();
774    
775 gezelter 419 }
776 gezelter 907
777     painCave.isEventLoop = 0;
778 tim 837
779     if( worldRank == 0 ) finalOut.close();
780 mmeineke 377 #endif // is_mpi
781     }
782 mmeineke 440
783    
784    
785     #ifdef IS_MPI
786    
787     // a couple of functions to let us escape the write loop
788    
789 gezelter 907 void dWrite::DieDieDie( void ){
790 tim 837
791 mmeineke 440 MPI_Finalize();
792     exit (0);
793     }
794    
795     #endif //is_mpi